PeriDyno 1.0.0
Loading...
Searching...
No Matches
HexahedronSet.cu
Go to the documentation of this file.
1#include "HexahedronSet.h"
2#include "QuadSet.h"
3#include <fstream>
4#include <iostream>
5#include <sstream>
6
7#include <thrust/sort.h>
8
9namespace dyno
10{
11 template<typename TDataType>
12 HexahedronSet<TDataType>::HexahedronSet()
13 : QuadSet<TDataType>()
14 {
15
16 }
17
18 template<typename TDataType>
19 HexahedronSet<TDataType>::~HexahedronSet()
20 {
21 }
22
23 template<typename TDataType>
24 void HexahedronSet<TDataType>::setHexahedrons(std::vector<Hexahedron>& hexahedrons)
25 {
26 std::vector<Quad> quads;
27
28 m_hexahedrons.resize(hexahedrons.size());
29 m_hexahedrons.assign(hexahedrons);
30
31 this->updateQuads();
32 }
33
34 template<typename TDataType>
35 void HexahedronSet<TDataType>::setHexahedrons(DArray<Hexahedron>& hexahedrons)
36 {
37 if (hexahedrons.size() != m_hexahedrons.size())
38 {
39 m_hexahedrons.resize(hexahedrons.size());
40 }
41
42 m_hexahedrons.assign(hexahedrons);
43
44 this->updateQuads();
45 }
46
47 template<typename Hexahedron>
48 __global__ void HS_CountHexs(
49 DArray<uint> counter,
50 DArray<Hexahedron> hexs)
51 {
52 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
53 if (tId >= hexs.size()) return;
54
55 Hexahedron t = hexs[tId];
56
57 atomicAdd(&counter[t[0]], 1);
58 atomicAdd(&counter[t[1]], 1);
59 atomicAdd(&counter[t[2]], 1);
60 atomicAdd(&counter[t[3]], 1);
61 atomicAdd(&counter[t[4]], 1);
62 atomicAdd(&counter[t[5]], 1);
63 atomicAdd(&counter[t[6]], 1);
64 atomicAdd(&counter[t[7]], 1);
65 }
66
67 template<typename Hexahedron>
68 __global__ void HS_SetupHexIds(
69 DArrayList<int> hexIds,
70 DArray<Hexahedron> hexs)
71 {
72 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
73 if (tId >= hexs.size()) return;
74
75 Hexahedron t = hexs[tId];
76
77 hexIds[t[0]].atomicInsert(tId);
78 hexIds[t[1]].atomicInsert(tId);
79 hexIds[t[2]].atomicInsert(tId);
80 hexIds[t[3]].atomicInsert(tId);
81 hexIds[t[4]].atomicInsert(tId);
82 hexIds[t[5]].atomicInsert(tId);
83 hexIds[t[6]].atomicInsert(tId);
84 hexIds[t[7]].atomicInsert(tId);
85 }
86
87 template<typename TDataType>
88 DArrayList<int>& HexahedronSet<TDataType>::getVer2Hex()
89 {
90 DArray<uint> counter;
91 counter.resize(this->mCoords.size());
92 counter.reset();
93
94 cuExecute(m_hexahedrons.size(),
95 HS_CountHexs,
96 counter,
97 m_hexahedrons);
98
99 m_ver2Hex.resize(counter);
100
101 counter.reset();
102 cuExecute(m_hexahedrons.size(),
103 HS_SetupHexIds,
104 m_ver2Hex,
105 m_hexahedrons);
106
107 counter.clear();
108
109 return m_ver2Hex;
110 }
111
112 template<typename TDataType>
113 void HexahedronSet<TDataType>::getVolume(DArray<Real>& volume)
114 {
115
116 }
117
118 template<typename QKey, typename Hexahedron>
119 __global__ void HS_SetupKeys(
120 DArray<QKey> keys,
121 DArray<int> ids,
122 DArray<Hexahedron> hexs)
123 {
124 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
125 if (tId >= hexs.size()) return;
126
127 Hexahedron hex = hexs[tId];
128 keys[6 * tId] = QKey(hex[0], hex[1], hex[2], hex[3]);
129 keys[6 * tId + 1] = QKey(hex[4], hex[5], hex[6], hex[7]);
130 keys[6 * tId + 2] = QKey(hex[0], hex[1], hex[5], hex[4]);
131 keys[6 * tId + 3] = QKey(hex[1], hex[2], hex[6], hex[5]);
132 keys[6 * tId + 4] = QKey(hex[3], hex[2], hex[6], hex[7]);
133 keys[6 * tId + 5] = QKey(hex[0], hex[3], hex[7], hex[4]);
134
135 ids[6 * tId] = tId;
136 ids[6 * tId + 1] = tId;
137 ids[6 * tId + 2] = tId;
138 ids[6 * tId + 3] = tId;
139 ids[6 * tId + 4] = tId;
140 ids[6 * tId + 5] = tId;
141 }
142
143 template<typename QKey>
144 __global__ void HS_CountQuadNumber(
145 DArray<int> counter,
146 DArray<QKey> keys)
147 {
148 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
149 if (tId >= keys.size()) return;
150
151 if (tId == 0 || keys[tId] != keys[tId - 1])
152 counter[tId] = 1;
153 else
154 counter[tId] = 0;
155 }
156
157 template<typename Quad, typename Quad2Hex, typename QKey>
158 __global__ void HS_SetupQuads(
159 DArray<Quad> quads,
160 DArray<Quad2Hex> quad2Hex,
161 DArray<QKey> keys,
162 DArray<int> counter,
163 DArray<int> hexIds)
164 {
165 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
166 if (tId >= keys.size()) return;
167
168 int shift = counter[tId];
169 if (tId == 0 || keys[tId] != keys[tId - 1])
170 {
171 QKey key = keys[tId];
172 quads[shift] = Quad(key[0], key[1], key[2], key[3]);
173
174 Quad2Hex q2H(EMPTY, EMPTY);
175 q2H[0] = hexIds[tId];
176
177 if (tId + 1 < keys.size() && keys[tId + 1] == key)
178 q2H[1] = hexIds[tId + 1];
179
180 quad2Hex[shift] = q2H;
181
182 }
183 }
184
185 template<typename QKey>
186 void printTKey(DArray<QKey> keys, int maxLength) {
187 CArray<QKey> h_keys;
188 h_keys.resize(keys.size());
189 h_keys.assign(keys);
190
191 int psize = min((int)h_keys.size(), maxLength);
192 for (int i = 0; i < psize; i++)
193 {
194 printf("%d: %d %d %d %d \n", i, h_keys[i][0], h_keys[i][1], h_keys[i][2], h_keys[i][3]);
195 }
196
197 h_keys.clear();
198 }
199
200 /*void printCount(DArray<int> keys, int maxLength) {
201 CArray<int> h_keys;
202 h_keys.resize(keys.size());
203 h_keys.assign(keys);
204
205 int psize = minimum((int)h_keys.size(), maxLength);
206 for (int i = 0; i < psize; i++)
207 {
208 printf("%d: %d \n", i, h_keys[i]);
209 }
210
211 h_keys.clear();
212 }*/
213
214 template<typename TDataType>
215 void HexahedronSet<TDataType>::updateQuads()
216 {
217 uint hexSize = m_hexahedrons.size();
218
219 DArray<QKey> keys;
220 DArray<int> hexIds;
221
222 keys.resize(6 * hexSize);
223 hexIds.resize(6 * hexSize);
224
225 cuExecute(hexSize,
226 HS_SetupKeys,
227 keys,
228 hexIds,
229 m_hexahedrons);
230
231 thrust::sort_by_key(thrust::device, keys.begin(), keys.begin() + keys.size(), hexIds.begin());
232
233 DArray<int> counter;
234 counter.resize(6 * hexSize);
235
236 cuExecute(keys.size(),
237 HS_CountQuadNumber,
238 counter,
239 keys);
240
241 int quadNum = thrust::reduce(thrust::device, counter.begin(), counter.begin() + counter.size());
242 thrust::exclusive_scan(thrust::device, counter.begin(), counter.begin() + counter.size(), counter.begin());
243
244 quad2Hex.resize(quadNum);
245
246 auto& pQuad = this->getQuads();
247 pQuad.resize(quadNum);
248 cuExecute(keys.size(),
249 HS_SetupQuads,
250 pQuad,
251 quad2Hex,
252 keys,
253 counter,
254 hexIds);
255
256 counter.clear();
257 hexIds.clear();
258 keys.clear();
259
260
261// this->updateTriangles();
262// this->updateEdges();
263 }
264
265
266 template<typename TDataType>
267 void HexahedronSet<TDataType>::copyFrom(HexahedronSet<TDataType> hexSet)
268 {
269 m_hexahedrons.resize(hexSet.m_hexahedrons.size());
270 m_hexahedrons.assign(hexSet.m_hexahedrons);
271
272 quad2Hex.resize(hexSet.quad2Hex.size());
273 quad2Hex.assign(hexSet.quad2Hex);
274
275 m_ver2Hex.assign(hexSet.m_ver2Hex);
276
277 QuadSet<TDataType>::copyFrom(hexSet);
278 }
279
280 DEFINE_CLASS(HexahedronSet);
281}