TTK
Loading...
Searching...
No Matches
RegularGridTriangulation.cpp
Go to the documentation of this file.
2#include <numeric>
3#include <string>
4
6 : dimensionality_{-1} {
7 setDebugMsgPrefix("RegularGridTriangulation");
8}
9
11 const SimplexId v0, const SimplexId v1) const {
12 // loop over v0 edges to find the one between v0 and v1
13 const auto nEdges = this->getVertexEdgeNumberInternal(v0);
14 for(SimplexId i = 0; i < nEdges; ++i) {
15 SimplexId e{};
16 std::array<SimplexId, 2> eVerts{};
17 this->getVertexEdgeInternal(v0, i, e);
18 this->getEdgeVertexInternal(e, 0, eVerts[0]);
19 this->getEdgeVertexInternal(e, 1, eVerts[1]);
20 if((v0 == eVerts[0] && v1 == eVerts[1])
21 || (v0 == eVerts[1] && v1 == eVerts[0])) {
22 return e;
23 }
24 }
25
26 return -1;
27}
28
30 std::array<SimplexId, 3> &verts) const {
31 std::sort(verts.begin(), verts.end());
32
33 // loop over verts[0] triangles to find the one shared by all 3
34 const auto nTriangles = this->getVertexTriangleNumberInternal(verts[0]);
35 for(SimplexId i = 0; i < nTriangles; ++i) {
36 SimplexId t{};
37 std::array<SimplexId, 3> tVerts{};
38 this->getVertexTriangleInternal(verts[0], i, t);
39 this->getTriangleVertexInternal(t, 0, tVerts[0]);
40 this->getTriangleVertexInternal(t, 1, tVerts[1]);
41 this->getTriangleVertexInternal(t, 2, tVerts[2]);
42 std::sort(tVerts.begin(), tVerts.end());
43 if(tVerts == verts) {
44 return t;
45 }
46 }
47
48 return -1;
49}
50
51#ifdef TTK_ENABLE_MPI
52
53int ttk::RegularGridTriangulation::getVertexRankInternal(
54 const SimplexId lvid) const {
55
56 if(this->vertexGhost_[lvid] == 0) {
57 return ttk::MPIrank_;
58 }
59
60#ifndef TTK_ENABLE_KAMIKAZE
61 if(this->neighborRanks_.empty()) {
62 this->printErr("Empty neighborsRanks_!");
63 return -1;
64 }
65#endif // TTK_ENABLE_KAMIKAZE
66 const auto p{this->getVertGlobalCoords(lvid)};
67 for(const auto neigh : this->neighborRanks_) {
68 const auto &bbox{this->neighborVertexBBoxes_[neigh]};
69 if(p[0] >= bbox[0] && p[0] <= bbox[1] && p[1] >= bbox[2] && p[1] <= bbox[3]
70 && p[2] >= bbox[4] && p[2] <= bbox[5]) {
71 return neigh;
72 }
73 }
74 return -1;
75}
76
77ttk::SimplexId ttk::RegularGridTriangulation::getVertexGlobalIdInternal(
78 const SimplexId lvid) const {
79 if(!ttk::isRunningWithMPI()) {
80 return lvid;
81 }
82
83#ifndef TTK_ENABLE_KAMIKAZE
84 if(lvid > this->TTK_TRIANGULATION_INTERNAL(getNumberOfVertices)() - 1
85 || lvid < 0) {
86 return -1;
87 }
88 if(this->metaGrid_ == nullptr) {
89 return -1;
90 }
91#endif // TTK_ENABLE_KAMIKAZE
92
93 const auto p{this->getVertGlobalCoords(lvid)};
94 const auto &dims{this->metaGrid_->getGridDimensions()};
95
96 // global coordinates to identifier (inverse of vertexToPosition)
97 return p[0] + p[1] * dims[0] + p[2] * dims[0] * dims[1];
98}
99
100ttk::SimplexId ttk::RegularGridTriangulation::getVertexLocalIdInternal(
101 const SimplexId gvid) const {
102 if(!ttk::isRunningWithMPI()) {
103 return gvid;
104 }
105
106#ifndef TTK_ENABLE_KAMIKAZE
107 if(this->metaGrid_ == nullptr) {
108 return -1;
109 }
110 if(gvid
111 > this->metaGrid_->TTK_TRIANGULATION_INTERNAL(getNumberOfVertices)() - 1
112 || gvid < 0) {
113 return -1;
114 }
115#endif // TTK_ENABLE_KAMIKAZE
116
117 const auto p{this->getVertLocalCoords(gvid)};
118 const auto &dims{this->getGridDimensions()};
119
120 if(p[0] < 0 || p[1] < 0 || p[2] < 0 || p[0] > dims[0] - 1
121 || p[1] > dims[1] - 1 || p[2] > dims[2] - 1) {
122 return -1;
123 }
124
125 // local coordinates to identifier (inverse of vertexToPosition)
126 return p[0] + p[1] * dims[0] + p[2] * dims[0] * dims[1];
127}
128
129ttk::SimplexId ttk::RegularGridTriangulation::getCellGlobalIdInternal(
130 const SimplexId lcid) const {
131 if(!ttk::isRunningWithMPI()) {
132 return lcid;
133 }
134
135#ifndef TTK_ENABLE_KAMIKAZE
136 if(lcid > this->TTK_TRIANGULATION_INTERNAL(getNumberOfCells)() - 1
137 || lcid < 0) {
138 return -1;
139 }
140 if(this->metaGrid_ == nullptr) {
141 return -1;
142 }
143#endif // TTK_ENABLE_KAMIKAZE
144
145 // local cube coordinates
146 std::array<SimplexId, 3> p{};
147 if(this->dimensionality_ == 3) {
148 this->tetrahedronToPosition(lcid, p.data());
149 } else if(this->dimensionality_ == 2) {
150 this->triangleToPosition2d(lcid, p.data());
151 }
152
153 // global cube coordinates
154 p[0] += this->localGridOffset_[0];
155 p[1] += this->localGridOffset_[1];
156 p[2] += this->localGridOffset_[2];
157
158 const auto &dims{this->metaGrid_->getGridDimensions()};
159
160 // global coordinates to identifier (inverse of tetrahedronToPosition)
161 const auto globCubeId{p[0] + p[1] * (dims[0] - 1)
162 + p[2] * (dims[0] - 1) * (dims[1] - 1)};
163
164 const auto nCellsPerCube{this->dimensionality_ == 3 ? 6 : 2};
165 return globCubeId * nCellsPerCube + lcid % nCellsPerCube;
166}
167
168ttk::SimplexId ttk::RegularGridTriangulation::getCellLocalIdInternal(
169 const SimplexId gcid) const {
170 if(!ttk::isRunningWithMPI()) {
171 return gcid;
172 }
173
174#ifndef TTK_ENABLE_KAMIKAZE
175 if(this->metaGrid_ == nullptr) {
176 return -1;
177 }
178 if(gcid > this->metaGrid_->TTK_TRIANGULATION_INTERNAL(getNumberOfCells)() - 1
179 || gcid < 0) {
180 return -1;
181 }
182#endif // TTK_ENABLE_KAMIKAZE
183
184 // global cube coordinates
185 std::array<SimplexId, 3> p{};
186 if(this->dimensionality_ == 3) {
187 this->metaGrid_->tetrahedronToPosition(gcid, p.data());
188 } else if(this->dimensionality_ == 2) {
189 this->metaGrid_->triangleToPosition2d(gcid, p.data());
190 }
191
192 // local cube coordinates
193 p[0] -= this->localGridOffset_[0];
194 p[1] -= this->localGridOffset_[1];
195 p[2] -= this->localGridOffset_[2];
196
197 const auto &dims{this->getGridDimensions()};
198
199 // local coordinates to identifier (inverse of tetrahedronToPosition)
200 const auto locCubeId{p[0] + p[1] * (dims[0] - 1)
201 + p[2] * (dims[0] - 1) * (dims[1] - 1)};
202
203 const auto nCellsPerCube{this->dimensionality_ == 3 ? 6 : 2};
204 return locCubeId * nCellsPerCube + gcid % nCellsPerCube;
205}
206
207ttk::SimplexId ttk::RegularGridTriangulation::getEdgeGlobalIdInternal(
208 const SimplexId leid) const {
209 if(!ttk::isRunningWithMPI()) {
210 return leid;
211 }
212
213#ifndef TTK_ENABLE_KAMIKAZE
214 if(leid > this->getNumberOfEdgesInternal() - 1 || leid < 0) {
215 return -1;
216 }
217 if(this->metaGrid_ == nullptr) {
218 return -1;
219 }
220#endif // TTK_ENABLE_KAMIKAZE
221
222 if(this->dimensionality_ == 1) {
223 return this->getCellGlobalIdInternal(leid);
224 }
225
226 // local vertices ids
227 SimplexId lv0{}, lv1{};
228 this->getEdgeVertexInternal(leid, 0, lv0);
229 this->getEdgeVertexInternal(leid, 1, lv1);
230
231 // global vertices ids
232 const auto gv0 = this->getVertexGlobalId(lv0);
233 const auto gv1 = this->getVertexGlobalId(lv1);
234 if(gv0 == -1 || gv1 == -1) {
235 return -1;
236 }
237
238 return this->metaGrid_->findEdgeFromVertices(gv0, gv1);
239}
240
241ttk::SimplexId ttk::RegularGridTriangulation::getEdgeLocalIdInternal(
242 const SimplexId geid) const {
243 if(!ttk::isRunningWithMPI()) {
244 return geid;
245 }
246
247#ifndef TTK_ENABLE_KAMIKAZE
248 if(this->metaGrid_ == nullptr) {
249 return -1;
250 }
251 if(geid > this->metaGrid_->getNumberOfEdgesInternal() - 1 || geid < 0) {
252 return -1;
253 }
254#endif // TTK_ENABLE_KAMIKAZE
255
256 if(this->dimensionality_ == 1) {
257 return this->getCellLocalIdInternal(geid);
258 }
259
260 // global vertices ids
261 SimplexId gv0{}, gv1{};
262 this->metaGrid_->getEdgeVertexInternal(geid, 0, gv0);
263 this->metaGrid_->getEdgeVertexInternal(geid, 1, gv1);
264
265 // local vertices ids
266 const auto lv0 = this->getVertexLocalId(gv0);
267 const auto lv1 = this->getVertexLocalId(gv1);
268 if(lv0 == -1 || lv1 == -1) {
269 return -1;
270 }
271
272 return this->findEdgeFromVertices(lv0, lv1);
273}
274
275ttk::SimplexId ttk::RegularGridTriangulation::getTriangleGlobalIdInternal(
276 const SimplexId ltid) const {
277 if(!ttk::isRunningWithMPI()) {
278 return ltid;
279 }
280
281#ifndef TTK_ENABLE_KAMIKAZE
282 if(ltid > this->getNumberOfTrianglesInternal() - 1 || ltid < 0) {
283 return -1;
284 }
285 if(this->metaGrid_ == nullptr) {
286 return -1;
287 }
288#endif // TTK_ENABLE_KAMIKAZE
289
290 if(this->dimensionality_ == 2) {
291 return this->getCellGlobalIdInternal(ltid);
292 }
293
294 // local vertices ids
295 SimplexId lv0{}, lv1{}, lv2{};
296 this->getTriangleVertexInternal(ltid, 0, lv0);
297 this->getTriangleVertexInternal(ltid, 1, lv1);
298 this->getTriangleVertexInternal(ltid, 2, lv2);
299
300 // global vertices ids
301 std::array<SimplexId, 3> globVerts{
302 this->getVertexGlobalId(lv0),
303 this->getVertexGlobalId(lv1),
304 this->getVertexGlobalId(lv2),
305 };
306 for(const auto gv : globVerts) {
307 if(gv == -1) {
308 return -1;
309 }
310 }
311
312 return this->metaGrid_->findTriangleFromVertices(globVerts);
313}
314
315ttk::SimplexId ttk::RegularGridTriangulation::getTriangleLocalIdInternal(
316 const SimplexId gtid) const {
317 if(!ttk::isRunningWithMPI()) {
318 return gtid;
319 }
320
321#ifndef TTK_ENABLE_KAMIKAZE
322 if(this->metaGrid_ == nullptr) {
323 return -1;
324 }
325 if(gtid > this->metaGrid_->getNumberOfTrianglesInternal() - 1 || gtid < 0) {
326 return -1;
327 }
328#endif // TTK_ENABLE_KAMIKAZE
329
330 if(this->dimensionality_ == 2) {
331 return this->getCellGlobalIdInternal(gtid);
332 }
333
334 // local vertices ids
335 SimplexId gv0{}, gv1{}, gv2{};
336 this->metaGrid_->getTriangleVertexInternal(gtid, 0, gv0);
337 this->metaGrid_->getTriangleVertexInternal(gtid, 1, gv1);
338 this->metaGrid_->getTriangleVertexInternal(gtid, 2, gv2);
339
340 // global vertices ids
341 std::array<SimplexId, 3> locVerts{
342 this->getVertexLocalId(gv0),
343 this->getVertexLocalId(gv1),
344 this->getVertexLocalId(gv2),
345 };
346 for(const auto lv : locVerts) {
347 if(lv == -1) {
348 return -1;
349 }
350 }
351
352 return this->findTriangleFromVertices(locVerts);
353}
354
355int ttk::RegularGridTriangulation::preconditionDistributedVertices() {
356 if(this->hasPreconditionedDistributedVertices_) {
357 return 0;
358 }
359 if(!isRunningWithMPI()) {
360 return -1;
361 }
362 if(this->vertexGhost_ == nullptr) {
363 if(ttk::isRunningWithMPI()) {
364 this->printErr("Missing vertex ghost array!");
365 }
366 return -3;
367 }
368 // number of local vertices (with ghost vertices...)
369 const auto nLocVertices{this->getNumberOfVertices()};
370 this->neighborVertexBBoxes_.resize(ttk::MPIsize_);
371 auto &localBBox{this->neighborVertexBBoxes_[ttk::MPIrank_]};
372 // "good" starting values?
373 ttk::SimplexId localBBox_x_min{this->localGridOffset_[0]
374 + this->dimensions_[0]},
375 localBBox_y_min{this->localGridOffset_[1] + this->dimensions_[1]},
376 localBBox_z_min{this->localGridOffset_[2] + this->dimensions_[2]};
377 ttk::SimplexId localBBox_x_max{this->localGridOffset_[0]},
378 localBBox_y_max{this->localGridOffset_[1]},
379 localBBox_z_max{this->localGridOffset_[2]};
380#ifdef TTK_ENABLE_OPENMP
381#pragma omp parallel for reduction( \
382 min \
383 : localBBox_x_min, localBBox_y_min, localBBox_z_min) \
384 reduction(max \
385 : localBBox_x_max, localBBox_y_max, localBBox_z_max)
386#endif
387 for(SimplexId lvid = 0; lvid < nLocVertices; ++lvid) {
388 // only keep non-ghost vertices
389 if(this->vertexGhost_[lvid] != 0) {
390 continue;
391 }
392 // local vertex coordinates
393 std::array<SimplexId, 3> p{};
394 p = this->getVertGlobalCoords(lvid);
395
396 if(p[0] < localBBox_x_min) {
397 localBBox_x_min = p[0];
398 }
399 if(p[0] > localBBox_x_max) {
400 localBBox_x_max = p[0];
401 }
402 if(p[1] < localBBox_y_min) {
403 localBBox_y_min = p[1];
404 }
405 if(p[1] > localBBox_y_max) {
406 localBBox_y_max = p[1];
407 }
408 if(p[2] < localBBox_z_min) {
409 localBBox_z_min = p[2];
410 }
411 if(p[2] > localBBox_z_max) {
412 localBBox_z_max = p[2];
413 }
414 }
415
416 localBBox = {
417 localBBox_x_min, localBBox_x_max, localBBox_y_min,
418 localBBox_y_max, localBBox_z_min, localBBox_z_max,
419 };
420
421 for(size_t i = 0; i < this->neighborRanks_.size(); ++i) {
422 const auto neigh{this->neighborRanks_[i]};
423 MPI_Sendrecv(this->neighborVertexBBoxes_[ttk::MPIrank_].data(), 6,
424 ttk::getMPIType(SimplexId{}), neigh, ttk::MPIrank_,
425 this->neighborVertexBBoxes_[neigh].data(), 6,
426 ttk::getMPIType(SimplexId{}), neigh, neigh, ttk::MPIcomm_,
427 MPI_STATUS_IGNORE);
428 }
429
430 this->hasPreconditionedDistributedVertices_ = true;
431
432 return 0;
433}
434
435int ttk::RegularGridTriangulation::preconditionExchangeGhostCells() {
436
437 if(this->hasPreconditionedExchangeGhostCells_) {
438 return 0;
439 }
440 if(!ttk::hasInitializedMPI()) {
441 return -1;
442 }
443
444 // number of local cells (with ghost cells...)
445 const auto nLocCells{this->getNumberOfCells()};
446
447 int cellRank = 0;
448 for(LongSimplexId lcid = 0; lcid < nLocCells; ++lcid) {
449 cellRank = this->getCellRankInternal(lcid);
450 if(cellRank != ttk::MPIrank_) {
451 // store ghost cell global ids (per rank)
452 this->ghostCellsPerOwner_[cellRank].emplace_back(
453 this->getCellGlobalIdInternal(lcid));
454 }
455 }
456
457 // for each rank, store the global id of local cells that are ghost cells of
458 // other ranks.
459 const auto MIT{ttk::getMPIType(ttk::SimplexId{})};
460 this->remoteGhostCells_.resize(ttk::MPIsize_);
461 // number of owned cells that are ghost cells of other ranks
462 std::vector<size_t> nOwnedGhostCellsPerRank(ttk::MPIsize_);
463
464 for(const auto neigh : this->neighborRanks_) {
465 // 1. send to neigh number of ghost cells owned by neigh
466 const auto nCells{this->ghostCellsPerOwner_[neigh].size()};
467 MPI_Sendrecv(&nCells, 1, ttk::getMPIType(nCells), neigh, ttk::MPIrank_,
468 &nOwnedGhostCellsPerRank[neigh], 1, ttk::getMPIType(nCells),
469 neigh, neigh, ttk::MPIcomm_, MPI_STATUS_IGNORE);
470 this->remoteGhostCells_[neigh].resize(nOwnedGhostCellsPerRank[neigh]);
471
472 // 2. send to neigh list of ghost cells owned by neigh
473 MPI_Sendrecv(this->ghostCellsPerOwner_[neigh].data(),
474 this->ghostCellsPerOwner_[neigh].size(), MIT, neigh,
475 ttk::MPIrank_, this->remoteGhostCells_[neigh].data(),
476 this->remoteGhostCells_[neigh].size(), MIT, neigh, neigh,
477 ttk::MPIcomm_, MPI_STATUS_IGNORE);
478 }
479
480 this->hasPreconditionedExchangeGhostCells_ = true;
481 return 0;
482}
483
484int ttk::RegularGridTriangulation::preconditionExchangeGhostVertices() {
485
486 if(this->hasPreconditionedExchangeGhostVertices_) {
487 return 0;
488 }
489 if((!ttk::hasInitializedMPI()) || (!ttk::isRunningWithMPI())) {
490 return -1;
491 }
492
493 this->ghostVerticesPerOwner_.resize(ttk::MPIsize_);
494
495 // number of local vertices (with ghost vertices...)
496 const auto nLocVertices{this->getNumberOfVertices()};
497
498 for(LongSimplexId lvid = 0; lvid < nLocVertices; ++lvid) {
499 if(this->getVertexRankInternal(lvid) != ttk::MPIrank_) {
500 // store ghost cell global ids (per rank)
501 this->ghostVerticesPerOwner_[this->getVertexRankInternal(lvid)]
502 .emplace_back(this->getVertexGlobalIdInternal(lvid));
503 }
504 }
505
506 // for each rank, store the global id of local cells that are ghost cells of
507 // other ranks.
508 const auto MIT{ttk::getMPIType(ttk::SimplexId{})};
509 this->remoteGhostVertices_.resize(ttk::MPIsize_);
510 // number of owned cells that are ghost cells of other ranks
511 std::vector<size_t> nOwnedGhostVerticesPerRank(ttk::MPIsize_);
512
513 for(const auto neigh : this->neighborRanks_) {
514 // 1. send to neigh number of ghost cells owned by neigh
515 const auto nVerts{this->ghostVerticesPerOwner_[neigh].size()};
516 MPI_Sendrecv(&nVerts, 1, ttk::getMPIType(nVerts), neigh, ttk::MPIrank_,
517 &nOwnedGhostVerticesPerRank[neigh], 1, ttk::getMPIType(nVerts),
518 neigh, neigh, ttk::MPIcomm_, MPI_STATUS_IGNORE);
519 this->remoteGhostVertices_[neigh].resize(nOwnedGhostVerticesPerRank[neigh]);
520
521 // 2. send to neigh list of ghost cells owned by neigh
522 MPI_Sendrecv(this->ghostVerticesPerOwner_[neigh].data(),
523 this->ghostVerticesPerOwner_[neigh].size(), MIT, neigh,
524 ttk::MPIrank_, this->remoteGhostVertices_[neigh].data(),
525 this->remoteGhostVertices_[neigh].size(), MIT, neigh, neigh,
526 ttk::MPIcomm_, MPI_STATUS_IGNORE);
527 }
528
529 this->hasPreconditionedExchangeGhostVertices_ = true;
530 return 0;
531}
532#endif
#define TTK_TRIANGULATION_INTERNAL(NAME)
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
SimplexId findEdgeFromVertices(const SimplexId v0, const SimplexId v1) const
SimplexId findTriangleFromVertices(std::array< SimplexId, 3 > &verts) const
COMMON_EXPORTS int MPIsize_
Definition BaseClass.cpp:10
COMMON_EXPORTS int MPIrank_
Definition BaseClass.cpp:9
long long int LongSimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:15
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22