TTK
Loading...
Searching...
No Matches
CompactTriangulation.h
Go to the documentation of this file.
1
28
29#pragma once
30
31// base code includes
33#include <CellArray.h>
34#include <FlatJaggedArray.h>
35#include <algorithm>
36#include <boost/functional/hash.hpp>
37#include <boost/unordered_map.hpp>
38#include <boost/unordered_set.hpp>
39#include <list>
40
41namespace ttk {
42
44 private:
45 /* components */
46 SimplexId nid;
47 std::vector<std::array<SimplexId, 2>> internalEdgeList_;
48 std::vector<std::array<SimplexId, 3>> internalTriangleList_;
49 boost::unordered_map<std::array<SimplexId, 2>, SimplexId> internalEdgeMap_;
50 boost::unordered_map<std::array<SimplexId, 2>, SimplexId> externalEdgeMap_;
51 boost::unordered_map<std::array<SimplexId, 3>, SimplexId>
52 internalTriangleMap_;
53 boost::unordered_map<std::array<SimplexId, 3>, SimplexId>
54 externalTriangleMap_;
55 /* boundary cells */
56 std::vector<bool> boundaryVertices_;
57 std::vector<bool> boundaryEdges_;
58 std::vector<bool> boundaryTriangles_;
59 /* vertex relationships */
60 FlatJaggedArray vertexEdges_;
61 FlatJaggedArray vertexLinks_;
62 FlatJaggedArray vertexNeighbors_;
63 FlatJaggedArray vertexStars_;
64 FlatJaggedArray vertexTriangles_;
65 /* edge relationships */
66 // edgeVertex relation can be extracted from internal edge list
67 FlatJaggedArray edgeLinks_;
68 FlatJaggedArray edgeStars_;
69 FlatJaggedArray edgeTriangles_;
70 /* triangle relationships */
71 // triangleVertex relation can be extracted from internal triangle list
72 std::vector<std::array<SimplexId, 3>> triangleEdges_;
73 FlatJaggedArray triangleLinks_;
74 FlatJaggedArray triangleStars_;
75 /* cell relationships */
76 std::vector<std::array<SimplexId, 6>> tetraEdges_;
77 FlatJaggedArray cellNeighbors_;
78 std::vector<std::array<SimplexId, 4>> tetraTriangles_;
79
80 public:
81 ImplicitCluster() = default;
82 ImplicitCluster(SimplexId id) : nid(id) {
83 }
84 ~ImplicitCluster() = default;
85
87 };
88
90
91 // different id types for compact triangulation
92 enum class SIMPLEX_ID { EDGE_ID = 1, TRIANGLE_ID = 2 };
93
94 public:
101
105 inline int setInputPoints(const SimplexId &pointNumber,
106 const void *pointSet,
107 const int *indexArray,
108 const bool &doublePrecision = false) {
109
110 if(vertexNumber_)
111 clear();
112
113 vertexNumber_ = pointNumber;
114 pointSet_ = pointSet;
115 vertexIndices_ = indexArray;
116 doublePrecision_ = doublePrecision;
117
118 return 0;
119 }
120
125#ifdef TTK_CELL_ARRAY_NEW
126 // Layout with connectivity + offset array (new)
127 inline int setInputCells(const SimplexId &cellNumber,
128 const LongSimplexId *connectivity,
129 const LongSimplexId *offset) {
130
131 // Cell Check
132 {
133 if(cellNumber > 0) {
134 const auto &cellDimension = offset[1] - offset[0] - 1;
135
136 if(cellDimension < 0 || cellDimension > 3) {
137 this->printErr("Unable to create triangulation for cells of "
138 "dimension 4 or higher ("
139 + std::to_string(cellDimension) + ").");
140 return -1;
141 }
142
143 bool error = false;
144
145#ifdef TTK_ENABLE_OPENMP
146#pragma omp parallel for num_threads(this->threadNumber_)
147#endif
148 for(SimplexId i = 0; i < cellNumber; i++) {
149 if(offset[i + 1] - offset[i] - 1 != cellDimension) {
150#ifdef TTK_ENABLE_OPENMP
151#pragma omp atomic write
152#endif // TTK_ENABLE_OPENMP
153 error = true;
154 }
155 }
156
157 if(error) {
158 this->printErr("Unable to create triangulation for "
159 "inhomogeneous\ncell dimensions.");
160 return -2;
161 }
162 }
163 }
164
165 if(cellNumber_)
166 clear();
167
168 cellNumber_ = cellNumber;
169 std::vector<SimplexId> vertexMap(vertexNumber_);
170 reorderVertices(vertexMap);
171 reorderCells(vertexMap, cellNumber, connectivity, offset);
173 = std::make_shared<CellArray>(connectivity, offset, cellNumber);
174
175 // ASSUME Regular Mesh Here to compute dimension!
176 if(cellNumber) {
177 if(cellArray_->getCellVertexNumber(0) == 3) {
178 maxCellDim_ = 2;
180 } else {
181 maxCellDim_ = 3;
182 }
183 }
184
185 return 0;
186 }
187#else
188 // Flat layout with a single array (legacy & default one)
189 inline int setInputCells(const SimplexId &cellNumber,
190 const LongSimplexId *cellArray) {
191 if(cellNumber_)
192 clear();
193
194 cellNumber_ = cellNumber;
195
196 if(cellNumber) {
197 // assume regular mesh here to compute dimension
198 maxCellDim_ = cellArray[0] - 1;
199 std::vector<SimplexId> vertexMap(vertexNumber_);
200 reorderVertices(vertexMap);
201 reorderCells(vertexMap, cellArray);
202 cellArray_ = std::make_shared<CellArray>(
203 cellArray, cellNumber, cellArray[0] - 1);
204 }
205
206 return 0;
207 }
208#endif
209
213 int reorderVertices(std::vector<SimplexId> &vertexMap);
214
218 int reorderCells(const std::vector<SimplexId> &vertexMap,
219 const SimplexId &cellNumber,
220 const LongSimplexId *connectivity,
221 const LongSimplexId *offset);
222
226 int reorderCells(const std::vector<SimplexId> &vertexMap,
227 const LongSimplexId *cellArray);
228
229 inline int getCellEdgeInternal(const SimplexId &cellId,
230 const int &localEdgeId,
231 SimplexId &edgeId) const override {
232
233#ifndef TTK_ENABLE_KAMIKAZE
234 if((cellId < 0) || (cellId >= cellNumber_)) {
235 edgeId = -1;
236 return 0;
237 }
238 if(localEdgeId < 0) {
239 edgeId = -2;
240 return 0;
241 }
242#endif
243
244 SimplexId nid = vertexIndices_[cellArray_->getCellVertex(cellId, 0)];
245 SimplexId localCellId = cellId - cellIntervals_[nid - 1] - 1;
246 ImplicitCluster *exnode = searchCache(nid);
247 if(exnode->tetraEdges_.empty()) {
248 getClusterTetraEdges(exnode);
249 }
250
251 if(localEdgeId >= (int)(exnode->tetraEdges_)[localCellId].size()) {
252 edgeId = -2;
253 } else {
254 edgeId = (exnode->tetraEdges_)[localCellId][localEdgeId];
255 }
256 return 0;
257 }
258
259 inline SimplexId
260 getCellEdgeNumberInternal(const SimplexId &cellId) const override {
261
262#ifndef TTK_ENABLE_KAMIKAZE
263 if((cellId < 0) || (cellId >= cellNumber_))
264 return -1;
265#endif
266
267 (void)cellId;
268 return (maxCellDim_ + 1) * maxCellDim_ / 2;
269 }
270
271 inline const std::vector<std::vector<SimplexId>> *
273 if(cellEdgeVector_.empty()) {
275 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
276 ImplicitCluster *exnode = searchCache(nid);
277 if(exnode->tetraEdges_.empty()) {
278 getClusterTetraEdges(exnode);
279 }
280 for(size_t i = 0; i < exnode->tetraEdges_.size(); i++) {
281 cellEdgeVector_.emplace_back(exnode->tetraEdges_.at(i).begin(),
282 exnode->tetraEdges_.at(i).end());
283 }
284 }
285 }
286 return &cellEdgeVector_;
287 }
288
290 const SimplexId &cellId,
291 const int &localNeighborId,
292 SimplexId &neighborId) const override {
293
294#ifndef TTK_ENABLE_KAMIKAZE
295 if((cellId < 0) || (cellId >= cellNumber_)) {
296 neighborId = -1;
297 return 0;
298 }
299 if(localNeighborId < 0) {
300 neighborId = -2;
301 return 0;
302 }
303#endif
304
305 SimplexId nid = vertexIndices_[cellArray_->getCellVertex(cellId, 0)];
306 SimplexId localCellId = cellId - cellIntervals_[nid - 1] - 1;
307 ImplicitCluster *exnode = searchCache(nid);
308 if(exnode->cellNeighbors_.empty()) {
310 }
311
312 if(localNeighborId >= exnode->cellNeighbors_.size(localCellId)) {
313 neighborId = -2;
314 } else {
315 neighborId = exnode->cellNeighbors_.get(localCellId, localNeighborId);
316 }
317 return 0;
318 }
319
321 const SimplexId &cellId) const override {
322
323#ifndef TTK_ENABLE_KAMIKAZE
324 if((cellId < 0) || (cellId >= cellNumber_))
325 return -1;
326#endif
327
328 SimplexId nid = vertexIndices_[cellArray_->getCellVertex(cellId, 0)];
329 SimplexId localCellId = cellId - cellIntervals_[nid - 1] - 1;
330 ImplicitCluster *exnode = searchCache(nid);
331 if(exnode->cellNeighbors_.empty()) {
333 }
334 return exnode->cellNeighbors_.size(localCellId);
335 }
336
337 inline const std::vector<std::vector<SimplexId>> *
340 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
341 std::vector<std::vector<SimplexId>> localCellNeighbors;
342 ImplicitCluster *exnode = searchCache(nid);
343 if(exnode->cellNeighbors_.empty()) {
345 }
346 exnode->cellNeighbors_.copyTo(localCellNeighbors);
348 localCellNeighbors.begin(),
349 localCellNeighbors.end());
350 }
351 return &cellNeighborList_;
352 }
353
354 inline int getCellTriangleInternal(const SimplexId &cellId,
355 const int &localTriangleId,
356 SimplexId &triangleId) const override {
357
358#ifndef TTK_ENABLE_KAMIKAZE
359 if((cellId < 0) || (cellId >= cellNumber_)) {
360 triangleId = -1;
361 return 0;
362 }
363 if((localTriangleId < 0)
364 || (localTriangleId >= getCellTriangleNumber(cellId))) {
365 triangleId = -2;
366 return 0;
367 }
368#endif
369
370 SimplexId nid = vertexIndices_[cellArray_->getCellVertex(cellId, 0)];
371 SimplexId localCellId = cellId - cellIntervals_[nid - 1] - 1;
372 ImplicitCluster *exnode = searchCache(nid);
373 if(exnode->tetraTriangles_.empty()) {
375 }
376 triangleId = (exnode->tetraTriangles_)[localCellId][localTriangleId];
377 return 0;
378 }
379
380 inline SimplexId
381 getCellTriangleNumberInternal(const SimplexId &cellId) const override {
382
383#ifndef TTK_ENABLE_KAMIKAZE
384 if((cellId < 0) || (cellId >= cellNumber_))
385 return -1;
386#endif
387
388 (void)cellId;
389 return (maxCellDim_ + 1) * maxCellDim_ * (maxCellDim_ - 1) / 6;
390 }
391
392 inline const std::vector<std::vector<SimplexId>> *
394 if(cellTriangleVector_.empty()) {
396 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
397 ImplicitCluster *exnode = searchCache(nid);
398 if(exnode->tetraTriangles_.empty()) {
400 }
401 for(size_t i = 0; i < exnode->tetraTriangles_.size(); i++) {
402 cellTriangleVector_.emplace_back(
403 exnode->tetraTriangles_.at(i).begin(),
404 exnode->tetraTriangles_.at(i).end());
405 }
406 }
407 }
408 return &cellTriangleVector_;
409 }
410
412 const SimplexId &cellId,
413 const int &localVertexId,
414 SimplexId &vertexId) const override {
415
416#ifndef TTK_ENABLE_KAMIKAZE
417 if((cellId < 0) || (cellId >= cellNumber_)) {
418 vertexId = -1;
419 return 0;
420 }
421 if((localVertexId < 0)
422 || (localVertexId >= cellArray_->getCellVertexNumber(cellId))) {
423 vertexId = -2;
424 return 0;
425 }
426#endif
427
428 vertexId = cellArray_->getCellVertex(cellId, localVertexId);
429 return 0;
430 }
431
433 const SimplexId &cellId) const override {
434
435#ifndef TTK_ENABLE_KAMIKAZE
436 if((cellId < 0) || (cellId >= cellNumber_))
437 return -1;
438#endif
439
440 return cellArray_->getCellVertexNumber(cellId);
441 }
442
444 return maxCellDim_;
445 }
446
447 inline const std::vector<std::array<SimplexId, 2>> *
449 edgeList_.reserve(edgeIntervals_.back() + 1);
450 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
451 ImplicitCluster *exnode = searchCache(nid);
452 if(exnode->internalEdgeList_.empty())
453 buildInternalEdgeMap(exnode, true, false);
454 edgeList_.insert(edgeList_.end(), exnode->internalEdgeList_.begin(),
455 exnode->internalEdgeList_.end());
456 }
457 return &edgeList_;
458 }
459
461 const SimplexId &edgeId,
462 const int &localLinkId,
463 SimplexId &linkId) const override {
464
465#ifndef TTK_ENABLE_KAMIKAZE
466 if((edgeId < 0) || (edgeId > edgeIntervals_.back())) {
467 linkId = -1;
468 return 0;
469 }
470 if(localLinkId < 0) {
471 linkId = -2;
472 return 0;
473 }
474#endif
475
476 SimplexId nid = findNodeIndex(edgeId, SIMPLEX_ID::EDGE_ID);
477 SimplexId localEdgeId = edgeId - edgeIntervals_[nid - 1] - 1;
478 ImplicitCluster *exnode = searchCache(nid);
479 if(exnode->edgeLinks_.empty()) {
480 getClusterEdgeLinks(exnode);
481 }
482
483 if(localLinkId >= exnode->edgeLinks_.size(localEdgeId)) {
484 linkId = -2;
485 } else {
486 linkId = exnode->edgeLinks_.get(localEdgeId, localLinkId);
487 }
488 return 0;
489 }
490
492 const SimplexId &edgeId) const override {
493
494#ifndef TTK_ENABLE_KAMIKAZE
495 if((edgeId < 0) || (edgeId > edgeIntervals_.back()))
496 return -1;
497#endif
498
499 SimplexId nid = findNodeIndex(edgeId, SIMPLEX_ID::EDGE_ID);
500 SimplexId localEdgeId = edgeId - edgeIntervals_[nid - 1] - 1;
501 ImplicitCluster *exnode = searchCache(nid);
502 if(exnode->edgeLinks_.empty()) {
503 getClusterEdgeLinks(exnode);
504 }
505 return exnode->edgeLinks_.size(localEdgeId);
506 }
507
508 inline const std::vector<std::vector<SimplexId>> *
510 edgeLinkList_.reserve(edgeIntervals_.back() + 1);
511 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
512 std::vector<std::vector<SimplexId>> localEdgeLinks;
513 ImplicitCluster *exnode = searchCache(nid);
514 if(exnode->edgeLinks_.empty()) {
515 getClusterEdgeLinks(exnode);
516 }
517 exnode->edgeLinks_.copyTo(localEdgeLinks);
518 edgeLinkList_.insert(
519 edgeLinkList_.end(), localEdgeLinks.begin(), localEdgeLinks.end());
520 }
521 return &edgeLinkList_;
522 }
523
525 const SimplexId &edgeId,
526 const int &localStarId,
527 SimplexId &starId) const override {
528
529#ifndef TTK_ENABLE_KAMIKAZE
530 if((edgeId < 0) || (edgeId > edgeIntervals_.back())) {
531 starId = -1;
532 return 0;
533 }
534 if(localStarId < 0) {
535 starId = -2;
536 return 0;
537 }
538#endif
539
540 SimplexId nid = findNodeIndex(edgeId, SIMPLEX_ID::EDGE_ID);
541 SimplexId localEdgeId = edgeId - edgeIntervals_[nid - 1] - 1;
542 ImplicitCluster *exnode = searchCache(nid);
543 if(exnode->edgeStars_.empty()) {
544 getClusterEdgeStars(exnode);
545 }
546
547 if(localStarId >= exnode->edgeStars_.size(localEdgeId)) {
548 starId = -2;
549 } else {
550 starId = exnode->edgeStars_.get(localEdgeId, localStarId);
551 }
552 return 0;
553 }
554
556 const SimplexId &edgeId) const override {
557
558#ifndef TTK_ENABLE_KAMIKAZE
559 if((edgeId < 0) || (edgeId > edgeIntervals_.back()))
560 return -1;
561#endif
562 SimplexId nid = findNodeIndex(edgeId, SIMPLEX_ID::EDGE_ID);
563 ImplicitCluster *exnode = searchCache(nid);
564 SimplexId localEdgeId = edgeId - edgeIntervals_[nid - 1] - 1;
565 if(exnode->edgeStars_.empty()) {
566 getClusterEdgeStars(exnode);
567 }
568 return exnode->edgeStars_.size(localEdgeId);
569 }
570
571 inline const std::vector<std::vector<SimplexId>> *
573 edgeStarList_.reserve(edgeIntervals_.back() + 1);
574 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
575 std::vector<std::vector<SimplexId>> localEdgeStars;
576 ImplicitCluster *exnode = searchCache(nid);
577 if(exnode->edgeStars_.empty()) {
578 getClusterEdgeStars(exnode);
579 }
580 exnode->edgeStars_.copyTo(localEdgeStars);
581 edgeStarList_.insert(
582 edgeStarList_.end(), localEdgeStars.begin(), localEdgeStars.end());
583 }
584 return &edgeStarList_;
585 }
586
587 inline int getEdgeTriangleInternal(const SimplexId &edgeId,
588 const int &localTriangleId,
589 SimplexId &triangleId) const override {
590
591#ifndef TTK_ENABLE_KAMIKAZE
592 if((edgeId < 0) || (edgeId > (SimplexId)edgeIntervals_.back())) {
593 triangleId = -1;
594 return 0;
595 }
596 if(localTriangleId < 0) {
597 triangleId = -2;
598 return 0;
599 }
600#endif
601
602 SimplexId nid = findNodeIndex(edgeId, SIMPLEX_ID::EDGE_ID);
603 SimplexId localEdgeId = edgeId - edgeIntervals_[nid - 1] - 1;
604 ImplicitCluster *exnode = searchCache(nid);
605 if(exnode->edgeTriangles_.empty()) {
607 }
608
609 if(localTriangleId >= exnode->edgeTriangles_.size(localEdgeId)) {
610 triangleId = -2;
611 } else {
612 triangleId = exnode->edgeTriangles_.get(localEdgeId, localTriangleId);
613 }
614 return 0;
615 }
616
617 inline SimplexId
618 getEdgeTriangleNumberInternal(const SimplexId &edgeId) const override {
619
620#ifndef TTK_ENABLE_KAMIKAZE
621 if((edgeId < 0) || (edgeId > (SimplexId)edgeIntervals_.back()))
622 return -1;
623#endif
624
625 SimplexId nid = findNodeIndex(edgeId, SIMPLEX_ID::EDGE_ID);
626 SimplexId localEdgeId = edgeId - edgeIntervals_[nid - 1] - 1;
627 ImplicitCluster *exnode = searchCache(nid);
628 if(exnode->edgeTriangles_.empty()) {
630 }
631 return exnode->edgeTriangles_.size(localEdgeId);
632 }
633
634 inline const std::vector<std::vector<SimplexId>> *
636 edgeTriangleList_.reserve(edgeIntervals_.back() + 1);
637 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
638 std::vector<std::vector<SimplexId>> localEdgeTriangles;
639 ImplicitCluster *exnode = searchCache(nid);
640 if(exnode->edgeTriangles_.empty()) {
642 }
643 exnode->edgeTriangles_.copyTo(localEdgeTriangles);
645 localEdgeTriangles.begin(),
646 localEdgeTriangles.end());
647 }
648 return &edgeTriangleList_;
649 }
650
651 inline int getEdgeVertexInternal(const SimplexId &edgeId,
652 const int &localVertexId,
653 SimplexId &vertexId) const override {
654
655#ifndef TTK_ENABLE_KAMIKAZE
656 if((edgeId < 0) || (edgeId > (SimplexId)edgeIntervals_.back())) {
657 vertexId = -1;
658 return 0;
659 }
660 if((localVertexId != 0) && (localVertexId != 1)) {
661 vertexId = -2;
662 return 0;
663 }
664#endif
665
666 SimplexId nid = findNodeIndex(edgeId, SIMPLEX_ID::EDGE_ID);
667 SimplexId localEdgeId = edgeId - edgeIntervals_[nid - 1] - 1;
668 ImplicitCluster *exnode = searchCache(nid);
669 if(exnode->internalEdgeList_.empty()) {
670 buildInternalEdgeMap(exnode, true, false);
671 }
672
673 if(localVertexId) {
674 vertexId = exnode->internalEdgeList_.at(localEdgeId)[1];
675 } else {
676 vertexId = exnode->internalEdgeList_.at(localEdgeId)[0];
677 }
678 return 0;
679 }
680
681 inline SimplexId
683 return cellNumber_;
684 }
685
686 inline SimplexId getNumberOfEdgesInternal() const override {
687
688#ifndef TTK_ENABLE_KAMIKAZE
689 if(!edgeIntervals_.size())
690 return -1;
691#endif
692
693 return edgeIntervals_.back() + 1;
694 }
695
696 inline SimplexId getNumberOfTrianglesInternal() const override {
697
698#ifndef TTK_ENABLE_KAMIKAZE
699 if(!triangleIntervals_.size())
700 return -1;
701#endif
702
703 return triangleIntervals_.back() + 1;
704 }
705
706 inline SimplexId
708 return vertexNumber_;
709 }
710
711 inline const std::vector<std::array<SimplexId, 3>> *
713 // if it is a triangle mesh
714 if(getDimensionality() == 2) {
715 triangleList_.resize(cellNumber_, std::array<SimplexId, 3>());
716 for(SimplexId cid = 0; cid < cellNumber_; cid++) {
717 triangleList_[cid][0] = cellArray_->getCellVertex(cid, 0);
718 triangleList_[cid][1] = cellArray_->getCellVertex(cid, 1);
719 triangleList_[cid][2] = cellArray_->getCellVertex(cid, 2);
720 }
721 } else {
722 triangleList_.reserve(triangleIntervals_.back() + 1);
723 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
724 ImplicitCluster *exnode = searchCache(nid);
725 if(exnode->internalTriangleList_.empty()) {
726 buildInternalTriangleMap(exnode, true, false);
727 }
728 triangleList_.insert(triangleList_.end(),
729 exnode->internalTriangleList_.begin(),
730 exnode->internalTriangleList_.end());
731 }
732 }
733 return &triangleList_;
734 }
735
736 inline int getTriangleEdgeInternal(const SimplexId &triangleId,
737 const int &localEdgeId,
738 SimplexId &edgeId) const override {
739
740#ifndef TTK_ENABLE_KAMIKAZE
741 if(triangleId < 0 || triangleId > triangleIntervals_.back()) {
742 edgeId = -1;
743 return 0;
744 }
745 if((localEdgeId < 0) || (localEdgeId > 2)) {
746 edgeId = -2;
747 return 0;
748 }
749#endif
750
751 SimplexId nid = findNodeIndex(triangleId, SIMPLEX_ID::TRIANGLE_ID);
752 SimplexId localTriangleId = triangleId - triangleIntervals_[nid - 1] - 1;
753 ImplicitCluster *exnode = searchCache(nid);
754 if(exnode->triangleEdges_.empty()) {
756 }
757 edgeId = (exnode->triangleEdges_)[localTriangleId][localEdgeId];
758 return 0;
759 }
760
762 const SimplexId &triangleId) const override {
763#ifndef TTK_ENABLE_KAMIKAZE
764 if((triangleId < 0) || (triangleId > triangleIntervals_.back()))
765 return -1;
766#endif
767
768 (void)triangleId;
769 return 3;
770 }
771
772 inline const std::vector<std::vector<SimplexId>> *
774 if(triangleEdgeVector_.empty()) {
775 triangleEdgeVector_.reserve(triangleIntervals_.size() + 1);
776 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
777 ImplicitCluster *exnode = searchCache(nid);
778 if(exnode->triangleEdges_.empty()) {
780 }
781 for(size_t i = 0; i < exnode->triangleEdges_.size(); i++) {
782 triangleEdgeVector_.emplace_back(
783 exnode->triangleEdges_.at(i).begin(),
784 exnode->triangleEdges_.at(i).end());
785 }
786 }
787 }
788 return &triangleEdgeVector_;
789 }
790
792 const SimplexId &triangleId,
793 const int &localLinkId,
794 SimplexId &linkId) const override {
795
796#ifndef TTK_ENABLE_KAMIKAZE
797 if((triangleId < 0) || (triangleId > triangleIntervals_.back())) {
798 linkId = -1;
799 return 0;
800 }
801 if(localLinkId < 0) {
802 linkId = -2;
803 return 0;
804 }
805#endif
806
807 SimplexId nid = findNodeIndex(triangleId, SIMPLEX_ID::TRIANGLE_ID);
808 SimplexId localTriangleId = triangleId - triangleIntervals_[nid - 1] - 1;
809 ImplicitCluster *exnode = searchCache(nid);
810 if(exnode->triangleLinks_.empty()) {
812 }
813
814 if(localLinkId >= exnode->triangleLinks_.size(localTriangleId)) {
815 linkId = -2;
816 } else {
817 linkId = exnode->triangleLinks_.get(localTriangleId, localLinkId);
818 }
819 return 0;
820 }
821
823 const SimplexId &triangleId) const override {
824
825#ifndef TTK_ENABLE_KAMIKAZE
826 if((triangleId < 0) || (triangleId > triangleIntervals_.back()))
827 return -1;
828#endif
829
830 SimplexId nid = findNodeIndex(triangleId, SIMPLEX_ID::TRIANGLE_ID);
831 SimplexId localTriangleId = triangleId - triangleIntervals_[nid - 1] - 1;
832 ImplicitCluster *exnode = searchCache(nid);
833 if(exnode->triangleLinks_.empty()) {
835 }
836 return exnode->triangleLinks_.size(localTriangleId);
837 }
838
839 inline const std::vector<std::vector<SimplexId>> *
841 triangleLinkList_.reserve(triangleIntervals_.back() + 1);
842 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
843 std::vector<std::vector<SimplexId>> localTriangleLinks;
844 ImplicitCluster *exnode = searchCache(nid);
845 if(exnode->triangleLinks_.empty()) {
847 }
848 exnode->triangleLinks_.copyTo(localTriangleLinks);
850 localTriangleLinks.begin(),
851 localTriangleLinks.end());
852 }
853 return &triangleLinkList_;
854 }
855
857 const SimplexId &triangleId,
858 const int &localStarId,
859 SimplexId &starId) const override {
860
861#ifndef TTK_ENABLE_KAMIKAZE
862 if((triangleId < 0) || !triangleIntervals_.size()
863 || (triangleId > triangleIntervals_.back())) {
864 starId = -1;
865 return 0;
866 }
867 if(localStarId < 0) {
868 starId = -2;
869 return 0;
870 }
871#endif
872
873 SimplexId nid = findNodeIndex(triangleId, SIMPLEX_ID::TRIANGLE_ID);
874 SimplexId localTriangleId = triangleId - triangleIntervals_[nid - 1] - 1;
875 ImplicitCluster *exnode = searchCache(nid);
876 if(exnode->triangleStars_.empty()) {
878 }
879
880 if(localStarId >= exnode->triangleStars_.size(localTriangleId)) {
881 starId = -2;
882 } else {
883 starId = exnode->triangleStars_.get(localTriangleId, localStarId);
884 }
885 return 0;
886 }
887
889 const SimplexId &triangleId) const override {
890
891#ifndef TTK_ENABLE_KAMIKAZE
892 if((triangleId < 0) || !triangleIntervals_.size()
893 || (triangleId > triangleIntervals_.back()))
894 return -1;
895#endif
896
897 SimplexId nid = findNodeIndex(triangleId, SIMPLEX_ID::TRIANGLE_ID);
898 SimplexId localTriangleId = triangleId - triangleIntervals_[nid - 1] - 1;
899 ImplicitCluster *exnode = searchCache(nid);
900 if(exnode->triangleStars_.empty()) {
902 }
903 return exnode->triangleStars_.size(localTriangleId);
904 }
905
906 inline const std::vector<std::vector<SimplexId>> *
908 triangleStarList_.reserve(triangleIntervals_.back() + 1);
909 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
910 std::vector<std::vector<SimplexId>> localTriangleStars;
911 ImplicitCluster *exnode = searchCache(nid);
912 if(exnode->triangleStars_.empty()) {
914 }
916 localTriangleStars.begin(),
917 localTriangleStars.end());
918 }
919 return &triangleStarList_;
920 }
921
922 inline int getTriangleVertexInternal(const SimplexId &triangleId,
923 const int &localVertexId,
924 SimplexId &vertexId) const override {
925
926#ifndef TTK_ENABLE_KAMIKAZE
927 if((triangleId < 0) || (triangleId > triangleIntervals_.back())) {
928 vertexId = -1;
929 return 0;
930 }
931 if((localVertexId < 0) || (localVertexId > 2)) {
932 vertexId = -2;
933 return 0;
934 }
935#endif
936
937 SimplexId nid = findNodeIndex(triangleId, SIMPLEX_ID::TRIANGLE_ID);
938 SimplexId localTriangleId = triangleId - triangleIntervals_[nid - 1] - 1;
939 ImplicitCluster *exnode = searchCache(nid);
940 if(exnode->internalTriangleList_.empty()) {
941 buildInternalTriangleMap(exnode, true, false);
942 }
943 vertexId
944 = exnode->internalTriangleList_.at(localTriangleId)[localVertexId];
945 return 0;
946 }
947
948 inline int getVertexEdgeInternal(const SimplexId &vertexId,
949 const int &localEdgeId,
950 SimplexId &edgeId) const override {
951
952#ifndef TTK_ENABLE_KAMIKAZE
953 if((vertexId < 0) || (vertexId >= vertexNumber_)) {
954 edgeId = -1;
955 return 0;
956 }
957 if(localEdgeId < 0) {
958 edgeId = -1;
959 return 0;
960 }
961#endif
962
963 SimplexId nid = vertexIndices_[vertexId];
964 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
965 ImplicitCluster *exnode = searchCache(nid);
966 if(exnode->vertexEdges_.empty()) {
967 getClusterVertexEdges(exnode);
968 }
969 if(localEdgeId >= exnode->vertexEdges_.size(localVertexId)) {
970 edgeId = -2;
971 } else {
972 edgeId = exnode->vertexEdges_.get(localVertexId, localEdgeId);
973 }
974 return 0;
975 }
976
977 inline SimplexId
978 getVertexEdgeNumberInternal(const SimplexId &vertexId) const override {
979
980#ifndef TTK_ENABLE_KAMIKAZE
981 if((vertexId < 0) || (vertexId >= vertexNumber_))
982 return -1;
983#endif
984
985 SimplexId nid = vertexIndices_[vertexId];
986 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
987 ImplicitCluster *exnode = searchCache(nid);
988 if(exnode->vertexEdges_.empty()) {
989 getClusterVertexEdges(exnode);
990 }
991 return exnode->vertexEdges_.size(localVertexId);
992 }
993
994 inline const std::vector<std::vector<SimplexId>> *
997 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
998 std::vector<std::vector<SimplexId>> localVertexEdges;
999 ImplicitCluster *exnode = searchCache(nid);
1000 if(exnode->vertexEdges_.empty()) {
1001 getClusterVertexEdges(exnode);
1002 }
1003 exnode->vertexEdges_.copyTo(localVertexEdges);
1004 vertexEdgeList_.insert(vertexEdgeList_.end(), localVertexEdges.begin(),
1005 localVertexEdges.end());
1006 }
1007
1008 return &vertexEdgeList_;
1009 }
1010
1012 const SimplexId &vertexId,
1013 const int &localLinkId,
1014 SimplexId &linkId) const override {
1015
1016#ifndef TTK_ENABLE_KAMIKAZE
1017 if((vertexId < 0) || (vertexId > vertexIntervals_.back())) {
1018 linkId = -1;
1019 return 0;
1020 }
1021 if(localLinkId < 0) {
1022 linkId = -2;
1023 return 0;
1024 }
1025#endif
1026
1027 SimplexId nid = vertexIndices_[vertexId];
1028 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
1029 ImplicitCluster *exnode = searchCache(nid);
1030 if(exnode->vertexLinks_.empty()) {
1031 getClusterVertexLinks(exnode);
1032 }
1033 if(localLinkId >= exnode->vertexLinks_.size(localVertexId)) {
1034 linkId = -2;
1035 } else {
1036 linkId = exnode->vertexLinks_.get(localVertexId, localLinkId);
1037 }
1038 return 0;
1039 }
1040
1042 const SimplexId &vertexId) const override {
1043
1044#ifndef TTK_ENABLE_KAMIKAZE
1045 if((vertexId < 0) || (vertexId > vertexIntervals_.back()))
1046 return -1;
1047#endif
1048
1049 SimplexId nid = vertexIndices_[vertexId];
1050 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
1051 ImplicitCluster *exnode = searchCache(nid);
1052 if(exnode->vertexLinks_.empty()) {
1053 getClusterVertexLinks(exnode);
1054 }
1055 return exnode->vertexLinks_.size(localVertexId);
1056 }
1057
1058 inline const std::vector<std::vector<SimplexId>> *
1060 vertexLinkList_.reserve(vertexIntervals_.back() + 1);
1061 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
1062 std::vector<std::vector<SimplexId>> localVertexLinks;
1063 ImplicitCluster *exnode = searchCache(nid);
1064 if(exnode->vertexLinks_.empty()) {
1065 getClusterVertexLinks(exnode);
1066 }
1067 exnode->vertexLinks_.copyTo(localVertexLinks);
1068 vertexLinkList_.insert(vertexLinkList_.end(), localVertexLinks.begin(),
1069 localVertexLinks.end());
1070 }
1071 return &vertexLinkList_;
1072 }
1073
1075 const SimplexId &vertexId,
1076 const int &localNeighborId,
1077 SimplexId &neighborId) const override {
1078#ifndef TTK_ENABLE_KAMIKAZE
1079 if((vertexId < 0) || (vertexId >= vertexNumber_)) {
1080 neighborId = -1;
1081 return 0;
1082 }
1083 if(localNeighborId < 0) {
1084 neighborId = -2;
1085 return 0;
1086 }
1087#endif
1088
1089 SimplexId nid = vertexIndices_[vertexId];
1090 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
1091 ImplicitCluster *exnode = searchCache(nid);
1092 if(exnode == nullptr) {
1093 return -1;
1094 }
1095 if(exnode->vertexNeighbors_.empty()) {
1097 }
1098 if(localNeighborId >= exnode->vertexNeighbors_.size(localVertexId)) {
1099 neighborId = -2;
1100 } else {
1101 neighborId
1102 = exnode->vertexNeighbors_.get(localVertexId, localNeighborId);
1103 }
1104 return 0;
1105 }
1106
1108 const SimplexId &vertexId) const override {
1109
1110#ifndef TTK_ENABLE_KAMIKAZE
1111 if((vertexId < 0) || (vertexId >= vertexNumber_))
1112 return -1;
1113#endif
1114
1115 SimplexId nid = vertexIndices_[vertexId];
1116 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
1117 ImplicitCluster *exnode = searchCache(nid);
1118 if(exnode->vertexNeighbors_.empty()) {
1120 }
1121 return exnode->vertexNeighbors_.size(localVertexId);
1122 }
1123
1124 inline const std::vector<std::vector<SimplexId>> *
1127 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
1128 std::vector<std::vector<SimplexId>> localVertexNeighbors;
1129 ImplicitCluster *exnode = searchCache(nid);
1130 if(exnode->vertexNeighbors_.empty()) {
1132 }
1133 exnode->vertexNeighbors_.copyTo(localVertexNeighbors);
1135 localVertexNeighbors.begin(),
1136 localVertexNeighbors.end());
1137 }
1138 return &vertexNeighborList_;
1139 }
1140
1142 const SimplexId &vertexId, float &x, float &y, float &z) const override {
1143
1144#ifndef TTK_ENABLE_KAMIKAZE
1145 if((vertexId < 0) || (vertexId >= vertexNumber_))
1146 return -1;
1147#endif
1148
1149 if(doublePrecision_) {
1150 x = ((const double *)pointSet_)[3 * vertexId];
1151 y = ((const double *)pointSet_)[3 * vertexId + 1];
1152 z = ((const double *)pointSet_)[3 * vertexId + 2];
1153 } else {
1154 x = ((const float *)pointSet_)[3 * vertexId];
1155 y = ((const float *)pointSet_)[3 * vertexId + 1];
1156 z = ((const float *)pointSet_)[3 * vertexId + 2];
1157 }
1158
1159 return 0;
1160 }
1161
1163 const SimplexId &vertexId,
1164 const int &localStarId,
1165 SimplexId &starId) const override {
1166
1167#ifndef TTK_ENABLE_KAMIKAZE
1168 if((vertexId < 0) || (vertexId >= vertexNumber_)) {
1169 starId = -1;
1170 return 0;
1171 }
1172 if(localStarId < 0) {
1173 starId = -2;
1174 return 0;
1175 }
1176#endif
1177
1178 SimplexId nid = vertexIndices_[vertexId];
1179 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
1180 ImplicitCluster *exnode = searchCache(nid);
1181 if(exnode->vertexStars_.empty()) {
1182 getClusterVertexStars(exnode);
1183 }
1184 if(localStarId >= exnode->vertexStars_.size(localVertexId)) {
1185 starId = -2;
1186 } else {
1187 starId = exnode->vertexStars_.get(localVertexId, localStarId);
1188 }
1189 return 0;
1190 }
1191
1193 const SimplexId &vertexId) const override {
1194
1195#ifndef TTK_ENABLE_KAMIKAZE
1196 if((vertexId < 0) || (vertexId >= vertexNumber_))
1197 return -1;
1198#endif
1199
1200 SimplexId nid = vertexIndices_[vertexId];
1201 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
1202 ImplicitCluster *exnode = searchCache(nid);
1203 if(exnode->vertexStars_.empty()) {
1204 getClusterVertexStars(exnode);
1205 }
1206 return exnode->vertexStars_.size(localVertexId);
1207 }
1208
1209 inline const std::vector<std::vector<SimplexId>> *
1212 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
1213 std::vector<std::vector<SimplexId>> localVertexStars;
1214 ImplicitCluster *exnode = searchCache(nid);
1215 if(exnode->vertexStars_.empty()) {
1216 getClusterVertexStars(exnode);
1217 }
1218 exnode->vertexStars_.copyTo(localVertexStars);
1219 vertexStarList_.insert(vertexStarList_.end(), localVertexStars.begin(),
1220 localVertexStars.end());
1221 }
1222 return &vertexStarList_;
1223 }
1224
1225 inline int getVertexTriangleInternal(const SimplexId &vertexId,
1226 const int &localTriangleId,
1227 SimplexId &triangleId) const override {
1228
1229#ifndef TTK_ENABLE_KAMIKAZE
1230 if((vertexId < 0) || (vertexId >= vertexNumber_)) {
1231 triangleId = -1;
1232 return 0;
1233 }
1234 if(localTriangleId < 0) {
1235 triangleId = -2;
1236 return 0;
1237 }
1238#endif
1239
1240 SimplexId nid = vertexIndices_[vertexId];
1241 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
1242 ImplicitCluster *exnode = searchCache(nid);
1243 if(exnode->vertexTriangles_.empty()) {
1245 }
1246 if(localTriangleId >= exnode->vertexTriangles_.size(localVertexId)) {
1247 triangleId = -2;
1248 } else {
1249 triangleId
1250 = exnode->vertexTriangles_.get(localVertexId, localTriangleId);
1251 }
1252 return 0;
1253 }
1254
1256 const SimplexId &vertexId) const override {
1257
1258#ifndef TTK_ENABLE_KAMIKAZE
1259 if((vertexId < 0) || (vertexId >= vertexNumber_))
1260 return -1;
1261#endif
1262
1263 SimplexId nid = vertexIndices_[vertexId];
1264 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
1265 ImplicitCluster *exnode = searchCache(nid);
1266 if(exnode->vertexTriangles_.empty()) {
1268 }
1269 return exnode->vertexTriangles_.size(localVertexId);
1270 }
1271
1272 inline const std::vector<std::vector<SimplexId>> *
1275 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
1276 std::vector<std::vector<SimplexId>> localVertexTriangles;
1277 ImplicitCluster *exnode = searchCache(nid);
1278 if(exnode->vertexTriangles_.empty()) {
1280 }
1281 exnode->vertexTriangles_.copyTo(localVertexTriangles);
1283 localVertexTriangles.begin(),
1284 localVertexTriangles.end());
1285 }
1286 return &vertexTriangleList_;
1287 }
1288
1290 const SimplexId &edgeId) const override {
1291#ifndef TTK_ENABLE_KAMIKAZE
1292 if((edgeId < 0) || (edgeId > edgeIntervals_.back()))
1293 return false;
1294#endif
1295 SimplexId nid = findNodeIndex(edgeId, SIMPLEX_ID::EDGE_ID);
1296 SimplexId localedgeId = edgeId - edgeIntervals_[nid - 1] - 1;
1297 ImplicitCluster *exnode = searchCache(nid);
1298 getBoundaryCells(exnode, 1);
1299 return (exnode->boundaryEdges_)[localedgeId];
1300 }
1301
1302 bool isEmpty() const override {
1303 return !vertexNumber_;
1304 }
1305
1307 const SimplexId &triangleId) const override {
1308 if(getDimensionality() == 2)
1309 return false;
1310
1311#ifndef TTK_ENABLE_KAMIKAZE
1312 if((triangleId < 0) || (triangleId > triangleIntervals_.back()))
1313 return false;
1314#endif
1315 SimplexId nid = findNodeIndex(triangleId, SIMPLEX_ID::TRIANGLE_ID);
1316 SimplexId localtriangleId = triangleId - triangleIntervals_[nid - 1] - 1;
1317 ImplicitCluster *exnode = searchCache(nid);
1318 getBoundaryCells(exnode);
1319 return (exnode->boundaryTriangles_)[localtriangleId];
1320 }
1321
1323 const SimplexId &vertexId) const override {
1324#ifndef TTK_ENABLE_KAMIKAZE
1325 if((vertexId < 0) || (vertexId >= vertexNumber_))
1326 return false;
1327#endif
1328 SimplexId nid = vertexIndices_[vertexId];
1329 SimplexId localVertexId = vertexId - vertexIntervals_[nid - 1] - 1;
1330 ImplicitCluster *exnode = searchCache(nid);
1331 getBoundaryCells(exnode, 0);
1332 return (exnode->boundaryVertices_)[localVertexId];
1333 }
1334
1336 return 0;
1337 }
1338
1340 return 0;
1341 }
1342
1344 if(getDimensionality() == 2) {
1346 } else if(getDimensionality() == 3) {
1348 }
1349 return 0;
1350 }
1351
1352 inline int preconditionCellEdgesInternal() override {
1353 return 0;
1354 }
1355
1357 return 0;
1358 }
1359
1361 return 0;
1362 }
1363
1364 inline int preconditionEdgesInternal() override {
1365
1366#ifndef TTK_ENABLE_KAMIKAZE
1367 if(vertexNumber_ <= 0)
1368 return -1;
1369 if(cellNumber_ <= 0)
1370 return -2;
1371 if(!cellArray_)
1372 return -3;
1373 if(nodeNumber_ <= 0)
1374 return -4;
1375#endif
1376
1377 if(edgeIntervals_.empty()) {
1378 Timer t;
1379 edgeIntervals_.resize(nodeNumber_ + 1);
1380 edgeIntervals_[0] = -1;
1381 std::vector<SimplexId> edgeCount(nodeNumber_ + 1);
1382
1383#ifdef TTK_ENABLE_OPENMP
1384#pragma omp parallel for num_threads(threadNumber_)
1385#endif
1386 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
1387 edgeCount[nid] = countInternalEdges(nid);
1388 }
1389
1390 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
1391 edgeIntervals_[nid] = edgeIntervals_[nid - 1] + edgeCount[nid];
1392 }
1393
1394 this->printMsg("Edges preconditioned in "
1395 + std::to_string(t.getElapsedTime()) + " s.");
1396 }
1397
1398 return 0;
1399 }
1400
1401 inline int preconditionEdgeLinksInternal() override {
1402 return 0;
1403 }
1404
1405 inline int preconditionEdgeStarsInternal() override {
1406 return 0;
1407 }
1408
1410 return 0;
1411 }
1412
1413 inline int preconditionTrianglesInternal() override {
1414#ifndef TTK_ENABLE_KAMIKAZE
1415 if(vertexNumber_ <= 0)
1416 return -1;
1417 if(cellNumber_ <= 0)
1418 return -2;
1419 if(!cellArray_)
1420 return -3;
1421 if(nodeNumber_ <= 0)
1422 return -4;
1423#endif
1424
1425 // build triangle interval list
1426 if(triangleIntervals_.empty()) {
1427 Timer t;
1428 triangleIntervals_.resize(nodeNumber_ + 1);
1429 triangleIntervals_[0] = -1;
1430 std::vector<SimplexId> triangleCount(nodeNumber_ + 1);
1431#ifdef TTK_ENABLE_OPENMP
1432#pragma omp parallel for num_threads(threadNumber_)
1433#endif
1434 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
1435 triangleCount[nid] = countInternalTriangles(nid);
1436 }
1437
1438 for(SimplexId nid = 1; nid <= nodeNumber_; nid++) {
1440 = triangleIntervals_[nid - 1] + triangleCount[nid];
1441 }
1442
1443 this->printMsg("Triangles preconditioned in "
1444 + std::to_string(t.getElapsedTime()) + " s.");
1445 }
1446
1447 return 0;
1448 }
1449
1451 return 0;
1452 }
1453
1455 return 0;
1456 }
1457
1459 return 0;
1460 }
1461
1462 inline int preconditionVertexEdgesInternal() override {
1463 return 0;
1464 }
1465
1466 inline int preconditionVertexLinksInternal() override {
1467 if(getDimensionality() == 2) {
1469 } else if(getDimensionality() == 3) {
1471 } else {
1472 this->printErr("Unsupported dimension for vertex link precondition");
1473 return -1;
1474 }
1475 return 0;
1476 }
1477
1479 return 0;
1480 }
1481
1482 inline int preconditionVertexStarsInternal() override {
1483
1484#ifndef TTK_ENABLE_KAMIKAZE
1485 if(!cellArray_)
1486 return -1;
1487#endif
1488
1489 return 0;
1490 }
1491
1493 return 0;
1494 }
1495
1499 inline void initCache(const float ratio = 0.2) {
1500 cacheSize_ = nodeNumber_ * ratio + 1;
1501 caches_.resize(threadNumber_);
1502 cacheMaps_.resize(threadNumber_);
1503 for(int i = 0; i < threadNumber_; i++) {
1504 caches_[i].clear();
1505 cacheMaps_[i].clear();
1506 }
1507 this->printMsg("Initializing cache: " + std::to_string(cacheSize_));
1508 }
1509
1513 inline size_t getCacheSize() {
1514 return cacheSize_;
1515 }
1516
1521 inline int resetCache(int option) {
1522 for(int i = 0; i < threadNumber_; i++) {
1523 caches_[i].clear();
1524 cacheMaps_[i].clear();
1525 }
1526 if(option) {
1527 caches_.resize(1);
1528 cacheMaps_.resize(1);
1530 } else {
1531 caches_.resize(threadNumber_);
1532 cacheMaps_.resize(threadNumber_);
1534 }
1535 return 0;
1536 }
1537
1538 protected:
1539 inline int clear() {
1540 vertexIntervals_.clear();
1541 edgeIntervals_.clear();
1542 triangleIntervals_.clear();
1543 cellIntervals_.clear();
1544 externalCells_.clear();
1545 cacheMaps_.clear();
1546 cacheMaps_.clear();
1548 }
1549
1553 inline SimplexId findNodeIndex(SimplexId id, SIMPLEX_ID idType) const {
1554 const std::vector<SimplexId> *intervals = nullptr;
1555 // determine which vector to search
1556 if(idType == SIMPLEX_ID::EDGE_ID) {
1557 intervals = &edgeIntervals_;
1558 } else if(idType == SIMPLEX_ID::TRIANGLE_ID) {
1559 intervals = &triangleIntervals_;
1560 } else {
1561 return -1;
1562 }
1563
1564 std::vector<SimplexId>::const_iterator low
1565 = lower_bound(intervals->begin(), intervals->end(), id);
1566 return (low - intervals->begin());
1567 }
1568
1573 const SimplexId reservedId = 0) const {
1574 ThreadId threadId = 0;
1575#ifdef TTK_ENABLE_OPENMP
1576 threadId = omp_get_thread_num();
1577#endif
1578
1579 // cannot find the expanded node in the cache
1580 if(cacheMaps_[threadId].find(nodeId) == cacheMaps_[threadId].end()) {
1581 if(caches_[threadId].size() >= cacheSize_) {
1582 if(caches_[threadId].back().nid == reservedId) {
1583 return nullptr;
1584 }
1585 cacheMaps_[threadId].erase(caches_[threadId].back().nid);
1586 caches_[threadId].pop_back();
1587 }
1588 caches_[threadId].push_front(ImplicitCluster(nodeId));
1589 cacheMaps_[threadId][nodeId] = caches_[threadId].begin();
1590 }
1591 return &(*cacheMaps_[threadId][nodeId]);
1592 }
1593
1597 int buildInternalEdgeMap(ImplicitCluster *const nodePtr,
1598 bool computeInternalEdgeList,
1599 bool computeInternalEdgeMap) const;
1603 int buildExternalEdgeMap(ImplicitCluster *const nodePtr) const;
1604
1608 int buildInternalTriangleMap(ImplicitCluster *const nodePtr,
1609 bool computeInternalTriangleList,
1610 bool computeInternalTriangleMap) const;
1611
1615 int buildExternalTriangleMap(ImplicitCluster *const nodePtr) const;
1616
1621
1625 int countInternalTriangles(SimplexId nodeId) const;
1626
1630 int getClusterCellNeighbors(ImplicitCluster *const nodePtr) const;
1631
1635 int getClusterCellTriangles(ImplicitCluster *const nodePtr) const;
1636
1640 int getClusterEdgeLinks(ImplicitCluster *const nodePtr) const;
1641
1645 int getClusterEdgeStars(ImplicitCluster *const nodePtr) const;
1646
1650 int getClusterEdgeTriangles(ImplicitCluster *const nodePtr) const;
1651
1656 int getClusterTetraEdges(ImplicitCluster *const nodePtr) const;
1657
1661 int getClusterTriangleEdges(ImplicitCluster *const nodePtr) const;
1662
1666 int getClusterTriangleLinks(ImplicitCluster *const nodePtr) const;
1667
1671 int getClusterTriangleStars(ImplicitCluster *const nodePtr) const;
1672
1676 int getClusterVertexEdges(ImplicitCluster *const nodePtr) const;
1677
1681 int getClusterVertexLinks(ImplicitCluster *const nodePtr) const;
1682
1686 int getClusterVertexNeighbors(ImplicitCluster *const nodePtr) const;
1687
1692 int getClusterVertexStars(ImplicitCluster *const nodePtr) const;
1693
1697 int getClusterVertexTriangles(ImplicitCluster *const nodePtr) const;
1698
1702 int getBoundaryCells(ImplicitCluster *const nodePtr,
1703 const SimplexId dim = 2) const;
1704
1711 const void *pointSet_;
1712 const int *vertexIndices_;
1713 std::vector<SimplexId> vertexIntervals_;
1714 std::vector<SimplexId> edgeIntervals_;
1715 std::vector<SimplexId> triangleIntervals_;
1716 std::vector<SimplexId> cellIntervals_;
1717 std::shared_ptr<CellArray> cellArray_;
1718 std::vector<std::vector<SimplexId>> externalCells_;
1719
1720 // Cache system
1722 mutable std::vector<std::list<ImplicitCluster>> caches_;
1723 mutable std::vector<
1724 boost::unordered_map<SimplexId, std::list<ImplicitCluster>::iterator>>
1726 };
1727} // namespace ttk
#define TTK_TRIANGULATION_INTERNAL(NAME)
long long int idType
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< std::vector< SimplexId > > triangleLinkList_
std::vector< std::vector< SimplexId > > vertexNeighborList_
std::vector< std::vector< SimplexId > > vertexLinkList_
std::vector< std::vector< SimplexId > > edgeLinkList_
std::vector< std::vector< SimplexId > > vertexEdgeList_
std::vector< std::vector< SimplexId > > triangleEdgeVector_
std::vector< std::vector< SimplexId > > triangleStarList_
virtual SimplexId getCellTriangleNumber(const SimplexId &cellId) const
std::vector< std::vector< SimplexId > > cellNeighborList_
std::vector< std::array< SimplexId, 2 > > edgeList_
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::vector< SimplexId > > edgeTriangleList_
int threadNumber_
Definition: BaseClass.h:95
CompactTriangulation is a class implemented based on the TopoCluster data structure,...
int setInputCells(const SimplexId &cellNumber, const LongSimplexId *cellArray)
SimplexId findNodeIndex(SimplexId id, SIMPLEX_ID idType) const
void initCache(const float ratio=0.2)
int preconditionCellEdgesInternal() override
int preconditionVertexLinksInternal() override
int TTK_TRIANGULATION_INTERNAL() getCellVertex(const SimplexId &cellId, const int &localVertexId, SimplexId &vertexId) const override
int preconditionBoundaryTrianglesInternal() override
int preconditionEdgeLinksInternal() override
int TTK_TRIANGULATION_INTERNAL() getTriangleLink(const SimplexId &triangleId, const int &localLinkId, SimplexId &linkId) const override
const std::vector< std::vector< SimplexId > > * getCellEdgesInternal() override
SimplexId getVertexEdgeNumberInternal(const SimplexId &vertexId) const override
int getClusterEdgeStars(ImplicitCluster *const nodePtr) const
int getClusterVertexTriangles(ImplicitCluster *const nodePtr) const
int preconditionBoundaryVerticesInternal() override
SimplexId TTK_TRIANGULATION_INTERNAL() getNumberOfCells() const override
~CompactTriangulation() override
SimplexId getCellTriangleNumberInternal(const SimplexId &cellId) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getVertexLinks() override
int getTriangleVertexInternal(const SimplexId &triangleId, const int &localVertexId, SimplexId &vertexId) const override
SimplexId getCellEdgeNumberInternal(const SimplexId &cellId) const override
SimplexId getEdgeTriangleNumberInternal(const SimplexId &edgeId) const override
int getClusterVertexStars(ImplicitCluster *const nodePtr) const
int getClusterTriangleStars(ImplicitCluster *const nodePtr) const
CompactTriangulation & operator=(const CompactTriangulation &rhs)
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getCellNeighbors() override
int TTK_TRIANGULATION_INTERNAL() getVertexLink(const SimplexId &vertexId, const int &localLinkId, SimplexId &linkId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getCellNeighborNumber(const SimplexId &cellId) const override
SimplexId getVertexTriangleNumberInternal(const SimplexId &vertexId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getCellVertexNumber(const SimplexId &cellId) const override
int TTK_TRIANGULATION_INTERNAL() getEdgeLink(const SimplexId &edgeId, const int &localLinkId, SimplexId &linkId) const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getTriangleStars() override
const std::vector< std::array< SimplexId, 2 > > *TTK_TRIANGULATION_INTERNAL() getEdges() override
const std::vector< std::vector< SimplexId > > * getCellTrianglesInternal() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getVertexNeighbors() override
SimplexId TTK_TRIANGULATION_INTERNAL() getVertexLinkNumber(const SimplexId &vertexId) const override
int preconditionCellTrianglesInternal() override
std::vector< std::vector< SimplexId > > externalCells_
int TTK_TRIANGULATION_INTERNAL() getVertexNeighbor(const SimplexId &vertexId, const int &localNeighborId, SimplexId &neighborId) const override
int getEdgeVertexInternal(const SimplexId &edgeId, const int &localVertexId, SimplexId &vertexId) const override
int getVertexTriangleInternal(const SimplexId &vertexId, const int &localTriangleId, SimplexId &triangleId) const override
std::vector< SimplexId > cellIntervals_
int countInternalTriangles(SimplexId nodeId) const
int preconditionTriangleEdgesInternal() override
int getClusterTriangleEdges(ImplicitCluster *const nodePtr) const
CompactTriangulation & operator=(CompactTriangulation &&rhs)=default
int getClusterCellTriangles(ImplicitCluster *const nodePtr) const
SimplexId TTK_TRIANGULATION_INTERNAL() getNumberOfVertices() const override
int getBoundaryCells(ImplicitCluster *const nodePtr, const SimplexId dim=2) const
std::vector< boost::unordered_map< SimplexId, std::list< ImplicitCluster >::iterator > > cacheMaps_
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getVertexStars() override
int buildInternalTriangleMap(ImplicitCluster *const nodePtr, bool computeInternalTriangleList, bool computeInternalTriangleMap) const
int getClusterVertexNeighbors(ImplicitCluster *const nodePtr) const
int getClusterVertexEdges(ImplicitCluster *const nodePtr) const
const std::vector< std::vector< SimplexId > > * getTriangleEdgesInternal() override
SimplexId TTK_TRIANGULATION_INTERNAL() getEdgeLinkNumber(const SimplexId &edgeId) const override
std::vector< SimplexId > vertexIntervals_
int preconditionTrianglesInternal() override
int getCellTriangleInternal(const SimplexId &cellId, const int &localTriangleId, SimplexId &triangleId) const override
CompactTriangulation(CompactTriangulation &&rhs)=default
int getTriangleEdgeInternal(const SimplexId &triangleId, const int &localEdgeId, SimplexId &edgeId) const override
int preconditionEdgeStarsInternal() override
int TTK_TRIANGULATION_INTERNAL() getVertexStar(const SimplexId &vertexId, const int &localStarId, SimplexId &starId) const override
const std::vector< std::vector< SimplexId > > * getVertexTrianglesInternal() override
int getClusterEdgeLinks(ImplicitCluster *const nodePtr) const
int preconditionVertexStarsInternal() override
bool TTK_TRIANGULATION_INTERNAL() isTriangleOnBoundary(const SimplexId &triangleId) const override
int buildInternalEdgeMap(ImplicitCluster *const nodePtr, bool computeInternalEdgeList, bool computeInternalEdgeMap) const
int getClusterVertexLinks(ImplicitCluster *const nodePtr) const
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getEdgeLinks() override
int preconditionBoundaryEdgesInternal() override
int getVertexEdgeInternal(const SimplexId &vertexId, const int &localEdgeId, SimplexId &edgeId) const override
int preconditionCellNeighborsInternal() override
const std::vector< std::vector< SimplexId > > * getVertexEdgesInternal() override
std::vector< SimplexId > edgeIntervals_
int getEdgeTriangleInternal(const SimplexId &edgeId, const int &localTriangleId, SimplexId &triangleId) const override
int reorderVertices(std::vector< SimplexId > &vertexMap)
SimplexId TTK_TRIANGULATION_INTERNAL() getTriangleLinkNumber(const SimplexId &triangleId) const override
bool TTK_TRIANGULATION_INTERNAL() isEdgeOnBoundary(const SimplexId &edgeId) const override
int TTK_TRIANGULATION_INTERNAL() getCellNeighbor(const SimplexId &cellId, const int &localNeighborId, SimplexId &neighborId) const override
int preconditionTriangleLinksInternal() override
SimplexId TTK_TRIANGULATION_INTERNAL() getTriangleStarNumber(const SimplexId &triangleId) const override
int reorderCells(const std::vector< SimplexId > &vertexMap, const SimplexId &cellNumber, const LongSimplexId *connectivity, const LongSimplexId *offset)
SimplexId countInternalEdges(SimplexId nodeId) const
int getClusterEdgeTriangles(ImplicitCluster *const nodePtr) const
int preconditionTriangleStarsInternal() override
std::shared_ptr< CellArray > cellArray_
int preconditionVertexTrianglesInternal() override
SimplexId getTriangleEdgeNumberInternal(const SimplexId &triangleId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getVertexStarNumber(const SimplexId &vertexId) const override
std::vector< SimplexId > triangleIntervals_
int buildExternalEdgeMap(ImplicitCluster *const nodePtr) const
const std::vector< std::array< SimplexId, 3 > > *TTK_TRIANGULATION_INTERNAL() getTriangles() override
int setInputPoints(const SimplexId &pointNumber, const void *pointSet, const int *indexArray, const bool &doublePrecision=false)
int getCellEdgeInternal(const SimplexId &cellId, const int &localEdgeId, SimplexId &edgeId) const override
SimplexId TTK_TRIANGULATION_INTERNAL() getEdgeStarNumber(const SimplexId &edgeId) const override
int preconditionEdgeTrianglesInternal() override
int getClusterCellNeighbors(ImplicitCluster *const nodePtr) const
int TTK_TRIANGULATION_INTERNAL() getVertexPoint(const SimplexId &vertexId, float &x, float &y, float &z) const override
int preconditionVertexEdgesInternal() override
SimplexId getNumberOfTrianglesInternal() const override
int preconditionVertexNeighborsInternal() override
bool isEmpty() const override
SimplexId getNumberOfEdgesInternal() const override
int getClusterTriangleLinks(ImplicitCluster *const nodePtr) const
std::vector< std::list< ImplicitCluster > > caches_
ImplicitCluster * searchCache(const SimplexId &nodeId, const SimplexId reservedId=0) const
SimplexId TTK_TRIANGULATION_INTERNAL() getVertexNeighborNumber(const SimplexId &vertexId) const override
int TTK_TRIANGULATION_INTERNAL() getDimensionality() const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getTriangleLinks() override
int buildExternalTriangleMap(ImplicitCluster *const nodePtr) const
int TTK_TRIANGULATION_INTERNAL() getTriangleStar(const SimplexId &triangleId, const int &localStarId, SimplexId &starId) const override
const std::vector< std::vector< SimplexId > > * getEdgeTrianglesInternal() override
int getClusterTetraEdges(ImplicitCluster *const nodePtr) const
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL() getEdgeStars() override
bool TTK_TRIANGULATION_INTERNAL() isVertexOnBoundary(const SimplexId &vertexId) const override
int TTK_TRIANGULATION_INTERNAL() getEdgeStar(const SimplexId &edgeId, const int &localStarId, SimplexId &starId) const override
int printMsg(const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
Definition: Debug.h:118
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition: Debug.h:149
Replacement for std::vector<std::vector<SimplexId>>
void copyTo(std::vector< std::vector< SimplexId > > &dst, int threadNumber=1) const
Copy buffers to a std::vector<std::vector<SimplexId>>
SimplexId get(SimplexId id, SimplexId local) const
Returns the data inside the sub-vectors.
SimplexId size(SimplexId id) const
Get the size of a particular sub-vector.
bool empty() const
If the underlying buffers are empty.
ImplicitCluster()=default
~ImplicitCluster()=default
double getElapsedTime()
Definition: Timer.h:15
The Topology ToolKit.
long long int LongSimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:15
int ThreadId
Identifier type for threads (i.e. with OpenMP).
Definition: DataTypes.h:26
int SimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:22