TTK
Loading...
Searching...
No Matches
MorseSmaleComplex.h
Go to the documentation of this file.
1
84
85#pragma once
86
87// base code includes
88#include <DiscreteGradient.h>
90#include <Triangulation.h>
91
92#include <queue>
93
94namespace ttk {
95 class MorseSmaleComplex : public virtual Debug {
96 public:
98
101 std::vector<std::array<float, 3>> points_{};
102 std::vector<char> cellDimensions_{};
103 std::vector<SimplexId> cellIds_{};
104 std::vector<char> isOnBoundary_{};
105 std::vector<SimplexId> PLVertexIdentifiers_{};
106 std::vector<SimplexId> manifoldSize_{};
107 void clear() {
108 *this = {};
109 };
110 };
111
114 struct {
116 std::vector<float> points_{};
117 std::vector<char> smoothingMask_{};
118 std::vector<char> cellDimensions_{};
119 std::vector<ttk::SimplexId> cellIds_{};
120 } pt{}; // point data arrays
121 struct {
123 std::vector<ttk::SimplexId> connectivity_{};
124 std::vector<ttk::SimplexId> sourceIds_{};
125 std::vector<ttk::SimplexId> destinationIds_{};
126 std::vector<ttk::SimplexId> separatrixIds_{};
127 std::vector<char> separatrixTypes_{};
128 std::vector<char> isOnBoundary_{};
129 std::vector<SimplexId> sepFuncMaxId_{};
130 std::vector<SimplexId> sepFuncMinId_{};
131 } cl{}; // cell data arrays
132 void clear() {
133 *this = {};
134 };
135 };
136
139 struct {
141 std::vector<float> points_{};
142 } pt{}; // point data arrays
143 struct {
145 std::vector<ttk::SimplexId> offsets_{};
146 std::vector<ttk::SimplexId> connectivity_{};
147 std::vector<ttk::SimplexId> sourceIds_{};
148 std::vector<ttk::SimplexId> separatrixIds_{};
149 std::vector<char> separatrixTypes_{};
150 std::vector<char> isOnBoundary_{};
151 std::vector<SimplexId> sepFuncMaxId_{};
152 std::vector<SimplexId> sepFuncMinId_{};
153 } cl{}; // cell data arrays
154 void clear() {
155 *this = {};
156 };
157 };
158
165
172 template <typename dataType, typename triangulationType>
173 inline int execute(OutputCriticalPoints &outCP,
174 Output1Separatrices &outSeps1,
175 Output2Separatrices &outSeps2,
176 OutputManifold &outManifold,
177 const dataType *const scalars,
178 const size_t scalarsMTime,
179 const SimplexId *const offsets,
180 const triangulationType &triangulation,
181 const unsigned int &seed = 0);
182
187 inline void setComputeCriticalPoints(const bool state) {
188 this->ComputeCriticalPoints = state;
189 }
190
194 inline void setComputeSeparatrices1(const bool doAscending,
195 const bool doDescending,
196 const bool doSaddleConnectors) {
197 this->ComputeAscendingSeparatrices1 = doAscending;
198 this->ComputeDescendingSeparatrices1 = doDescending;
199 this->ComputeSaddleConnectors = doSaddleConnectors;
200 }
201
205 inline void setComputeSeparatrices2(const bool doAscending,
206 const bool doDescending) {
207 this->ComputeAscendingSeparatrices2 = doAscending;
208 this->ComputeDescendingSeparatrices2 = doDescending;
209 }
210
214 inline void setComputeSegmentation(const bool doAscending,
215 const bool doDescending,
216 const bool doMorseSmale) {
217 this->ComputeAscendingSegmentation = doAscending;
218 this->ComputeDescendingSegmentation = doDescending;
219 this->ComputeFinalSegmentation = doMorseSmale;
220 }
221
227 inline void setReturnSaddleConnectors(const bool state) {
229 }
230
235
237 const DiscreteGradient::BACKEND selectedBackend) {
238 DiscreteGradientBackend = selectedBackend;
239 }
240
246 inline void
250
256 this->discreteGradient_.preconditionTriangulation(data);
257 data->preconditionCellEdges();
259 }
260
261 protected:
266 struct Separatrix {
276 std::vector<dcg::Cell> geometry_;
277 };
278
283 template <typename triangulationType>
284 int
285 getDescendingSeparatrices1(const std::vector<SimplexId> &saddles,
286 std::vector<Separatrix> &separatrices,
287 const triangulationType &triangulation) const;
288
293 template <typename triangulationType>
294 int getAscendingSeparatrices1(const std::vector<SimplexId> &saddles,
295 std::vector<Separatrix> &separatrices,
296 const triangulationType &triangulation) const;
297
302 template <typename triangulationType>
303 int getSaddleConnectors(const std::vector<SimplexId> &saddles2,
304 std::vector<Separatrix> &separatrices,
305 const triangulationType &triangulation) const;
306
310 template <typename triangulationType>
312 const std::vector<Separatrix> &separatrices,
313 const SimplexId *const offsets,
314 const triangulationType &triangulation) const;
315
320 template <typename triangulationType>
322 const std::vector<SimplexId> &saddles2,
323 std::vector<Separatrix> &separatrices,
324 std::vector<std::vector<SimplexId>> &separatricesSaddles,
325 const triangulationType &triangulation) const;
326
331 template <typename triangulationType>
333 Output2Separatrices &outSeps2,
334 const std::vector<Separatrix> &separatrices,
335 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
336 const SimplexId *const offsets,
337 const triangulationType &triangulation) const;
338
344 template <typename triangulationType>
345 int getDualPolygon(const SimplexId edgeId,
346 SimplexId *const polygon,
347 const size_t polSize,
348 const triangulationType &triangulation) const;
349
353 template <typename triangulationType>
354 int sortDualPolygonVertices(SimplexId *const polygon,
355 const size_t polSize,
356 const triangulationType &triangulation) const;
357
362 template <typename triangulationType>
364 const std::vector<SimplexId> &saddles1,
365 std::vector<Separatrix> &separatrices,
366 std::vector<std::vector<SimplexId>> &separatricesSaddles,
367 const triangulationType &triangulation) const;
368
373 template <typename triangulationType>
375 Output2Separatrices &outSeps2,
376 const std::vector<Separatrix> &separatrices,
377 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
378 const SimplexId *const offsets,
379 const triangulationType &triangulation) const;
380
385 std::vector<std::vector<Separatrix>> &separatrices) const;
386
390 template <typename triangulationType>
391 int setAscendingSegmentation(const std::vector<SimplexId> &maxima,
392 SimplexId *const morseSmaleManifold,
393 const triangulationType &triangulation) const;
394
398 template <typename triangulationType>
399 int setDescendingSegmentation(const std::vector<SimplexId> &minima,
400 SimplexId *const morseSmaleManifold,
401 const triangulationType &triangulation) const;
402
407 template <typename triangulationType>
408 int setFinalSegmentation(const SimplexId numberOfMaxima,
409 const SimplexId *const ascendingManifold,
410 const SimplexId *const descendingManifold,
411 SimplexId *const morseSmaleManifold,
412 const triangulationType &triangulation) const;
413
414 template <typename dataType, typename triangulationType>
415 int returnSaddleConnectors(const double persistenceThreshold,
416 const dataType *const scalars,
417 const SimplexId *const offsets,
418 const triangulationType &triangulation);
419
420 template <typename dataType, typename triangulationType>
421 int returnSaddleConnectorsStochastic(const triangulationType &triangulation,
422 const double &persistenceThreshold,
423 const bool &detectCriticalPoints
424 = true);
425
427
437
443 };
444} // namespace ttk
445
446// ---------------- //
447// Execute method //
448// ---------------- //
449
450template <typename dataType, typename triangulationType>
452 Output1Separatrices &outSeps1,
453 Output2Separatrices &outSeps2,
454 OutputManifold &outManifold,
455 const dataType *const scalars,
456 const size_t scalarsMTime,
457 const SimplexId *const offsets,
458 const triangulationType &triangulation,
459 const unsigned int &seed) {
460#ifndef TTK_ENABLE_KAMIKAZE
461 if(scalars == nullptr) {
462 this->printErr("Input scalar field pointer is null.");
463 return -1;
464 }
465
466 if(offsets == nullptr) {
467 this->printErr("Input offset field pointer is null.");
468 return -1;
469 }
470#endif
471 Timer t;
472
473 outCP.clear();
474 outSeps1.clear();
475 outSeps2.clear();
476 const auto dim = triangulation.getDimensionality();
477
478 this->discreteGradient_.setThreadNumber(threadNumber_);
479 this->discreteGradient_.setDebugLevel(debugLevel_);
480 this->discreteGradient_.setInputScalarField(scalars, scalarsMTime);
481 this->discreteGradient_.setInputOffsets(offsets);
482 this->discreteGradient_.setBackend(this->DiscreteGradientBackend);
484 == DiscreteGradient::BACKEND::STOCHASTIC_BACKEND) {
485 this->discreteGradient_.setSeed(seed);
486 }
487 this->discreteGradient_.buildGradient<dataType, triangulationType>(
488 triangulation, this->ReturnSaddleConnectors);
489 if(this->ReturnSaddleConnectors) {
490 auto persistenceThreshold{this->SaddleConnectorsPersistenceThreshold};
491 if(!this->ThresholdIsAbsolute) {
492 const auto nVerts{triangulation.getNumberOfVertices()};
493 // global extrema are (generally) faster computed on offsets
494 // than on scalar field
495 const auto pair{std::minmax_element(offsets, offsets + nVerts)};
496 // global extrema vertex ids
497 const auto globmin = std::distance(offsets, pair.first);
498 const auto globmax = std::distance(offsets, pair.second);
499 persistenceThreshold *= (scalars[globmax] - scalars[globmin]);
500 this->printMsg("Absolute saddle connectors persistence threshold is "
501 + std::to_string(persistenceThreshold),
503 }
504
506 == DiscreteGradient::BACKEND::CLASSIC_BACKEND)
508 persistenceThreshold, scalars, offsets, triangulation);
510 == DiscreteGradient::BACKEND::STOCHASTIC_BACKEND) {
512 triangulation, persistenceThreshold, true);
513 }
514 }
515
516 std::array<std::vector<SimplexId>, 4> criticalPoints{};
517 {
518 Timer tm{};
519 discreteGradient_.getCriticalPoints(criticalPoints, triangulation);
520 this->printMsg(" Critical points extracted", 1.0, tm.getElapsedTime(),
521 this->threadNumber_, debug::LineMode::NEW,
523 }
524
525 std::vector<std::vector<Separatrix>> separatrices1{};
526
527 // 1-separatrices
528 Timer tm1sep{};
529
530 if(dim > 1 && ComputeDescendingSeparatrices1) {
531 Timer tmp;
532 separatrices1.emplace_back();
533
535 criticalPoints[1], separatrices1.back(), triangulation);
536
537 this->printMsg(" Descending 1-separatrices computed", 1.0,
538 tmp.getElapsedTime(), this->threadNumber_,
540 }
541
542 if(dim > 1 && ComputeAscendingSeparatrices1) {
543 Timer tmp;
544 separatrices1.emplace_back();
545
547 criticalPoints[dim - 1], separatrices1.back(), triangulation);
548
549 this->printMsg(" Ascending 1-separatrices computed", 1.0,
550 tmp.getElapsedTime(), this->threadNumber_,
552 }
553
554 // saddle-connectors
555 if(dim == 3 && ComputeSaddleConnectors) {
556 Timer tmp;
557 separatrices1.emplace_back();
558
559 getSaddleConnectors(criticalPoints[2], separatrices1.back(), triangulation);
560
561 this->printMsg(" Saddle connectors computed", 1.0, tmp.getElapsedTime(),
562 this->threadNumber_, debug::LineMode::NEW,
564 }
565
566 if(dim > 1
569 Timer tmp{};
570
571 this->flattenSeparatricesVectors(separatrices1);
572 setSeparatrices1(outSeps1, separatrices1[0], offsets, triangulation);
573
574 this->printMsg(" 1-separatrices set", 1.0, tmp.getElapsedTime(),
575 this->threadNumber_, debug::LineMode::NEW,
577
578 this->printMsg("1-separatrices computed", 1.0, tm1sep.getElapsedTime(),
579 this->threadNumber_);
580 }
581
582 // 2-separatrices
583 Timer tm2sep{};
584
585 if(dim == 3 && ComputeDescendingSeparatrices2) {
586 Timer tmp;
587 std::vector<Separatrix> separatrices;
588 std::vector<std::vector<SimplexId>> separatricesSaddles;
590 criticalPoints[2], separatrices, separatricesSaddles, triangulation);
592 outSeps2, separatrices, separatricesSaddles, offsets, triangulation);
593
594 this->printMsg(" Descending 2-separatrices computed", 1.0,
595 tmp.getElapsedTime(), this->threadNumber_,
597 }
598
599 if(dim == 3 && ComputeAscendingSeparatrices2) {
600 Timer tmp;
601 std::vector<Separatrix> separatrices;
602 std::vector<std::vector<SimplexId>> separatricesSaddles;
604 criticalPoints[1], separatrices, separatricesSaddles, triangulation);
606 outSeps2, separatrices, separatricesSaddles, offsets, triangulation);
607
608 this->printMsg(" Ascending 2-separatrices computed", 1.0,
609 tmp.getElapsedTime(), this->threadNumber_,
611 }
612
615 this->printMsg("2-separatrices computed", 1.0, tm2sep.getElapsedTime(),
616 this->threadNumber_);
617 }
618
620 Timer tmp;
621
624 criticalPoints[dim], outManifold.ascending_, triangulation);
625 }
628 criticalPoints[0], outManifold.descending_, triangulation);
629 }
632 setFinalSegmentation(criticalPoints[dim].size(), outManifold.ascending_,
633 outManifold.descending_, outManifold.morseSmale_,
634 triangulation);
635 }
636
637 this->printMsg(
638 "Segmentation computed", 1.0, tmp.getElapsedTime(), this->threadNumber_);
639 }
640
642 discreteGradient_.setCriticalPoints(
643 criticalPoints, outCP.points_, outCP.cellDimensions_, outCP.cellIds_,
644 outCP.isOnBoundary_, outCP.PLVertexIdentifiers_, triangulation);
645
647 discreteGradient_.setManifoldSize(criticalPoints, outManifold.ascending_,
648 outManifold.descending_,
649 outCP.manifoldSize_);
650 }
651 }
652
653 this->printMsg("Data-set ("
654 + std::to_string(triangulation.getNumberOfVertices())
655 + " points) processed",
656 1.0, t.getElapsedTime(), this->threadNumber_);
657
658 return 0;
659}
660
661// ---------------- //
662// 1-Separatrices //
663// ---------------- //
664
665template <typename triangulationType>
667 const std::vector<SimplexId> &saddles,
668 std::vector<Separatrix> &separatrices,
669 const triangulationType &triangulation) const {
670
671 const SimplexId numberOfSaddles = saddles.size();
672
673 // only 2 descending separatrices per 1-saddle
674 const SimplexId numberOfSeparatrices = 2 * numberOfSaddles;
675 separatrices.resize(numberOfSeparatrices);
676
677 // apriori: by default construction, the separatrices are not valid
678#ifdef TTK_ENABLE_OPENMP
679#pragma omp parallel for num_threads(threadNumber_)
680#endif
681 for(SimplexId i = 0; i < numberOfSaddles; ++i) {
682 const Cell saddle{1, saddles[i]};
683
684 // add descending vpaths
685 {
686 const Cell &saddle1 = saddle;
687
688 for(int j = 0; j < 2; ++j) {
689 SimplexId vertexId;
690 triangulation.getEdgeVertex(saddle1.id_, j, vertexId);
691
692 std::vector<Cell> vpath;
693 vpath.push_back(saddle1);
694 discreteGradient_.getDescendingPath(
695 Cell(0, vertexId), vpath, triangulation);
696
697 const Cell &lastCell = vpath.back();
698 if(lastCell.dim_ == 0 and discreteGradient_.isCellCritical(lastCell)) {
699 separatrices[2 * i + j].source_ = saddle;
700 separatrices[2 * i + j].destination_ = lastCell;
701 separatrices[2 * i + j].geometry_ = std::move(vpath);
702 }
703 }
704 }
705 }
706
707 return 0;
708}
709
710template <typename triangulationType>
712 const std::vector<SimplexId> &saddles,
713 std::vector<Separatrix> &separatrices,
714 const triangulationType &triangulation) const {
715
716 if(saddles.empty())
717 return 0;
718
719 const auto dim{triangulation.getDimensionality()};
720
721 // Triangulation method pointers for 3D
722 auto getFaceStarNumber = &triangulationType::getTriangleStarNumber;
723 auto getFaceStar = &triangulationType::getTriangleStar;
724 if(dim == 2) {
725 // Triangulation method pointers for 2D
726 getFaceStarNumber = &triangulationType::getEdgeStarNumber;
727 getFaceStar = &triangulationType::getEdgeStar;
728 }
729
730 const SimplexId numberOfSaddles = saddles.size();
731
732 std::vector<std::vector<Separatrix>> sepsPerSaddle(numberOfSaddles);
733
734 // apriori: by default construction, the separatrices are not valid
735#ifdef TTK_ENABLE_OPENMP
736#pragma omp parallel for num_threads(threadNumber_)
737#endif
738 for(SimplexId i = 0; i < numberOfSaddles; ++i) {
739 const Cell saddle{dim - 1, saddles[i]};
740
741 // add ascending vpaths
742 const auto starNumber{(triangulation.*getFaceStarNumber)(saddle.id_)};
743 for(SimplexId j = 0; j < starNumber; ++j) {
744
745 SimplexId sId{};
746 (triangulation.*getFaceStar)(saddle.id_, j, sId);
747
748 std::vector<Cell> vpath{saddle};
749 discreteGradient_.getAscendingPath(Cell(dim, sId), vpath, triangulation);
750
751 const Cell &lastCell = vpath.back();
752 if(lastCell.dim_ == dim and discreteGradient_.isCellCritical(lastCell)) {
753 sepsPerSaddle[i].emplace_back();
754 sepsPerSaddle[i].back().source_ = saddle;
755 sepsPerSaddle[i].back().destination_ = lastCell;
756 sepsPerSaddle[i].back().geometry_ = std::move(vpath);
757 }
758 }
759 }
760
761 this->flattenSeparatricesVectors(sepsPerSaddle);
762
763 separatrices = std::move(sepsPerSaddle[0]);
764
765 return 0;
766}
767
768template <typename triangulationType>
770 const std::vector<SimplexId> &saddles2,
771 std::vector<Separatrix> &separatrices,
772 const triangulationType &triangulation) const {
773
774 if(saddles2.empty())
775 return 0;
776
777 const auto nTriangles = triangulation.getNumberOfTriangles();
778 // visited triangles (one vector per thread)
779 std::vector<bool> isVisited(nTriangles, false);
780 std::vector<SimplexId> visitedTriangles{};
781
782 using Vpath = std::vector<Cell>;
783
784 const auto dim{triangulation.getDimensionality()};
785
786 std::vector<std::vector<Separatrix>> sepsByThread(saddles2.size());
787 std::vector<SimplexId> saddles1{};
788
789#ifdef TTK_ENABLE_OPENMP
790#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
791 firstprivate(isVisited, visitedTriangles, saddles1)
792#endif // TTK_ENABLE_OPENMP
793 for(size_t i = 0; i < saddles2.size(); ++i) {
794 const Cell s2{dim - 1, saddles2[i]};
795
796 VisitedMask mask{isVisited, visitedTriangles};
797 discreteGradient_.getDescendingWall(
798 s2, mask, triangulation, nullptr, &saddles1);
799
800 for(const auto saddle1Id : saddles1) {
801 const Cell s1{1, saddle1Id};
802
803 Vpath vpath;
804 const bool isMultiConnected
805 = discreteGradient_.getAscendingPathThroughWall(
806 s1, s2, isVisited, &vpath, triangulation);
807
808 if(vpath.empty()) {
809 // safety, should be unreachable
810 continue;
811 }
812 const auto &last = vpath.back();
813
814 if(!isMultiConnected && last.dim_ == s2.dim_ && last.id_ == s2.id_) {
815 sepsByThread[i].emplace_back();
816 sepsByThread[i].back().source_ = s1;
817 sepsByThread[i].back().destination_ = s2;
818 sepsByThread[i].back().geometry_ = std::move(vpath);
819 }
820 }
821 }
822
823 this->flattenSeparatricesVectors(sepsByThread);
824
825 separatrices = std::move(sepsByThread[0]);
826
827 return 0;
828}
829
830template <typename triangulationType>
832 Output1Separatrices &outSeps1,
833 const std::vector<Separatrix> &separatrices,
834 const SimplexId *const offsets,
835 const triangulationType &triangulation) const {
836
837 auto &separatrixFunctionMaxima = outSeps1.cl.sepFuncMaxId_;
838 auto &separatrixFunctionMinima = outSeps1.cl.sepFuncMinId_;
839
840 // max existing separatrix id + 1 or 0
841 const SimplexId separatrixId
842 = !outSeps1.cl.separatrixIds_.empty()
843 ? *std::max_element(outSeps1.cl.separatrixIds_.begin(),
844 outSeps1.cl.separatrixIds_.end())
845 + 1
846 : 0;
847
848 // total number of separatrices points
849 auto npoints{static_cast<size_t>(outSeps1.pt.numberOfPoints_)};
850 // total number of separatrices cells
851 auto ncells{static_cast<size_t>(outSeps1.cl.numberOfCells_)};
852 // points beginning id for each separatrix geometry
853 std::vector<size_t> geomPointsBegId{npoints};
854 // cells beginning id for each separatrix geometry
855 std::vector<size_t> geomCellsBegId{ncells};
856
857 // count total number of points and cells, flatten geometryId loops
858 for(size_t i = 0; i < separatrices.size(); ++i) {
859 const auto &sep = separatrices[i];
860 const auto sepSize = sep.geometry_.size();
861 npoints += sepSize;
862 ncells += sepSize - 1;
863 geomPointsBegId.emplace_back(npoints);
864 geomCellsBegId.emplace_back(ncells);
865 }
866
867 const int dimensionality = triangulation.getCellVertexNumber(0) - 1;
868
869 // resize arrays
870 outSeps1.pt.points_.resize(3 * npoints);
871 auto &points = outSeps1.pt.points_;
872 outSeps1.cl.connectivity_.resize(2 * ncells);
873 auto &cellsConn = outSeps1.cl.connectivity_;
874 outSeps1.pt.smoothingMask_.resize(npoints);
875 outSeps1.pt.cellDimensions_.resize(npoints);
876 outSeps1.pt.cellIds_.resize(npoints);
877 outSeps1.cl.sourceIds_.resize(ncells);
878 outSeps1.cl.destinationIds_.resize(ncells);
879 outSeps1.cl.separatrixIds_.resize(ncells);
880 outSeps1.cl.separatrixTypes_.resize(ncells);
881 separatrixFunctionMaxima.resize(separatrixId + separatrices.size());
882 separatrixFunctionMinima.resize(separatrixId + separatrices.size());
883 outSeps1.cl.isOnBoundary_.resize(ncells);
884
885#ifdef TTK_ENABLE_OPENMP
886#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic)
887#endif // TTK_ENABLE_OPENMP
888 for(size_t i = 0; i < separatrices.size(); ++i) {
889 const auto &sep = separatrices[i];
890 const auto &sepGeom = sep.geometry_;
891 const auto sepId = separatrixId + i;
892 // saddle (asc/desc sep) or saddle1 (saddle connector)
893 const dcg::Cell &src = sep.source_;
894 // extremum (asc/desc sep) or saddle2 (saddle connector)
895 const dcg::Cell &dst = sep.destination_;
896
897 // get separatrix type
898 const auto saddleConnector
899 = dimensionality == 3 && src.dim_ == 1 && dst.dim_ == 2;
900 const char sepType
901 = saddleConnector ? 1 : std::min(dst.dim_, dimensionality - 1);
902
903 // compute separatrix function diff
904 const auto vertsOrder = [offsets](const SimplexId a, const SimplexId b) {
905 return offsets[a] < offsets[b];
906 };
907 const std::array<SimplexId, 2> gVerts{
908 discreteGradient_.getCellGreaterVertex(src, triangulation),
909 discreteGradient_.getCellGreaterVertex(dst, triangulation)};
910 const auto sepFuncMax
911 = *std::max_element(gVerts.begin(), gVerts.end(), vertsOrder);
912 const std::array<SimplexId, 2> lVerts{
913 discreteGradient_.getCellLowerVertex(src, triangulation),
914 discreteGradient_.getCellLowerVertex(dst, triangulation)};
915 const auto sepFuncMin
916 = *std::min_element(lVerts.begin(), lVerts.end(), vertsOrder);
917 separatrixFunctionMaxima[sepId] = sepFuncMax;
918 separatrixFunctionMinima[sepId] = sepFuncMin;
919
920 // get boundary condition
921 const auto onBoundary
922 = static_cast<char>(discreteGradient_.isBoundary(src, triangulation))
923 + static_cast<char>(discreteGradient_.isBoundary(dst, triangulation));
924
925 for(size_t j = 0; j < sepGeom.size(); ++j) {
926 const auto &cell = sepGeom[j];
927 std::array<float, 3> pt{};
928 triangulation.getCellIncenter(cell.id_, cell.dim_, pt.data());
929
930 // index of current point in point data arrays
931 const auto k = geomPointsBegId[i] + j;
932
933 points[3 * k + 0] = pt[0];
934 points[3 * k + 1] = pt[1];
935 points[3 * k + 2] = pt[2];
936
937 outSeps1.pt.smoothingMask_[k]
938 = (j == 0 || j == sepGeom.size() - 1) ? 0 : 1;
939 outSeps1.pt.cellDimensions_[k] = cell.dim_;
940 outSeps1.pt.cellIds_[k] = cell.id_;
941
942 // skip filling cell data for first geometry point
943 if(j == 0)
944 continue;
945
946 // index of current cell in cell data arrays
947 const auto l = geomCellsBegId[i] + j - 1;
948
949 cellsConn[2 * l + 0] = k - 1;
950 cellsConn[2 * l + 1] = k;
951
952 outSeps1.cl.sourceIds_[l] = src.id_;
953 outSeps1.cl.destinationIds_[l] = dst.id_;
954 outSeps1.cl.separatrixIds_[l] = sepId;
955 outSeps1.cl.separatrixTypes_[l] = sepType;
956 outSeps1.cl.isOnBoundary_[l] = onBoundary;
957 }
958 }
959
960 // update pointers
961 outSeps1.pt.numberOfPoints_ = npoints;
962 outSeps1.cl.numberOfCells_ = ncells;
963
964 return 0;
965}
966
967// ---------------- //
968// 2-Separatrices //
969// ---------------- //
970
971template <typename triangulationType>
973 const std::vector<SimplexId> &saddles1,
974 std::vector<Separatrix> &separatrices,
975 std::vector<std::vector<SimplexId>> &separatricesSaddles,
976 const triangulationType &triangulation) const {
977 const Cell emptyCell;
978
979 const SimplexId numberOfSaddles = saddles1.size();
980
981 // estimation of the number of separatrices, apriori : numberOfWalls =
982 // numberOfSaddles
983 const SimplexId numberOfSeparatrices = numberOfSaddles;
984 separatrices.resize(numberOfSeparatrices);
985 separatricesSaddles.resize(numberOfSeparatrices);
986
987 const auto nEdges = triangulation.getNumberOfEdges();
988 std::vector<bool> isVisited(nEdges, false);
989 std::vector<SimplexId> visitedEdges{};
990
991 // apriori: by default construction, the separatrices are not valid
992#ifdef TTK_ENABLE_OPENMP
993#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
994 firstprivate(isVisited, visitedEdges)
995#endif // TTK_ENABLE_OPENMP
996 for(SimplexId i = 0; i < numberOfSaddles; ++i) {
997 const Cell saddle1{1, saddles1[i]};
998
999 std::vector<Cell> wall;
1000 VisitedMask mask{isVisited, visitedEdges};
1001 discreteGradient_.getAscendingWall(
1002 saddle1, mask, triangulation, &wall, &separatricesSaddles[i]);
1003
1004 separatrices[i].source_ = saddle1;
1005 separatrices[i].destination_ = emptyCell;
1006 separatrices[i].geometry_ = std::move(wall);
1007 }
1008
1009 return 0;
1010}
1011
1012template <typename triangulationType>
1014 const std::vector<SimplexId> &saddles2,
1015 std::vector<Separatrix> &separatrices,
1016 std::vector<std::vector<SimplexId>> &separatricesSaddles,
1017 const triangulationType &triangulation) const {
1018 const Cell emptyCell;
1019
1020 const SimplexId numberOfSaddles = saddles2.size();
1021
1022 // estimation of the number of separatrices, apriori : numberOfWalls =
1023 // numberOfSaddles
1024 const SimplexId numberOfSeparatrices = numberOfSaddles;
1025 separatrices.resize(numberOfSeparatrices);
1026 separatricesSaddles.resize(numberOfSeparatrices);
1027
1028 const auto nTriangles = triangulation.getNumberOfTriangles();
1029 std::vector<bool> isVisited(nTriangles, false);
1030 std::vector<SimplexId> visitedTriangles{};
1031
1032 const auto dim{triangulation.getDimensionality()};
1033
1034 // apriori: by default construction, the separatrices are not valid
1035#ifdef TTK_ENABLE_OPENMP
1036#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
1037 firstprivate(isVisited, visitedTriangles)
1038#endif // TTK_ENABLE_OPENMP
1039 for(SimplexId i = 0; i < numberOfSaddles; ++i) {
1040 const Cell saddle2{dim - 1, saddles2[i]};
1041
1042 std::vector<Cell> wall;
1043 VisitedMask mask{isVisited, visitedTriangles};
1044 discreteGradient_.getDescendingWall(
1045 saddle2, mask, triangulation, &wall, &separatricesSaddles[i]);
1046
1047 separatrices[i].source_ = saddle2;
1048 separatrices[i].destination_ = emptyCell;
1049 separatrices[i].geometry_ = std::move(wall);
1050 }
1051
1052 return 0;
1053}
1054
1055template <typename triangulationType>
1057 const SimplexId edgeId,
1058 SimplexId *const polygon,
1059 const size_t polSize,
1060 const triangulationType &triangulation) const {
1061
1062 for(size_t i = 0; i < polSize; ++i) {
1063 SimplexId starId;
1064 triangulation.getEdgeStar(edgeId, i, starId);
1065 polygon[i] = starId;
1066 }
1067
1068 return 0;
1069}
1070
1071template <typename triangulationType>
1073 SimplexId *const polygon,
1074 const size_t polSize,
1075 const triangulationType &triangulation) const {
1076
1077 for(size_t i = 1; i < polSize; ++i) {
1078
1079 // find polygon[i - 1] neighboring tetra in polygon[i..]
1080 bool isFound = false;
1081 size_t j = i;
1082 for(; j < polSize; ++j) {
1083 // check if current is the neighbor
1084 for(SimplexId k = 0;
1085 k < triangulation.getCellNeighborNumber(polygon[i - 1]); ++k) {
1086 SimplexId neighborId{};
1087 triangulation.getCellNeighbor(polygon[i - 1], k, neighborId);
1088 if(neighborId == polygon[j]) {
1089 isFound = true;
1090 break;
1091 }
1092 }
1093 if(isFound)
1094 break;
1095 }
1096
1097 // place polygon[j] next to polygon[i - 1]
1098 if(isFound) {
1099 std::swap(polygon[j], polygon[i]);
1100 }
1101 }
1102
1103 return 0;
1104}
1105
1106template <typename triangulationType>
1108 Output2Separatrices &outSeps2,
1109 const std::vector<Separatrix> &separatrices,
1110 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
1111 const SimplexId *const offsets,
1112 const triangulationType &triangulation) const {
1113
1114 auto &separatrixFunctionMaxima = outSeps2.cl.sepFuncMaxId_;
1115 auto &separatrixFunctionMinima = outSeps2.cl.sepFuncMinId_;
1116
1117 // max existing separatrix id + 1 or 0 if no previous separatrices
1118 const SimplexId separatrixId
1119 = !outSeps2.cl.separatrixIds_.empty()
1120 ? *std::max_element(outSeps2.cl.separatrixIds_.begin(),
1121 outSeps2.cl.separatrixIds_.end())
1122 + 1
1123 : 0;
1124
1125 // total number of separatrices points
1126 auto npoints{static_cast<size_t>(outSeps2.pt.numberOfPoints_)};
1127 // total number of separatrices cells
1128 auto ncells{static_cast<size_t>(outSeps2.cl.numberOfCells_)};
1129 // old number of separatrices cells
1130 const auto noldcells{ncells};
1131 // index of last vertex of last old cell + 1
1132 const auto firstCellId{outSeps2.cl.connectivity_.size()};
1133 // cells beginning id for each separatrix geometry
1134 std::vector<size_t> geomCellsBegId{ncells};
1135
1136 // count total number of cells, flatten geometryId loops
1137 for(size_t i = 0; i < separatrices.size(); ++i) {
1138 const auto &sep = separatrices[i];
1139 ncells += sep.geometry_.size();
1140 geomCellsBegId.emplace_back(ncells);
1141 }
1142
1143 // store the separatrices info (one per separatrix)
1144 std::vector<SimplexId> sepSourceIds(separatrices.size());
1145 std::vector<SimplexId> sepIds(separatrices.size());
1146 std::vector<char> sepOnBoundary(separatrices.size());
1147 separatrixFunctionMaxima.resize(separatrixId + separatrices.size());
1148 separatrixFunctionMinima.resize(separatrixId + separatrices.size());
1149 // store the polygonal cells tetras SimplexId
1150 std::vector<SimplexId> polygonNTetras(ncells - noldcells);
1151 std::vector<SimplexId> polygonEdgeIds(ncells - noldcells);
1152 std::vector<SimplexId> polygonSepInfosIds(ncells - noldcells);
1153
1154#ifdef TTK_ENABLE_OPENMP
1155#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic)
1156#endif // TTK_ENABLE_OPENMP
1157 for(size_t i = 0; i < separatrices.size(); ++i) {
1158 const auto &sep = separatrices[i];
1159 const auto &sepGeom = sep.geometry_;
1160 const auto &sepSaddles = separatricesSaddles[i];
1161 const auto sepId = separatrixId + i;
1162 const dcg::Cell &src = sep.source_; // saddle1
1163
1164 // compute separatrix function diff
1165 const auto sepFuncMin
1166 = discreteGradient_.getCellLowerVertex(src, triangulation);
1167 SimplexId sepFuncMax{};
1168 if(!sepSaddles.empty()) {
1169 // find minimum vertex on the critical triangles of the 2-separatrix
1170 const auto maxId = *std::max_element(
1171 sepSaddles.begin(), sepSaddles.end(),
1172 [&triangulation, offsets, this](const SimplexId a, const SimplexId b) {
1173 return offsets[discreteGradient_.getCellGreaterVertex(
1174 Cell{2, a}, triangulation)]
1175 < offsets[discreteGradient_.getCellGreaterVertex(
1176 Cell{2, b}, triangulation)];
1177 });
1178 sepFuncMax
1179 = discreteGradient_.getCellGreaterVertex(Cell{2, maxId}, triangulation);
1180 } else {
1181 // find maximum vertex by iterating over all the edges in the
1182 // 2-separatrix
1183 const auto maxId = *std::max_element(
1184 sepGeom.begin(), sepGeom.end(),
1185 [&triangulation, offsets, this](const Cell &a, const Cell &b) {
1186 return offsets[discreteGradient_.getCellGreaterVertex(
1187 a, triangulation)]
1188 < offsets[discreteGradient_.getCellGreaterVertex(
1189 b, triangulation)];
1190 });
1191 sepFuncMax = discreteGradient_.getCellGreaterVertex(maxId, triangulation);
1192 }
1193
1194 // get boundary condition
1195 const char onBoundary
1196 = (sepSaddles.empty()
1197 ? 0
1198 : std::count_if(sepSaddles.begin(), sepSaddles.end(),
1199 [&triangulation](const SimplexId a) {
1200 return triangulation.isTriangleOnBoundary(a);
1201 }))
1202 + triangulation.isEdgeOnBoundary(src.id_);
1203
1204 sepIds[i] = sepId;
1205 sepSourceIds[i] = src.id_;
1206 separatrixFunctionMaxima[sepId] = sepFuncMax;
1207 separatrixFunctionMinima[sepId] = sepFuncMin;
1208 sepOnBoundary[i] = onBoundary;
1209
1210 for(size_t j = 0; j < sepGeom.size(); ++j) {
1211 const auto &cell = sepGeom[j];
1212 // index of current cell in cell data arrays
1213 const auto k = geomCellsBegId[i] + j - noldcells;
1214
1215 polygonNTetras[k] = triangulation.getEdgeStarNumber(cell.id_);
1216
1217 if(polygonNTetras[k] > 2) {
1218 polygonEdgeIds[k] = cell.id_;
1219 polygonSepInfosIds[k] = i;
1220 }
1221 }
1222 }
1223
1224 // indices of valid polygon tetras
1225 std::vector<SimplexId> validTetraIds{};
1226 validTetraIds.reserve(polygonNTetras.size());
1227
1228 for(size_t i = 0; i < polygonNTetras.size(); ++i) {
1229 if(polygonNTetras[i] > 2) {
1230 validTetraIds.emplace_back(i);
1231 }
1232 }
1233
1234 // count number of valid new cells and new points
1235 size_t nnewpoints{};
1236 std::vector<SimplexId> pointsPerCell(validTetraIds.size() + 1);
1237 for(size_t i = 0; i < validTetraIds.size(); ++i) {
1238 nnewpoints += polygonNTetras[validTetraIds[i]];
1239 pointsPerCell[i + 1] = nnewpoints;
1240 }
1241
1242 // resize connectivity array
1243 outSeps2.cl.connectivity_.resize(firstCellId + nnewpoints);
1244 auto cellsConn = &outSeps2.cl.connectivity_[firstCellId];
1245 // copy of cell connectivity array (for removing duplicates vertices)
1246 std::vector<SimplexId> cellVertsIds(nnewpoints);
1247
1248#ifdef TTK_ENABLE_OPENMP
1249#pragma omp parallel for num_threads(threadNumber_)
1250#endif // TTK_ENABLE_OPENMP
1251 for(size_t i = 0; i < validTetraIds.size(); ++i) {
1252 const auto k = validTetraIds[i];
1253
1254 // get tetras in edge star
1255 getDualPolygon(polygonEdgeIds[k], &cellVertsIds[pointsPerCell[i]],
1256 polygonNTetras[k], triangulation);
1257 // sort tetras (in-place)
1258 sortDualPolygonVertices(
1259 &cellVertsIds[pointsPerCell[i]], polygonNTetras[k], triangulation);
1260
1261 for(SimplexId j = 0; j < polygonNTetras[k]; ++j) {
1262 cellsConn[pointsPerCell[i] + j] = cellVertsIds[pointsPerCell[i] + j];
1263 }
1264 }
1265
1266 TTK_PSORT(this->threadNumber_, cellVertsIds.begin(), cellVertsIds.end());
1267 const auto last = std::unique(cellVertsIds.begin(), cellVertsIds.end());
1268 cellVertsIds.erase(last, cellVertsIds.end());
1269
1270 // vertex Id to index in points array
1271 std::vector<SimplexId> vertId2PointsId(triangulation.getNumberOfCells());
1272
1273 const auto noldpoints{npoints};
1274 npoints += cellVertsIds.size();
1275 ncells = noldcells + validTetraIds.size();
1276
1277 // resize arrays
1278 outSeps2.pt.points_.resize(3 * npoints);
1279 auto points = &outSeps2.pt.points_[3 * noldpoints];
1280 outSeps2.cl.offsets_.resize(ncells + 1);
1281 outSeps2.cl.offsets_[0] = 0;
1282 auto cellsOff = &outSeps2.cl.offsets_[noldcells];
1283 outSeps2.cl.sourceIds_.resize(ncells);
1284 outSeps2.cl.separatrixIds_.resize(ncells);
1285 outSeps2.cl.separatrixTypes_.resize(ncells);
1286 outSeps2.cl.isOnBoundary_.resize(ncells);
1287
1288#ifdef TTK_ENABLE_OPENMP
1289#pragma omp parallel for num_threads(threadNumber_)
1290#endif // TTK_ENABLE_OPENMP
1291 for(size_t i = 0; i < cellVertsIds.size(); ++i) {
1292 // vertex 3D coords
1293 triangulation.getTetraIncenter(cellVertsIds[i], &points[3 * i]);
1294 // vertex index in cellVertsIds array (do not forget offset)
1295 vertId2PointsId[cellVertsIds[i]] = i + noldpoints;
1296 }
1297
1298#ifdef TTK_ENABLE_OPENMP
1299#pragma omp parallel for num_threads(threadNumber_)
1300#endif // TTK_ENABLE_OPENMP
1301 for(size_t i = 0; i < validTetraIds.size(); ++i) {
1302 const auto m = validTetraIds[i];
1303 const auto k = pointsPerCell[i];
1304 for(SimplexId j = 0; j < polygonNTetras[m]; ++j) {
1305 cellsConn[k + j] = vertId2PointsId[cellsConn[k + j]];
1306 }
1307 const auto l = i + noldcells;
1308 const auto n = polygonSepInfosIds[m];
1309 outSeps2.cl.sourceIds_[l] = sepSourceIds[n];
1310 outSeps2.cl.separatrixIds_[l] = sepIds[n];
1311 outSeps2.cl.separatrixTypes_[l] = 1;
1312 outSeps2.cl.isOnBoundary_[l] = sepOnBoundary[n];
1313 }
1314
1315 for(size_t i = 0; i < validTetraIds.size(); ++i) {
1316 // fill offsets sequentially (due to iteration dependencies)
1317 cellsOff[i + 1] = cellsOff[i] + polygonNTetras[validTetraIds[i]];
1318 }
1319
1320 outSeps2.pt.numberOfPoints_ = npoints;
1321 outSeps2.cl.numberOfCells_ = ncells;
1322
1323 return 0;
1324}
1325
1326template <typename triangulationType>
1328 Output2Separatrices &outSeps2,
1329 const std::vector<Separatrix> &separatrices,
1330 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
1331 const SimplexId *const offsets,
1332 const triangulationType &triangulation) const {
1333
1334 auto &separatrixFunctionMaxima = outSeps2.cl.sepFuncMaxId_;
1335 auto &separatrixFunctionMinima = outSeps2.cl.sepFuncMinId_;
1336
1337 // max existing separatrix id + 1 or 0 if no previous separatrices
1338 const SimplexId separatrixId
1339 = !outSeps2.cl.separatrixIds_.empty()
1340 ? *std::max_element(outSeps2.cl.separatrixIds_.begin(),
1341 outSeps2.cl.separatrixIds_.end())
1342 + 1
1343 : 0;
1344
1345 // total number of separatrices points
1346 auto npoints{static_cast<size_t>(outSeps2.pt.numberOfPoints_)};
1347 // total number of separatrices cells
1348 auto ncells{static_cast<size_t>(outSeps2.cl.numberOfCells_)};
1349 // old number of separatrices cells
1350 const auto noldcells{ncells};
1351 // index of last vertex of last old cell + 1
1352 const auto firstCellId{outSeps2.cl.connectivity_.size()};
1353
1354 // cells beginning id for each separatrix geometry
1355 std::vector<size_t> geomCellsBegId{ncells};
1356
1357 // count total number of cells, flatten geometryId loops
1358 for(size_t i = 0; i < separatrices.size(); ++i) {
1359 const auto &sep = separatrices[i];
1360 ncells += sep.geometry_.size();
1361 geomCellsBegId.emplace_back(ncells);
1362 }
1363
1364 // resize arrays
1365 outSeps2.cl.offsets_.resize(ncells + 1);
1366 outSeps2.cl.offsets_[0] = 0;
1367 outSeps2.cl.connectivity_.resize(
1368 firstCellId + 3 * (ncells - noldcells)); // triangles cells
1369 auto cellsOff = &outSeps2.cl.offsets_[noldcells];
1370 auto cellsConn = &outSeps2.cl.connectivity_[firstCellId];
1371 outSeps2.cl.sourceIds_.resize(ncells);
1372 outSeps2.cl.separatrixIds_.resize(ncells);
1373 outSeps2.cl.separatrixTypes_.resize(ncells);
1374 separatrixFunctionMaxima.resize(separatrixId + separatrices.size());
1375 separatrixFunctionMinima.resize(separatrixId + separatrices.size());
1376 outSeps2.cl.isOnBoundary_.resize(ncells);
1377
1378 // store the cells/triangles vertices vertexId
1379 std::vector<SimplexId> cellVertsIds(3 * (ncells - noldcells));
1380
1381#ifdef TTK_ENABLE_OPENMP
1382#pragma omp parallel for num_threads(threadNumber_)
1383#endif // TTK_ENABLE_OPENMP
1384 for(size_t i = 0; i < separatrices.size(); ++i) {
1385 const auto &sep = separatrices[i];
1386 const auto &sepGeom = sep.geometry_;
1387 const auto &sepSaddles = separatricesSaddles[i];
1388 const auto sepId = separatrixId + i;
1389 const dcg::Cell &src = sep.source_; // saddle2
1390 const char sepType = 2;
1391
1392 // compute separatrix function diff
1393 const auto sepFuncMax
1394 = discreteGradient_.getCellGreaterVertex(src, triangulation);
1395 SimplexId sepFuncMin{};
1396 if(!sepSaddles.empty()) {
1397 // find minimum vertex on the critical edges of the 2-separatrix
1398 const auto minId = *std::min_element(
1399 sepSaddles.begin(), sepSaddles.end(),
1400 [&triangulation, offsets, this](const SimplexId a, const SimplexId b) {
1401 return offsets[discreteGradient_.getCellLowerVertex(
1402 Cell{1, a}, triangulation)]
1403 < offsets[discreteGradient_.getCellLowerVertex(
1404 Cell{1, b}, triangulation)];
1405 });
1406 sepFuncMin
1407 = discreteGradient_.getCellLowerVertex(Cell{1, minId}, triangulation);
1408 } else {
1409 // find minimum vertex by iterating over all the triangles in the
1410 // 2-separatrix
1411 const auto minId = *std::min_element(
1412 sepGeom.begin(), sepGeom.end(),
1413 [&triangulation, offsets, this](const Cell &a, const Cell &b) {
1414 return offsets[discreteGradient_.getCellLowerVertex(a, triangulation)]
1415 < offsets[discreteGradient_.getCellLowerVertex(
1416 b, triangulation)];
1417 });
1418 sepFuncMin = discreteGradient_.getCellLowerVertex(minId, triangulation);
1419 }
1420 separatrixFunctionMaxima[sepId] = sepFuncMax;
1421 separatrixFunctionMinima[sepId] = sepFuncMin;
1422
1423 // get boundary condition
1424 const char onBoundary
1425 = (sepSaddles.empty()
1426 ? 0
1427 : std::count_if(sepSaddles.begin(), sepSaddles.end(),
1428 [&triangulation](const SimplexId a) {
1429 return triangulation.isEdgeOnBoundary(a);
1430 }))
1431 + triangulation.isTriangleOnBoundary(src.id_);
1432
1433 for(size_t j = 0; j < sepGeom.size(); ++j) {
1434 const auto &cell = sepGeom[j];
1435
1436 // first store the SimplexId of the cell/triangle vertices
1437 SimplexId v0{}, v1{}, v2{};
1438 triangulation.getTriangleVertex(cell.id_, 0, v0);
1439 triangulation.getTriangleVertex(cell.id_, 1, v1);
1440 triangulation.getTriangleVertex(cell.id_, 2, v2);
1441
1442 // index of current cell in cell data arrays
1443 const auto l = geomCellsBegId[i] + j;
1444 // index of current cell among all new cells
1445 const auto m = l - noldcells;
1446
1447 cellsConn[3 * m + 0] = v0;
1448 cellsConn[3 * m + 1] = v1;
1449 cellsConn[3 * m + 2] = v2;
1450 cellVertsIds[3 * m + 0] = v0;
1451 cellVertsIds[3 * m + 1] = v1;
1452 cellVertsIds[3 * m + 2] = v2;
1453
1454 outSeps2.cl.sourceIds_[l] = src.id_;
1455 outSeps2.cl.separatrixIds_[l] = sepId;
1456 outSeps2.cl.separatrixTypes_[l] = sepType;
1457 outSeps2.cl.isOnBoundary_[l] = onBoundary;
1458 }
1459 }
1460
1461 // reduce the cell vertices ids
1462 // (cells are triangles sharing two vertices)
1463 TTK_PSORT(this->threadNumber_, cellVertsIds.begin(), cellVertsIds.end());
1464 const auto last = std::unique(cellVertsIds.begin(), cellVertsIds.end());
1465 cellVertsIds.erase(last, cellVertsIds.end());
1466
1467 // vertex Id to index in points array
1468 std::vector<size_t> vertId2PointsId(triangulation.getNumberOfVertices());
1469
1470 const auto noldpoints{npoints};
1471 npoints += cellVertsIds.size();
1472 outSeps2.pt.points_.resize(3 * npoints);
1473 auto points = &outSeps2.pt.points_[3 * noldpoints];
1474
1475#ifdef TTK_ENABLE_OPENMP
1476#pragma omp parallel for num_threads(threadNumber_)
1477#endif // TTK_ENABLE_OPENMP
1478 for(size_t i = 0; i < cellVertsIds.size(); ++i) {
1479 // vertex 3D coords
1480 triangulation.getVertexPoint(
1481 cellVertsIds[i], points[3 * i + 0], points[3 * i + 1], points[3 * i + 2]);
1482 // vertex index in cellVertsIds array (do not forget offset)
1483 vertId2PointsId[cellVertsIds[i]] = i + noldpoints;
1484 }
1485
1486 const auto lastOffset = noldcells == 0 ? 0 : cellsOff[-1];
1487
1488#ifdef TTK_ENABLE_OPENMP
1489#pragma omp parallel for num_threads(threadNumber_)
1490#endif // TTK_ENABLE_OPENMP
1491 for(size_t i = 0; i < ncells - noldcells; ++i) {
1492 cellsOff[i] = 3 * i + lastOffset;
1493 cellsConn[3 * i + 0] = vertId2PointsId[cellsConn[3 * i + 0]];
1494 cellsConn[3 * i + 1] = vertId2PointsId[cellsConn[3 * i + 1]];
1495 cellsConn[3 * i + 2] = vertId2PointsId[cellsConn[3 * i + 2]];
1496 }
1497
1498 cellsOff[ncells - noldcells] = cellsOff[ncells - noldcells - 1] + 3;
1499
1500 outSeps2.pt.numberOfPoints_ = npoints;
1501 outSeps2.cl.numberOfCells_ = ncells;
1502
1503 return 0;
1504}
1505
1506// ---------------- //
1507// Segmentation //
1508// ---------------- //
1509
1510template <typename triangulationType>
1512 const std::vector<SimplexId> &maxima,
1513 SimplexId *const morseSmaleManifold,
1514 const triangulationType &triangulation) const {
1515
1516 if(morseSmaleManifold == nullptr) {
1517 this->printErr("Could not compute ascending segmentation");
1518 return 1;
1519 }
1520
1521 Timer tm{};
1522
1523 const auto nVerts{triangulation.getNumberOfVertices()};
1524 std::fill(morseSmaleManifold, morseSmaleManifold + nVerts, -1);
1525 if(maxima.empty()) {
1526 // shortcut for elevation
1527 return 0;
1528 }
1529
1530 const auto nCells{triangulation.getNumberOfCells()};
1531 std::vector<SimplexId> morseSmaleManifoldOnCells(nCells, -1);
1532
1533 size_t nMax{};
1534 for(const auto &id : maxima) {
1535 // mark the maxima
1536 morseSmaleManifoldOnCells[id] = nMax++;
1537 }
1538
1539 const auto dim{triangulation.getDimensionality()};
1540
1541 // Triangulation method pointers for 3D
1542 auto getFaceStarNumber = &triangulationType::getTriangleStarNumber;
1543 auto getFaceStar = &triangulationType::getTriangleStar;
1544 if(dim == 2) {
1545 // Triangulation method pointers for 2D
1546 getFaceStarNumber = &triangulationType::getEdgeStarNumber;
1547 getFaceStar = &triangulationType::getEdgeStar;
1548 } else if(dim == 1) {
1549 // Triangulation method pointers for 1D
1550 getFaceStarNumber = &triangulationType::getVertexStarNumber;
1551 getFaceStar = &triangulationType::getVertexStar;
1552 }
1553
1554 // cells visited during the propagation alongside one integral line
1555 std::vector<SimplexId> visited{};
1556 // all marked cells
1557 std::vector<uint8_t> isMarked(nCells, 0);
1558
1559#ifdef TTK_ENABLE_OPENMP
1560#pragma omp parallel for num_threads(threadNumber_) firstprivate(visited)
1561#endif // TTK_ENABLE_OPENMP
1562 for(SimplexId i = 0; i < nCells; ++i) {
1563 if(isMarked[i] == 1) {
1564 continue;
1565 }
1566 visited.clear();
1567 auto curr{i};
1568 while(morseSmaleManifoldOnCells[curr] == -1) {
1569 if(isMarked[curr] == 1) {
1570 break;
1571 }
1572 // follow a V-path till an already marked cell is reached
1573 const auto paired{this->discreteGradient_.getPairedCell(
1574 Cell{dim, curr}, triangulation, true)};
1575 SimplexId next{curr};
1576 const auto nStars{(triangulation.*getFaceStarNumber)(paired)};
1577 for(SimplexId j = 0; j < nStars; ++j) {
1578 (triangulation.*getFaceStar)(paired, j, next);
1579 // get the first cell != curr (what of non-manifold datasets?)
1580 if(next != curr) {
1581 break;
1582 }
1583 }
1584 visited.emplace_back(curr);
1585 if(next == curr) {
1586 // on the boundary?
1587 break;
1588 }
1589 curr = next;
1590 }
1591 for(const auto el : visited) {
1592 morseSmaleManifoldOnCells[el] = morseSmaleManifoldOnCells[curr];
1593 isMarked[el] = 1;
1594 }
1595 }
1596
1597#ifdef TTK_ENABLE_OPENMP
1598#pragma omp parallel for num_threads(threadNumber_)
1599#endif // TTK_ENABLE_OPENMP
1600 for(SimplexId i = 0; i < nVerts; ++i) {
1601 if(triangulation.getVertexStarNumber(i) < 1) {
1602 // handle non-manifold datasets?
1603 continue;
1604 }
1605 SimplexId starId{};
1606 triangulation.getVertexStar(i, 0, starId);
1607 // put segmentation infos from cells to points
1608 morseSmaleManifold[i] = morseSmaleManifoldOnCells[starId];
1609 }
1610
1611 this->printMsg(" Ascending segmentation computed", 1.0, tm.getElapsedTime(),
1612 this->threadNumber_, debug::LineMode::NEW,
1614
1615 return 0;
1616}
1617
1618template <typename triangulationType>
1620 const std::vector<SimplexId> &minima,
1621 SimplexId *const morseSmaleManifold,
1622 const triangulationType &triangulation) const {
1623
1624 if(morseSmaleManifold == nullptr) {
1625 this->printErr("Could not compute descending segmentation");
1626 return 1;
1627 }
1628
1629 Timer tm{};
1630
1631 const auto nVerts{triangulation.getNumberOfVertices()};
1632
1633 if(minima.size() == 1) {
1634 // shortcut for elevation
1635 std::fill(morseSmaleManifold, morseSmaleManifold + nVerts, 0);
1636 return 0;
1637 }
1638
1639 std::fill(morseSmaleManifold, morseSmaleManifold + nVerts, -1);
1640
1641 size_t nMin{};
1642 for(const auto &cp : minima) {
1643 // mark the minima
1644 morseSmaleManifold[cp] = nMin++;
1645 }
1646
1647 std::vector<SimplexId> visited{};
1648
1649#ifdef TTK_ENABLE_OPENMP
1650#pragma omp parallel for num_threads(threadNumber_) firstprivate(visited)
1651#endif // TTK_ENABLE_OPENMP
1652 for(SimplexId i = 0; i < nVerts; ++i) {
1653 if(morseSmaleManifold[i] != -1) {
1654 continue;
1655 }
1656 visited.clear();
1657 auto curr{i};
1658 while(morseSmaleManifold[curr] == -1) {
1659 // follow a V-path till an already marked vertex is reached
1660 const auto pairedEdge{
1661 this->discreteGradient_.getPairedCell(Cell{0, curr}, triangulation)};
1662 SimplexId next{};
1663 triangulation.getEdgeVertex(pairedEdge, 0, next);
1664 if(next == curr) {
1665 triangulation.getEdgeVertex(pairedEdge, 1, next);
1666 }
1667 visited.emplace_back(curr);
1668 curr = next;
1669 }
1670 for(const auto el : visited) {
1671 morseSmaleManifold[el] = morseSmaleManifold[curr];
1672 }
1673 }
1674
1675 this->printMsg(" Descending segmentation computed", 1.0, tm.getElapsedTime(),
1676 this->threadNumber_, debug::LineMode::NEW,
1678
1679 return 0;
1680}
1681
1682template <typename triangulationType>
1684 const SimplexId numberOfMaxima,
1685 const SimplexId *const ascendingManifold,
1686 const SimplexId *const descendingManifold,
1687 SimplexId *const morseSmaleManifold,
1688 const triangulationType &triangulation) const {
1689
1690 if(ascendingManifold == nullptr || descendingManifold == nullptr
1691 || morseSmaleManifold == nullptr) {
1692 this->printErr("Could not compute final segmentation");
1693 return 1;
1694 }
1695
1696 Timer tm{};
1697
1698 const size_t nVerts = triangulation.getNumberOfVertices();
1699
1700 // associate a unique "sparse region id" to each (ascending, descending) pair
1701
1702#ifdef TTK_ENABLE_OPENMP
1703#pragma omp parallel for num_threads(threadNumber_)
1704#endif // TTK_ENABLE_OPENMP
1705 for(size_t i = 0; i < nVerts; ++i) {
1706 const auto d = ascendingManifold[i];
1707 const auto a = descendingManifold[i];
1708 if(a == -1 || d == -1) {
1709 morseSmaleManifold[i] = -1;
1710 } else {
1711 morseSmaleManifold[i] = a * numberOfMaxima + d;
1712 }
1713 }
1714
1715 // store the "sparse region ids" by copying the morseSmaleManifold output
1716 std::vector<SimplexId> sparseRegionIds(
1717 morseSmaleManifold, morseSmaleManifold + nVerts);
1718
1719 // get unique "sparse region ids"
1720 TTK_PSORT(
1721 this->threadNumber_, sparseRegionIds.begin(), sparseRegionIds.end());
1722 const auto last = std::unique(sparseRegionIds.begin(), sparseRegionIds.end());
1723 sparseRegionIds.erase(last, sparseRegionIds.end());
1724
1725 // "sparse region id" -> "dense region id"
1726 std::map<SimplexId, size_t> sparseToDenseRegionId{};
1727
1728 for(size_t i = 0; i < sparseRegionIds.size(); ++i) {
1729 sparseToDenseRegionId[sparseRegionIds[i]] = i;
1730 }
1731
1732 // update region id on all vertices: "sparse id" -> "dense id"
1733
1734#ifdef TTK_ENABLE_OPENMP
1735#pragma omp parallel for num_threads(threadNumber_)
1736#endif // TTK_ENABLE_OPENMP
1737 for(size_t i = 0; i < nVerts; ++i) {
1738 morseSmaleManifold[i] = sparseToDenseRegionId[morseSmaleManifold[i]];
1739 }
1740
1741 this->printMsg(" Final segmentation computed", 1.0, tm.getElapsedTime(),
1742 this->threadNumber_, debug::LineMode::NEW,
1744
1745 return 0;
1746}
1747template <typename dataType, typename triangulationType>
1749 const double persistenceThreshold,
1750 const dataType *const scalars,
1751 const SimplexId *const offsets,
1752 const triangulationType &triangulation) {
1753
1754 Timer tm{};
1755
1756 const auto dim{triangulation.getDimensionality()};
1757 if(dim != 3) {
1758 this->printWrn("Can't return saddle connectors without a 3D dataset");
1759 return 0;
1760 }
1761
1762 // compute saddle-saddle Persistence Pairs from DiscreteMorseSandwich
1764 dms.setThreadNumber(this->threadNumber_);
1765 dms.setDebugLevel(this->debugLevel_);
1766 dms.setGradient(std::move(this->discreteGradient_));
1767
1768 using PersPairType = DiscreteMorseSandwich::PersistencePair;
1769
1770 std::vector<PersPairType> dms_pairs{};
1771 dms.computePersistencePairs(dms_pairs, offsets, triangulation, false, true);
1772 this->discreteGradient_ = dms.getGradient();
1773 // reset gradient pointer to local storage
1774 this->discreteGradient_.setLocalGradient();
1775
1776 const auto getPersistence
1777 = [this, &triangulation, scalars](const PersPairType &p) {
1778 return this->discreteGradient_.getPersistence(
1779 Cell{2, p.death}, Cell{1, p.birth}, scalars, triangulation);
1780 };
1781
1782 // saddle-saddle pairs should be in one continuous block inside dms_pairs
1783 const auto firstSadSadPair{std::distance(
1784 dms_pairs.begin(),
1785 std::find_if(dms_pairs.begin(), dms_pairs.end(),
1786 [](const auto &pair) { return pair.type == 1; }))};
1787
1788 // DMS pairs output are already sorted by 2-saddles in ascending order
1789
1790 std::vector<bool> isVisited(triangulation.getNumberOfTriangles(), false);
1791 std::vector<SimplexId> visitedTriangles{};
1792
1793 size_t nReturned{};
1794
1795 // Sort pairs to process by persistence
1796 std::vector<std::tuple<size_t, dataType>> pairs;
1797 for(size_t i = firstSadSadPair; i < dms_pairs.size(); ++i) {
1798 const auto &pair{dms_pairs[i]};
1799 pairs.emplace_back(std::make_tuple(i, getPersistence(pair)));
1800 }
1801 const auto comparePersistence
1802 = [](const std::tuple<size_t, dataType> &pair1,
1803 const std::tuple<size_t, dataType> &pair2) {
1804 return std::get<1>(pair1) < std::get<1>(pair2);
1805 };
1806 TTK_PSORT(
1807 this->threadNumber_, pairs.begin(), pairs.end(), comparePersistence);
1808
1809 std::vector<std::tuple<dataType, SimplexId, SimplexId>> skippedPairsPers;
1810
1811 // Process pairs
1812 for(const auto &pairTup : pairs) {
1813 const auto &pairIndex = std::get<0>(pairTup);
1814 const auto &pair{dms_pairs[pairIndex]};
1815 const auto &pairPersistence = std::get<1>(pairTup);
1816
1817 if(pair.type != 1 || pairPersistence > persistenceThreshold) {
1818 continue;
1819 }
1820
1821 const Cell birth{1, pair.birth};
1822 const Cell death{2, pair.death};
1823
1824 // 1. get the 2-saddle wall
1825 VisitedMask mask{isVisited, visitedTriangles};
1826 this->discreteGradient_.getDescendingWall(death, mask, triangulation);
1827 // 2. get the saddle connector
1828 std::vector<Cell> vpath{};
1829 bool disableForkReversal = not ForceLoopFreeGradient;
1830 this->discreteGradient_.getAscendingPathThroughWall(
1831 birth, death, isVisited, &vpath, triangulation, disableForkReversal);
1832 // 3. reverse the gradient on the saddle connector path
1833 if(vpath.back() == death) {
1834 this->discreteGradient_.reverseAscendingPathOnWall(vpath, triangulation);
1835
1836 // 3.1 Detect cycle
1837 bool cycle = false;
1839 cycle
1840 = this->discreteGradient_.detectGradientCycle(birth, triangulation);
1841 } else {
1842 std::vector<bool> isVisitedUpdated(
1843 triangulation.getNumberOfTriangles(), false);
1844 std::vector<SimplexId> visitedTrianglesUpdated{};
1845 VisitedMask maskUpdated{isVisitedUpdated, visitedTrianglesUpdated};
1846 discreteGradient_.getDescendingWall(death, maskUpdated, triangulation);
1847 discreteGradient_.getAscendingPathThroughWall(
1848 birth, death, isVisitedUpdated, nullptr, triangulation, false, true,
1849 &cycle);
1850 }
1851 if(cycle) {
1852 this->discreteGradient_.reverseAscendingPathOnWall(
1853 vpath, triangulation, true);
1854 this->printMsg("Could not return saddle connector " + birth.to_string()
1855 + " -> " + death.to_string()
1856 + " without creating a cycle.",
1858 if(this->debugLevel_ == (int)debug::Priority::DETAIL)
1859 skippedPairsPers.emplace_back(
1860 std::make_tuple(pairPersistence, pair.birth, pair.death));
1861 } else
1862 nReturned++;
1863
1864 } else {
1865 this->printMsg("Could not return saddle connector " + birth.to_string()
1866 + " -> " + death.to_string(),
1868 if(this->debugLevel_ == (int)debug::Priority::DETAIL)
1869 skippedPairsPers.emplace_back(
1870 std::make_tuple(pairPersistence, pair.birth, pair.death));
1871 }
1872 }
1873
1874 if(this->debugLevel_ == (int)debug::Priority::DETAIL) {
1875 TTK_PSORT(
1876 this->threadNumber_, skippedPairsPers.begin(), skippedPairsPers.end());
1877 for(unsigned int i = 0; i < skippedPairsPers.size(); ++i)
1878 this->printMsg(std::to_string(i) + " "
1879 + std::to_string(std::get<1>(skippedPairsPers[i])) + " "
1880 + std::to_string(std::get<2>(skippedPairsPers[i])) + " "
1881 + std::to_string(std::get<0>(skippedPairsPers[i])));
1882 }
1883
1884 this->printMsg("Returned " + std::to_string(nReturned) + " saddle connectors",
1885 1.0, tm.getElapsedTime(), this->threadNumber_);
1886
1887 return 0;
1888}
1889
1890template <typename dataType, typename triangulationType>
1892 const triangulationType &triangulation,
1893 const double &persistenceThreshold,
1894 const bool &detectCriticalPoints) {
1895
1896 std::vector<std::pair<SimplexId, char>> criticalPoints{};
1897 this->discreteGradient_.setSaddleConnectorsPersistenceThreshold(
1898 persistenceThreshold);
1899
1900 if(detectCriticalPoints) {
1901
1902 std::vector<Cell> criticalCells{};
1903 this->discreteGradient_.getCriticalPoints(criticalCells, triangulation);
1904
1905 criticalPoints.resize(criticalCells.size());
1906
1907 // iterate over cells to get points (max vertex) and type
1908 for(size_t i = 0; i < criticalCells.size(); ++i) {
1909 const auto &c = criticalCells[i];
1910 criticalPoints[i]
1911 = {this->discreteGradient_.getCellGreaterVertex(c, triangulation),
1912 static_cast<char>(
1913 this->discreteGradient_.criticalTypeFromCellDimension(c.dim_))};
1914 }
1915
1916 // print number of critical cells
1917 {
1918 // foreach dimension
1919 const int numberOfDimensions
1920 = this->discreteGradient_.getNumberOfDimensions();
1921 std::vector<SimplexId> nDMTCriticalPoints(numberOfDimensions, 0);
1922 for(const auto &c : criticalCells) {
1923 ++nDMTCriticalPoints[c.dim_];
1924 }
1925
1926 std::vector<SimplexId> nPLInteriorCriticalPoints(numberOfDimensions, 0);
1927 for(const auto &cp : criticalPoints) {
1928 if(!triangulation.isVertexOnBoundary(cp.first)) {
1929 ++nPLInteriorCriticalPoints[cp.second];
1930 }
1931 }
1932
1933 std::vector<std::vector<std::string>> rows(numberOfDimensions);
1934 for(int i = 0; i < numberOfDimensions; ++i) {
1935 rows[i] = std::vector<std::string>{
1936 "#" + std::to_string(i) + "-cell(s)",
1937 std::to_string(nDMTCriticalPoints[i]) + " (with "
1938 + std::to_string(nPLInteriorCriticalPoints[i]) + " interior PL)"};
1939 }
1940 this->printMsg(rows);
1941 }
1942 }
1943
1944 Timer t;
1945
1946 const bool allowBoundary = true;
1947
1948 std::vector<char> isPL;
1949 this->discreteGradient_.getCriticalPointMap(criticalPoints, isPL);
1950
1951 if(triangulation.getDimensionality() == 3) {
1952 this->discreteGradient_.filterSaddleConnectors<dataType>(
1953 allowBoundary, triangulation);
1954 }
1955
1956 this->printMsg(
1957 "Gradient reversed", 1.0, t.getElapsedTime(), this->threadNumber_);
1958
1959 return 0;
1960}
#define TTK_PSORT(NTHREADS,...)
Parallel sort macro.
Definition OpenMP.h:46
AbstractTriangulation is an interface class that defines an interface for efficient traversal methods...
virtual int setThreadNumber(const int threadNumber)
Definition BaseClass.h:80
int debugLevel_
Definition Debug.h:379
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
TTK DiscreteMorseSandwich processing package.
int execute(OutputCriticalPoints &outCP, Output1Separatrices &outSeps1, Output2Separatrices &outSeps2, OutputManifold &outManifold, const dataType *const scalars, const size_t scalarsMTime, const SimplexId *const offsets, const triangulationType &triangulation, const unsigned int &seed=0)
int setAscendingSegmentation(const std::vector< SimplexId > &maxima, SimplexId *const morseSmaleManifold, const triangulationType &triangulation) const
int setFinalSegmentation(const SimplexId numberOfMaxima, const SimplexId *const ascendingManifold, const SimplexId *const descendingManifold, SimplexId *const morseSmaleManifold, const triangulationType &triangulation) const
void flattenSeparatricesVectors(std::vector< std::vector< Separatrix > > &separatrices) const
Flatten the vectors of vectors into their first component.
int getDualPolygon(const SimplexId edgeId, SimplexId *const polygon, const size_t polSize, const triangulationType &triangulation) const
void setComputeCriticalPoints(const bool state)
int returnSaddleConnectorsStochastic(const triangulationType &triangulation, const double &persistenceThreshold, const bool &detectCriticalPoints=true)
int getAscendingSeparatrices2(const std::vector< SimplexId > &saddles1, std::vector< Separatrix > &separatrices, std::vector< std::vector< SimplexId > > &separatricesSaddles, const triangulationType &triangulation) const
int getSaddleConnectors(const std::vector< SimplexId > &saddles2, std::vector< Separatrix > &separatrices, const triangulationType &triangulation) const
int getDescendingSeparatrices1(const std::vector< SimplexId > &saddles, std::vector< Separatrix > &separatrices, const triangulationType &triangulation) const
dcg::DiscreteGradient discreteGradient_
int returnSaddleConnectors(const double persistenceThreshold, const dataType *const scalars, const SimplexId *const offsets, const triangulationType &triangulation)
int setAscendingSeparatrices2(Output2Separatrices &outSeps2, const std::vector< Separatrix > &separatrices, const std::vector< std::vector< SimplexId > > &separatricesSaddles, const SimplexId *const offsets, const triangulationType &triangulation) const
int setDescendingSegmentation(const std::vector< SimplexId > &minima, SimplexId *const morseSmaleManifold, const triangulationType &triangulation) const
int setSeparatrices1(Output1Separatrices &outSeps1, const std::vector< Separatrix > &separatrices, const SimplexId *const offsets, const triangulationType &triangulation) const
void setComputeSegmentation(const bool doAscending, const bool doDescending, const bool doMorseSmale)
void setSaddleConnectorsPersistenceThreshold(const double threshold)
int sortDualPolygonVertices(SimplexId *const polygon, const size_t polSize, const triangulationType &triangulation) const
void setComputeSeparatrices2(const bool doAscending, const bool doDescending)
void preconditionTriangulation(AbstractTriangulation *const data)
int setDescendingSeparatrices2(Output2Separatrices &outSeps2, const std::vector< Separatrix > &separatrices, const std::vector< std::vector< SimplexId > > &separatricesSaddles, const SimplexId *const offsets, const triangulationType &triangulation) const
int getDescendingSeparatrices2(const std::vector< SimplexId > &saddles2, std::vector< Separatrix > &separatrices, std::vector< std::vector< SimplexId > > &separatricesSaddles, const triangulationType &triangulation) const
void setComputeSeparatrices1(const bool doAscending, const bool doDescending, const bool doSaddleConnectors)
void setReturnSaddleConnectors(const bool state)
void setDiscreteGradientBackend(const DiscreteGradient::BACKEND selectedBackend)
int getAscendingSeparatrices1(const std::vector< SimplexId > &saddles, std::vector< Separatrix > &separatrices, const triangulationType &triangulation) const
DiscreteGradient::BACKEND DiscreteGradientBackend
double getElapsedTime()
Definition Timer.h:15
TTK discreteGradient processing package.
TTK base package defining the standard types.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
Persistence pair struct as exported by DiscreteGradient.
1-Separatrices point and cell data arrays
std::vector< ttk::SimplexId > destinationIds_
std::vector< ttk::SimplexId > separatrixIds_
struct ttk::MorseSmaleComplex::Output1Separatrices::@346010247245230204106122273277326155246155345264 cl
std::vector< ttk::SimplexId > connectivity_
struct ttk::MorseSmaleComplex::Output1Separatrices::@336365226070351077156040303020054073121141066106 pt
2-Separatrices point and cell data arrays
struct ttk::MorseSmaleComplex::Output2Separatrices::@005131327202314304030346343376373171222130207232 cl
struct ttk::MorseSmaleComplex::Output2Separatrices::@103373205041225224020060067072033253227034337377 pt
std::vector< ttk::SimplexId > connectivity_
std::vector< ttk::SimplexId > separatrixIds_
std::vector< std::array< float, 3 > > points_
Pointers to pre-allocated segmentation point data arrays.
Auto-cleaning re-usable graph propagations data structure.
Definition VisitedMask.h:27
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)