PeriDyno 1.0.0
Loading...
Searching...
No Matches
DiscreteElementsToTriangleSet.cu
Go to the documentation of this file.
1#include "DiscreteElementsToTriangleSet.h"
2
3namespace dyno
4{
5 typedef typename ::dyno::TOrientedBox3D<Real> Box3D;
6
7 template<typename TDataType>
8 DiscreteElementsToTriangleSet<TDataType>::DiscreteElementsToTriangleSet()
9 : TopologyMapping()
10 {
11 mStandardSphere.loadObjFile(getAssetPath() + "standard/standard_icosahedron.obj");
12 mStandardCapsule.loadObjFile(getAssetPath() + "standard/standard_capsule.obj");
13 }
14
15 template<typename Triangle>
16 __global__ void SetupCubeInstances(
17 DArray<Vec3f> vertices,
18 DArray<Triangle> indices,
19 DArray<Box3D> boxes,
20 uint pointOffset,
21 uint indexOffset,
22 uint cubeOffset)
23 {
24 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
25 if (tId >= boxes.size()) return;
26
27 int idx = tId;
28 Box3D box = boxes[idx];
29
30 Vec3f hx = box.u * box.extent[0];
31 Vec3f hy = box.v * box.extent[1];
32 Vec3f hz = box.w * box.extent[2];
33
34 Vec3f hyz = hy + hz;
35 Vec3f hxy = hx + hy;
36 Vec3f hxz = hx + hz;
37
38 Vec3f c = box.center;
39
40 Vec3f v0 = c - hx - hyz;
41 Vec3f v1 = c + hx - hyz;
42 Vec3f v2 = c + hxz - hy;
43 Vec3f v3 = c - hxy + hz;
44
45 Vec3f v4 = c - hxz + hy;
46 Vec3f v5 = c + hxy - hz;
47 Vec3f v6 = c + hx + hyz;
48 Vec3f v7 = c - hx + hyz;
49
50 vertices[pointOffset + idx * 8] = v0;
51 vertices[pointOffset + idx * 8 + 1] = v1;
52 vertices[pointOffset + idx * 8 + 2] = v2;
53 vertices[pointOffset + idx * 8 + 3] = v3;
54 vertices[pointOffset + idx * 8 + 4] = v4;
55 vertices[pointOffset + idx * 8 + 5] = v5;
56 vertices[pointOffset + idx * 8 + 6] = v6;
57 vertices[pointOffset + idx * 8 + 7] = v7;
58
59 uint offset = idx * 8 + pointOffset;
60
61 indices[indexOffset + idx * 12] = Triangle(offset + 0, offset + 1, offset + 2);
62 indices[indexOffset + idx * 12 + 1] = Triangle(offset + 0, offset + 2, offset + 3);
63
64 indices[indexOffset + idx * 12 + 2] = Triangle(offset + 0, offset + 4, offset + 5);
65 indices[indexOffset + idx * 12 + 3] = Triangle(offset + 0, offset + 5, offset + 1);
66
67 indices[indexOffset + idx * 12 + 4] = Triangle(offset + 4, offset + 7, offset + 6);
68 indices[indexOffset + idx * 12 + 5] = Triangle(offset + 4, offset + 6, offset + 5);
69
70 indices[indexOffset + idx * 12 + 6] = Triangle(offset + 1, offset + 5, offset + 6);
71 indices[indexOffset + idx * 12 + 7] = Triangle(offset + 1, offset + 6, offset + 2);
72
73 indices[indexOffset + idx * 12 + 8] = Triangle(offset + 2, offset + 6, offset + 7);
74 indices[indexOffset + idx * 12 + 9] = Triangle(offset + 2, offset + 7, offset + 3);
75
76 indices[indexOffset + idx * 12 + 10] = Triangle(offset + 0, offset + 3, offset + 7);
77 indices[indexOffset + idx * 12 + 11] = Triangle(offset + 0, offset + 7, offset + 4);
78 }
79
80 template<typename Triangle>
81 __global__ void SetupTetInstances(
82 DArray<Vec3f> vertices,
83 DArray<Triangle> indices,
84 DArray<Tet3D> tets,
85 uint pointOffset,
86 uint indexOffset,
87 uint tetOffset)
88 {
89 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
90 if (tId >= tets.size()) return;
91
92 int idx = tId;
93 Tet3D tet = tets[idx];
94
95 Vec3f v0 = tet.v[0];
96 Vec3f v1 = tet.v[1];
97 Vec3f v2 = tet.v[2];
98 Vec3f v3 = tet.v[3];
99
100 vertices[pointOffset + idx * 4] = v0;
101 vertices[pointOffset + idx * 4 + 1] = v1;
102 vertices[pointOffset + idx * 4 + 2] = v2;
103 vertices[pointOffset + idx * 4 + 3] = v3;
104
105 uint offset = idx * 4 + pointOffset;
106
107 indices[indexOffset + idx * 4] = Triangle(offset + 0, offset + 1, offset + 2);
108 indices[indexOffset + idx * 4 + 1] = Triangle(offset + 0, offset + 1, offset + 3);
109 indices[indexOffset + idx * 4 + 2] = Triangle(offset + 1, offset + 2, offset + 3);
110 indices[indexOffset + idx * 4 + 3] = Triangle(offset + 0, offset + 2, offset + 3);
111 }
112
113 __global__ void SetupVerticesForSphereInstances(
114 DArray<Vec3f> vertices,
115 DArray<Vec3f> sphereVertices,
116 DArray<Sphere3D> sphereInstances,
117 uint pointOffset,
118 uint sphereOffset)
119 {
120 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
121 if (tId >= sphereInstances.size() * sphereVertices.size()) return;
122
123 uint instanceId = tId / sphereVertices.size();
124 uint vertexId = tId % sphereVertices.size();
125
126 Sphere3D sphere = sphereInstances[instanceId];
127
128 Vec3f v = sphereVertices[vertexId];
129 vertices[pointOffset + tId] = sphere.center + sphere.radius * sphere.rotation.rotate(v);
130 }
131
132 template<typename Triangle>
133 __global__ void SetupIndicesForSphereInstances(
134 DArray<Triangle> indices,
135 DArray<Triangle> sphereIndices,
136 DArray<Sphere3D> sphereInstances,
137 uint vertexSize, //vertex size of the instance sphere
138 uint indexOffset)
139 {
140 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
141 if (tId >= sphereInstances.size() * sphereIndices.size()) return;
142
143 uint instanceId = tId / sphereIndices.size();
144 uint indexId = tId % sphereIndices.size();
145
146 int vertexOffset = indexOffset + instanceId * vertexSize;
147
148 Triangle tIndex = sphereIndices[indexId];
149 indices[indexOffset + tId] = Triangle(tIndex[0] + vertexOffset, tIndex[1] + vertexOffset, tIndex[2] + vertexOffset);
150 }
151
152 __global__ void SetupVerticesForCapsuleInstances(
153 DArray<Vec3f> vertices,
154 DArray<Vec3f> capsuleVertices,
155 DArray<Capsule3D> capsuleInstances,
156 uint pointOffset,
157 uint capsuleOffset)
158 {
159 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
160 if (tId >= capsuleInstances.size() * capsuleVertices.size()) return;
161
162 uint instanceId = tId / capsuleVertices.size();
163 uint vertexId = tId % capsuleVertices.size();
164
165 Capsule3D capsule = capsuleInstances[instanceId];
166 float r = capsule.radius;
167 float h = capsule.halfLength;
168 auto rot = capsule.rotation.toMatrix3x3();
169 Vec3f center = capsule.center;
170
171 Vec3f v = capsuleVertices[vertexId];
172 Vec3f orignZ = Vec3f(0, 1, 0);
173 Vec3f newZ = Vec3f(0, h, 0);
174
175 if (v.y >= 1)
176 {
177 vertices[pointOffset + tId] = rot * ((v - orignZ) * r + newZ) + center;
178 }
179 else if (v.y <= -1)
180 {
181 vertices[pointOffset + tId] = rot * ((v + orignZ) * r - newZ) + center;
182 }
183 else
184 {
185 vertices[pointOffset + tId] = rot * (v * Vec3f(r, h, r)) + center;
186 }
187 }
188
189 template<typename Triangle>
190 __global__ void SetupIndicesForCapsuleInstances(
191 DArray<Triangle> indices,
192 DArray<Triangle> capsuleIndices,
193 DArray<Capsule3D> capsuleInstances,
194 uint vertexSize, //vertex size of the instance sphere
195 uint vertexOffset,
196 uint indexOffset)
197 {
198 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
199 if (tId >= capsuleInstances.size() * capsuleIndices.size()) return;
200
201 uint instanceId = tId / capsuleIndices.size();
202 uint indexId = tId % capsuleIndices.size();
203
204 vertexOffset += instanceId * vertexSize;
205
206 Triangle tIndex = capsuleIndices[indexId];
207 indices[indexOffset + tId] = Triangle(tIndex[0] + vertexOffset, tIndex[1] + vertexOffset, tIndex[2] + vertexOffset);
208 }
209
210 template<typename TDataType>
211 bool DiscreteElementsToTriangleSet<TDataType>::apply()
212 {
213 if (this->outTriangleSet()->isEmpty())
214 {
215 this->outTriangleSet()->allocate();
216 }
217
218 auto inTopo = this->inDiscreteElements()->constDataPtr();
219
220 DArray<Box3D>& boxInGlobal = inTopo->boxesInGlobal();
221 DArray<Sphere3D>& sphereInGlobal = inTopo->spheresInGlobal();
222 DArray<Tet3D>& tetInGlobal = inTopo->tetsInGlobal();
223 DArray<Capsule3D>& capsuleInGlobal = inTopo->capsulesInGlobal();
224
225 ElementOffset elementOffset = inTopo->calculateElementOffset();
226
227 int numOfSpheres = sphereInGlobal.size();
228 int numofCaps = capsuleInGlobal.size();
229 int numOfBoxes = boxInGlobal.size();
230 int numOfTets = tetInGlobal.size();
231
232 auto triSet = this->outTriangleSet()->getDataPtr();
233
234 auto& vertices = triSet->getPoints();
235 auto& indices = triSet->getTriangles();
236
237 auto& sphereVertices = mStandardSphere.getPoints();
238 auto& sphereIndices = mStandardSphere.getTriangles();
239
240 auto& capsuleVertices = mStandardCapsule.getPoints();
241 auto& capsuleIndices = mStandardCapsule.getTriangles();
242
243 int numOfVertices = 8 * numOfBoxes + 4 * numOfTets + sphereVertices.size() * numOfSpheres + capsuleVertices.size() * numofCaps;
244 int numOfTriangles = 12 * numOfBoxes + 4 * numOfTets + sphereIndices.size() * numOfSpheres + capsuleIndices.size() * numofCaps;
245
246 vertices.resize(numOfVertices);
247 indices.resize(numOfTriangles);
248
249 uint vertexOffset = 0;
250 uint indexOffset = 0;
251
252 //Setup spheres
253 cuExecute(numOfSpheres * sphereVertices.size(),
254 SetupVerticesForSphereInstances,
255 vertices,
256 sphereVertices,
257 sphereInGlobal,
258 vertexOffset,
259 elementOffset.sphereIndex());
260
261 cuExecute(numOfSpheres * sphereIndices.size(),
262 SetupIndicesForSphereInstances,
263 indices,
264 sphereIndices,
265 sphereInGlobal,
266 sphereVertices.size(),
267 indexOffset);
268
269 vertexOffset += numOfSpheres * sphereVertices.size();
270 indexOffset += numOfSpheres * sphereIndices.size();
271
272 //Setup boxes
273 cuExecute(numOfBoxes,
274 SetupCubeInstances,
275 vertices,
276 indices,
277 boxInGlobal,
278 vertexOffset,
279 indexOffset,
280 elementOffset.boxIndex());
281
282 vertexOffset += numOfBoxes * 8;
283 indexOffset += numOfBoxes * 12;
284
285 //Setup tets
286 cuExecute(numOfTets,
287 SetupTetInstances,
288 vertices,
289 indices,
290 tetInGlobal,
291 vertexOffset,
292 indexOffset,
293 elementOffset.tetIndex());
294
295 vertexOffset += numOfTets * 4;
296 indexOffset += numOfTets * 4;
297
298 cuExecute(numofCaps * capsuleVertices.size(),
299 SetupVerticesForCapsuleInstances,
300 vertices,
301 capsuleVertices,
302 capsuleInGlobal,
303 vertexOffset,
304 elementOffset.capsuleIndex());
305
306 cuExecute(numofCaps * capsuleIndices.size(),
307 SetupIndicesForCapsuleInstances,
308 indices,
309 capsuleIndices,
310 capsuleInGlobal,
311 capsuleVertices.size(),
312 vertexOffset,
313 indexOffset);
314
315 vertexOffset += numofCaps * capsuleVertices.size();
316 indexOffset += numofCaps * capsuleIndices.size();
317
318 this->outTriangleSet()->getDataPtr()->update();
319
320 return true;
321 }
322
323 DEFINE_CLASS(DiscreteElementsToTriangleSet);
324}