PeriDyno 1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
ParticleIntegrator.cu
Go to the documentation of this file.
1#include <cuda_runtime.h>
2#include "ParticleIntegrator.h"
3#include "Node.h"
4#include "SceneGraphFactory.h"
5
6namespace dyno
7{
8 IMPLEMENT_TCLASS(ParticleIntegrator, TDataType)
9
10 template<typename TDataType>
11 ParticleIntegrator<TDataType>::ParticleIntegrator()
12 : ComputeModule()
13 {
14 this->inAttribute()->tagOptional(true);
15 }
16
17 template<typename Real, typename Coord>
18 __global__ void K_UpdateVelocity(
19 DArray<Coord> vel,
20 Coord gravity,
21 Real dt)
22 {
23 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
24 if (pId >= vel.size()) return;
25
26 vel[pId] += dt * (gravity);
27 }
28
29
30
31
32 template<typename Real, typename Coord>
33 __global__ void K_UpdateVelocity(
34 DArray<Coord> vel,
35 DArray<Attribute> atts,
36 Coord gravity,
37 Real dt)
38 {
39 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
40 if (pId >= vel.size()) return;
41
42 Attribute att = atts[pId];
43
44 if (att.isDynamic())
45 {
46 vel[pId] += dt * (gravity);
47 }
48 }
49
50
51
52
53 template<typename TDataType>
54 bool ParticleIntegrator<TDataType>::updateVelocity()
55 {
56 Real dt = this->inTimeStep()->getData();
57
58 auto scn = dyno::SceneGraphFactory::instance()->active();
59 Coord gravity = scn->getGravity();
60
61 int total_num = this->inPosition()->size();
62
63 if (this->inAttribute()->isEmpty())
64 {
65 cuExecute(total_num,
66 K_UpdateVelocity,
67 this->inVelocity()->getData(),
68 gravity,
69 dt);
70 }
71 else
72 {
73 cuExecute(total_num,
74 K_UpdateVelocity,
75 this->inVelocity()->getData(),
76 this->inAttribute()->getData(),
77 gravity,
78 dt);
79 }
80
81
82 return true;
83 }
84
85 template<typename Real, typename Coord>
86 __global__ void K_UpdatePosition(
87 DArray<Coord> pos,
88 DArray<Coord> vel,
89 Real dt)
90 {
91 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
92 if (pId >= pos.size()) return;
93
94 pos[pId] += dt * vel[pId];
95 }
96
97 template<typename Real, typename Coord>
98 __global__ void K_UpdatePosition(
99 DArray<Coord> pos,
100 DArray<Coord> vel,
101 DArray<Attribute> att,
102 Real dt)
103 {
104 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
105 if (pId >= pos.size()) return;
106
107 Attribute att_i = att[pId];
108
109 if (!att_i.isFixed())
110 {
111 pos[pId] += dt * vel[pId];
112 }
113 }
114
115 template<typename TDataType>
116 bool ParticleIntegrator<TDataType>::updatePosition()
117 {
118 Real dt = this->inTimeStep()->getData();
119
120 int total_num = this->inPosition()->getDataPtr()->size();
121
122
123 if (this->inAttribute()->isEmpty())
124 {
125 cuExecute(total_num,
126 K_UpdatePosition,
127 this->inPosition()->getData(),
128 this->inVelocity()->getData(),
129 dt);
130 }
131 else
132 {
133 cuExecute(total_num,
134 K_UpdatePosition,
135 this->inPosition()->getData(),
136 this->inVelocity()->getData(),
137 this->inAttribute()->getData(),
138 dt);
139 }
140
141
142 return true;
143 }
144
145 template<typename TDataType>
146 void ParticleIntegrator<TDataType>::compute()
147 {
148 updatePosition();
149 updateVelocity();
150 }
151
152 DEFINE_CLASS(ParticleIntegrator);
153}