TTK
Loading...
Searching...
No Matches
ExplicitTriangulation.cpp
Go to the documentation of this file.
2
3#include <OneSkeleton.h>
4#include <ThreeSkeleton.h>
5#include <TwoSkeleton.h>
6#include <ZeroSkeleton.h>
7
8#include <cstring>
9#include <numeric>
10
11using namespace ttk;
12
14
15 setDebugMsgPrefix("ExplicitTriangulation");
16
17 clear();
18}
19
21
23 vertexNumber_ = 0;
24 cellNumber_ = 0;
25 doublePrecision_ = false;
26
27 printMsg("Triangulation cleared.", debug::Priority::DETAIL);
28
30}
31
32size_t ExplicitTriangulation::footprint(size_t size) const {
33
34 const auto printArrayFootprint
35 = [this](const FlatJaggedArray &array, const std::string &name) {
36 if(!array.empty() && !name.empty()) {
37 this->printMsg(name + std::string{": "}
38 + std::to_string(array.footprint()) + " bytes");
39 }
40 return array.footprint();
41 };
42
43 size += printArrayFootprint(vertexNeighborData_, "vertexNeighborData_");
44 size += printArrayFootprint(cellNeighborData_, "cellNeighborData_");
45 size += printArrayFootprint(vertexEdgeData_, "vertexEdgeData_");
46 size += printArrayFootprint(vertexTriangleData_, "vertexTriangleData_");
47 size += printArrayFootprint(edgeTriangleData_, "edgeTriangleData_");
48 size += printArrayFootprint(vertexStarData_, "vertexStarData_");
49 size += printArrayFootprint(edgeStarData_, "edgeStarData_");
50 size += printArrayFootprint(triangleStarData_, "triangleStarData_");
51 size += printArrayFootprint(vertexLinkData_, "vertexLinkData_");
52 size += printArrayFootprint(edgeLinkData_, "edgeLinkData_");
53 size += printArrayFootprint(triangleLinkData_, "triangleLinkData_");
54
56}
57
59
60 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
61#ifdef TTK_ENABLE_MPI
62 if(!(ttk::isRunningWithMPI()))
63#endif
64 this->printErr("Empty dataset, precondition skipped");
65 return 1;
66 }
67
68 if((!boundaryEdges_.empty()) && (boundaryEdges_.size() == edgeList_.size())) {
69 return 0;
70 }
71
72 Timer tm{};
73
75 boundaryEdges_.resize(edgeList_.size(), false);
76
77 if(getDimensionality() == 2) {
79 for(SimplexId i = 0; i < (SimplexId)edgeStarData_.size(); i++) {
80 if(edgeStarData_[i].size() == 1) {
81 boundaryEdges_[i] = true;
82 }
83 }
84 } else if(getDimensionality() == 3) {
87
88 for(size_t i = 0; i < triangleStarData_.size(); i++) {
89 if(triangleStarData_[i].size() == 1) {
90 for(int j = 0; j < 3; j++) {
92 }
93 }
94 }
95 } else {
96 // unsupported dimension
97 printErr("Unsupported dimension for edge boundary precondition");
98 return -1;
99 }
100
101#ifdef TTK_ENABLE_MPI
102
103 if(ttk::isRunningWithMPI()) {
104 ttk::SimplexId edgeNumber = edgeList_.size();
105 std::vector<unsigned char> charBoundary(edgeNumber, false);
106 for(int i = 0; i < edgeNumber; ++i) {
107 charBoundary[i] = boundaryEdges_[i] ? '1' : '0';
108 }
109 ttk::exchangeGhostDataWithoutTriangulation(
110 charBoundary.data(),
111 [this](const SimplexId a) { return this->edgeRankArray_[a]; },
112 [this](const SimplexId a) { return this->edgeLidToGid_[a]; },
113 [this](const SimplexId a) { return this->edgeGidToLid_[a]; }, edgeNumber,
114 ttk::MPIcomm_, this->getNeighborRanks());
115 for(int i = 0; i < edgeNumber; ++i) {
116 boundaryEdges_[i] = (charBoundary[i] == '1');
117 }
118 }
119
120#endif // TTK_ENABLE_MPI
121
122 this->printMsg("Extracted boundary edges", 1.0, tm.getElapsedTime(), 1);
123
124 return 0;
125}
126
127#ifdef TTK_ENABLE_MPI
128
129int ExplicitTriangulation::preconditionVertexRankArray() {
130 if(vertexRankArray_.size() == 0) {
131 this->vertexRankArray_.resize(this->vertexNumber_, 0);
132 if(ttk::isRunningWithMPI()) {
133 ttk::produceRankArray(this->vertexRankArray_, this->vertGid_,
134 this->vertexGhost_, this->vertexNumber_,
135 this->boundingBox_.data(), this->neighborRanks_,
136 this->neighborsToId_);
137 ttk::preconditionNeighborsUsingRankArray<ttk::SimplexId>(
138 this->neighborRanks_, this->neighborsToId_,
139 [this](const ttk::SimplexId a) { return this->getVertexRank(a); },
140 this->vertexNumber_, ttk::MPIcomm_);
141 }
142 }
143 return 0;
144}
145
146int ExplicitTriangulation::preconditionCellRankArray() {
147 if(cellRankArray_.size() == 0) {
148 this->cellRankArray_.resize(this->cellNumber_, 0);
149 if(ttk::isRunningWithMPI()) {
150 ttk::produceRankArray(this->cellRankArray_, this->cellGid_,
151 this->cellGhost_, this->cellNumber_,
152 this->boundingBox_.data(), this->neighborRanks_,
153 this->neighborsToId_);
154 }
155 }
156 return 0;
157}
158
159int ExplicitTriangulation::preconditionEdgeRankArray() {
160 ttk::SimplexId edgeNumber = this->getNumberOfEdgesInternal();
161 edgeRankArray_.resize(edgeNumber, 0);
162 if(ttk::isRunningWithMPI()) {
163 ttk::SimplexId min_id;
164 for(ttk::SimplexId id = 0; id < edgeNumber; id++) {
165 this->TTK_TRIANGULATION_INTERNAL(getEdgeStar)(id, 0, min_id);
166 const auto nStar{this->TTK_TRIANGULATION_INTERNAL(getEdgeStarNumber)(id)};
167 for(SimplexId i = 1; i < nStar; ++i) {
168 SimplexId sid{-1};
169 this->TTK_TRIANGULATION_INTERNAL(getEdgeStar)(id, i, sid);
170 // rule: an edge is owned by the cell in its star with the
171 // lowest rank
172 if(this->cellRankArray_[sid] < this->cellRankArray_[min_id]) {
173 min_id = sid;
174 }
175 }
176 edgeRankArray_[id] = cellRankArray_[min_id];
177 }
178 }
179 return 0;
180}
181
182int ExplicitTriangulation::preconditionTriangleRankArray() {
183 ttk::SimplexId triangleNumber = this->getNumberOfTrianglesInternal();
184 triangleRankArray_.resize(triangleNumber, 0);
185 if(ttk::isRunningWithMPI()) {
186 ttk::SimplexId min_id;
187 for(ttk::SimplexId id = 0; id < triangleNumber; id++) {
188 this->TTK_TRIANGULATION_INTERNAL(getTriangleStar)(id, 0, min_id);
189 const auto nStar{
191 for(SimplexId i = 1; i < nStar; ++i) {
192 SimplexId sid{-1};
194 // rule: an triangle is owned by the cell in its star with the
195 // lowest rank
196 if(this->cellRankArray_[sid] < this->cellRankArray_[min_id]) {
197 min_id = sid;
198 }
199 }
200 triangleRankArray_[id] = cellRankArray_[min_id];
201 }
202 }
203 return 0;
204}
205
206#endif // TTK_ENABLE_MPI
207
209
210 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
211#ifdef TTK_ENABLE_MPI
212 if(!(ttk::isRunningWithMPI()))
213#endif
214 this->printErr("Empty dataset, precondition skipped");
215 return 1;
216 }
217
218 if(getDimensionality() == 2)
219 return 0;
220
221 Timer tm{};
222
223 if((!boundaryTriangles_.empty())
224 && (boundaryTriangles_.size() == triangleList_.size())) {
225 return 0;
226 }
227
229 boundaryTriangles_.resize(triangleList_.size(), false);
230
231 if(getDimensionality() == 3) {
233
234 for(size_t i = 0; i < triangleStarData_.size(); i++) {
235 if(triangleStarData_[i].size() == 1) {
236 boundaryTriangles_[i] = true;
237 }
238 }
239 } else {
240 // unsupported dimension
241 printErr("Unsupported dimension for triangle boundary precondition");
242 return -1;
243 }
244
245#ifdef TTK_ENABLE_MPI
246
247 if(ttk::isRunningWithMPI()) {
248 ttk::SimplexId triangleNumber = triangleList_.size();
249 std::vector<unsigned char> charBoundary(triangleNumber, false);
250 for(int i = 0; i < triangleNumber; ++i) {
251 charBoundary[i] = boundaryTriangles_[i] ? '1' : '0';
252 }
253 ttk::exchangeGhostDataWithoutTriangulation(
254 charBoundary.data(),
255 [this](const SimplexId a) { return this->triangleRankArray_[a]; },
256 [this](const SimplexId a) { return this->triangleLidToGid_[a]; },
257 [this](const SimplexId a) { return this->triangleGidToLid_[a]; },
258 triangleNumber, ttk::MPIcomm_, this->getNeighborRanks());
259 for(int i = 0; i < triangleNumber; ++i) {
260 boundaryTriangles_[i] = (charBoundary[i] == '1');
261 }
262 }
263
264#endif // TTK_ENABLE_MPI
265
266 this->printMsg("Extracted boundary triangles", 1.0, tm.getElapsedTime(), 1);
267
268 return 0;
269}
270
272
273 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
274#ifdef TTK_ENABLE_MPI
275 if(!(ttk::isRunningWithMPI()))
276#endif
277 this->printErr("Empty dataset, precondition skipped");
278 return 1;
279 }
280
281 if((!boundaryVertices_.empty())
282 && ((SimplexId)boundaryVertices_.size() == vertexNumber_))
283 return 0;
284
285 Timer tm{};
286
287 boundaryVertices_.resize(vertexNumber_, false);
288
289 // create the list of boundary elements
290 // create their star
291 // look for singletons
292 if(getDimensionality() == 1) {
294 for(size_t i = 0; i < vertexStarData_.size(); i++) {
295 if(vertexStarData_[i].size() == 1) {
296 boundaryVertices_[i] = true;
297 }
298 }
299 } else if(getDimensionality() == 2) {
302
303 for(SimplexId i = 0; i < (SimplexId)edgeStarData_.size(); i++) {
304 if(edgeStarData_[i].size() == 1) {
305 boundaryVertices_[edgeList_[i][0]] = true;
306 boundaryVertices_[edgeList_[i][1]] = true;
307 }
308 }
309 } else if(getDimensionality() == 3) {
312
313 for(size_t i = 0; i < triangleStarData_.size(); i++) {
314 if(triangleStarData_[i].size() == 1) {
315 boundaryVertices_[triangleList_[i][0]] = true;
316 boundaryVertices_[triangleList_[i][1]] = true;
317 boundaryVertices_[triangleList_[i][2]] = true;
318 }
319 }
320 } else {
321 // unsupported dimension
322 printErr("Unsupported dimension for vertex boundary precondition");
323 return -1;
324 }
325
326#ifdef TTK_ENABLE_MPI
327
328 if(ttk::isRunningWithMPI()) {
329 this->preconditionDistributedVertices();
330 std::vector<unsigned char> charBoundary(vertexNumber_, false);
331 for(int i = 0; i < vertexNumber_; ++i) {
332 charBoundary[i] = boundaryVertices_[i] ? '1' : '0';
333 }
334 ttk::exchangeGhostVertices<unsigned char, ExplicitTriangulation>(
335 charBoundary.data(), this, ttk::MPIcomm_);
336
337 for(int i = 0; i < vertexNumber_; ++i) {
338 if(vertexRankArray_[i] != ttk::MPIrank_) {
339 boundaryVertices_[i] = (charBoundary[i] == '1');
340 }
341 }
342 }
343
344#endif // TTK_ENABLE_MPI
345
346 this->printMsg("Extracted boundary vertices", 1.0, tm.getElapsedTime(), 1);
347
348 return 0;
349}
350
352
353 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
354#ifdef TTK_ENABLE_MPI
355 if(!(ttk::isRunningWithMPI()))
356#endif
357 this->printErr("Empty dataset, precondition skipped");
358 return 1;
359 }
360
361 if((tetraEdgeList_.empty() && getDimensionality() == 3)
362 || (triangleEdgeList_.empty() && getDimensionality() == 2)) {
364 }
365
366 return 0;
367}
368
370
371 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
372#ifdef TTK_ENABLE_MPI
373 if(!(ttk::isRunningWithMPI()))
374#endif
375 this->printErr("Empty dataset, precondition skipped");
376 return 1;
377 }
378
379 if(cellNeighborData_.empty()) {
380 if(getDimensionality() == 3) {
381 ThreeSkeleton threeSkeleton;
382 threeSkeleton.setWrapper(this);
384 vertexNumber_, *cellArray_, cellNeighborData_, &triangleStarData_);
385 } else if(getDimensionality() == 2) {
387 TwoSkeleton twoSkeleton;
388 twoSkeleton.setWrapper(this);
389 twoSkeleton.buildCellNeighborsFromEdges(
390 *cellArray_, cellNeighborData_, edgeStarData_);
391 }
392 }
393
394 return 0;
395}
396
398
399 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
400#ifdef TTK_ENABLE_MPI
401 if(!(ttk::isRunningWithMPI()))
402#endif
403 this->printErr("Empty dataset, precondition skipped");
404 return 1;
405 }
406
407 if(tetraTriangleList_.empty()) {
408
409 TwoSkeleton twoSkeleton;
410 twoSkeleton.setWrapper(this);
411
412 if(triangleList_.size()) {
413 // we already computed this guy, let's just get the cell triangles
414 if(!triangleStarData_.empty()) {
415 return twoSkeleton.buildTriangleList(
416 vertexNumber_, *cellArray_, nullptr, nullptr, &tetraTriangleList_);
417 } else {
418 // let's compute the triangle star while we're at it...
419 // it's just a tiny overhead.
420 return twoSkeleton.buildTriangleList(vertexNumber_, *cellArray_,
421 nullptr, &triangleStarData_,
423 }
424 } else {
425 // we have not computed this guy, let's do it while we're at it
426 if(!triangleStarData_.empty()) {
427 return twoSkeleton.buildTriangleList(vertexNumber_, *cellArray_,
428 &triangleList_, nullptr,
430 } else {
431 // let's compute the triangle star while we're at it...
432 // it's just a tiny overhead.
433 return twoSkeleton.buildTriangleList(vertexNumber_, *cellArray_,
434 &triangleList_, &triangleStarData_,
436 }
437 }
438 }
439
440 return 0;
441}
442
444
445 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
446#ifdef TTK_ENABLE_MPI
447 if(!(ttk::isRunningWithMPI()))
448#endif
449 this->printErr("Empty dataset, precondition skipped");
450 return 1;
451 }
452
453 if(edgeList_.empty()) {
454 OneSkeleton oneSkeleton;
455 oneSkeleton.setWrapper(this);
456 // also computes edgeStar and triangleEdge / tetraEdge lists for free...
457 int ret{};
458 if(getDimensionality() == 1) {
459 std::vector<std::array<SimplexId, 1>> tmp{};
460 return oneSkeleton.buildEdgeList<1>(
461 vertexNumber_, *cellArray_, edgeList_, edgeStarData_, tmp);
462 } else if(getDimensionality() == 2) {
463 ret = oneSkeleton.buildEdgeList(vertexNumber_, *cellArray_, edgeList_,
464 edgeStarData_, triangleEdgeList_);
465 } else if(getDimensionality() == 3) {
466 ret = oneSkeleton.buildEdgeList(
467 vertexNumber_, *cellArray_, edgeList_, edgeStarData_, tetraEdgeList_);
468 }
469
470 if(ret != 0) {
471 return ret;
472 }
473
474#ifdef TTK_ENABLE_MPI
475 if(this->getDimensionality() == 2 || this->getDimensionality() == 3) {
476 return this->preconditionDistributedEdges();
477 }
478#endif // TTK_ENABLE_MPI
479 }
480
481 return 0;
482}
483
485
486 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
487#ifdef TTK_ENABLE_MPI
488 if(!(ttk::isRunningWithMPI()))
489#endif
490 this->printErr("Empty dataset, precondition skipped");
491 return 1;
492 }
493
494 if(edgeLinkData_.empty()) {
495
496 if(getDimensionality() == 2) {
499
500 OneSkeleton oneSkeleton;
501 oneSkeleton.setWrapper(this);
502 return oneSkeleton.buildEdgeLinks(
503 edgeList_, edgeStarData_, *cellArray_, edgeLinkData_);
504 } else if(getDimensionality() == 3) {
508
509 OneSkeleton oneSkeleton;
510 oneSkeleton.setWrapper(this);
511 return oneSkeleton.buildEdgeLinks(
512 edgeList_, edgeStarData_, tetraEdgeList_, edgeLinkData_);
513 } else {
514 // unsupported dimension
515 printErr("Unsupported dimension for edge link precondition");
516 return -1;
517 }
518 }
519
520 return 0;
521}
522
524
525 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
526#ifdef TTK_ENABLE_MPI
527 if(!(ttk::isRunningWithMPI()))
528#endif
529 this->printErr("Empty dataset, precondition skipped");
530 return 1;
531 }
532
533 if(edgeStarData_.empty()) {
535 }
536 return 0;
537}
538
540
541 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
542#ifdef TTK_ENABLE_MPI
543 if(!(ttk::isRunningWithMPI()))
544#endif
545 this->printErr("Empty dataset, precondition skipped");
546 return 1;
547 }
548
549 if(edgeTriangleData_.empty()) {
552
553 TwoSkeleton twoSkeleton;
554 twoSkeleton.setWrapper(this);
555 return twoSkeleton.buildEdgeTriangles(vertexNumber_, *cellArray_,
556 edgeTriangleData_, edgeList_,
558 }
559
560 return 0;
561}
562
564
565 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
566#ifdef TTK_ENABLE_MPI
567 if(!(ttk::isRunningWithMPI()))
568#endif
569 this->printErr("Empty dataset, precondition skipped");
570 return 1;
571 }
572
573 if(triangleList_.empty()) {
574
575 TwoSkeleton twoSkeleton;
576 twoSkeleton.setWrapper(this);
577
578 twoSkeleton.buildTriangleList(vertexNumber_, *cellArray_, &triangleList_,
579 &triangleStarData_, &tetraTriangleList_);
580
581#ifdef TTK_ENABLE_MPI
582 this->preconditionDistributedTriangles();
583#endif // TTK_ENABLE_MPI
584 }
585
586 return 0;
587}
588
590
591 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
592#ifdef TTK_ENABLE_MPI
593 if(!(ttk::isRunningWithMPI()))
594#endif
595 this->printErr("Empty dataset, precondition skipped");
596 return 1;
597 }
598
599 if(triangleEdgeList_.empty()) {
601
602 // WARNING
603 // here triangleStarList and cellTriangleList will be computed (for
604 // free) although they are not requireed to get the edgeTriangleList.
605 // if memory usage is an issue, please change these pointers by nullptr.
606
607 TwoSkeleton twoSkeleton;
608 twoSkeleton.setWrapper(this);
609
610 return twoSkeleton.buildTriangleEdgeList(
611 vertexNumber_, *cellArray_, triangleEdgeList_, edgeList_,
612 &vertexEdgeData_, &triangleList_, &triangleStarData_,
614 }
615
616 return 0;
617}
618
620
621 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
622#ifdef TTK_ENABLE_MPI
623 if(!(ttk::isRunningWithMPI()))
624#endif
625 this->printErr("Empty dataset, precondition skipped");
626 return 1;
627 }
628
629 if(triangleLinkData_.empty()) {
630
632
633 TwoSkeleton twoSkeleton;
634 twoSkeleton.setWrapper(this);
635 return twoSkeleton.buildTriangleLinks(
636 triangleList_, triangleStarData_, *cellArray_, triangleLinkData_);
637 }
638
639 return 0;
640}
641
643
644 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
645#ifdef TTK_ENABLE_MPI
646 if(!(ttk::isRunningWithMPI()))
647#endif
648 this->printErr("Empty dataset, precondition skipped");
649 return 1;
650 }
651
652 if(triangleStarData_.empty()) {
653
654 TwoSkeleton twoSkeleton;
655 twoSkeleton.setWrapper(this);
656 return twoSkeleton.buildTriangleList(
657 vertexNumber_, *cellArray_, &triangleList_, &triangleStarData_);
658 }
659
660 return 0;
661}
662
664
665 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
666#ifdef TTK_ENABLE_MPI
667 if(!(ttk::isRunningWithMPI()))
668#endif
669 this->printErr("Empty dataset, precondition skipped");
670 return 1;
671 }
672
673 if((SimplexId)vertexEdgeData_.size() != vertexNumber_) {
674 ZeroSkeleton zeroSkeleton;
675
676 if(edgeList_.empty()) {
678 }
679
680 zeroSkeleton.setWrapper(this);
681 return zeroSkeleton.buildVertexEdges(
682 vertexNumber_, edgeList_, vertexEdgeData_);
683 }
684 return 0;
685}
686
688
689 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
690#ifdef TTK_ENABLE_MPI
691 if(!(ttk::isRunningWithMPI()))
692#endif
693 this->printErr("Empty dataset, precondition skipped");
694 return 1;
695 }
696
697 if((SimplexId)vertexLinkData_.size() != vertexNumber_) {
698
699 if(getDimensionality() == 2) {
702
703 ZeroSkeleton zeroSkeleton;
704 zeroSkeleton.setWrapper(this);
705 return zeroSkeleton.buildVertexLinks(
706 vertexStarData_, triangleEdgeList_, edgeList_, vertexLinkData_);
707 } else if(getDimensionality() == 3) {
710
711 ZeroSkeleton zeroSkeleton;
712 zeroSkeleton.setWrapper(this);
713 return zeroSkeleton.buildVertexLinks(
714 vertexStarData_, tetraTriangleList_, triangleList_, vertexLinkData_);
715 } else {
716 // unsupported dimension
717 printErr("Unsupported dimension for vertex link precondition");
718 return -1;
719 }
720 }
721 return 0;
722}
723
725
726 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
727#ifdef TTK_ENABLE_MPI
728 if(!(ttk::isRunningWithMPI()))
729#endif
730 this->printErr("Empty dataset, precondition skipped");
731 return 1;
732 }
733
734 if((SimplexId)vertexNeighborData_.size() != vertexNumber_) {
736 ZeroSkeleton zeroSkeleton;
737 zeroSkeleton.setWrapper(this);
738 return zeroSkeleton.buildVertexNeighbors(
739 vertexNumber_, vertexNeighborData_, edgeList_);
740 }
741 return 0;
742}
743
745
746 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
747#ifdef TTK_ENABLE_MPI
748 if(!(ttk::isRunningWithMPI()))
749#endif
750 this->printErr("Empty dataset, precondition skipped");
751 return 1;
752 }
753
754 if((SimplexId)vertexStarData_.size() != vertexNumber_) {
755 ZeroSkeleton zeroSkeleton;
756 zeroSkeleton.setWrapper(this);
757
758 return zeroSkeleton.buildVertexStars(
759 vertexNumber_, *cellArray_, vertexStarData_);
760 }
761 return 0;
762}
763
765
766 if(this->cellArray_ == nullptr || this->vertexNumber_ == 0) {
767#ifdef TTK_ENABLE_MPI
768 if(!(ttk::isRunningWithMPI()))
769#endif
770 this->printErr("Empty dataset, precondition skipped");
771 return 1;
772 }
773
774 if((SimplexId)vertexTriangleData_.size() != vertexNumber_) {
775
777
778 TwoSkeleton twoSkeleton;
779 twoSkeleton.setWrapper(this);
780
781 twoSkeleton.buildVertexTriangles(
782 vertexNumber_, triangleList_, vertexTriangleData_);
783 }
784
785 return 0;
786}
787
789
790 // quick check by numbering (d-1)-simplices star
791 FlatJaggedArray *simplexStar{};
792
793 if(this->getDimensionality() == 3) {
794 this->preconditionTriangleStarsInternal();
795 simplexStar = &this->triangleStarData_;
796 } else if(this->getDimensionality() == 2) {
797 this->preconditionEdgeStarsInternal();
798 simplexStar = &this->edgeStarData_;
799 } else if(this->getDimensionality() == 1) {
800 this->preconditionVertexStarsInternal();
801 simplexStar = &this->vertexStarData_;
802 }
803
804 if(simplexStar == nullptr) {
805 return 0;
806 }
807
808 for(const auto star : *simplexStar) {
809 if(star.size() < 1 || star.size() > 2) {
810 this->isManifold_ = false;
811 this->printWrn("Non manifold data-set detected");
812 break;
813 }
814 }
815
816 return 0;
817}
818
819#ifdef TTK_ENABLE_MPI
820
821int ExplicitTriangulation::preconditionDistributedCells() {
822 if(this->hasPreconditionedDistributedCells_) {
823 return 0;
824 }
825 if(!ttk::hasInitializedMPI()) {
826 return -1;
827 }
828 if(this->cellGid_ == nullptr) {
829 this->printErr("Missing global cell identifiers array!");
830 return -2;
831 }
832
833 Timer tm{};
834
835 this->preconditionCellRankArray();
836
837 // number of local cells (with ghost cells...)
838 const auto nLocCells{this->getNumberOfCells()};
839
840 // global cell id -> local cell id (reverse of this->cellGid_)
841 this->cellGidToLid_.reserve(nLocCells);
842 for(LongSimplexId lcid = 0; lcid < nLocCells; ++lcid) {
843 this->cellGidToLid_[this->cellGid_[lcid]] = lcid;
844 }
845
846 this->preconditionExchangeGhostCells();
847
848 this->preconditionDistributedCellRanges();
849
850 this->hasPreconditionedDistributedCells_ = true;
851
852 if(ttk::MPIrank_ == 0) {
853 this->printMsg("Domain contains "
854 + std::to_string(this->gatheredCellRanges_.back().end + 1)
855 + " cells",
856 1.0, tm.getElapsedTime(), this->threadNumber_);
857 }
858
859 return 0;
860}
861
862int ttk::ExplicitTriangulation::preconditionExchangeGhostCells() {
863
864 if(this->hasPreconditionedExchangeGhostCells_) {
865 return 0;
866 }
867 if(!ttk::hasInitializedMPI()) {
868 return -1;
869 }
870 if(this->cellGid_ == nullptr) {
871 this->printErr("Missing global cell identifiers array!");
872 return -2;
873 }
874
875 this->ghostCellsPerOwner_.resize(ttk::MPIsize_);
876
877 // number of local cells (with ghost cells...)
878 const auto nLocCells{this->getNumberOfCells()};
879
880 for(LongSimplexId lcid = 0; lcid < nLocCells; ++lcid) {
881 if(this->cellRankArray_[lcid] != ttk::MPIrank_) {
882 // store ghost cell global ids (per rank)
883 this->ghostCellsPerOwner_[this->cellRankArray_[lcid]].emplace_back(
884 this->cellGid_[lcid]);
885 }
886 }
887
888 // for each rank, store the global id of local cells that are ghost cells of
889 // other ranks.
890 const auto MIT{ttk::getMPIType(ttk::SimplexId{})};
891 this->remoteGhostCells_.resize(ttk::MPIsize_);
892 // number of owned cells that are ghost cells of other ranks
893 std::vector<size_t> nOwnedGhostCellsPerRank(ttk::MPIsize_);
894
895 for(const auto neigh : this->neighborRanks_) {
896 // 1. send to neigh number of ghost cells owned by neigh
897 const auto nCells{this->ghostCellsPerOwner_[neigh].size()};
898 MPI_Sendrecv(&nCells, 1, ttk::getMPIType(nCells), neigh, ttk::MPIrank_,
899 &nOwnedGhostCellsPerRank[neigh], 1, ttk::getMPIType(nCells),
900 neigh, neigh, ttk::MPIcomm_, MPI_STATUS_IGNORE);
901 this->remoteGhostCells_[neigh].resize(nOwnedGhostCellsPerRank[neigh]);
902
903 // 2. send to neigh list of ghost cells owned by neigh
904 MPI_Sendrecv(this->ghostCellsPerOwner_[neigh].data(),
905 this->ghostCellsPerOwner_[neigh].size(), MIT, neigh,
906 ttk::MPIrank_, this->remoteGhostCells_[neigh].data(),
907 this->remoteGhostCells_[neigh].size(), MIT, neigh, neigh,
908 ttk::MPIcomm_, MPI_STATUS_IGNORE);
909 }
910
911 this->hasPreconditionedExchangeGhostCells_ = true;
912 return 0;
913}
914
915int ttk::ExplicitTriangulation::preconditionDistributedCellRanges() {
916
917 // 1. store all local cells owned by current rank by global id
918
919 std::vector<SimplexId> localCellIds{};
920 localCellIds.reserve(this->getNumberOfCells());
921 for(SimplexId i = 0; i < this->getNumberOfCells(); ++i) {
922 if(this->cellRankArray_[i] == ttk::MPIrank_) {
923 localCellIds.emplace_back(i);
924 }
925 }
926
927 TTK_PSORT(this->threadNumber_, localCellIds.begin(), localCellIds.end(),
928 [this](const SimplexId a, const SimplexId b) {
929 return this->cellGid_[a] < this->cellGid_[b];
930 });
931
932 // 2. determine ranges of contiguous cell global ids
933
934 size_t begRange{};
935 while(begRange < localCellIds.size()) {
936 size_t endRange{begRange + 1};
937
938 if(begRange < localCellIds.size() - 1) {
939 for(size_t j = begRange + 1; j < localCellIds.size(); ++j) {
940 if(this->cellGid_[localCellIds[j]]
941 > this->cellGid_[localCellIds[j - 1]] + 1) {
942 endRange = j;
943 break;
944 }
945 }
946 if(endRange == begRange + 1
947 && this->cellGid_[localCellIds[endRange]]
948 == this->cellGid_[localCellIds[endRange - 1]] + 1) {
949 endRange = localCellIds.size();
950 }
951 }
952
953 const size_t gbeg = this->cellGid_[localCellIds[begRange]];
954 const size_t gend = this->cellGid_[localCellIds[endRange - 1]];
955 const auto nRanges{this->localCellRanges_.size()};
956
957 // inclusive range
958 this->localCellRanges_.emplace_back(
959 CellRange{nRanges, gbeg, gend, static_cast<size_t>(ttk::MPIrank_)});
960
961 begRange = endRange;
962 }
963
964 // 3. send to rank 0 the vector of ranges so it can compute range offsets
965
966 if(ttk::MPIrank_ == 0) {
967 this->nRangesPerRank_.resize(ttk::MPIsize_);
968 }
969
970 const int rangeSize = this->localCellRanges_.size();
971 MPI_Gather(&rangeSize, 1, ttk::getMPIType(rangeSize),
972 this->nRangesPerRank_.data(), 1, ttk::getMPIType(rangeSize), 0,
973 ttk::MPIcomm_);
974
975 std::vector<int> displacements{};
976
977 if(ttk::MPIrank_ == 0) {
978 const auto nRanges{std::accumulate(
979 this->nRangesPerRank_.begin(), this->nRangesPerRank_.end(), 0)};
980 this->gatheredCellRanges_.resize(nRanges);
981 displacements.resize(this->nRangesPerRank_.size());
982
983 for(size_t i = 0; i < this->nRangesPerRank_.size() - 1; ++i) {
984 displacements[i + 1] = displacements[i] + this->nRangesPerRank_[i];
985 }
986 }
987
988 auto cellRangeDT{CellRange::getMPIType()};
989 MPI_Type_commit(&cellRangeDT);
990
991 MPI_Gatherv(this->localCellRanges_.data(), this->localCellRanges_.size(),
992 cellRangeDT, this->gatheredCellRanges_.data(),
993 this->nRangesPerRank_.data(), displacements.data(), cellRangeDT,
994 0, ttk::MPIcomm_);
995
996 MPI_Type_free(&cellRangeDT);
997
998 // 4. sort range vector on rank 0
999
1000 if(ttk::MPIrank_ == 0) {
1001 TTK_PSORT(
1002 this->threadNumber_, this->gatheredCellRanges_.begin(),
1003 this->gatheredCellRanges_.end(),
1004 [](const CellRange &a, const CellRange &b) { return a.begin < b.begin; });
1005 }
1006
1007 return 0;
1008}
1009
1010size_t ttk::ExplicitTriangulation::computeCellRangeOffsets(
1011 std::vector<size_t> &nSimplicesPerRange) const {
1012
1013 // 1. send to rank 0 number of edges per cell range
1014
1015 std::vector<std::vector<size_t>> nSimplicesPerRangePerRank{};
1016
1017 if(ttk::MPIrank_ == 0) {
1018 nSimplicesPerRangePerRank.resize(this->nRangesPerRank_.size());
1019 for(int i = 0; i < ttk::MPIsize_; ++i) {
1020 nSimplicesPerRangePerRank[i].resize(this->nRangesPerRank_[i]);
1021 }
1022 for(int i = 0; i < ttk::MPIsize_; ++i) {
1023 if(i == 0) {
1024 continue;
1025 }
1026 // receive src content from other ranks
1027 MPI_Recv(nSimplicesPerRangePerRank[i].data(),
1028 nSimplicesPerRangePerRank[i].size(), ttk::getMPIType(size_t{}),
1029 i, MPI_ANY_TAG, ttk::MPIcomm_, MPI_STATUS_IGNORE);
1030 }
1031 std::swap(nSimplicesPerRangePerRank[0], nSimplicesPerRange);
1032 } else {
1033 MPI_Send(nSimplicesPerRange.data(), nSimplicesPerRange.size(),
1034 ttk::getMPIType(size_t{}), 0, 0, ttk::MPIcomm_);
1035 }
1036
1037 // 2. compute range offsets on rank 0
1038
1039 size_t nSimplices{};
1040 if(ttk::MPIrank_ == 0) {
1041
1042 for(const auto &range : this->gatheredCellRanges_) {
1043 auto &pSum{nSimplicesPerRangePerRank[range.rank][range.id]};
1044 std::swap(pSum, nSimplices);
1045 nSimplices += pSum;
1046 }
1047 }
1048
1049 // 3. send back range offsets to other ranks
1050
1051 if(ttk::MPIrank_ == 0) {
1052 for(int i = 1; i < ttk::MPIsize_; ++i) {
1053 MPI_Send(nSimplicesPerRangePerRank[i].data(),
1054 nSimplicesPerRangePerRank[i].size(), ttk::getMPIType(size_t{}),
1055 i, 0, ttk::MPIcomm_);
1056 }
1057 std::swap(nSimplicesPerRange, nSimplicesPerRangePerRank[0]);
1058 } else {
1059 MPI_Recv(nSimplicesPerRange.data(), nSimplicesPerRange.size(),
1060 ttk::getMPIType(size_t{}), MPI_ANY_TAG, 0, ttk::MPIcomm_,
1061 MPI_STATUS_IGNORE);
1062 }
1063
1064 return nSimplices;
1065}
1066
1067template <typename Func0, typename Func1, typename Func2>
1068int ttk::ExplicitTriangulation::exchangeDistributedInternal(
1069 const Func0 &getGlobalSimplexId,
1070 const Func1 &storeGlobalSimplexId,
1071 const Func2 &iterCond,
1072 const int nSimplicesPerCell) {
1073
1074 // per neighbor, owned ghost cell simplex global ids to transfer back
1075 std::vector<std::vector<SimplexId>> globalIdPerOwnedGhostCell(ttk::MPIsize_);
1076 // per neighbor, non-owned ghost cell simplex global ids to transfer back
1077 std::vector<std::vector<SimplexId>> globalIdPerLocalGhostCell(ttk::MPIsize_);
1078
1079 const auto MIT{ttk::getMPIType(ttk::SimplexId{})};
1080
1081 // make sure that all simplices are correctly labelled: for a given
1082 // rank, a simplex can be owned by a ghost cell from a neighboring
1083 // rank but in reality can be owned by another ghost cell in a third
1084 // rank
1085 bool doIter{true};
1086
1087 while(doIter) {
1088
1089 doIter = false;
1090
1091 // 3. for each list of ghost cell, accumulate the global simplex id
1092 for(const auto neigh : this->neighborRanks_) {
1093 // sending side
1094 globalIdPerOwnedGhostCell[neigh].resize(
1095 nSimplicesPerCell * this->remoteGhostCells_[neigh].size());
1096 for(size_t i = 0; i < this->remoteGhostCells_[neigh].size(); ++i) {
1097 const auto lcid{this->cellGidToLid_[this->remoteGhostCells_[neigh][i]]};
1098 for(int j = 0; j < nSimplicesPerCell; ++j) {
1099 globalIdPerOwnedGhostCell[neigh][nSimplicesPerCell * i + j]
1100 = getGlobalSimplexId(lcid, j);
1101 }
1102 }
1103 // receiving side
1104 globalIdPerLocalGhostCell[neigh].resize(
1105 nSimplicesPerCell * this->ghostCellsPerOwner_[neigh].size());
1106
1107 // 4. transfer back global simplex ids
1108 MPI_Sendrecv(globalIdPerOwnedGhostCell[neigh].data(),
1109 globalIdPerOwnedGhostCell[neigh].size(), MIT, neigh,
1110 ttk::MPIrank_, globalIdPerLocalGhostCell[neigh].data(),
1111 globalIdPerLocalGhostCell[neigh].size(), MIT, neigh, neigh,
1112 ttk::MPIcomm_, MPI_STATUS_IGNORE);
1113 }
1114
1115 // 5. extend local <-> global simplex ids mappings
1116 for(const auto neigh : this->neighborRanks_) {
1117 for(size_t i = 0; i < this->ghostCellsPerOwner_[neigh].size(); ++i) {
1118 const auto gcid{this->ghostCellsPerOwner_[neigh][i]};
1119 const auto lcid{this->cellGidToLid_[gcid]};
1120 for(int j = 0; j < nSimplicesPerCell; ++j) {
1121 const auto geid{
1122 globalIdPerLocalGhostCell[neigh][nSimplicesPerCell * i + j]};
1123 storeGlobalSimplexId(lcid, geid, j);
1124 }
1125 }
1126 }
1127
1128 // do an additional transmission if there still is some locally
1129 // non-labelled simplices
1130 int doNextIter{0};
1131 if(iterCond()) {
1132 doNextIter = 1;
1133 doIter = true;
1134 }
1135 for(int i = 0; i < ttk::MPIsize_; ++i) {
1136 if(doIter) {
1137 // reset doNextIter (might have been erased by the MPI_Bcast)
1138 doNextIter = 1;
1139 }
1140 MPI_Bcast(&doNextIter, 1, ttk::getMPIType(doNextIter), i, ttk::MPIcomm_);
1141 doIter |= (doNextIter == 1);
1142 }
1143
1144 if(doIter && ttk::MPIrank_ == 0) {
1145 this->printMsg("Re-sending global ids to neighbors...");
1146 }
1147 }
1148
1149 return 0;
1150}
1151
1152int ExplicitTriangulation::preconditionDistributedEdges() {
1153 if(this->hasPreconditionedDistributedEdges_) {
1154 return 0;
1155 }
1156 if(!ttk::hasInitializedMPI()) {
1157 return -1;
1158 }
1159 if(this->cellGid_ == nullptr) {
1160 this->printErr("Missing global cell identifiers array!");
1161 return -2;
1162 }
1163
1164 if(this->getDimensionality() != 2 && this->getDimensionality() != 3) {
1165 return -3;
1166 }
1167
1168 if(this->getDimensionality() == 2) {
1170 }
1171
1172 Timer tm{};
1173
1174 this->preconditionDistributedCells();
1175
1176 // allocate memory
1177 this->edgeLidToGid_.resize(this->getNumberOfEdgesInternal(), -1);
1178 this->edgeGidToLid_.reserve(this->getNumberOfEdgesInternal());
1179
1180 // 1. for every range of local cells, number the edges locally
1181
1182 std::vector<SimplexId> edgeLidToRangeId(this->getNumberOfEdgesInternal(), -1);
1183 std::vector<size_t> nEdgesPerRange(this->localCellRanges_.size());
1184
1185 const auto edgeAlreadyProcessed = [this](const SimplexId leid,
1186 const SimplexId lcid) {
1187 const auto nStar{this->TTK_TRIANGULATION_INTERNAL(getEdgeStarNumber)(leid)};
1188 for(SimplexId i = 0; i < nStar; ++i) {
1189 SimplexId sid{-1};
1190 this->TTK_TRIANGULATION_INTERNAL(getEdgeStar)(leid, i, sid);
1191 if(sid == -1 || sid == lcid) {
1192 continue;
1193 }
1194 // rule: an edge is owned by the cell in its star with the
1195 // lowest rank
1196 if(this->cellRankArray_[sid] < this->cellRankArray_[lcid]) {
1197 return true;
1198 break;
1199 }
1200 }
1201 return false;
1202 };
1203
1204 const auto countCellEdges
1205 = [this, &edgeAlreadyProcessed](const SimplexId lcid,
1206 std::vector<SimplexId> &edgeGid,
1207 std::vector<SimplexId> &edgeRangeId,
1208 const size_t rangeId, size_t &edgeCount) {
1209 SimplexId nEdges{};
1210 if(this->maxCellDim_ == 3) {
1211 nEdges = this->getCellEdgeNumberInternal(lcid);
1212 } else if(this->maxCellDim_ == 2) {
1213 nEdges = this->getTriangleEdgeNumberInternal(lcid);
1214 }
1215 for(SimplexId k = 0; k < nEdges; ++k) {
1216 SimplexId leid{-1};
1217 if(this->maxCellDim_ == 3) {
1218 this->getCellEdgeInternal(lcid, k, leid);
1219 } else if(this->maxCellDim_ == 2) {
1220 this->getTriangleEdge(lcid, k, leid);
1221 }
1222 const auto alreadyProcessed = edgeAlreadyProcessed(leid, lcid);
1223 if(!alreadyProcessed) {
1224 edgeGid[leid] = edgeCount;
1225 edgeRangeId[leid] = rangeId;
1226 edgeCount++;
1227 }
1228 }
1229 };
1230
1231#ifdef TTK_ENABLE_OPENMP
1232#pragma omp parallel for num_threads(this->threadNumber_)
1233#endif // TTK_ENABLE_OPENMP
1234 for(size_t i = 0; i < this->localCellRanges_.size(); ++i) {
1235 auto &range{this->localCellRanges_[i]};
1236 range.id = i;
1237 for(size_t j = range.begin; j <= range.end; ++j) {
1238 // local cell id
1239 const auto lcid{this->cellGidToLid_[j]};
1240 countCellEdges(lcid, this->edgeLidToGid_, edgeLidToRangeId, range.id,
1241 nEdgesPerRange[i]);
1242 }
1243 }
1244
1245 // 2. compute range offset on rank 0
1246 const auto nEdges = this->computeCellRangeOffsets(nEdgesPerRange);
1247
1248 // 3. locally edit the edge global id with range offsets
1249
1250 for(SimplexId leid = 0; leid < this->getNumberOfEdgesInternal(); ++leid) {
1251 if(this->edgeLidToGid_[leid] == -1) {
1252 // not owned by a cell of this rank
1253 continue;
1254 }
1255 const auto geid{this->edgeLidToGid_[leid]
1256 + nEdgesPerRange[edgeLidToRangeId[leid]]};
1257 this->edgeLidToGid_[leid] = geid;
1258 this->edgeGidToLid_[geid] = leid;
1259 }
1260
1261 // 4. exchange global ids between ghost cells
1262
1263 const auto nEdgesPerCell{this->getDimensionality() == 3 ? 6 : 3};
1264 this->exchangeDistributedInternal(
1265 [this](const SimplexId lcid, const int j) {
1266 SimplexId leid{};
1267 this->getCellEdgeInternal(lcid, j, leid);
1268 return this->edgeLidToGid_[leid];
1269 },
1270 [this](const SimplexId lcid, const SimplexId geid, const int j) {
1271 SimplexId leid{};
1272 this->getCellEdgeInternal(lcid, j, leid);
1273 if(this->edgeLidToGid_[leid] == -1 && geid != -1) {
1274 this->edgeLidToGid_[leid] = geid;
1275 this->edgeGidToLid_[geid] = leid;
1276 }
1277 },
1278 [this]() {
1279 return std::count(
1280 this->edgeLidToGid_.begin(), this->edgeLidToGid_.end(), -1)
1281 > 0;
1282 },
1283 nEdgesPerCell);
1284
1285 this->preconditionEdgeRankArray();
1286
1287 if(MPIrank_ == 0) {
1288 this->printMsg("Domain contains " + std::to_string(nEdges) + " edges", 1.0,
1289 tm.getElapsedTime(), 1);
1290 }
1291
1292 this->hasPreconditionedDistributedEdges_ = true;
1293
1294 return 0;
1295}
1296
1297int ExplicitTriangulation::preconditionDistributedTriangles() {
1298 if(this->hasPreconditionedDistributedTriangles_) {
1299 return 0;
1300 }
1301 if(!ttk::hasInitializedMPI()) {
1302 return -1;
1303 }
1304 if(this->cellGid_ == nullptr) {
1305 this->printErr("Missing global cell identifiers array!");
1306 return -2;
1307 }
1308
1309 if(this->getDimensionality() != 3) {
1310 return -3;
1311 }
1312
1313 Timer tm{};
1314
1315 this->preconditionDistributedCells();
1316
1317 // allocate memory
1318 this->triangleLidToGid_.resize(this->getNumberOfTrianglesInternal(), -1);
1319 this->triangleGidToLid_.reserve(this->getNumberOfTrianglesInternal());
1320
1321 // 1. for every range of local cells, number the edges locally
1322
1323 std::vector<SimplexId> triangleLidToRangeId(
1324 this->getNumberOfTrianglesInternal(), -1);
1325 std::vector<size_t> nTrianglesPerRange(this->localCellRanges_.size());
1326
1327 const auto triangleAlreadyProcessed
1328 = [this](const SimplexId leid, const SimplexId lcid) {
1329 const auto nStar{
1331 for(SimplexId i = 0; i < nStar; ++i) {
1332 SimplexId sid{-1};
1333 this->TTK_TRIANGULATION_INTERNAL(getTriangleStar)(leid, i, sid);
1334 if(sid == -1 || sid == lcid) {
1335 continue;
1336 }
1337 // rule: an triangle is owned by the cell in its star with the
1338 // lowest rank
1339 if(this->cellRankArray_[sid] < this->cellRankArray_[lcid]) {
1340 return true;
1341 break;
1342 }
1343 }
1344 return false;
1345 };
1346
1347 const auto countCellTriangles
1348 = [this, &triangleAlreadyProcessed](
1349 const SimplexId lcid, std::vector<SimplexId> &triangleGid,
1350 std::vector<SimplexId> &triangleRangeId, const size_t rangeId,
1351 size_t &triangleCount) {
1352 const auto nTriangles{this->getCellTriangleNumberInternal(lcid)};
1353 for(SimplexId k = 0; k < nTriangles; ++k) {
1354 SimplexId leid{-1};
1355 this->getCellTriangleInternal(lcid, k, leid);
1356 const auto alreadyProcessed = triangleAlreadyProcessed(leid, lcid);
1357 if(!alreadyProcessed) {
1358 triangleGid[leid] = triangleCount;
1359 triangleRangeId[leid] = rangeId;
1360 triangleCount++;
1361 }
1362 }
1363 };
1364
1365#ifdef TTK_ENABLE_OPENMP
1366#pragma omp parallel for num_threads(this->threadNumber_)
1367#endif // TTK_ENABLE_OPENMP
1368 for(size_t i = 0; i < this->localCellRanges_.size(); ++i) {
1369 auto &range{this->localCellRanges_[i]};
1370 range.id = i;
1371 for(size_t j = range.begin; j <= range.end; ++j) {
1372 // local cell id
1373 const auto lcid{this->cellGidToLid_[j]};
1374 countCellTriangles(lcid, this->triangleLidToGid_, triangleLidToRangeId,
1375 range.id, nTrianglesPerRange[i]);
1376 }
1377 }
1378
1379 // 2. compute range offset on rank 0
1380 const auto nTriangles = this->computeCellRangeOffsets(nTrianglesPerRange);
1381
1382 // 3. locally edit the triangle global id with range offsets
1383
1384 for(SimplexId leid = 0; leid < this->getNumberOfTrianglesInternal(); ++leid) {
1385 if(this->triangleLidToGid_[leid] == -1) {
1386 // not owned by a cell of this rank
1387 continue;
1388 }
1389 const auto geid{this->triangleLidToGid_[leid]
1390 + nTrianglesPerRange[triangleLidToRangeId[leid]]};
1391 this->triangleLidToGid_[leid] = geid;
1392 this->triangleGidToLid_[geid] = leid;
1393 }
1394
1395 // 4. exchange global ids between ghost cells
1396
1397 const auto nTrianglesPerCell{4};
1398 this->exchangeDistributedInternal(
1399 [this](const SimplexId lcid, const int j) {
1400 SimplexId ltid{};
1401 this->getCellTriangleInternal(lcid, j, ltid);
1402 return this->triangleLidToGid_[ltid];
1403 },
1404 [this](const SimplexId lcid, const SimplexId gtid, const int j) {
1405 SimplexId ltid{};
1406 this->getCellTriangleInternal(lcid, j, ltid);
1407 if(this->triangleLidToGid_[ltid] == -1 && gtid != -1) {
1408 this->triangleLidToGid_[ltid] = gtid;
1409 this->triangleGidToLid_[gtid] = ltid;
1410 }
1411 },
1412 [this]() {
1413 return std::count(this->triangleLidToGid_.begin(),
1414 this->triangleLidToGid_.end(), -1)
1415 > 0;
1416 },
1417 nTrianglesPerCell);
1418
1419 this->preconditionTriangleRankArray();
1420
1421 if(MPIrank_ == 0) {
1422 this->printMsg(
1423 "Domain contains " + std::to_string(nTriangles) + " triangles", 1.0,
1424 tm.getElapsedTime(), 1);
1425 }
1426
1427 this->hasPreconditionedDistributedTriangles_ = true;
1428
1429 return 0;
1430}
1431
1432int ExplicitTriangulation::preconditionDistributedVertices() {
1433 if(this->hasPreconditionedDistributedVertices_) {
1434 return 0;
1435 }
1436 if(!hasInitializedMPI()) {
1437 return -1;
1438 }
1439 if(this->vertGid_ == nullptr) {
1440 this->printErr("Missing global vertex identifiers array!");
1441 return -2;
1442 }
1443
1444 this->hasPreconditionedDistributedVertices_ = true;
1445 this->preconditionVertexRankArray();
1446
1447 // number of local vertices (with ghost vertices...)
1448 const auto nLocVertices{this->getNumberOfVertices()};
1449
1450 // global vertex id -> local vertex id (reverse of this->vertGid_)
1451 this->vertexGidToLid_.reserve(nLocVertices);
1452 for(LongSimplexId lvid = 0; lvid < nLocVertices; ++lvid) {
1453 this->vertexGidToLid_[this->vertGid_[lvid]] = lvid;
1454 }
1455
1456 this->preconditionExchangeGhostVertices();
1457
1458 return 0;
1459}
1460
1461int ttk::ExplicitTriangulation::preconditionExchangeGhostVertices() {
1462
1463 if(this->hasPreconditionedExchangeGhostVertices_) {
1464 return 0;
1465 }
1466 if(!ttk::hasInitializedMPI()) {
1467 return -1;
1468 }
1469 if(this->vertGid_ == nullptr) {
1470 this->printErr("Missing global vertex identifiers array!");
1471 return -2;
1472 }
1473
1474 // number of local vertices (with ghost vertices...)
1475 const auto nLocVertices{this->getNumberOfVertices()};
1476
1477 this->ghostVerticesPerOwner_.resize(ttk::MPIsize_);
1478
1479 for(LongSimplexId lvid = 0; lvid < nLocVertices; ++lvid) {
1480 if(this->vertexRankArray_[lvid] != ttk::MPIrank_) {
1481 // store ghost cell global ids (per rank)
1482 this->ghostVerticesPerOwner_[this->vertexRankArray_[lvid]].emplace_back(
1483 this->vertGid_[lvid]);
1484 }
1485 }
1486
1487 // for each rank, store the global id of local cells that are ghost cells of
1488 // other ranks.
1489 const auto MIT{ttk::getMPIType(ttk::SimplexId{})};
1490 this->remoteGhostVertices_.resize(ttk::MPIsize_);
1491 // number of owned cells that are ghost cells of other ranks
1492 std::vector<size_t> nOwnedGhostVerticesPerRank(ttk::MPIsize_);
1493
1494 for(const auto neigh : this->neighborRanks_) {
1495 // 1. send to neigh number of ghost cells owned by neigh
1496 const auto nVerts{this->ghostVerticesPerOwner_[neigh].size()};
1497 MPI_Sendrecv(&nVerts, 1, ttk::getMPIType(nVerts), neigh, ttk::MPIrank_,
1498 &nOwnedGhostVerticesPerRank[neigh], 1, ttk::getMPIType(nVerts),
1499 neigh, neigh, ttk::MPIcomm_, MPI_STATUS_IGNORE);
1500 this->remoteGhostVertices_[neigh].resize(nOwnedGhostVerticesPerRank[neigh]);
1501
1502 // 2. send to neigh list of ghost cells owned by neigh
1503 MPI_Sendrecv(this->ghostVerticesPerOwner_[neigh].data(),
1504 this->ghostVerticesPerOwner_[neigh].size(), MIT, neigh,
1505 ttk::MPIrank_, this->remoteGhostVertices_[neigh].data(),
1506 this->remoteGhostVertices_[neigh].size(), MIT, neigh, neigh,
1507 ttk::MPIcomm_, MPI_STATUS_IGNORE);
1508 }
1509
1510 this->hasPreconditionedExchangeGhostVertices_ = true;
1511
1512 return 0;
1513}
1514
1515#endif // TTK_ENABLE_MPI
1516
1517template <typename T>
1518void writeBin(std::ofstream &stream, const T var) {
1519 stream.write(reinterpret_cast<const char *>(&var), sizeof(var));
1520}
1521
1522template <typename T>
1523void writeBinArray(std::ofstream &stream,
1524 const T *const buff,
1525 const size_t size) {
1526 stream.write(reinterpret_cast<const char *>(buff), size * sizeof(T));
1527}
1528
1529// initialize static member variables
1530const char *ExplicitTriangulation::magicBytes_ = "TTKTriangulationFileFormat";
1531const unsigned long ExplicitTriangulation::formatVersion_ = 1;
1532
1533int ExplicitTriangulation::writeToFile(std::ofstream &stream) const {
1534
1535 // 1. magic bytes (char *)
1536 stream.write(ttk::ExplicitTriangulation::magicBytes_,
1537 std::strlen(ttk::ExplicitTriangulation::magicBytes_));
1538 // 2. format version (unsigned long)
1539 writeBin(stream, ttk::ExplicitTriangulation::formatVersion_);
1540 // 3. dimensionality (int)
1541 const auto dim = this->getDimensionality();
1542 writeBin(stream, dim);
1543 // 4. number of vertices (SimplexId)
1544 const auto nVerts = this->getNumberOfVertices();
1545 writeBin(stream, nVerts);
1546 // 5. number of edges (SimplexId)
1547 const auto edgesNumber = [this, dim]() -> SimplexId {
1548 if(dim == 1) {
1549 return this->getNumberOfCells();
1550 } else if(dim > 1) {
1551 return this->edgeList_.size();
1552 }
1553 return 0;
1554 };
1555 const auto nEdges = edgesNumber();
1556 writeBin(stream, nEdges);
1557 // 6. number of triangles (SimplexId, 0 in 1D)
1558 const auto trianglesNumber = [this, dim]() -> SimplexId {
1559 if(dim == 2) {
1560 return this->getNumberOfCells();
1561 } else if(dim == 3) {
1562 return this->triangleList_.size();
1563 }
1564 return 0;
1565 };
1566 const auto nTriangles = trianglesNumber();
1567 writeBin(stream, nTriangles);
1568 // 7. number of tetrahedron (SimplexId, 0 in 2D)
1569 const auto nTetras = dim > 2 ? this->getNumberOfCells() : 0;
1570 writeBin(stream, nTetras);
1571
1572 // only write buffers oustside this->cellArray_ (cellVertex, vertexCoords),
1573 // those ones will be provided by VTK
1574
1575 // fixed-size arrays (in AbstractTriangulation.h)
1576
1577#define WRITE_FIXED(ARRAY) \
1578 if(ARRAY.empty()) { \
1579 writeBin(stream, char{0}); \
1580 } else { \
1581 writeBin(stream, char{1}); \
1582 writeBinArray(stream, ARRAY.data(), ARRAY.size()); \
1583 }
1584
1585 // 8. edgeList (SimplexId array)
1586 WRITE_FIXED(this->edgeList_);
1587 // 9. triangleList (SimplexId array)
1589 // 10. triangleEdgeList (SimplexId array)
1591 // 11. tetraEdgeList (SimplexId array)
1593 // 12. tetraTriangleList (SimplexId array)
1595
1596 // variable-size arrays (FlatJaggedArray in ExplicitTriangulation.h)
1597
1598#define WRITE_GUARD(ARRAY) \
1599 if(ARRAY.empty()) { \
1600 writeBin(stream, char{0}); \
1601 return; \
1602 } else { \
1603 writeBin(stream, char{1}); \
1604 }
1605
1606 const auto write_variable = [&stream](const FlatJaggedArray &arr) {
1607 // empty array guard
1608 WRITE_GUARD(arr);
1609 writeBinArray(stream, arr.offset_ptr(), arr.size() + 1);
1610 writeBinArray(stream, arr.get_ptr(0, 0), arr.dataSize());
1611 };
1612
1613 // 13. vertexNeighbors (SimplexId array, offsets then data)
1614 write_variable(this->vertexNeighborData_);
1615 // 14. cellNeighbors (SimplexId array, offsets then data)
1616 write_variable(this->cellNeighborData_);
1617 // 15. vertexEdges (SimplexId array, offsets then data)
1618 write_variable(this->vertexEdgeData_);
1619 // 16. vertexTriangles (SimplexId array, offsets then data)
1620 write_variable(this->vertexTriangleData_);
1621 // 17. edgeTriangles (SimplexId array, offsets then data)
1622 write_variable(this->edgeTriangleData_);
1623 // 18. vertexStars (SimplexId array, offsets then data)
1624 write_variable(this->vertexStarData_);
1625 // 19. edgeStars (SimplexId array, offsets then data)
1626 write_variable(this->edgeStarData_);
1627 // 20. triangleStars (SimplexId array, offsets then data)
1628 write_variable(this->triangleStarData_);
1629 // 21. vertexLinks (SimplexId array, offsets then data)
1630 write_variable(this->vertexLinkData_);
1631 // 22. edgeLinks (SimplexId array, offsets then data)
1632 write_variable(this->edgeLinkData_);
1633 // 23. triangleLinks (SimplexId array, offsets then data)
1634 write_variable(this->triangleLinkData_);
1635
1636 const auto write_bool = [&stream](const std::vector<bool> &arr) {
1637 // empty array guard
1638 WRITE_GUARD(arr);
1639 for(size_t i = 0; i < arr.size(); ++i) {
1640 const auto b = static_cast<char>(arr[i]);
1641 writeBin(stream, b);
1642 }
1643 };
1644
1645 // 24. boundary vertices (bool array)
1646 write_bool(this->boundaryVertices_);
1647 // 25. boundary edges (bool array)
1648 write_bool(this->boundaryEdges_);
1649 // 26. boundary triangles (bool array)
1650 write_bool(this->boundaryTriangles_);
1651
1652 return 0;
1653}
1654
1655int ExplicitTriangulation::writeToFileASCII(std::ofstream &stream) const {
1656 // 1. magic bytes
1657 stream << ttk::ExplicitTriangulation::magicBytes_ << '\n';
1658 // 2. format version
1659 stream << ttk::ExplicitTriangulation::formatVersion_ + 1 << '\n';
1660 // 3. dimensionality
1661 const auto dim = this->getDimensionality();
1662 stream << "dim " << dim << '\n';
1663 // 4. -> 7. number of simplices
1664 const auto nVerts = this->getNumberOfVertices();
1665 const auto edgesNumber = [this, dim]() -> SimplexId {
1666 if(dim == 1) {
1667 return this->getNumberOfCells();
1668 } else if(dim > 1) {
1669 return this->edgeList_.size();
1670 }
1671 return 0;
1672 };
1673 const auto nEdges = edgesNumber();
1674 const auto trianglesNumber = [this, dim]() -> SimplexId {
1675 if(dim == 2) {
1676 return this->getNumberOfCells();
1677 } else if(dim == 3) {
1678 return this->triangleList_.size();
1679 }
1680 return 0;
1681 };
1682 const auto nTriangles = trianglesNumber();
1683 const auto nTetras = dim > 2 ? this->getNumberOfCells() : 0;
1684 stream << nVerts << ' ' << nEdges << ' ' << nTriangles << ' ' << nTetras
1685 << '\n';
1686
1687#define WRITE_ASCII(ARRAY) \
1688 stream << #ARRAY << '\n'; \
1689 for(const auto &slice : ARRAY) { \
1690 for(size_t i = 0; i < slice.size(); ++i) { \
1691 if(i > 0) { \
1692 stream << ' '; \
1693 } \
1694 stream << slice[i]; \
1695 } \
1696 stream << '\n'; \
1697 }
1698
1699 const auto writeCellArray = [this, &stream]() {
1700 stream << "this->cellArray_\n";
1701 for(SimplexId i = 0; i < this->cellNumber_; ++i) {
1702 for(SimplexId j = 0; j < this->cellArray_->getCellVertexNumber(i); ++j) {
1703 stream << this->cellArray_->getCellVertex(i, j) << ' ';
1704 }
1705 stream << '\n';
1706 }
1707 };
1708
1709 // ?. cellArray_
1710 writeCellArray();
1711
1712 // 8. -> 12. fixed-size arrays
1713 WRITE_ASCII(this->edgeList_);
1717
1718 // 13. -> 23. variable-size arrays
1719 WRITE_ASCII(this->vertexNeighborData_);
1720 WRITE_ASCII(this->cellNeighborData_);
1721 WRITE_ASCII(this->vertexEdgeData_);
1722 WRITE_ASCII(this->vertexTriangleData_);
1723 WRITE_ASCII(this->edgeTriangleData_);
1724 WRITE_ASCII(this->vertexStarData_);
1725 WRITE_ASCII(this->edgeStarData_);
1726 WRITE_ASCII(this->triangleStarData_);
1727 WRITE_ASCII(this->vertexLinkData_);
1728 WRITE_ASCII(this->edgeLinkData_);
1729 WRITE_ASCII(this->triangleLinkData_);
1730
1731#define WRITE_BOOL(ARRAY) \
1732 stream << #ARRAY << '\n'; \
1733 for(const auto el : ARRAY) { \
1734 stream << el << '\n'; \
1735 }
1736
1737 // 24. -> 26. boolean arrays
1741
1742 return 0;
1743}
1744
1745template <typename T>
1746void readBin(std::ifstream &stream, T &res) {
1747 stream.read(reinterpret_cast<char *>(&res), sizeof(res));
1748}
1749
1750template <typename T>
1751void readBinArray(std::ifstream &stream, T *const res, const size_t size) {
1752 stream.read(reinterpret_cast<char *>(res), size * sizeof(T));
1753}
1754
1755int ExplicitTriangulation::readFromFile(std::ifstream &stream) {
1756
1757 // 1. magic bytes (char *)
1758 const auto magicBytesLen
1759 = std::strlen(ttk::ExplicitTriangulation::magicBytes_);
1760 std::vector<char> mBytes(magicBytesLen + 1);
1761 stream.read(mBytes.data(), magicBytesLen);
1762 const auto hasMagicBytes
1763 = std::strcmp(mBytes.data(), ttk::ExplicitTriangulation::magicBytes_) == 0;
1764 if(!hasMagicBytes) {
1765 this->printErr("Could not find magic bytes in input files!");
1766 this->printErr("Aborting...");
1767 return 0;
1768 }
1769 // 2. format version (unsigned long)
1770 unsigned long version{};
1771 readBin(stream, version);
1772 if(version != ttk::ExplicitTriangulation::formatVersion_) {
1773 this->printWrn("File format version (" + std::to_string(version)
1774 + ") and software version ("
1775 + std::to_string(ttk::ExplicitTriangulation::formatVersion_)
1776 + ") are different!");
1777 }
1778
1779 int dim{};
1780 SimplexId nVerts{}, nEdges{}, nTriangles{}, nTetras{};
1781
1782 // 3. dimensionality (int)
1783 readBin(stream, dim);
1784 // 4. number of vertices (SimplexId)
1785 readBin(stream, nVerts);
1786 // 5. number of edges (SimplexId)
1787 readBin(stream, nEdges);
1788 // 6. number of triangles (SimplexId, 0 in 1D)
1789 readBin(stream, nTriangles);
1790 // 7. number of tetrahedron (SimplexId, 0 in 2D)
1791 readBin(stream, nTetras);
1792
1793 if(dim != this->getDimensionality()) {
1794 this->printErr("Incorrect dimension!");
1795 return 0;
1796 }
1797 if(nVerts != this->getNumberOfVertices()) {
1798 this->printErr("Incorrect number of vertices!");
1799 return 0;
1800 }
1801 if((dim == 2 && nTriangles != this->getNumberOfCells())
1802 || (dim == 3 && nTetras != this->getNumberOfCells())) {
1803 this->printErr("Incorrect number of cells!");
1804 return 0;
1805 }
1806
1807 // fixed-size arrays (in AbstractTriangulation.h)
1808
1809 const auto read_guard = [&stream]() {
1810 char g{};
1811 readBin(stream, g);
1812 return g == 0;
1813 };
1814
1815#define READ_FIXED(ARRAY, N_ITEMS) \
1816 if(!read_guard()) { \
1817 ARRAY.resize(N_ITEMS); \
1818 readBinArray(stream, ARRAY.data(), N_ITEMS); \
1819 }
1820
1821 // 8. edgeList (SimplexId array)
1822 READ_FIXED(this->edgeList_, nEdges);
1823 // 9. triangleList (SimplexId array)
1824 READ_FIXED(this->triangleList_, nTriangles);
1825 // 10. triangleEdgeList (SimplexId array)
1826 READ_FIXED(this->triangleEdgeList_, nTriangles);
1827 // 11. tetraEdgeList (SimplexId array)
1828 READ_FIXED(this->tetraEdgeList_, nTetras);
1829 // 12. tetraTriangleList (SimplexId array)
1830 READ_FIXED(this->tetraTriangleList_, nTetras);
1831
1832 // variable-size arrays (FlagJaggedArrays in ExplicitTriangulation.h)
1833
1834 const auto read_variable
1835 = [&stream, &read_guard](FlatJaggedArray &arr, const SimplexId n_items) {
1836 // empty array guard
1837 if(read_guard()) {
1838 return;
1839 }
1840 std::vector<SimplexId> offsets{}, data{};
1841 offsets.resize(n_items + 1);
1842 readBinArray(stream, offsets.data(), offsets.size());
1843 data.resize(offsets.back());
1844 readBinArray(stream, data.data(), data.size());
1845 arr.setData(std::move(data), std::move(offsets));
1846 };
1847
1848 // 13. vertexNeighbors (SimplexId array, offsets then data)
1849 read_variable(this->vertexNeighborData_, nVerts);
1850 // 14. cellNeighbors (SimplexId array, offsets then data)
1851 read_variable(this->cellNeighborData_, this->getNumberOfCells());
1852 // 15. vertexEdges (SimplexId array, offsets then data)
1853 read_variable(this->vertexEdgeData_, nVerts);
1854 // 16. vertexTriangles (SimplexId array, offsets then data)
1855 read_variable(this->vertexTriangleData_, nVerts);
1856 // 17. edgeTriangles (SimplexId array, offsets then data)
1857 read_variable(this->edgeTriangleData_, nEdges);
1858 // 18. vertexStars (SimplexId array, offsets then data)
1859 read_variable(this->vertexStarData_, nVerts);
1860 // 19. edgeStars (SimplexId array, offsets then data)
1861 read_variable(this->edgeStarData_, nEdges);
1862 // 20. triangleStars (SimplexId array, offsets then data)
1863 read_variable(this->triangleStarData_, nTriangles);
1864 // 21. vertexLinks (SimplexId array, offsets then data)
1865 read_variable(this->vertexLinkData_, nVerts);
1866 // 22. edgeLinks (SimplexId array, offsets then data)
1867 read_variable(this->edgeLinkData_, nEdges);
1868 // 23. triangleLinks (SimplexId array, offsets then data)
1869 read_variable(this->triangleLinkData_, nTriangles);
1870
1871 const auto read_bool
1872 = [&stream, &read_guard](std::vector<bool> &arr, const SimplexId n_items) {
1873 // empty array guard
1874 if(read_guard()) {
1875 return;
1876 }
1877 // resize vector
1878 arr.resize(n_items);
1879 for(SimplexId i = 0; i < n_items; ++i) {
1880 char b{};
1881 stream.read(&b, sizeof(b));
1882 arr[i] = static_cast<bool>(b);
1883 }
1884 };
1885
1886 // 24. boundary vertices (bool array)
1887 read_bool(this->boundaryVertices_, nVerts);
1888 // 25. boundary edges (bool array)
1889 read_bool(this->boundaryEdges_, nEdges);
1890 // 26. boundary triangles (bool array)
1891 read_bool(this->boundaryTriangles_, nTriangles);
1892
1893 return 0;
1894}
#define TTK_TRIANGULATION_INTERNAL(NAME)
void readBinArray(std::ifstream &stream, T *const res, const size_t size)
void writeBin(std::ofstream &stream, const T var)
void readBin(std::ifstream &stream, T &res)
#define READ_FIXED(ARRAY, N_ITEMS)
#define WRITE_GUARD(ARRAY)
#define WRITE_FIXED(ARRAY)
void writeBinArray(std::ofstream &stream, const T *const buff, const size_t size)
#define WRITE_BOOL(ARRAY)
#define WRITE_ASCII(ARRAY)
#define TTK_PSORT(NTHREADS,...)
Parallel sort macro.
Definition OpenMP.h:46
std::vector< std::array< SimplexId, 3 > > triangleList_
virtual int getTriangleEdge(const SimplexId &triangleId, const int &localEdgeId, SimplexId &edgeId) const
std::vector< std::array< SimplexId, 4 > > tetraTriangleList_
std::vector< std::array< SimplexId, 2 > > edgeList_
std::vector< std::array< SimplexId, 6 > > tetraEdgeList_
size_t footprint(size_t size=0) const
std::vector< std::array< SimplexId, 3 > > triangleEdgeList_
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
int setWrapper(const Wrapper *wrapper) override
Definition Debug.cpp:155
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
SimplexId getTriangleEdgeNumberInternal(const SimplexId &triangleId) const override
size_t footprint(size_t size=0) const
SimplexId getNumberOfTrianglesInternal() const override
SimplexId getCellEdgeNumberInternal(const SimplexId &cellId) const override
int writeToFile(std::ofstream &stream) const
Write internal state to disk.
int preconditionBoundaryTrianglesInternal() override
SimplexId TTK_TRIANGULATION_INTERNAL() getNumberOfCells() const override
int TTK_TRIANGULATION_INTERNAL() getTriangleStar(const SimplexId &triangleId, const int &localStarId, SimplexId &starId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getTriangleStarNumber(const SimplexId &triangleId) const override
int getCellEdgeInternal(const SimplexId &cellId, const int &localEdgeId, SimplexId &edgeId) const override
int preconditionBoundaryVerticesInternal() override
SimplexId getNumberOfEdgesInternal() const override
SimplexId TTK_TRIANGULATION_INTERNAL() getEdgeStarNumber(const SimplexId &edgeId) const override
int TTK_TRIANGULATION_INTERNAL() getDimensionality() const override
int TTK_TRIANGULATION_INTERNAL() getEdgeStar(const SimplexId &edgeId, const int &localStarId, SimplexId &starId) const override
int writeToFileASCII(std::ofstream &stream) const
Write internal state to disk using an ASCII format.
SimplexId getCellTriangleNumberInternal(const SimplexId &cellId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getNumberOfVertices() const override
int getCellTriangleInternal(const SimplexId &cellId, const int &localTriangleId, SimplexId &triangleId) const override
int readFromFile(std::ifstream &stream)
Read from disk into internal state.
Replacement for std::vector<std::vector<SimplexId>>
size_t footprint() const
Computes the memory footprint of the array.
SimplexId size(SimplexId id) const
Get the size of a particular sub-vector.
bool empty() const
If the underlying buffers are empty.
OneSkeleton processing package.
Definition OneSkeleton.h:25
int buildEdgeLinks(const std::vector< std::array< SimplexId, 2 > > &edgeList, const FlatJaggedArray &edgeStars, const CellArray &cellArray, FlatJaggedArray &edgeLinks) const
int buildEdgeList(const SimplexId &vertexNumber, const CellArray &cellArray, std::vector< std::array< SimplexId, 2 > > &edgeList, FlatJaggedArray &edgeStars, std::vector< std::array< SimplexId, n > > &cellEdgeList) const
ThreeSkeleton processing package.
int buildCellNeighborsFromTriangles(const SimplexId &vertexNumber, const CellArray &cellArray, FlatJaggedArray &cellNeighbors, FlatJaggedArray *triangleStars=nullptr) const
TwoSkeleton processing package.
Definition TwoSkeleton.h:24
int buildTriangleEdgeList(const SimplexId &vertexNumber, const CellArray &cellArray, std::vector< std::array< SimplexId, 3 > > &triangleEdgeList, const std::vector< std::array< SimplexId, 2 > > &edgeList, FlatJaggedArray *vertexEdgeList=nullptr, std::vector< std::array< SimplexId, 3 > > *triangleList=nullptr, FlatJaggedArray *triangleStarList=nullptr, std::vector< std::array< SimplexId, 4 > > *cellTriangleList=nullptr) const
int buildVertexTriangles(const SimplexId &vertexNumber, const std::vector< std::array< SimplexId, 3 > > &triangleList, FlatJaggedArray &vertexTriangles) const
int buildCellNeighborsFromEdges(const CellArray &cellArray, FlatJaggedArray &cellNeighbors, const FlatJaggedArray &edgeStars) const
int buildTriangleLinks(const std::vector< std::array< SimplexId, 3 > > &triangleList, const FlatJaggedArray &triangleStars, const CellArray &cellArray, FlatJaggedArray &triangleLinks) const
int buildTriangleList(const SimplexId &vertexNumber, const CellArray &cellArray, std::vector< std::array< SimplexId, 3 > > *triangleList=nullptr, FlatJaggedArray *triangleStars=nullptr, std::vector< std::array< SimplexId, 4 > > *cellTriangleList=nullptr) const
int buildEdgeTriangles(const SimplexId &vertexNumber, const CellArray &cellArray, FlatJaggedArray &edgeTriangleList, const std::vector< std::array< SimplexId, 2 > > &edgeList, std::vector< std::array< SimplexId, 3 > > *triangleEdgeList=nullptr) const
ZeroSkeleton processing package.
int buildVertexEdges(const SimplexId &vertexNumber, const std::vector< std::array< SimplexId, 2 > > &edgeList, FlatJaggedArray &vertexEdges) const
int buildVertexStars(const SimplexId &vertexNumber, const CellArray &cellArray, FlatJaggedArray &vertexStars) const
int buildVertexLinks(const FlatJaggedArray &vertexStars, const std::vector< std::array< SimplexId, 3 > > &cellEdges, const std::vector< std::array< SimplexId, 2 > > &edgeList, FlatJaggedArray &vertexLinks) const
int buildVertexNeighbors(const SimplexId &vertexNumber, FlatJaggedArray &vertexNeighbors, const std::vector< std::array< SimplexId, 2 > > &edgeList) const
The Topology ToolKit.
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
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)