PeriDyno 1.0.0
Loading...
Searching...
No Matches
NeighborTriangleQuery.cu
Go to the documentation of this file.
1#include "NeighborTriangleQuery.h"
2
3#include "Topology/SparseOctree.h"
4
5#include "Collision/CollisionDetectionBroadPhase.h"
6
7namespace dyno
8{
9 IMPLEMENT_TCLASS(NeighborTriangleQuery, TDataType)
10
11 template<typename TDataType>
12 NeighborTriangleQuery<TDataType>::NeighborTriangleQuery()
13 : ComputeModule()
14 {
15 mBroadPhaseCD = std::make_shared<CollisionDetectionBroadPhase<TDataType>>();
16 }
17
18 template<typename TDataType>
19 NeighborTriangleQuery<TDataType>::~NeighborTriangleQuery()
20 {
21 }
22
23 template<typename Real, typename Coord>
24 __global__ void NTQ_SetupAABB(
25 DArray<AABB> boundingBox,
26 DArray<Coord> position,
27 Real radius)
28 {
29 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
30 if (pId >= position.size()) return;
31
32 AABB box;
33 Coord p = position[pId];
34 box.v0 = p - radius;
35 box.v1 = p + radius;
36
37 boundingBox[pId] = box;
38 }
39
40 template<typename Coord>
41 __global__ void NTQ_SetupAABB(
42 DArray<AABB> boundingBox,
43 DArray<Coord> vertex,
44 DArray<TopologyModule::Triangle> tIndex)
45 {
46 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
47 if (tId >= boundingBox.size()) return;
48
49 AABB box;
50 TopologyModule::Triangle index = tIndex[tId];
51
52 Coord v0 = vertex[index[0]];
53 Coord v1 = vertex[index[1]];
54 Coord v2 = vertex[index[2]];
55
56 box.v0 = minimum(v0, minimum(v1, v2));
57 box.v1 = maximum(v0, maximum(v1, v2));
58
59 boundingBox[tId] = box;
60 }
61
62 template<typename Coord>
63 __global__ void NTQ_Narrow_Count(
64 DArrayList<int> nbr,
65 DArray<Coord> position,
66 DArray<Coord> vertex,
67 DArray<TopologyModule::Triangle> triangles,
68 DArray<uint> count,
69 Real radius)
70 {
71 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
72 if (tId >= position.size()) return;
73 uint cnt = 0;
74
75 List<int>& nbrIds_i = nbr[tId];
76 int nbSize = nbrIds_i.size();
77 //printf("CNT_OKKKK nbSize: %d\n", nbSize);
78 for (int ne = 0; ne < nbSize; ne++)
79 {
80 int j = nbrIds_i[ne];
81 Point3D point(position[tId]);
82 Triangle3D t3d_n(vertex[triangles[j][0]], vertex[triangles[j][1]], vertex[triangles[j][2]]);
83 Real p_dis_t = abs(point.distance(t3d_n));
84 if (p_dis_t< radius && p_dis_t > EPSILON)
85 {
86 cnt++;
87 }
88 }
89 count[tId] = cnt;
90 }
91
92 template<typename Coord>
93 __global__ void NTQ_Narrow_Set(
94 DArrayList<int> nbr,
95 DArrayList<int> nbr_out,
96 DArray<Coord> position,
97 DArray<Coord> vertex,
98 DArray<TopologyModule::Triangle> triangles,
99 Real radius)
100 {
101 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
102 if (tId >= position.size()) return;
103 int cnt = 0;
104
105 List<int>& nbrIds_i = nbr[tId];
106 int nbSize = nbrIds_i.size(); //int nbSize = nbr.getNeighborSize(tId);
107 List<int>& nbrOutIds_i = nbr_out[tId];
108 for (int ne = 0; ne < nbSize; ne++)
109 {
110 int j = nbrIds_i[ne];
111
112 Point3D point(position[tId]);
113 Triangle3D t3d_n(vertex[triangles[j][0]], vertex[triangles[j][1]], vertex[triangles[j][2]]);
114 Real proj_dist = abs(point.distance(t3d_n));
115
116 if (proj_dist < radius && proj_dist > EPSILON)
117 {
118 nbrOutIds_i.insert(j);
119 }
120 }
121 }
122
123 template<typename TDataType>
124 void NeighborTriangleQuery<TDataType>::compute()
125 {
126 int pNum = this->inPosition()->size();
127 if (pNum == 0) return;
128
129 if (mQueryAABB.size() != pNum) {
130 mQueryAABB.resize(pNum);
131 }
132
133 auto ts = this->inTriangleSet()->constDataPtr();
134 auto& triVertex = ts->getPoints();
135 auto& triIndex = ts->getTriangles();
136
137 int tNum = triIndex.size();
138 if (tNum == 0) return;
139 if (mQueriedAABB.size() != tNum) {
140 mQueriedAABB.resize(tNum);
141 }
142
143 cuExecute(pNum,
144 NTQ_SetupAABB,
145 mQueryAABB,
146 this->inPosition()->constData(),
147 this->inRadius()->getValue());
148
149 cuExecute(tNum,
150 NTQ_SetupAABB,
151 mQueriedAABB,
152 triVertex,
153 triIndex);
154
155 Real radius = this->inRadius()->getValue();
156
157 mBroadPhaseCD->varGridSizeLimit()->setValue(2 * radius);
158 mBroadPhaseCD->inSource()->assign(mQueryAABB);
159
160 if (this->inTriangleSet()->isModified()) {
161 mBroadPhaseCD->inTarget()->assign(mQueriedAABB);
162 }
163
164
165 auto type = this->varSpatial()->getDataPtr()->currentKey();
166
167 switch (type)
168 {
169 case Spatial::BVH:
170 mBroadPhaseCD->varAccelerationStructure()->setCurrentKey(CollisionDetectionBroadPhase<TDataType>::BVH);
171 break;
172 case Spatial::OCTREE:
173 mBroadPhaseCD->varAccelerationStructure()->setCurrentKey(CollisionDetectionBroadPhase<TDataType>::Octree);
174 break;
175 default:
176 break;
177 }
178
179 mBroadPhaseCD->update();
180
181 auto& nbr = mBroadPhaseCD->outContactList()->getData();
182
183 if (this->outNeighborIds()->size() != pNum)
184 {
185 this->outNeighborIds()->allocate();
186 DArray<uint> nbrNum;
187 nbrNum.resize(pNum);
188 nbrNum.reset();
189 auto& nbrIds = this->outNeighborIds()->getData();
190 nbrIds.resize(nbrNum);
191 nbrNum.clear();
192
193 //this->outNeighborIds()->getData().resize(p_num);
194 }
195 auto& nbrIds = this->outNeighborIds()->getData();
196
197 //new
198 DArray<uint> nbrNum;
199 nbrNum.resize(pNum);
200 nbrNum.reset();
201 int sum1 = nbr.elementSize();
202 if (sum1 > 0)
203 {
204 //printf("one %d %d\n", nbrNum.size(), p_num);
205 cuExecute(pNum,
206 NTQ_Narrow_Count,
207 nbr,
208 this->inPosition()->getData(),
209 triVertex,
210 triIndex,
211 nbrNum,
212 this->inRadius()->getData());
213
214 nbrIds.resize(nbrNum);
215
216 int sum = mReduce.accumulate(nbrNum.begin(), nbrNum.size());
217 if (sum > 0)
218 {
219 cuExecute(pNum,
220 NTQ_Narrow_Set,
221 nbr,
222 nbrIds,
223 this->inPosition()->getData(),
224 triVertex,
225 triIndex,
226 this->inRadius()->getData());
227
228 nbrNum.clear();
229 }
230 }
231 }
232
233 DEFINE_CLASS(NeighborTriangleQuery);
234}