PeriDyno 1.0.0
Loading...
Searching...
No Matches
HeightFieldToTriangleSet.cu
Go to the documentation of this file.
1#include "HeightFieldToTriangleSet.h"
2
3namespace dyno
4{
5 template<typename TDataType>
6 HeightFieldToTriangleSet<TDataType>::HeightFieldToTriangleSet()
7 : TopologyMapping()
8 {
9 }
10
11 template<typename Real, typename Coord>
12 __global__ void SetupVerticesForHeightField(
13 DArray<Coord> vertices,
14 DArray2D<Coord> displacement,
15 Coord origin,
16 Real h,
17 Coord translation,
18 Real scale)
19 {
20 int i = threadIdx.x + (blockIdx.x * blockDim.x);
21 int j = threadIdx.y + (blockIdx.y * blockDim.y);
22
23 if (i >= displacement.nx() || j >= displacement.ny()) return;
24
25 Coord di = displacement(i, j);
26 Coord v = Coord(origin.x + i * h + di.x, origin.y + di.y, origin.z + j * h + di.z);
27
28 vertices[i + j * displacement.nx()] = scale * v + translation;
29 }
30
31 template<typename Triangle>
32 __global__ void SetupTrianglesForHeightField(
33 DArray<Triangle> vertices,
34 uint nx,
35 uint ny)
36 {
37 int i = threadIdx.x + (blockIdx.x * blockDim.x);
38 int j = threadIdx.y + (blockIdx.y * blockDim.y);
39
40 if (i >= nx - 1 || j >= ny - 1) return;
41
42 uint v0 = i + j * nx;
43 uint v1 = i + 1 + j * nx;
44
45 uint v2 = i + (j + 1) * nx;
46 uint v3 = i + 1 + (j + 1) * nx;
47
48 Triangle t0(v1, v0, v2);
49 Triangle t1(v1, v2, v3);
50
51 uint offset = 2 * i + 2 * j * (nx - 1);
52 vertices[offset] = t0;
53 vertices[offset + 1] = t1;
54 }
55
56 template<typename TDataType>
57 bool HeightFieldToTriangleSet<TDataType>::apply()
58 {
59 if (this->outTriangleSet()->isEmpty())
60 {
61 this->outTriangleSet()->allocate();
62 }
63
64 auto heights = this->inHeightField()->getDataPtr();
65
66 auto triSet = this->outTriangleSet()->getDataPtr();
67
68 auto& vertices = triSet->getPoints();
69 auto& indices = triSet->getTriangles();
70
71 int numOfVertices = heights->width() * heights->height();
72 int numOfTriangles = 2 * (heights->width() - 1) * (heights->height() - 1);
73
74 vertices.resize(numOfVertices);
75 indices.resize(numOfTriangles);
76
77
78 auto& disp = heights->getDisplacement();
79
80 Real scale = this->varScale()->getData();
81 Coord translation = this->varTranslation()->getData();
82
83 uint2 dim;
84 dim.x = disp.nx();
85 dim.y = disp.ny();
86 cuExecute2D(dim,
87 SetupVerticesForHeightField,
88 vertices,
89 disp,
90 heights->getOrigin(),
91 heights->getGridSpacing(),
92 translation,
93 scale);
94
95 dim.x = disp.nx() - 1;
96 dim.y = disp.ny() - 1;
97 cuExecute2D(dim,
98 SetupTrianglesForHeightField,
99 indices,
100 heights->width(),
101 heights->height());
102
103 triSet->update();
104
105 return true;
106 }
107
108 DEFINE_CLASS(HeightFieldToTriangleSet);
109}