TTK
Loading...
Searching...
No Matches
ExplicitTriangulation.h
Go to the documentation of this file.
1
9
10#pragma once
11
12// base code includes
14#include <CellArray.h>
15#include <FlatJaggedArray.h>
16
17#include <memory>
18
19namespace ttk {
20
22
23 public:
25
27
32
33 int clear();
34
35 size_t footprint(size_t size = 0) const;
36
37 inline int getCellEdgeInternal(const SimplexId &cellId,
38 const int &localEdgeId,
39 SimplexId &edgeId) const override {
40
41#ifndef TTK_ENABLE_KAMIKAZE
42 if((cellId < 0) || (cellId >= (SimplexId)tetraEdgeList_.size()))
43 return -1;
44 if((localEdgeId < 0)
45 || (localEdgeId >= (SimplexId)tetraEdgeList_[cellId].size()))
46 return -2;
47#endif
48 edgeId = tetraEdgeList_[cellId][localEdgeId];
49 return 0;
50 }
51
52 inline SimplexId
53 getCellEdgeNumberInternal(const SimplexId &cellId) const override {
54#ifndef TTK_ENABLE_KAMIKAZE
55 if((cellId < 0) || (cellId >= (SimplexId)tetraEdgeList_.size()))
56 return -1;
57#endif
58 return tetraEdgeList_[cellId].size();
59 }
60
61 template <std::size_t N>
62 inline void
63 convertToVector(const std::vector<std::array<SimplexId, N>> &table,
64 std::vector<std::vector<SimplexId>> &vec) {
65 for(size_t i = 0; i < table.size(); ++i) {
66 vec[i] = {table[i].begin(), table[i].end()};
67 }
68 }
69
70 inline const std::vector<std::vector<SimplexId>> *
72
74 return &cellEdgeVector_;
75 }
76
78 const SimplexId &cellId,
79 const int &localNeighborId,
80 SimplexId &neighborId) const override {
81
82 neighborId = cellNeighborData_[cellId][localNeighborId];
83 return 0;
84 }
85
87 const SimplexId &cellId) const override {
88 return cellNeighborData_.size(cellId);
89 }
90
91 inline const std::vector<std::vector<SimplexId>> *
93 cellNeighborData_.copyTo(cellNeighborList_);
94 return &cellNeighborList_;
95 }
96
97 inline int getCellTriangleInternal(const SimplexId &cellId,
98 const int &localTriangleId,
99 SimplexId &triangleId) const override {
100
101#ifndef TTK_ENABLE_KAMIKAZE
102 if((cellId < 0) || (cellId >= (SimplexId)tetraTriangleList_.size()))
103 return -1;
104 if((localTriangleId < 0)
105 || (localTriangleId >= (SimplexId)tetraTriangleList_[cellId].size()))
106 return -2;
107#endif
108 triangleId = tetraTriangleList_[cellId][localTriangleId];
109
110 return 0;
111 }
112
113 inline SimplexId
114 getCellTriangleNumberInternal(const SimplexId &cellId) const override {
115
116#ifndef TTK_ENABLE_KAMIKAZE
117 if((cellId < 0) || (cellId >= (SimplexId)tetraTriangleList_.size()))
118 return -1;
119#endif
120
121 return tetraTriangleList_[cellId].size();
122 }
123
124 inline const std::vector<std::vector<SimplexId>> *
126
128 return &cellTriangleVector_;
129 }
130
132 const SimplexId &cellId,
133 const int &localVertexId,
134 SimplexId &vertexId) const override {
135#ifndef TTK_ENABLE_KAMIKAZE
136 if((!cellArray_) || (!cellNumber_))
137 return -1;
138#endif
139 vertexId = cellArray_->getCellVertex(cellId, localVertexId);
140 return 0;
141 }
142
144 const SimplexId &cellId) const override {
145#ifndef TTK_ENABLE_KAMIKAZE
146 if((!cellArray_) || (!cellNumber_))
147 return -1;
148#endif
149 return cellArray_->getCellVertexNumber(cellId);
150 }
151
153 return maxCellDim_;
154 }
155
156 inline const std::vector<std::array<SimplexId, 2>> *
158 return &edgeList_;
159 }
160
162 const SimplexId &edgeId,
163 const int &localLinkId,
164 SimplexId &linkId) const override {
165
166 linkId = edgeLinkData_[edgeId][localLinkId];
167 return 0;
168 }
169
171 const SimplexId &edgeId) const override {
172 return edgeLinkData_.size(edgeId);
173 }
174
175 inline const std::vector<std::vector<SimplexId>> *
177 edgeLinkData_.copyTo(edgeLinkList_);
178 return &edgeLinkList_;
179 }
180
182 const SimplexId &edgeId,
183 const int &localStarId,
184 SimplexId &starId) const override {
185
186 starId = edgeStarData_[edgeId][localStarId];
187 return 0;
188 }
189
191 const SimplexId &edgeId) const override {
192 return edgeStarData_.size(edgeId);
193 }
194
195 inline const std::vector<std::vector<SimplexId>> *
197 edgeStarData_.copyTo(edgeStarList_);
198 return &edgeStarList_;
199 }
200
201 inline int getEdgeTriangleInternal(const SimplexId &edgeId,
202 const int &localTriangleId,
203 SimplexId &triangleId) const override {
204 triangleId = edgeTriangleData_[edgeId][localTriangleId];
205 return 0;
206 }
207
208 inline SimplexId
209 getEdgeTriangleNumberInternal(const SimplexId &edgeId) const override {
210 return edgeTriangleData_.size(edgeId);
211 }
212
213 inline const std::vector<std::vector<SimplexId>> *
215 edgeTriangleData_.copyTo(edgeTriangleList_);
216 return &edgeTriangleList_;
217 }
218
219 inline int getEdgeVertexInternal(const SimplexId &edgeId,
220 const int &localVertexId,
221 SimplexId &vertexId) const override {
222#ifndef TTK_ENABLE_KAMIKAZE
223 if((edgeId < 0) || (edgeId >= (SimplexId)edgeList_.size()))
224 return -1;
225 if((localVertexId != 0) && (localVertexId != 1))
226 return -2;
227#endif
228 if(!localVertexId)
229 vertexId = edgeList_[edgeId][0];
230 else
231 vertexId = edgeList_[edgeId][1];
232 return 0;
233 }
234
235 inline SimplexId
237 return cellNumber_;
238 }
239
240 inline SimplexId getNumberOfEdgesInternal() const override {
241 return edgeList_.size();
242 }
243
244 inline SimplexId getNumberOfTrianglesInternal() const override {
245 return triangleList_.size();
246 }
247
248 inline SimplexId
250 return vertexNumber_;
251 }
252
253 inline const std::vector<std::array<SimplexId, 3>> *
255 return &triangleList_;
256 }
257
258 inline int getTriangleEdgeInternal(const SimplexId &triangleId,
259 const int &localEdgeId,
260 SimplexId &edgeId) const override {
261
262#ifndef TTK_ENABLE_KAMIKAZE
263 if((triangleId < 0)
264 || (triangleId >= (SimplexId)triangleEdgeList_.size()))
265 return -1;
266 if((localEdgeId < 0) || (localEdgeId > 2))
267 return -2;
268#endif
269
270 edgeId = triangleEdgeList_[triangleId][localEdgeId];
271
272 return 0;
273 }
274
276 const SimplexId &triangleId) const override {
277
278#ifndef TTK_ENABLE_KAMIKAZE
279 if((triangleId < 0)
280 || (triangleId >= (SimplexId)triangleEdgeList_.size()))
281 return -1;
282#endif
283
284 return triangleEdgeList_[triangleId].size();
285 }
286
287 inline const std::vector<std::vector<SimplexId>> *
289
291 return &triangleEdgeVector_;
292 }
293
295 const SimplexId &triangleId,
296 const int &localLinkId,
297 SimplexId &linkId) const override {
298
299 linkId = triangleLinkData_[triangleId][localLinkId];
300 return 0;
301 }
302
304 const SimplexId &triangleId) const override {
305 return triangleLinkData_[triangleId].size();
306 }
307
308 inline const std::vector<std::vector<SimplexId>> *
310 triangleLinkData_.copyTo(triangleLinkList_);
311 return &triangleLinkList_;
312 }
313
315 const SimplexId &triangleId,
316 const int &localStarId,
317 SimplexId &starId) const override {
318 starId = triangleStarData_[triangleId][localStarId];
319 return 0;
320 }
321
323 const SimplexId &triangleId) const override {
324 return triangleStarData_[triangleId].size();
325 }
326
327 inline const std::vector<std::vector<SimplexId>> *
329 triangleStarData_.copyTo(triangleStarList_);
330 return &triangleStarList_;
331 }
332
333 inline int getTriangleVertexInternal(const SimplexId &triangleId,
334 const int &localVertexId,
335 SimplexId &vertexId) const override {
336#ifndef TTK_ENABLE_KAMIKAZE
337 if((triangleId < 0) || (triangleId >= (SimplexId)triangleList_.size()))
338 return -1;
339 if((localVertexId < 0)
340 || (localVertexId >= (SimplexId)triangleList_[triangleId].size()))
341 return -2;
342#endif
343 vertexId = triangleList_[triangleId][localVertexId];
344 return 0;
345 }
346
347 inline int getVertexEdgeInternal(const SimplexId &vertexId,
348 const int &localEdgeId,
349 SimplexId &edgeId) const override {
350 edgeId = vertexEdgeData_[vertexId][localEdgeId];
351 return 0;
352 }
353
354 inline SimplexId
355 getVertexEdgeNumberInternal(const SimplexId &vertexId) const override {
356 return vertexEdgeData_[vertexId].size();
357 }
358
359 inline const std::vector<std::vector<SimplexId>> *
361 vertexEdgeData_.copyTo(vertexEdgeList_);
362 return &vertexEdgeList_;
363 }
364
366 const SimplexId &vertexId,
367 const int &localLinkId,
368 SimplexId &linkId) const override {
369
370 linkId = vertexLinkData_[vertexId][localLinkId];
371 return 0;
372 }
373
375 const SimplexId &vertexId) const override {
376 return vertexLinkData_[vertexId].size();
377 }
378
379 inline const std::vector<std::vector<SimplexId>> *
381 vertexLinkData_.copyTo(vertexLinkList_);
382 return &vertexLinkList_;
383 }
384
386 const SimplexId &vertexId,
387 const int &localNeighborId,
388 SimplexId &neighborId) const override {
389
390 neighborId = vertexNeighborData_[vertexId][localNeighborId];
391 return 0;
392 }
393
395 const SimplexId &vertexId) const override {
396 return vertexNeighborData_[vertexId].size();
397 }
398
399 inline const std::vector<std::vector<SimplexId>> *
401 vertexNeighborData_.copyTo(vertexNeighborList_);
402 return &vertexNeighborList_;
403 }
404
406 const SimplexId &vertexId, float &x, float &y, float &z) const override {
407#ifndef TTK_ENABLE_KAMIKAZE
408 if((vertexId < 0) || (vertexId >= vertexNumber_))
409 return -1;
410#endif
411
412 if(doublePrecision_) {
413 x = ((const double *)pointSet_)[3 * vertexId];
414 y = ((const double *)pointSet_)[3 * vertexId + 1];
415 z = ((const double *)pointSet_)[3 * vertexId + 2];
416 } else {
417 x = ((const float *)pointSet_)[3 * vertexId];
418 y = ((const float *)pointSet_)[3 * vertexId + 1];
419 z = ((const float *)pointSet_)[3 * vertexId + 2];
420 }
421
422 return 0;
423 }
424
426 const SimplexId &vertexId,
427 const int &localStarId,
428 SimplexId &starId) const override {
429 starId = vertexStarData_[vertexId][localStarId];
430 return 0;
431 }
432
434 const SimplexId &vertexId) const override {
435 return vertexStarData_[vertexId].size();
436 }
437
438 inline const std::vector<std::vector<SimplexId>> *
440 vertexStarData_.copyTo(vertexStarList_);
441 return &vertexStarList_;
442 }
443
444 inline int getVertexTriangleInternal(const SimplexId &vertexId,
445 const int &localTriangleId,
446 SimplexId &triangleId) const override {
447 triangleId = vertexTriangleData_[vertexId][localTriangleId];
448 return 0;
449 }
450
452 const SimplexId &vertexId) const override {
453 return vertexTriangleData_[vertexId].size();
454 }
455
456 inline const std::vector<std::vector<SimplexId>> *
458 vertexTriangleData_.copyTo(vertexTriangleList_);
459 return &vertexTriangleList_;
460 }
461
463 const SimplexId &edgeId) const override {
464#ifndef TTK_ENABLE_KAMIKAZE
465 if((edgeId < 0) || (edgeId >= (SimplexId)boundaryEdges_.size()))
466 return false;
467#endif
468 return boundaryEdges_[edgeId];
469 }
470
471 inline bool isEmpty() const override {
472 return !vertexNumber_;
473 }
474
476 const SimplexId &triangleId) const override {
477#ifndef TTK_ENABLE_KAMIKAZE
478 if((triangleId < 0)
479 || (triangleId >= (SimplexId)boundaryTriangles_.size()))
480 return false;
481#endif
482 return boundaryTriangles_[triangleId];
483 }
484
486 const SimplexId &vertexId) const override {
487#ifndef TTK_ENABLE_KAMIKAZE
488 if((vertexId < 0) || (vertexId >= (SimplexId)boundaryVertices_.size()))
489 return false;
490#endif
491 return boundaryVertices_[vertexId];
492 }
493
497
498 int preconditionCellEdgesInternal() override;
501
502 int preconditionEdgesInternal() override;
503 int preconditionEdgeLinksInternal() override;
504 int preconditionEdgeStarsInternal() override;
506
507 int preconditionTrianglesInternal() override;
511
512 int preconditionVertexEdgesInternal() override;
513 int preconditionVertexLinksInternal() override;
515 int preconditionVertexStarsInternal() override;
517
518 int preconditionManifoldInternal() override;
519
520#ifdef TTK_CELL_ARRAY_NEW
521 // Layout with connectivity + offset array (new)
522 inline int setInputCells(const SimplexId &cellNumber,
523 const LongSimplexId *connectivity,
524 const LongSimplexId *offset) {
525
526 // Cell Check
527 {
528 if(cellNumber > 0) {
529 const auto &cellDimension = offset[1] - offset[0] - 1;
530
531 if(cellDimension < 0 || cellDimension > 3) {
532 this->printErr("Unable to create triangulation for cells of "
533 "dimension 4 or higher ("
534 + std::to_string(cellDimension) + ").");
535 return -1;
536 }
537
538 bool error = false;
539
540#ifdef TTK_ENABLE_OPENMP
541#pragma omp parallel for num_threads(this->threadNumber_)
542#endif
543 for(SimplexId i = 0; i < cellNumber; i++) {
544 if(offset[i + 1] - offset[i] - 1 != cellDimension) {
545#ifdef TTK_ENABLE_OPENMP
546#pragma omp atomic write
547#endif // TTK_ENABLE_OPENMP
548 error = true;
549 }
550 }
551
552 if(error) {
553 this->printErr("Unable to create triangulation for "
554 "inhomogeneous\ncell dimensions.");
555 return -2;
556 }
557 }
558 }
559
560 if(cellNumber_)
561 clear();
562
563 cellNumber_ = cellNumber;
564
565 cellArray_
566 = std::make_shared<CellArray>(connectivity, offset, cellNumber);
567
568 // TODO: ASSUME Regular Mesh Here to compute dimension!
569 if(cellNumber) {
570 maxCellDim_ = cellArray_->getCellVertexNumber(0) - 1;
571 }
572 return 0;
573 }
574#else
575 // Flat layout with a single array (legacy & default one)
576 inline int setInputCells(const SimplexId &cellNumber,
577 const LongSimplexId *cellArray) {
578 if(cellNumber_)
579 clear();
580
581 cellNumber_ = cellNumber;
582
583 if(cellNumber) {
584 // assume regular mesh here to compute dimension
585 cellArray_ = std::make_shared<CellArray>(
586 cellArray, cellNumber, cellArray[0] - 1);
587 maxCellDim_ = cellArray[0] - 1;
588 }
589 return 0;
590 }
591#endif
592
593 inline int setInputPoints(const SimplexId &pointNumber,
594 const void *pointSet,
595 const bool &doublePrecision = false) {
596
597 if(vertexNumber_)
598 clear();
599
600 vertexNumber_ = pointNumber;
601 pointSet_ = pointSet;
602 doublePrecision_ = doublePrecision;
603 return 0;
604 }
605
611 int writeToFile(std::ofstream &stream) const;
615 int writeToFileASCII(std::ofstream &stream) const;
621 int readFromFile(std::ifstream &stream);
622
623#ifdef TTK_ENABLE_MPI
624
625 inline void setCellsGlobalIds(const LongSimplexId *const cellGid) {
626 this->cellGid_ = cellGid;
627 }
628 inline void setVertsGlobalIds(const LongSimplexId *array) {
629 this->vertGid_ = array;
630 }
631
632 inline SimplexId
633 getVertexGlobalIdInternal(const SimplexId lvid) const override {
634 return this->vertGid_[lvid];
635 }
636
637 inline SimplexId
638 getVertexLocalIdInternal(const SimplexId gvid) const override {
639 const auto it{this->vertexGidToLid_.find(gvid)};
640 if(it == this->vertexGidToLid_.end()) {
641 return -1;
642 }
643 return it->second;
644 }
645
646 inline SimplexId
647 getCellGlobalIdInternal(const SimplexId lcid) const override {
648 return this->cellGid_[lcid];
649 }
650
651 inline SimplexId
652 getCellLocalIdInternal(const SimplexId gcid) const override {
653 const auto it{this->cellGidToLid_.find(gcid)};
654#ifndef TTK_ENABLE_KAMIKAZE
655 if(it == this->cellGidToLid_.end()) {
656 return -1;
657 }
658#endif // TTK_ENABLE_KAMIKAZE
659 return it->second;
660 }
661
662 inline SimplexId
663 getEdgeGlobalIdInternal(const SimplexId leid) const override {
664 return this->edgeLidToGid_[leid];
665 }
666
667 inline SimplexId
668 getEdgeLocalIdInternal(const SimplexId geid) const override {
669 const auto it = this->edgeGidToLid_.find(geid);
670#ifndef TTK_ENABLE_KAMIKAZE
671 if(it == this->edgeGidToLid_.end()) {
672 return -1;
673 }
674#endif // TTK_ENABLE_KAMIKAZE
675 return it->second;
676 }
677
678 inline SimplexId
679 getTriangleGlobalIdInternal(const SimplexId ltid) const override {
680 return this->triangleLidToGid_[ltid];
681 }
682
683 inline SimplexId
684 getTriangleLocalIdInternal(const SimplexId gtid) const override {
685 const auto it = this->triangleGidToLid_.find(gtid);
686#ifndef TTK_ENABLE_KAMIKAZE
687 if(it == this->triangleGidToLid_.end()) {
688 return -1;
689 }
690#endif // TTK_ENABLE_KAMIKAZE
691 return it->second;
692 }
693
694 inline int setVertexRankArray(const int *rankArray) override {
695 vertexRankArray_.resize(vertexNumber_);
696 std::copy(rankArray, rankArray + vertexNumber_, vertexRankArray_.begin());
697 return 0;
698 }
699
700 inline int setCellRankArray(const int *rankArray) override {
701 cellRankArray_.resize(cellNumber_);
702 std::copy(rankArray, rankArray + cellNumber_, cellRankArray_.begin());
703 return 0;
704 }
705
706 inline int getVertexRankInternal(const SimplexId lvid) const override {
707 return this->vertexRankArray_[lvid];
708 }
709
710 inline std::unordered_map<SimplexId, SimplexId> &getVertexGlobalIdMap() {
711 return this->vertexGidToLid_;
712 }
713
714 inline void setBoundingBox(const double *const bBox) {
715 this->boundingBox_
716 = {bBox[0], bBox[1], bBox[2], bBox[3], bBox[4], bBox[5]};
717 }
718
719 protected:
720 template <typename Func0, typename Func1, typename Func2>
721 int exchangeDistributedInternal(const Func0 &getGlobalSimplexId,
722 const Func1 &storeGlobalSimplexId,
723 const Func2 &iterCond,
724 const int nSimplicesPerCell);
725
726 int preconditionDistributedCellRanges();
727 size_t
728 computeCellRangeOffsets(std::vector<size_t> &nSimplicesPerRange) const;
729
730 int preconditionDistributedCells() override;
731 int preconditionExchangeGhostCells() override;
732 int preconditionDistributedEdges() override;
733 int preconditionDistributedVertices() override;
734 int preconditionExchangeGhostVertices() override;
735 int preconditionDistributedTriangles() override;
736 int preconditionVertexRankArray();
737 int preconditionCellRankArray();
738 int preconditionEdgeRankArray() override;
739 int preconditionTriangleRankArray() override;
740
741 // range of (local) cells owned by the current rank that have
742 // contiguous global ids (to label edges & triangles)
743 struct CellRange {
744 // rank-local range id
745 size_t id;
746 // range beginning (global cell id)
747 size_t begin;
748 // range end (inclusive, global cell id)
749 size_t end;
750 // owner rank
751 size_t rank;
752
753 static inline MPI_Datatype getMPIType() {
754 MPI_Datatype res{};
755 const auto cellRangeSize = sizeof(CellRange) / sizeof(size_t);
756 MPI_Type_contiguous(cellRangeSize, ttk::getMPIType(size_t{}), &res);
757 return res;
758 }
759 };
760 // cell ranges per rank
761 std::vector<CellRange> localCellRanges_{};
762 // cell ranges from all ranks (gathered on rank 0)
763 std::vector<CellRange> gatheredCellRanges_{};
764 // number of CellRanges per rank
765 std::vector<int> nRangesPerRank_{};
766
767 // "GlobalCellIds" from "Generate Global Ids"
768 const LongSimplexId *cellGid_{};
769 // "GlobalPointIds" from "Generate Global Ids"
770 const LongSimplexId *vertGid_{};
771
772 // inverse of vertGid_
773 std::unordered_map<SimplexId, SimplexId> vertexGidToLid_{};
774 // inverse of cellGid_
775 std::unordered_map<SimplexId, SimplexId> cellGidToLid_{};
776
777 std::vector<SimplexId> edgeLidToGid_{};
778 std::unordered_map<SimplexId, SimplexId> edgeGidToLid_{};
779 std::vector<SimplexId> triangleLidToGid_{};
780 std::unordered_map<SimplexId, SimplexId> triangleGidToLid_{};
781
782 std::array<double, 6> boundingBox_{};
783
784 std::vector<int> vertexRankArray_{};
785 std::vector<int> cellRankArray_{};
786 std::vector<int> edgeRankArray_{};
787 std::vector<int> triangleRankArray_{};
788
789#endif // TTK_ENABLE_MPI
790
791 private:
792 bool doublePrecision_;
793 SimplexId cellNumber_, vertexNumber_;
794 const void *pointSet_;
795 int maxCellDim_;
796 std::shared_ptr<CellArray> cellArray_;
797
798 FlatJaggedArray vertexNeighborData_{};
799 FlatJaggedArray cellNeighborData_{};
800 FlatJaggedArray vertexEdgeData_{};
801 FlatJaggedArray vertexTriangleData_{};
802 FlatJaggedArray edgeTriangleData_{};
803 FlatJaggedArray vertexStarData_{};
804 FlatJaggedArray edgeStarData_{};
805 FlatJaggedArray triangleStarData_{};
806 FlatJaggedArray vertexLinkData_{};
807 FlatJaggedArray edgeLinkData_{};
808 FlatJaggedArray triangleLinkData_{};
809
810 // Char array that identifies the file format.
811 static const char *magicBytes_;
812 // Current version of the file format. To be incremented at every
813 // breaking change to keep backward compatibility.
814 static const unsigned long formatVersion_;
815 };
816} // namespace ttk
#define TTK_TRIANGULATION_INTERNAL(NAME)
AbstractTriangulation is an interface class that defines an interface for efficient traversal methods...
std::vector< std::vector< SimplexId > > vertexStarList_
std::vector< std::array< SimplexId, 3 > > triangleList_
std::vector< bool > boundaryTriangles_
std::vector< std::array< SimplexId, 4 > > tetraTriangleList_
std::vector< std::vector< SimplexId > > triangleLinkList_
std::vector< std::vector< SimplexId > > vertexNeighborList_
std::vector< std::vector< SimplexId > > vertexLinkList_
std::vector< std::vector< SimplexId > > edgeLinkList_
std::vector< bool > boundaryVertices_
std::vector< std::vector< SimplexId > > vertexEdgeList_
std::vector< std::vector< SimplexId > > triangleEdgeVector_
std::vector< std::vector< SimplexId > > triangleStarList_
std::vector< std::vector< SimplexId > > cellNeighborList_
std::vector< std::array< SimplexId, 2 > > edgeList_
std::vector< std::array< SimplexId, 6 > > tetraEdgeList_
std::vector< std::vector< SimplexId > > edgeStarList_
std::vector< std::vector< SimplexId > > vertexTriangleList_
std::vector< std::vector< SimplexId > > cellEdgeVector_
std::vector< std::vector< SimplexId > > cellTriangleVector_
std::vector< std::array< SimplexId, 3 > > triangleEdgeList_
std::vector< std::vector< SimplexId > > edgeTriangleList_
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition: Debug.h:149
ExplicitTriangulation is a class that provides time efficient traversal methods on triangulations of ...
SimplexId getTriangleEdgeNumberInternal(const SimplexId &triangleId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getVertexNeighborNumber(const SimplexId &vertexId) const override
int getVertexTriangleInternal(const SimplexId &vertexId, const int &localTriangleId, SimplexId &triangleId) const override
size_t footprint(size_t size=0) const
const std::vector< std::vector< SimplexId > > * getTriangleEdgesInternal() override
int TTK_TRIANGULATION_INTERNAL() getVertexPoint(const SimplexId &vertexId, float &x, float &y, float &z) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getVertexStars() override
int TTK_TRIANGULATION_INTERNAL() getVertexNeighbor(const SimplexId &vertexId, const int &localNeighborId, SimplexId &neighborId) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getEdgeStars() override
int TTK_TRIANGULATION_INTERNAL() getVertexStar(const SimplexId &vertexId, const int &localStarId, SimplexId &starId) const override
int preconditionTriangleStarsInternal() override
int TTK_TRIANGULATION_INTERNAL() getCellVertex(const SimplexId &cellId, const int &localVertexId, SimplexId &vertexId) const override
int getVertexEdgeInternal(const SimplexId &vertexId, const int &localEdgeId, SimplexId &edgeId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getTriangleLinkNumber(const SimplexId &triangleId) const override
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
const std::vector< std::vector< SimplexId > > * getVertexTrianglesInternal() override
int setInputCells(const SimplexId &cellNumber, const LongSimplexId *cellArray)
SimplexId TTK_TRIANGULATION_INTERNAL() getVertexLinkNumber(const SimplexId &vertexId) const override
int preconditionCellTrianglesInternal() override
int preconditionVertexNeighborsInternal() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getTriangleStars() 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
bool TTK_TRIANGULATION_INTERNAL() isVertexOnBoundary(const SimplexId &vertexId) const override
const std::vector< std::array< SimplexId, 3 > > *TTK_TRIANGULATION_INTERNAL() getTriangles() override
int preconditionTriangleLinksInternal() override
int getTriangleVertexInternal(const SimplexId &triangleId, const int &localVertexId, SimplexId &vertexId) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getVertexLinks() override
ExplicitTriangulation & operator=(ExplicitTriangulation &&)=default
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getEdgeLinks() override
SimplexId getVertexTriangleNumberInternal(const SimplexId &vertexId) const override
ExplicitTriangulation & operator=(const ExplicitTriangulation &)=default
const std::vector< std::array< SimplexId, 2 > > *TTK_TRIANGULATION_INTERNAL() getEdges() override
SimplexId TTK_TRIANGULATION_INTERNAL() getCellNeighborNumber(const SimplexId &cellId) const override
bool TTK_TRIANGULATION_INTERNAL() isTriangleOnBoundary(const SimplexId &triangleId) const override
int getEdgeVertexInternal(const SimplexId &edgeId, const int &localVertexId, SimplexId &vertexId) const override
const std::vector< std::vector< SimplexId > > * getCellEdgesInternal() override
int setInputPoints(const SimplexId &pointNumber, const void *pointSet, const bool &doublePrecision=false)
int preconditionBoundaryEdgesInternal() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getCellNeighbors() override
void convertToVector(const std::vector< std::array< SimplexId, N > > &table, std::vector< std::vector< SimplexId > > &vec)
int preconditionEdgeTrianglesInternal() override
int getCellEdgeInternal(const SimplexId &cellId, const int &localEdgeId, SimplexId &edgeId) const override
int preconditionBoundaryVerticesInternal() override
SimplexId getNumberOfEdgesInternal() const override
const std::vector< std::vector< SimplexId > > * getCellTrianglesInternal() override
int preconditionTriangleEdgesInternal() override
ExplicitTriangulation(const ExplicitTriangulation &)=default
SimplexId TTK_TRIANGULATION_INTERNAL() getEdgeStarNumber(const SimplexId &edgeId) const override
int TTK_TRIANGULATION_INTERNAL() getDimensionality() const override
const std::vector< std::vector< SimplexId > > * getEdgeTrianglesInternal() override
int TTK_TRIANGULATION_INTERNAL() getEdgeStar(const SimplexId &edgeId, const int &localStarId, SimplexId &starId) const override
int TTK_TRIANGULATION_INTERNAL() getEdgeLink(const SimplexId &edgeId, const int &localLinkId, SimplexId &linkId) const override
ExplicitTriangulation(ExplicitTriangulation &&)=default
int preconditionCellNeighborsInternal() override
const std::vector< std::vector< SimplexId > > * getVertexEdgesInternal() override
int getEdgeTriangleInternal(const SimplexId &edgeId, const int &localTriangleId, SimplexId &triangleId) const override
int TTK_TRIANGULATION_INTERNAL() getCellNeighbor(const SimplexId &cellId, const int &localNeighborId, SimplexId &neighborId) const override
bool TTK_TRIANGULATION_INTERNAL() isEdgeOnBoundary(const SimplexId &edgeId) const override
int writeToFileASCII(std::ofstream &stream) const
Write internal state to disk using an ASCII format.
SimplexId getEdgeTriangleNumberInternal(const SimplexId &edgeId) const override
SimplexId getCellTriangleNumberInternal(const SimplexId &cellId) const override
SimplexId getVertexEdgeNumberInternal(const SimplexId &vertexId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getNumberOfVertices() const override
SimplexId TTK_TRIANGULATION_INTERNAL() getEdgeLinkNumber(const SimplexId &edgeId) const override
int getTriangleEdgeInternal(const SimplexId &triangleId, const int &localEdgeId, SimplexId &edgeId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getVertexStarNumber(const SimplexId &vertexId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getCellVertexNumber(const SimplexId &cellId) const override
int getCellTriangleInternal(const SimplexId &cellId, const int &localTriangleId, SimplexId &triangleId) const override
int TTK_TRIANGULATION_INTERNAL() getVertexLink(const SimplexId &vertexId, const int &localLinkId, SimplexId &linkId) const override
int preconditionVertexTrianglesInternal() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getVertexNeighbors() override
int readFromFile(std::ifstream &stream)
Read from disk into internal state.
int TTK_TRIANGULATION_INTERNAL() getTriangleLink(const SimplexId &triangleId, const int &localLinkId, SimplexId &linkId) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getTriangleLinks() override
void copyTo(std::vector< std::vector< SimplexId > > &dst, int threadNumber=1) const
Copy buffers to a std::vector<std::vector<SimplexId>>
SimplexId size(SimplexId id) const
Get the size of a particular sub-vector.
The Topology ToolKit.
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