TTK
Loading...
Searching...
No Matches
MorseSmaleComplex.h
Go to the documentation of this file.
1
76
77#pragma once
78
79// base code includes
80#include <DiscreteGradient.h>
82#include <Triangulation.h>
83
84#include <queue>
85
86namespace ttk {
87 class MorseSmaleComplex : public virtual Debug {
88 public:
90
93 std::vector<std::array<float, 3>> points_{};
94 std::vector<char> cellDimensions_{};
95 std::vector<SimplexId> cellIds_{};
96 std::vector<char> isOnBoundary_{};
97 std::vector<SimplexId> PLVertexIdentifiers_{};
98 std::vector<SimplexId> manifoldSize_{};
99 void clear() {
100 *this = {};
101 };
102 };
103
106 struct {
108 std::vector<float> points_{};
109 std::vector<char> smoothingMask_{};
110 std::vector<char> cellDimensions_{};
111 std::vector<ttk::SimplexId> cellIds_{};
112 } pt{}; // point data arrays
113 struct {
115 std::vector<ttk::SimplexId> connectivity_{};
116 std::vector<ttk::SimplexId> sourceIds_{};
117 std::vector<ttk::SimplexId> destinationIds_{};
118 std::vector<ttk::SimplexId> separatrixIds_{};
119 std::vector<char> separatrixTypes_{};
120 std::vector<char> isOnBoundary_{};
121 std::vector<SimplexId> sepFuncMaxId_{};
122 std::vector<SimplexId> sepFuncMinId_{};
123 } cl{}; // cell data arrays
124 void clear() {
125 *this = {};
126 };
127 };
128
131 struct {
133 std::vector<float> points_{};
134 } pt{}; // point data arrays
135 struct {
137 std::vector<ttk::SimplexId> offsets_{};
138 std::vector<ttk::SimplexId> connectivity_{};
139 std::vector<ttk::SimplexId> sourceIds_{};
140 std::vector<ttk::SimplexId> separatrixIds_{};
141 std::vector<char> separatrixTypes_{};
142 std::vector<char> isOnBoundary_{};
143 std::vector<SimplexId> sepFuncMaxId_{};
144 std::vector<SimplexId> sepFuncMinId_{};
145 } cl{}; // cell data arrays
146 void clear() {
147 *this = {};
148 };
149 };
150
156 };
157
164 template <typename dataType, typename triangulationType>
165 inline int execute(OutputCriticalPoints &outCP,
166 Output1Separatrices &outSeps1,
167 Output2Separatrices &outSeps2,
168 OutputManifold &outManifold,
169 const dataType *const scalars,
170 const size_t scalarsMTime,
171 const SimplexId *const offsets,
172 const triangulationType &triangulation);
173
178 inline void setComputeCriticalPoints(const bool state) {
179 this->ComputeCriticalPoints = state;
180 }
185 inline void setComputeSeparatrices1(const bool doAscending,
186 const bool doDescending,
187 const bool doSaddleConnectors) {
188 this->ComputeAscendingSeparatrices1 = doAscending;
189 this->ComputeDescendingSeparatrices1 = doDescending;
190 this->ComputeSaddleConnectors = doSaddleConnectors;
191 }
196 inline void setComputeSeparatrices2(const bool doAscending,
197 const bool doDescending) {
198 this->ComputeAscendingSeparatrices2 = doAscending;
199 this->ComputeDescendingSeparatrices2 = doDescending;
200 }
205 inline void setComputeSegmentation(const bool doAscending,
206 const bool doDescending,
207 const bool doMorseSmale) {
208 this->ComputeAscendingSegmentation = doAscending;
209 this->ComputeDescendingSegmentation = doDescending;
210 this->ComputeFinalSegmentation = doMorseSmale;
211 }
212
218 inline void setReturnSaddleConnectors(const bool state) {
220 }
221
227 inline void
230 }
231
238 data->preconditionCellEdges();
240 }
241
242 protected:
247 struct Separatrix {
257 std::vector<dcg::Cell> geometry_;
258 };
259
264 template <typename triangulationType>
265 int
266 getDescendingSeparatrices1(const std::vector<SimplexId> &saddles,
267 std::vector<Separatrix> &separatrices,
268 const triangulationType &triangulation) const;
269
274 template <typename triangulationType>
275 int getAscendingSeparatrices1(const std::vector<SimplexId> &saddles,
276 std::vector<Separatrix> &separatrices,
277 const triangulationType &triangulation) const;
278
283 template <typename triangulationType>
284 int getSaddleConnectors(const std::vector<SimplexId> &saddles2,
285 std::vector<Separatrix> &separatrices,
286 const triangulationType &triangulation) const;
287
291 template <typename triangulationType>
293 const std::vector<Separatrix> &separatrices,
294 const SimplexId *const offsets,
295 const triangulationType &triangulation) const;
296
301 template <typename triangulationType>
303 const std::vector<SimplexId> &saddles2,
304 std::vector<Separatrix> &separatrices,
305 std::vector<std::vector<SimplexId>> &separatricesSaddles,
306 const triangulationType &triangulation) const;
307
312 template <typename triangulationType>
314 Output2Separatrices &outSeps2,
315 const std::vector<Separatrix> &separatrices,
316 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
317 const SimplexId *const offsets,
318 const triangulationType &triangulation) const;
319
325 template <typename triangulationType>
326 int getDualPolygon(const SimplexId edgeId,
327 SimplexId *const polygon,
328 const size_t polSize,
329 const triangulationType &triangulation) const;
330
334 template <typename triangulationType>
335 int sortDualPolygonVertices(SimplexId *const polygon,
336 const size_t polSize,
337 const triangulationType &triangulation) const;
338
343 template <typename triangulationType>
345 const std::vector<SimplexId> &saddles1,
346 std::vector<Separatrix> &separatrices,
347 std::vector<std::vector<SimplexId>> &separatricesSaddles,
348 const triangulationType &triangulation) const;
349
354 template <typename triangulationType>
356 Output2Separatrices &outSeps2,
357 const std::vector<Separatrix> &separatrices,
358 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
359 const SimplexId *const offsets,
360 const triangulationType &triangulation) const;
361
366 std::vector<std::vector<Separatrix>> &separatrices) const;
367
371 template <typename triangulationType>
372 int setAscendingSegmentation(const std::vector<SimplexId> &maxima,
373 SimplexId *const morseSmaleManifold,
374 const triangulationType &triangulation) const;
375
379 template <typename triangulationType>
380 int setDescendingSegmentation(const std::vector<SimplexId> &minima,
381 SimplexId *const morseSmaleManifold,
382 const triangulationType &triangulation) const;
383
388 template <typename triangulationType>
389 int setFinalSegmentation(const SimplexId numberOfMaxima,
390 const SimplexId *const ascendingManifold,
391 const SimplexId *const descendingManifold,
392 SimplexId *const morseSmaleManifold,
393 const triangulationType &triangulation) const;
394
395 template <typename dataType, typename triangulationType>
396 int returnSaddleConnectors(const double persistenceThreshold,
397 const dataType *const scalars,
398 const SimplexId *const offsets,
399 const triangulationType &triangulation);
400
402
412
416 };
417} // namespace ttk
418
419// ---------------- //
420// Execute method //
421// ---------------- //
422
423template <typename dataType, typename triangulationType>
425 Output1Separatrices &outSeps1,
426 Output2Separatrices &outSeps2,
427 OutputManifold &outManifold,
428 const dataType *const scalars,
429 const size_t scalarsMTime,
430 const SimplexId *const offsets,
431 const triangulationType &triangulation) {
432#ifndef TTK_ENABLE_KAMIKAZE
433 if(scalars == nullptr) {
434 this->printErr("Input scalar field pointer is null.");
435 return -1;
436 }
437
438 if(offsets == nullptr) {
439 this->printErr("Input offset field pointer is null.");
440 return -1;
441 }
442#endif
443 Timer t;
444
445 outCP.clear();
446 outSeps1.clear();
447 outSeps2.clear();
448 const auto dim = triangulation.getDimensionality();
449
452 this->discreteGradient_.setInputScalarField(scalars, scalarsMTime);
453 this->discreteGradient_.setInputOffsets(offsets);
455 triangulation, this->ReturnSaddleConnectors);
456
457 if(this->ReturnSaddleConnectors) {
458 auto persistenceThreshold{this->SaddleConnectorsPersistenceThreshold};
459 if(!this->ThresholdIsAbsolute) {
460 const auto nVerts{triangulation.getNumberOfVertices()};
461 // global extrema are (generally) faster computed on offsets
462 // than on scalar field
463 const auto pair{std::minmax_element(offsets, offsets + nVerts)};
464 // global extrema vertex ids
465 const auto globmin = std::distance(offsets, pair.first);
466 const auto globmax = std::distance(offsets, pair.second);
467 persistenceThreshold *= (scalars[globmax] - scalars[globmin]);
468 this->printMsg("Absolute saddle connectors persistence threshold is "
469 + std::to_string(persistenceThreshold),
471 }
472
474 persistenceThreshold, scalars, offsets, triangulation);
475 }
476
477 std::array<std::vector<SimplexId>, 4> criticalPoints{};
478 {
479 Timer tm{};
480 discreteGradient_.getCriticalPoints(criticalPoints, triangulation);
481 this->printMsg(" Critical points extracted", 1.0, tm.getElapsedTime(),
482 this->threadNumber_, debug::LineMode::NEW,
484 }
485
486 std::vector<std::vector<Separatrix>> separatrices1{};
487
488 // 1-separatrices
489 Timer tm1sep{};
490
491 if(dim > 1 && ComputeDescendingSeparatrices1) {
492 Timer tmp;
493 separatrices1.emplace_back();
494
496 criticalPoints[1], separatrices1.back(), triangulation);
497
498 this->printMsg(" Descending 1-separatrices computed", 1.0,
499 tmp.getElapsedTime(), this->threadNumber_,
501 }
502
503 if(dim > 1 && ComputeAscendingSeparatrices1) {
504 Timer tmp;
505 separatrices1.emplace_back();
506
508 criticalPoints[dim - 1], separatrices1.back(), triangulation);
509
510 this->printMsg(" Ascending 1-separatrices computed", 1.0,
511 tmp.getElapsedTime(), this->threadNumber_,
513 }
514
515 // saddle-connectors
516 if(dim == 3 && ComputeSaddleConnectors) {
517 Timer tmp;
518 separatrices1.emplace_back();
519
520 getSaddleConnectors(criticalPoints[2], separatrices1.back(), triangulation);
521
522 this->printMsg(" Saddle connectors computed", 1.0, tmp.getElapsedTime(),
523 this->threadNumber_, debug::LineMode::NEW,
525 }
526
527 if(dim > 1
530 Timer tmp{};
531
532 this->flattenSeparatricesVectors(separatrices1);
533 setSeparatrices1(outSeps1, separatrices1[0], offsets, triangulation);
534
535 this->printMsg(" 1-separatrices set", 1.0, tmp.getElapsedTime(),
536 this->threadNumber_, debug::LineMode::NEW,
538
539 this->printMsg("1-separatrices computed", 1.0, tm1sep.getElapsedTime(),
540 this->threadNumber_);
541 }
542
543 // 2-separatrices
544 Timer tm2sep{};
545
546 if(dim == 3 && ComputeDescendingSeparatrices2) {
547 Timer tmp;
548 std::vector<Separatrix> separatrices;
549 std::vector<std::vector<SimplexId>> separatricesSaddles;
551 criticalPoints[2], separatrices, separatricesSaddles, triangulation);
553 outSeps2, separatrices, separatricesSaddles, offsets, triangulation);
554
555 this->printMsg(" Descending 2-separatrices computed", 1.0,
556 tmp.getElapsedTime(), this->threadNumber_,
558 }
559
560 if(dim == 3 && ComputeAscendingSeparatrices2) {
561 Timer tmp;
562 std::vector<Separatrix> separatrices;
563 std::vector<std::vector<SimplexId>> separatricesSaddles;
565 criticalPoints[1], separatrices, separatricesSaddles, triangulation);
567 outSeps2, separatrices, separatricesSaddles, offsets, triangulation);
568
569 this->printMsg(" Ascending 2-separatrices computed", 1.0,
570 tmp.getElapsedTime(), this->threadNumber_,
572 }
573
576 this->printMsg("2-separatrices computed", 1.0, tm2sep.getElapsedTime(),
577 this->threadNumber_);
578 }
579
581 Timer tmp;
582
585 criticalPoints[dim], outManifold.ascending_, triangulation);
586 }
589 criticalPoints[0], outManifold.descending_, triangulation);
590 }
593 setFinalSegmentation(criticalPoints[dim].size(), outManifold.ascending_,
594 outManifold.descending_, outManifold.morseSmale_,
595 triangulation);
596 }
597
598 this->printMsg(
599 "Segmentation computed", 1.0, tmp.getElapsedTime(), this->threadNumber_);
600 }
601
604 criticalPoints, outCP.points_, outCP.cellDimensions_, outCP.cellIds_,
605 outCP.isOnBoundary_, outCP.PLVertexIdentifiers_, triangulation);
606
608 discreteGradient_.setManifoldSize(criticalPoints, outManifold.ascending_,
609 outManifold.descending_,
610 outCP.manifoldSize_);
611 }
612 }
613
614 this->printMsg("Data-set ("
615 + std::to_string(triangulation.getNumberOfVertices())
616 + " points) processed",
617 1.0, t.getElapsedTime(), this->threadNumber_);
618
619 return 0;
620}
621
622// ---------------- //
623// 1-Separatrices //
624// ---------------- //
625
626template <typename triangulationType>
628 const std::vector<SimplexId> &saddles,
629 std::vector<Separatrix> &separatrices,
630 const triangulationType &triangulation) const {
631
632 const SimplexId numberOfSaddles = saddles.size();
633
634 // only 2 descending separatrices per 1-saddle
635 const SimplexId numberOfSeparatrices = 2 * numberOfSaddles;
636 separatrices.resize(numberOfSeparatrices);
637
638 // apriori: by default construction, the separatrices are not valid
639#ifdef TTK_ENABLE_OPENMP
640#pragma omp parallel for num_threads(threadNumber_)
641#endif
642 for(SimplexId i = 0; i < numberOfSaddles; ++i) {
643 const Cell saddle{1, saddles[i]};
644
645 // add descending vpaths
646 {
647 const Cell &saddle1 = saddle;
648
649 for(int j = 0; j < 2; ++j) {
650 SimplexId vertexId;
651 triangulation.getEdgeVertex(saddle1.id_, j, vertexId);
652
653 std::vector<Cell> vpath;
654 vpath.push_back(saddle1);
655 discreteGradient_.getDescendingPath(
656 Cell(0, vertexId), vpath, triangulation);
657
658 const Cell &lastCell = vpath.back();
659 if(lastCell.dim_ == 0 and discreteGradient_.isCellCritical(lastCell)) {
660 separatrices[2 * i + j].source_ = saddle;
661 separatrices[2 * i + j].destination_ = lastCell;
662 separatrices[2 * i + j].geometry_ = std::move(vpath);
663 }
664 }
665 }
666 }
667
668 return 0;
669}
670
671template <typename triangulationType>
673 const std::vector<SimplexId> &saddles,
674 std::vector<Separatrix> &separatrices,
675 const triangulationType &triangulation) const {
676
677 const auto dim{triangulation.getDimensionality()};
678
679 // Triangulation method pointers for 3D
680 auto getFaceStarNumber = &triangulationType::getTriangleStarNumber;
681 auto getFaceStar = &triangulationType::getTriangleStar;
682 if(dim == 2) {
683 // Triangulation method pointers for 2D
684 getFaceStarNumber = &triangulationType::getEdgeStarNumber;
685 getFaceStar = &triangulationType::getEdgeStar;
686 }
687
688 const SimplexId numberOfSaddles = saddles.size();
689
690 std::vector<std::vector<Separatrix>> sepsPerSaddle(numberOfSaddles);
691
692 // apriori: by default construction, the separatrices are not valid
693#ifdef TTK_ENABLE_OPENMP
694#pragma omp parallel for num_threads(threadNumber_)
695#endif
696 for(SimplexId i = 0; i < numberOfSaddles; ++i) {
697 const Cell saddle{dim - 1, saddles[i]};
698
699 // add ascending vpaths
700 const auto starNumber{(triangulation.*getFaceStarNumber)(saddle.id_)};
701 for(SimplexId j = 0; j < starNumber; ++j) {
702
703 SimplexId sId{};
704 (triangulation.*getFaceStar)(saddle.id_, j, sId);
705
706 std::vector<Cell> vpath{saddle};
707 discreteGradient_.getAscendingPath(Cell(dim, sId), vpath, triangulation);
708
709 const Cell &lastCell = vpath.back();
710 if(lastCell.dim_ == dim and discreteGradient_.isCellCritical(lastCell)) {
711 sepsPerSaddle[i].emplace_back();
712 sepsPerSaddle[i].back().source_ = saddle;
713 sepsPerSaddle[i].back().destination_ = lastCell;
714 sepsPerSaddle[i].back().geometry_ = std::move(vpath);
715 }
716 }
717 }
718
719 this->flattenSeparatricesVectors(sepsPerSaddle);
720
721 separatrices = std::move(sepsPerSaddle[0]);
722
723 return 0;
724}
725
726template <typename triangulationType>
728 const std::vector<SimplexId> &saddles2,
729 std::vector<Separatrix> &separatrices,
730 const triangulationType &triangulation) const {
731
732 const auto nTriangles = triangulation.getNumberOfTriangles();
733 // visited triangles (one vector per thread)
734 std::vector<bool> isVisited(nTriangles, false);
735 std::vector<SimplexId> visitedTriangles{};
736
737 using Vpath = std::vector<Cell>;
738
739 const auto dim{triangulation.getDimensionality()};
740
741 std::vector<std::vector<Separatrix>> sepsByThread(saddles2.size());
742 std::vector<SimplexId> saddles1{};
743
744#ifdef TTK_ENABLE_OPENMP
745#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
746 firstprivate(isVisited, visitedTriangles, saddles1)
747#endif // TTK_ENABLE_OPENMP
748 for(size_t i = 0; i < saddles2.size(); ++i) {
749 const Cell s2{dim - 1, saddles2[i]};
750
751 VisitedMask mask{isVisited, visitedTriangles};
752 discreteGradient_.getDescendingWall(
753 s2, mask, triangulation, nullptr, &saddles1);
754
755 for(const auto saddle1Id : saddles1) {
756 const Cell s1{1, saddle1Id};
757
758 Vpath vpath;
759 const bool isMultiConnected
760 = discreteGradient_.getAscendingPathThroughWall(
761 s1, s2, isVisited, &vpath, triangulation);
762
763 if(vpath.empty()) {
764 // safety, should be unreachable
765 continue;
766 }
767 const auto &last = vpath.back();
768
769 if(!isMultiConnected && last.dim_ == s2.dim_ && last.id_ == s2.id_) {
770 sepsByThread[i].emplace_back();
771 sepsByThread[i].back().source_ = s1;
772 sepsByThread[i].back().destination_ = s2;
773 sepsByThread[i].back().geometry_ = std::move(vpath);
774 }
775 }
776 }
777
778 this->flattenSeparatricesVectors(sepsByThread);
779
780 separatrices = std::move(sepsByThread[0]);
781
782 return 0;
783}
784
785template <typename triangulationType>
787 Output1Separatrices &outSeps1,
788 const std::vector<Separatrix> &separatrices,
789 const SimplexId *const offsets,
790 const triangulationType &triangulation) const {
791
792 auto &separatrixFunctionMaxima = outSeps1.cl.sepFuncMaxId_;
793 auto &separatrixFunctionMinima = outSeps1.cl.sepFuncMinId_;
794
795 // max existing separatrix id + 1 or 0
796 const SimplexId separatrixId
797 = !outSeps1.cl.separatrixIds_.empty()
798 ? *std::max_element(outSeps1.cl.separatrixIds_.begin(),
799 outSeps1.cl.separatrixIds_.end())
800 + 1
801 : 0;
802
803 // total number of separatrices points
804 auto npoints{static_cast<size_t>(outSeps1.pt.numberOfPoints_)};
805 // total number of separatrices cells
806 auto ncells{static_cast<size_t>(outSeps1.cl.numberOfCells_)};
807 // points beginning id for each separatrix geometry
808 std::vector<size_t> geomPointsBegId{npoints};
809 // cells beginning id for each separatrix geometry
810 std::vector<size_t> geomCellsBegId{ncells};
811
812 // count total number of points and cells, flatten geometryId loops
813 for(size_t i = 0; i < separatrices.size(); ++i) {
814 const auto &sep = separatrices[i];
815 const auto sepSize = sep.geometry_.size();
816 npoints += sepSize;
817 ncells += sepSize - 1;
818 geomPointsBegId.emplace_back(npoints);
819 geomCellsBegId.emplace_back(ncells);
820 }
821
822 const int dimensionality = triangulation.getCellVertexNumber(0) - 1;
823
824 // resize arrays
825 outSeps1.pt.points_.resize(3 * npoints);
826 auto &points = outSeps1.pt.points_;
827 outSeps1.cl.connectivity_.resize(2 * ncells);
828 auto &cellsConn = outSeps1.cl.connectivity_;
829 outSeps1.pt.smoothingMask_.resize(npoints);
830 outSeps1.pt.cellDimensions_.resize(npoints);
831 outSeps1.pt.cellIds_.resize(npoints);
832 outSeps1.cl.sourceIds_.resize(ncells);
833 outSeps1.cl.destinationIds_.resize(ncells);
834 outSeps1.cl.separatrixIds_.resize(ncells);
835 outSeps1.cl.separatrixTypes_.resize(ncells);
836 separatrixFunctionMaxima.resize(separatrixId + separatrices.size());
837 separatrixFunctionMinima.resize(separatrixId + separatrices.size());
838 outSeps1.cl.isOnBoundary_.resize(ncells);
839
840#ifdef TTK_ENABLE_OPENMP
841#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic)
842#endif // TTK_ENABLE_OPENMP
843 for(size_t i = 0; i < separatrices.size(); ++i) {
844 const auto &sep = separatrices[i];
845 const auto &sepGeom = sep.geometry_;
846 const auto sepId = separatrixId + i;
847 // saddle (asc/desc sep) or saddle1 (saddle connector)
848 const dcg::Cell &src = sep.source_;
849 // extremum (asc/desc sep) or saddle2 (saddle connector)
850 const dcg::Cell &dst = sep.destination_;
851
852 // get separatrix type
853 const auto saddleConnector
854 = dimensionality == 3 && src.dim_ == 1 && dst.dim_ == 2;
855 const char sepType
856 = saddleConnector ? 1 : std::min(dst.dim_, dimensionality - 1);
857
858 // compute separatrix function diff
859 const auto vertsOrder = [offsets](const SimplexId a, const SimplexId b) {
860 return offsets[a] < offsets[b];
861 };
862 const std::array<SimplexId, 2> gVerts{
863 discreteGradient_.getCellGreaterVertex(src, triangulation),
864 discreteGradient_.getCellGreaterVertex(dst, triangulation)};
865 const auto sepFuncMax
866 = *std::max_element(gVerts.begin(), gVerts.end(), vertsOrder);
867 const std::array<SimplexId, 2> lVerts{
868 discreteGradient_.getCellLowerVertex(src, triangulation),
869 discreteGradient_.getCellLowerVertex(dst, triangulation)};
870 const auto sepFuncMin
871 = *std::min_element(lVerts.begin(), lVerts.end(), vertsOrder);
872 separatrixFunctionMaxima[sepId] = sepFuncMax;
873 separatrixFunctionMinima[sepId] = sepFuncMin;
874
875 // get boundary condition
876 const auto onBoundary
877 = static_cast<char>(discreteGradient_.isBoundary(src, triangulation))
878 + static_cast<char>(discreteGradient_.isBoundary(dst, triangulation));
879
880 for(size_t j = 0; j < sepGeom.size(); ++j) {
881 const auto &cell = sepGeom[j];
882 std::array<float, 3> pt{};
883 triangulation.getCellIncenter(cell.id_, cell.dim_, pt.data());
884
885 // index of current point in point data arrays
886 const auto k = geomPointsBegId[i] + j;
887
888 points[3 * k + 0] = pt[0];
889 points[3 * k + 1] = pt[1];
890 points[3 * k + 2] = pt[2];
891
892 outSeps1.pt.smoothingMask_[k]
893 = (j == 0 || j == sepGeom.size() - 1) ? 0 : 1;
894 outSeps1.pt.cellDimensions_[k] = cell.dim_;
895 outSeps1.pt.cellIds_[k] = cell.id_;
896
897 // skip filling cell data for first geometry point
898 if(j == 0)
899 continue;
900
901 // index of current cell in cell data arrays
902 const auto l = geomCellsBegId[i] + j - 1;
903
904 cellsConn[2 * l + 0] = k - 1;
905 cellsConn[2 * l + 1] = k;
906
907 outSeps1.cl.sourceIds_[l] = src.id_;
908 outSeps1.cl.destinationIds_[l] = dst.id_;
909 outSeps1.cl.separatrixIds_[l] = sepId;
910 outSeps1.cl.separatrixTypes_[l] = sepType;
911 outSeps1.cl.isOnBoundary_[l] = onBoundary;
912 }
913 }
914
915 // update pointers
916 outSeps1.pt.numberOfPoints_ = npoints;
917 outSeps1.cl.numberOfCells_ = ncells;
918
919 return 0;
920}
921
922// ---------------- //
923// 2-Separatrices //
924// ---------------- //
925
926template <typename triangulationType>
928 const std::vector<SimplexId> &saddles1,
929 std::vector<Separatrix> &separatrices,
930 std::vector<std::vector<SimplexId>> &separatricesSaddles,
931 const triangulationType &triangulation) const {
932 const Cell emptyCell;
933
934 const SimplexId numberOfSaddles = saddles1.size();
935
936 // estimation of the number of separatrices, apriori : numberOfWalls =
937 // numberOfSaddles
938 const SimplexId numberOfSeparatrices = numberOfSaddles;
939 separatrices.resize(numberOfSeparatrices);
940 separatricesSaddles.resize(numberOfSeparatrices);
941
942 const auto nEdges = triangulation.getNumberOfEdges();
943 std::vector<bool> isVisited(nEdges, false);
944 std::vector<SimplexId> visitedEdges{};
945
946 // apriori: by default construction, the separatrices are not valid
947#ifdef TTK_ENABLE_OPENMP
948#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
949 firstprivate(isVisited, visitedEdges)
950#endif // TTK_ENABLE_OPENMP
951 for(SimplexId i = 0; i < numberOfSaddles; ++i) {
952 const Cell saddle1{1, saddles1[i]};
953
954 std::vector<Cell> wall;
955 VisitedMask mask{isVisited, visitedEdges};
956 discreteGradient_.getAscendingWall(
957 saddle1, mask, triangulation, &wall, &separatricesSaddles[i]);
958
959 separatrices[i].source_ = saddle1;
960 separatrices[i].destination_ = emptyCell;
961 separatrices[i].geometry_ = std::move(wall);
962 }
963
964 return 0;
965}
966
967template <typename triangulationType>
969 const std::vector<SimplexId> &saddles2,
970 std::vector<Separatrix> &separatrices,
971 std::vector<std::vector<SimplexId>> &separatricesSaddles,
972 const triangulationType &triangulation) const {
973 const Cell emptyCell;
974
975 const SimplexId numberOfSaddles = saddles2.size();
976
977 // estimation of the number of separatrices, apriori : numberOfWalls =
978 // numberOfSaddles
979 const SimplexId numberOfSeparatrices = numberOfSaddles;
980 separatrices.resize(numberOfSeparatrices);
981 separatricesSaddles.resize(numberOfSeparatrices);
982
983 const auto nTriangles = triangulation.getNumberOfTriangles();
984 std::vector<bool> isVisited(nTriangles, false);
985 std::vector<SimplexId> visitedTriangles{};
986
987 const auto dim{triangulation.getDimensionality()};
988
989 // apriori: by default construction, the separatrices are not valid
990#ifdef TTK_ENABLE_OPENMP
991#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
992 firstprivate(isVisited, visitedTriangles)
993#endif // TTK_ENABLE_OPENMP
994 for(SimplexId i = 0; i < numberOfSaddles; ++i) {
995 const Cell saddle2{dim - 1, saddles2[i]};
996
997 std::vector<Cell> wall;
998 VisitedMask mask{isVisited, visitedTriangles};
999 discreteGradient_.getDescendingWall(
1000 saddle2, mask, triangulation, &wall, &separatricesSaddles[i]);
1001
1002 separatrices[i].source_ = saddle2;
1003 separatrices[i].destination_ = emptyCell;
1004 separatrices[i].geometry_ = std::move(wall);
1005 }
1006
1007 return 0;
1008}
1009
1010template <typename triangulationType>
1012 const SimplexId edgeId,
1013 SimplexId *const polygon,
1014 const size_t polSize,
1015 const triangulationType &triangulation) const {
1016
1017 for(size_t i = 0; i < polSize; ++i) {
1018 SimplexId starId;
1019 triangulation.getEdgeStar(edgeId, i, starId);
1020 polygon[i] = starId;
1021 }
1022
1023 return 0;
1024}
1025
1026template <typename triangulationType>
1028 SimplexId *const polygon,
1029 const size_t polSize,
1030 const triangulationType &triangulation) const {
1031
1032 for(size_t i = 1; i < polSize; ++i) {
1033
1034 // find polygon[i - 1] neighboring tetra in polygon[i..]
1035 bool isFound = false;
1036 size_t j = i;
1037 for(; j < polSize; ++j) {
1038 // check if current is the neighbor
1039 for(SimplexId k = 0;
1040 k < triangulation.getCellNeighborNumber(polygon[i - 1]); ++k) {
1041 SimplexId neighborId{};
1042 triangulation.getCellNeighbor(polygon[i - 1], k, neighborId);
1043 if(neighborId == polygon[j]) {
1044 isFound = true;
1045 break;
1046 }
1047 }
1048 if(isFound)
1049 break;
1050 }
1051
1052 // place polygon[j] next to polygon[i - 1]
1053 if(isFound) {
1054 std::swap(polygon[j], polygon[i]);
1055 }
1056 }
1057
1058 return 0;
1059}
1060
1061template <typename triangulationType>
1063 Output2Separatrices &outSeps2,
1064 const std::vector<Separatrix> &separatrices,
1065 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
1066 const SimplexId *const offsets,
1067 const triangulationType &triangulation) const {
1068
1069 auto &separatrixFunctionMaxima = outSeps2.cl.sepFuncMaxId_;
1070 auto &separatrixFunctionMinima = outSeps2.cl.sepFuncMinId_;
1071
1072 // max existing separatrix id + 1 or 0 if no previous separatrices
1073 const SimplexId separatrixId
1074 = !outSeps2.cl.separatrixIds_.empty()
1075 ? *std::max_element(outSeps2.cl.separatrixIds_.begin(),
1076 outSeps2.cl.separatrixIds_.end())
1077 + 1
1078 : 0;
1079
1080 // total number of separatrices points
1081 auto npoints{static_cast<size_t>(outSeps2.pt.numberOfPoints_)};
1082 // total number of separatrices cells
1083 auto ncells{static_cast<size_t>(outSeps2.cl.numberOfCells_)};
1084 // old number of separatrices cells
1085 const auto noldcells{ncells};
1086 // index of last vertex of last old cell + 1
1087 const auto firstCellId{outSeps2.cl.connectivity_.size()};
1088 // cells beginning id for each separatrix geometry
1089 std::vector<size_t> geomCellsBegId{ncells};
1090
1091 // count total number of cells, flatten geometryId loops
1092 for(size_t i = 0; i < separatrices.size(); ++i) {
1093 const auto &sep = separatrices[i];
1094 ncells += sep.geometry_.size();
1095 geomCellsBegId.emplace_back(ncells);
1096 }
1097
1098 // store the separatrices info (one per separatrix)
1099 std::vector<SimplexId> sepSourceIds(separatrices.size());
1100 std::vector<SimplexId> sepIds(separatrices.size());
1101 std::vector<char> sepOnBoundary(separatrices.size());
1102 separatrixFunctionMaxima.resize(separatrixId + separatrices.size());
1103 separatrixFunctionMinima.resize(separatrixId + separatrices.size());
1104 // store the polygonal cells tetras SimplexId
1105 std::vector<SimplexId> polygonNTetras(ncells - noldcells);
1106 std::vector<SimplexId> polygonEdgeIds(ncells - noldcells);
1107 std::vector<SimplexId> polygonSepInfosIds(ncells - noldcells);
1108
1109#ifdef TTK_ENABLE_OPENMP
1110#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic)
1111#endif // TTK_ENABLE_OPENMP
1112 for(size_t i = 0; i < separatrices.size(); ++i) {
1113 const auto &sep = separatrices[i];
1114 const auto &sepGeom = sep.geometry_;
1115 const auto &sepSaddles = separatricesSaddles[i];
1116 const auto sepId = separatrixId + i;
1117 const dcg::Cell &src = sep.source_; // saddle1
1118
1119 // compute separatrix function diff
1120 const auto sepFuncMin
1121 = discreteGradient_.getCellLowerVertex(src, triangulation);
1122 SimplexId sepFuncMax{};
1123 if(!sepSaddles.empty()) {
1124 // find minimum vertex on the critical triangles of the 2-separatrix
1125 const auto maxId = *std::max_element(
1126 sepSaddles.begin(), sepSaddles.end(),
1127 [&triangulation, offsets, this](const SimplexId a, const SimplexId b) {
1128 return offsets[discreteGradient_.getCellGreaterVertex(
1129 Cell{2, a}, triangulation)]
1130 < offsets[discreteGradient_.getCellGreaterVertex(
1131 Cell{2, b}, triangulation)];
1132 });
1133 sepFuncMax
1134 = discreteGradient_.getCellGreaterVertex(Cell{2, maxId}, triangulation);
1135 } else {
1136 // find maximum vertex by iterating over all the edges in the
1137 // 2-separatrix
1138 const auto maxId = *std::max_element(
1139 sepGeom.begin(), sepGeom.end(),
1140 [&triangulation, offsets, this](const Cell &a, const Cell &b) {
1141 return offsets[discreteGradient_.getCellGreaterVertex(
1142 a, triangulation)]
1143 < offsets[discreteGradient_.getCellGreaterVertex(
1144 b, triangulation)];
1145 });
1146 sepFuncMax = discreteGradient_.getCellGreaterVertex(maxId, triangulation);
1147 }
1148
1149 // get boundary condition
1150 const char onBoundary
1151 = (sepSaddles.empty()
1152 ? 0
1153 : std::count_if(sepSaddles.begin(), sepSaddles.end(),
1154 [&triangulation](const SimplexId a) {
1155 return triangulation.isTriangleOnBoundary(a);
1156 }))
1157 + triangulation.isEdgeOnBoundary(src.id_);
1158
1159 sepIds[i] = sepId;
1160 sepSourceIds[i] = src.id_;
1161 separatrixFunctionMaxima[sepId] = sepFuncMax;
1162 separatrixFunctionMinima[sepId] = sepFuncMin;
1163 sepOnBoundary[i] = onBoundary;
1164
1165 for(size_t j = 0; j < sepGeom.size(); ++j) {
1166 const auto &cell = sepGeom[j];
1167 // index of current cell in cell data arrays
1168 const auto k = geomCellsBegId[i] + j - noldcells;
1169
1170 polygonNTetras[k] = triangulation.getEdgeStarNumber(cell.id_);
1171
1172 if(polygonNTetras[k] > 2) {
1173 polygonEdgeIds[k] = cell.id_;
1174 polygonSepInfosIds[k] = i;
1175 }
1176 }
1177 }
1178
1179 // indices of valid polygon tetras
1180 std::vector<SimplexId> validTetraIds{};
1181 validTetraIds.reserve(polygonNTetras.size());
1182
1183 for(size_t i = 0; i < polygonNTetras.size(); ++i) {
1184 if(polygonNTetras[i] > 2) {
1185 validTetraIds.emplace_back(i);
1186 }
1187 }
1188
1189 // count number of valid new cells and new points
1190 size_t nnewpoints{};
1191 std::vector<SimplexId> pointsPerCell(validTetraIds.size() + 1);
1192 for(size_t i = 0; i < validTetraIds.size(); ++i) {
1193 nnewpoints += polygonNTetras[validTetraIds[i]];
1194 pointsPerCell[i + 1] = nnewpoints;
1195 }
1196
1197 // resize connectivity array
1198 outSeps2.cl.connectivity_.resize(firstCellId + nnewpoints);
1199 auto cellsConn = &outSeps2.cl.connectivity_[firstCellId];
1200 // copy of cell connectivity array (for removing duplicates vertices)
1201 std::vector<SimplexId> cellVertsIds(nnewpoints);
1202
1203#ifdef TTK_ENABLE_OPENMP
1204#pragma omp parallel for num_threads(threadNumber_)
1205#endif // TTK_ENABLE_OPENMP
1206 for(size_t i = 0; i < validTetraIds.size(); ++i) {
1207 const auto k = validTetraIds[i];
1208
1209 // get tetras in edge star
1210 getDualPolygon(polygonEdgeIds[k], &cellVertsIds[pointsPerCell[i]],
1211 polygonNTetras[k], triangulation);
1212 // sort tetras (in-place)
1213 sortDualPolygonVertices(
1214 &cellVertsIds[pointsPerCell[i]], polygonNTetras[k], triangulation);
1215
1216 for(SimplexId j = 0; j < polygonNTetras[k]; ++j) {
1217 cellsConn[pointsPerCell[i] + j] = cellVertsIds[pointsPerCell[i] + j];
1218 }
1219 }
1220
1221 TTK_PSORT(this->threadNumber_, cellVertsIds.begin(), cellVertsIds.end());
1222 const auto last = std::unique(cellVertsIds.begin(), cellVertsIds.end());
1223 cellVertsIds.erase(last, cellVertsIds.end());
1224
1225 // vertex Id to index in points array
1226 std::vector<SimplexId> vertId2PointsId(triangulation.getNumberOfCells());
1227
1228 const auto noldpoints{npoints};
1229 npoints += cellVertsIds.size();
1230 ncells = noldcells + validTetraIds.size();
1231
1232 // resize arrays
1233 outSeps2.pt.points_.resize(3 * npoints);
1234 auto points = &outSeps2.pt.points_[3 * noldpoints];
1235 outSeps2.cl.offsets_.resize(ncells + 1);
1236 outSeps2.cl.offsets_[0] = 0;
1237 auto cellsOff = &outSeps2.cl.offsets_[noldcells];
1238 outSeps2.cl.sourceIds_.resize(ncells);
1239 outSeps2.cl.separatrixIds_.resize(ncells);
1240 outSeps2.cl.separatrixTypes_.resize(ncells);
1241 outSeps2.cl.isOnBoundary_.resize(ncells);
1242
1243#ifdef TTK_ENABLE_OPENMP
1244#pragma omp parallel for num_threads(threadNumber_)
1245#endif // TTK_ENABLE_OPENMP
1246 for(size_t i = 0; i < cellVertsIds.size(); ++i) {
1247 // vertex 3D coords
1248 triangulation.getTetraIncenter(cellVertsIds[i], &points[3 * i]);
1249 // vertex index in cellVertsIds array (do not forget offset)
1250 vertId2PointsId[cellVertsIds[i]] = i + noldpoints;
1251 }
1252
1253#ifdef TTK_ENABLE_OPENMP
1254#pragma omp parallel for num_threads(threadNumber_)
1255#endif // TTK_ENABLE_OPENMP
1256 for(size_t i = 0; i < validTetraIds.size(); ++i) {
1257 const auto m = validTetraIds[i];
1258 const auto k = pointsPerCell[i];
1259 for(SimplexId j = 0; j < polygonNTetras[m]; ++j) {
1260 cellsConn[k + j] = vertId2PointsId[cellsConn[k + j]];
1261 }
1262 const auto l = i + noldcells;
1263 const auto n = polygonSepInfosIds[m];
1264 outSeps2.cl.sourceIds_[l] = sepSourceIds[n];
1265 outSeps2.cl.separatrixIds_[l] = sepIds[n];
1266 outSeps2.cl.separatrixTypes_[l] = 1;
1267 outSeps2.cl.isOnBoundary_[l] = sepOnBoundary[n];
1268 }
1269
1270 for(size_t i = 0; i < validTetraIds.size(); ++i) {
1271 // fill offsets sequentially (due to iteration dependencies)
1272 cellsOff[i + 1] = cellsOff[i] + polygonNTetras[validTetraIds[i]];
1273 }
1274
1275 outSeps2.pt.numberOfPoints_ = npoints;
1276 outSeps2.cl.numberOfCells_ = ncells;
1277
1278 return 0;
1279}
1280
1281template <typename triangulationType>
1283 Output2Separatrices &outSeps2,
1284 const std::vector<Separatrix> &separatrices,
1285 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
1286 const SimplexId *const offsets,
1287 const triangulationType &triangulation) const {
1288
1289 auto &separatrixFunctionMaxima = outSeps2.cl.sepFuncMaxId_;
1290 auto &separatrixFunctionMinima = outSeps2.cl.sepFuncMinId_;
1291
1292 // max existing separatrix id + 1 or 0 if no previous separatrices
1293 const SimplexId separatrixId
1294 = !outSeps2.cl.separatrixIds_.empty()
1295 ? *std::max_element(outSeps2.cl.separatrixIds_.begin(),
1296 outSeps2.cl.separatrixIds_.end())
1297 + 1
1298 : 0;
1299
1300 // total number of separatrices points
1301 auto npoints{static_cast<size_t>(outSeps2.pt.numberOfPoints_)};
1302 // total number of separatrices cells
1303 auto ncells{static_cast<size_t>(outSeps2.cl.numberOfCells_)};
1304 // old number of separatrices cells
1305 const auto noldcells{ncells};
1306 // index of last vertex of last old cell + 1
1307 const auto firstCellId{outSeps2.cl.connectivity_.size()};
1308
1309 // cells beginning id for each separatrix geometry
1310 std::vector<size_t> geomCellsBegId{ncells};
1311
1312 // count total number of cells, flatten geometryId loops
1313 for(size_t i = 0; i < separatrices.size(); ++i) {
1314 const auto &sep = separatrices[i];
1315 ncells += sep.geometry_.size();
1316 geomCellsBegId.emplace_back(ncells);
1317 }
1318
1319 // resize arrays
1320 outSeps2.cl.offsets_.resize(ncells + 1);
1321 outSeps2.cl.offsets_[0] = 0;
1322 outSeps2.cl.connectivity_.resize(
1323 firstCellId + 3 * (ncells - noldcells)); // triangles cells
1324 auto cellsOff = &outSeps2.cl.offsets_[noldcells];
1325 auto cellsConn = &outSeps2.cl.connectivity_[firstCellId];
1326 outSeps2.cl.sourceIds_.resize(ncells);
1327 outSeps2.cl.separatrixIds_.resize(ncells);
1328 outSeps2.cl.separatrixTypes_.resize(ncells);
1329 separatrixFunctionMaxima.resize(separatrixId + separatrices.size());
1330 separatrixFunctionMinima.resize(separatrixId + separatrices.size());
1331 outSeps2.cl.isOnBoundary_.resize(ncells);
1332
1333 // store the cells/triangles vertices vertexId
1334 std::vector<SimplexId> cellVertsIds(3 * (ncells - noldcells));
1335
1336#ifdef TTK_ENABLE_OPENMP
1337#pragma omp parallel for num_threads(threadNumber_)
1338#endif // TTK_ENABLE_OPENMP
1339 for(size_t i = 0; i < separatrices.size(); ++i) {
1340 const auto &sep = separatrices[i];
1341 const auto &sepGeom = sep.geometry_;
1342 const auto &sepSaddles = separatricesSaddles[i];
1343 const auto sepId = separatrixId + i;
1344 const dcg::Cell &src = sep.source_; // saddle2
1345 const char sepType = 2;
1346
1347 // compute separatrix function diff
1348 const auto sepFuncMax
1349 = discreteGradient_.getCellGreaterVertex(src, triangulation);
1350 SimplexId sepFuncMin{};
1351 if(!sepSaddles.empty()) {
1352 // find minimum vertex on the critical edges of the 2-separatrix
1353 const auto minId = *std::min_element(
1354 sepSaddles.begin(), sepSaddles.end(),
1355 [&triangulation, offsets, this](const SimplexId a, const SimplexId b) {
1356 return offsets[discreteGradient_.getCellLowerVertex(
1357 Cell{1, a}, triangulation)]
1358 < offsets[discreteGradient_.getCellLowerVertex(
1359 Cell{1, b}, triangulation)];
1360 });
1361 sepFuncMin
1362 = discreteGradient_.getCellLowerVertex(Cell{1, minId}, triangulation);
1363 } else {
1364 // find minimum vertex by iterating over all the triangles in the
1365 // 2-separatrix
1366 const auto minId = *std::min_element(
1367 sepGeom.begin(), sepGeom.end(),
1368 [&triangulation, offsets, this](const Cell &a, const Cell &b) {
1369 return offsets[discreteGradient_.getCellLowerVertex(a, triangulation)]
1370 < offsets[discreteGradient_.getCellLowerVertex(
1371 b, triangulation)];
1372 });
1373 sepFuncMin = discreteGradient_.getCellLowerVertex(minId, triangulation);
1374 }
1375 separatrixFunctionMaxima[sepId] = sepFuncMax;
1376 separatrixFunctionMinima[sepId] = sepFuncMin;
1377
1378 // get boundary condition
1379 const char onBoundary
1380 = (sepSaddles.empty()
1381 ? 0
1382 : std::count_if(sepSaddles.begin(), sepSaddles.end(),
1383 [&triangulation](const SimplexId a) {
1384 return triangulation.isEdgeOnBoundary(a);
1385 }))
1386 + triangulation.isTriangleOnBoundary(src.id_);
1387
1388 for(size_t j = 0; j < sepGeom.size(); ++j) {
1389 const auto &cell = sepGeom[j];
1390
1391 // first store the SimplexId of the cell/triangle vertices
1392 SimplexId v0{}, v1{}, v2{};
1393 triangulation.getTriangleVertex(cell.id_, 0, v0);
1394 triangulation.getTriangleVertex(cell.id_, 1, v1);
1395 triangulation.getTriangleVertex(cell.id_, 2, v2);
1396
1397 // index of current cell in cell data arrays
1398 const auto l = geomCellsBegId[i] + j;
1399 // index of current cell among all new cells
1400 const auto m = l - noldcells;
1401
1402 cellsConn[3 * m + 0] = v0;
1403 cellsConn[3 * m + 1] = v1;
1404 cellsConn[3 * m + 2] = v2;
1405 cellVertsIds[3 * m + 0] = v0;
1406 cellVertsIds[3 * m + 1] = v1;
1407 cellVertsIds[3 * m + 2] = v2;
1408
1409 outSeps2.cl.sourceIds_[l] = src.id_;
1410 outSeps2.cl.separatrixIds_[l] = sepId;
1411 outSeps2.cl.separatrixTypes_[l] = sepType;
1412 outSeps2.cl.isOnBoundary_[l] = onBoundary;
1413 }
1414 }
1415
1416 // reduce the cell vertices ids
1417 // (cells are triangles sharing two vertices)
1418 TTK_PSORT(this->threadNumber_, cellVertsIds.begin(), cellVertsIds.end());
1419 const auto last = std::unique(cellVertsIds.begin(), cellVertsIds.end());
1420 cellVertsIds.erase(last, cellVertsIds.end());
1421
1422 // vertex Id to index in points array
1423 std::vector<size_t> vertId2PointsId(triangulation.getNumberOfVertices());
1424
1425 const auto noldpoints{npoints};
1426 npoints += cellVertsIds.size();
1427 outSeps2.pt.points_.resize(3 * npoints);
1428 auto points = &outSeps2.pt.points_[3 * noldpoints];
1429
1430#ifdef TTK_ENABLE_OPENMP
1431#pragma omp parallel for num_threads(threadNumber_)
1432#endif // TTK_ENABLE_OPENMP
1433 for(size_t i = 0; i < cellVertsIds.size(); ++i) {
1434 // vertex 3D coords
1435 triangulation.getVertexPoint(
1436 cellVertsIds[i], points[3 * i + 0], points[3 * i + 1], points[3 * i + 2]);
1437 // vertex index in cellVertsIds array (do not forget offset)
1438 vertId2PointsId[cellVertsIds[i]] = i + noldpoints;
1439 }
1440
1441 const auto lastOffset = noldcells == 0 ? 0 : cellsOff[-1];
1442
1443#ifdef TTK_ENABLE_OPENMP
1444#pragma omp parallel for num_threads(threadNumber_)
1445#endif // TTK_ENABLE_OPENMP
1446 for(size_t i = 0; i < ncells - noldcells; ++i) {
1447 cellsOff[i] = 3 * i + lastOffset;
1448 cellsConn[3 * i + 0] = vertId2PointsId[cellsConn[3 * i + 0]];
1449 cellsConn[3 * i + 1] = vertId2PointsId[cellsConn[3 * i + 1]];
1450 cellsConn[3 * i + 2] = vertId2PointsId[cellsConn[3 * i + 2]];
1451 }
1452
1453 cellsOff[ncells - noldcells] = cellsOff[ncells - noldcells - 1] + 3;
1454
1455 outSeps2.pt.numberOfPoints_ = npoints;
1456 outSeps2.cl.numberOfCells_ = ncells;
1457
1458 return 0;
1459}
1460
1461// ---------------- //
1462// Segmentation //
1463// ---------------- //
1464
1465template <typename triangulationType>
1467 const std::vector<SimplexId> &maxima,
1468 SimplexId *const morseSmaleManifold,
1469 const triangulationType &triangulation) const {
1470
1471 if(morseSmaleManifold == nullptr) {
1472 this->printErr("Could not compute ascending segmentation");
1473 return 1;
1474 }
1475
1476 Timer tm{};
1477
1478 const auto nVerts{triangulation.getNumberOfVertices()};
1479 std::fill(morseSmaleManifold, morseSmaleManifold + nVerts, -1);
1480 if(maxima.empty()) {
1481 // shortcut for elevation
1482 return 0;
1483 }
1484
1485 const auto nCells{triangulation.getNumberOfCells()};
1486 std::vector<SimplexId> morseSmaleManifoldOnCells(nCells, -1);
1487
1488 size_t nMax{};
1489 for(const auto &id : maxima) {
1490 // mark the maxima
1491 morseSmaleManifoldOnCells[id] = nMax++;
1492 }
1493
1494 const auto dim{triangulation.getDimensionality()};
1495
1496 // Triangulation method pointers for 3D
1497 auto getFaceStarNumber = &triangulationType::getTriangleStarNumber;
1498 auto getFaceStar = &triangulationType::getTriangleStar;
1499 if(dim == 2) {
1500 // Triangulation method pointers for 2D
1501 getFaceStarNumber = &triangulationType::getEdgeStarNumber;
1502 getFaceStar = &triangulationType::getEdgeStar;
1503 } else if(dim == 1) {
1504 // Triangulation method pointers for 1D
1505 getFaceStarNumber = &triangulationType::getVertexStarNumber;
1506 getFaceStar = &triangulationType::getVertexStar;
1507 }
1508
1509 // cells visited during the propagation alongside one integral line
1510 std::vector<SimplexId> visited{};
1511 // all marked cells
1512 std::vector<uint8_t> isMarked(nCells, 0);
1513
1514#ifdef TTK_ENABLE_OPENMP
1515#pragma omp parallel for num_threads(threadNumber_) firstprivate(visited)
1516#endif // TTK_ENABLE_OPENMP
1517 for(SimplexId i = 0; i < nCells; ++i) {
1518 if(isMarked[i] == 1) {
1519 continue;
1520 }
1521 visited.clear();
1522 auto curr{i};
1523 while(morseSmaleManifoldOnCells[curr] == -1) {
1524 if(isMarked[curr] == 1) {
1525 break;
1526 }
1527 // follow a V-path till an already marked cell is reached
1528 const auto paired{this->discreteGradient_.getPairedCell(
1529 Cell{dim, curr}, triangulation, true)};
1530 SimplexId next{curr};
1531 const auto nStars{(triangulation.*getFaceStarNumber)(paired)};
1532 for(SimplexId j = 0; j < nStars; ++j) {
1533 (triangulation.*getFaceStar)(paired, j, next);
1534 // get the first cell != curr (what of non-manifold datasets?)
1535 if(next != curr) {
1536 break;
1537 }
1538 }
1539 visited.emplace_back(curr);
1540 if(next == curr) {
1541 // on the boundary?
1542 break;
1543 }
1544 curr = next;
1545 }
1546 for(const auto el : visited) {
1547 morseSmaleManifoldOnCells[el] = morseSmaleManifoldOnCells[curr];
1548 isMarked[el] = 1;
1549 }
1550 }
1551
1552#ifdef TTK_ENABLE_OPENMP
1553#pragma omp parallel for num_threads(threadNumber_)
1554#endif // TTK_ENABLE_OPENMP
1555 for(SimplexId i = 0; i < nVerts; ++i) {
1556 if(triangulation.getVertexStarNumber(i) < 1) {
1557 // handle non-manifold datasets?
1558 continue;
1559 }
1560 SimplexId starId{};
1561 triangulation.getVertexStar(i, 0, starId);
1562 // put segmentation infos from cells to points
1563 morseSmaleManifold[i] = morseSmaleManifoldOnCells[starId];
1564 }
1565
1566 this->printMsg(" Ascending segmentation computed", 1.0, tm.getElapsedTime(),
1567 this->threadNumber_, debug::LineMode::NEW,
1569
1570 return 0;
1571}
1572
1573template <typename triangulationType>
1575 const std::vector<SimplexId> &minima,
1576 SimplexId *const morseSmaleManifold,
1577 const triangulationType &triangulation) const {
1578
1579 if(morseSmaleManifold == nullptr) {
1580 this->printErr("Could not compute descending segmentation");
1581 return 1;
1582 }
1583
1584 Timer tm{};
1585
1586 const auto nVerts{triangulation.getNumberOfVertices()};
1587
1588 if(minima.size() == 1) {
1589 // shortcut for elevation
1590 std::fill(morseSmaleManifold, morseSmaleManifold + nVerts, 0);
1591 return 0;
1592 }
1593
1594 std::fill(morseSmaleManifold, morseSmaleManifold + nVerts, -1);
1595
1596 size_t nMin{};
1597 for(const auto &cp : minima) {
1598 // mark the minima
1599 morseSmaleManifold[cp] = nMin++;
1600 }
1601
1602 std::vector<SimplexId> visited{};
1603
1604#ifdef TTK_ENABLE_OPENMP
1605#pragma omp parallel for num_threads(threadNumber_) firstprivate(visited)
1606#endif // TTK_ENABLE_OPENMP
1607 for(SimplexId i = 0; i < nVerts; ++i) {
1608 if(morseSmaleManifold[i] != -1) {
1609 continue;
1610 }
1611 visited.clear();
1612 auto curr{i};
1613 while(morseSmaleManifold[curr] == -1) {
1614 // follow a V-path till an already marked vertex is reached
1615 const auto pairedEdge{
1616 this->discreteGradient_.getPairedCell(Cell{0, curr}, triangulation)};
1617 SimplexId next{};
1618 triangulation.getEdgeVertex(pairedEdge, 0, next);
1619 if(next == curr) {
1620 triangulation.getEdgeVertex(pairedEdge, 1, next);
1621 }
1622 visited.emplace_back(curr);
1623 curr = next;
1624 }
1625 for(const auto el : visited) {
1626 morseSmaleManifold[el] = morseSmaleManifold[curr];
1627 }
1628 }
1629
1630 this->printMsg(" Descending segmentation computed", 1.0, tm.getElapsedTime(),
1631 this->threadNumber_, debug::LineMode::NEW,
1633
1634 return 0;
1635}
1636
1637template <typename triangulationType>
1639 const SimplexId numberOfMaxima,
1640 const SimplexId *const ascendingManifold,
1641 const SimplexId *const descendingManifold,
1642 SimplexId *const morseSmaleManifold,
1643 const triangulationType &triangulation) const {
1644
1645 if(ascendingManifold == nullptr || descendingManifold == nullptr
1646 || morseSmaleManifold == nullptr) {
1647 this->printErr("Could not compute final segmentation");
1648 return 1;
1649 }
1650
1651 Timer tm{};
1652
1653 const size_t nVerts = triangulation.getNumberOfVertices();
1654
1655 // associate a unique "sparse region id" to each (ascending, descending) pair
1656
1657#ifdef TTK_ENABLE_OPENMP
1658#pragma omp parallel for num_threads(threadNumber_)
1659#endif // TTK_ENABLE_OPENMP
1660 for(size_t i = 0; i < nVerts; ++i) {
1661 const auto d = ascendingManifold[i];
1662 const auto a = descendingManifold[i];
1663 if(a == -1 || d == -1) {
1664 morseSmaleManifold[i] = -1;
1665 } else {
1666 morseSmaleManifold[i] = a * numberOfMaxima + d;
1667 }
1668 }
1669
1670 // store the "sparse region ids" by copying the morseSmaleManifold output
1671 std::vector<SimplexId> sparseRegionIds(
1672 morseSmaleManifold, morseSmaleManifold + nVerts);
1673
1674 // get unique "sparse region ids"
1675 TTK_PSORT(
1676 this->threadNumber_, sparseRegionIds.begin(), sparseRegionIds.end());
1677 const auto last = std::unique(sparseRegionIds.begin(), sparseRegionIds.end());
1678 sparseRegionIds.erase(last, sparseRegionIds.end());
1679
1680 // "sparse region id" -> "dense region id"
1681 std::map<SimplexId, size_t> sparseToDenseRegionId{};
1682
1683 for(size_t i = 0; i < sparseRegionIds.size(); ++i) {
1684 sparseToDenseRegionId[sparseRegionIds[i]] = i;
1685 }
1686
1687 // update region id on all vertices: "sparse id" -> "dense id"
1688
1689#ifdef TTK_ENABLE_OPENMP
1690#pragma omp parallel for num_threads(threadNumber_)
1691#endif // TTK_ENABLE_OPENMP
1692 for(size_t i = 0; i < nVerts; ++i) {
1693 morseSmaleManifold[i] = sparseToDenseRegionId[morseSmaleManifold[i]];
1694 }
1695
1696 this->printMsg(" Final segmentation computed", 1.0, tm.getElapsedTime(),
1697 this->threadNumber_, debug::LineMode::NEW,
1699
1700 return 0;
1701}
1702
1703template <typename dataType, typename triangulationType>
1705 const double persistenceThreshold,
1706 const dataType *const scalars,
1707 const SimplexId *const offsets,
1708 const triangulationType &triangulation) {
1709
1710 Timer tm{};
1711
1712 const auto dim{triangulation.getDimensionality()};
1713 if(dim != 3) {
1714 this->printWrn("Can't return saddle connectors without a 3D dataset");
1715 return 0;
1716 }
1717
1718 // compute saddle-saddle Persistence Pairs from DiscreteMorseSandwich
1720 dms.setThreadNumber(this->threadNumber_);
1721 dms.setDebugLevel(this->debugLevel_);
1722 dms.setGradient(std::move(this->discreteGradient_));
1723
1724 using PersPairType = DiscreteMorseSandwich::PersistencePair;
1725
1726 std::vector<PersPairType> dms_pairs{};
1727 dms.computePersistencePairs(dms_pairs, offsets, triangulation, false, true);
1728 this->discreteGradient_ = dms.getGradient();
1729 // reset gradient pointer to local storage
1730 this->discreteGradient_.setLocalGradient();
1731
1732 const auto getPersistence
1733 = [this, &triangulation, scalars](const PersPairType &p) {
1734 return this->discreteGradient_.getPersistence(
1735 Cell{2, p.death}, Cell{1, p.birth}, scalars, triangulation);
1736 };
1737
1738 // saddle-saddle pairs should be in one continuous block inside dms_pairs
1739 const auto firstSadSadPair{std::distance(
1740 dms_pairs.begin(),
1741 std::find_if(dms_pairs.begin(), dms_pairs.end(),
1742 [](const auto &pair) { return pair.type == 1; }))};
1743
1744 // DMS pairs output are already sorted by 2-saddles in ascending order
1745
1746 std::vector<bool> isVisited(triangulation.getNumberOfTriangles(), false);
1747 std::vector<SimplexId> visitedTriangles{};
1748
1749 size_t nReturned{};
1750
1751 // offset from firstSadSadPair should index s2Children and simplifyS2
1752 const auto &s2Children{dms.get2SaddlesChildren()};
1753 std::vector<bool> simplifyS2(s2Children.size(), true);
1754
1755 for(size_t i = firstSadSadPair; i < dms_pairs.size(); ++i) {
1756 const auto &pair{dms_pairs[i]};
1757
1758 if(pair.type != 1 || getPersistence(pair) > persistenceThreshold) {
1759 continue;
1760 }
1761
1762 const auto o{i - firstSadSadPair};
1763 const Cell birth{1, pair.birth};
1764 const Cell death{2, pair.death};
1765 bool skip{false};
1766 for(const auto child : s2Children[o]) {
1767 if(!simplifyS2[child]) {
1768 skip = true;
1769 break;
1770 }
1771 }
1772 if(skip) {
1773 this->printMsg("Skipping saddle connector " + birth.to_string() + " -> "
1774 + death.to_string(),
1776 simplifyS2[o] = false;
1777 continue;
1778 }
1779 // 1. get the 2-saddle wall
1780 VisitedMask mask{isVisited, visitedTriangles};
1781 this->discreteGradient_.getDescendingWall(death, mask, triangulation);
1782 // 2. get the saddle connector
1783 std::vector<Cell> vpath{};
1784 this->discreteGradient_.getAscendingPathThroughWall(
1785 birth, death, isVisited, &vpath, triangulation, true);
1786 // 3. reverse the gradient on the saddle connector path
1787 if(vpath.back() == death) {
1788 this->discreteGradient_.reverseAscendingPathOnWall(vpath, triangulation);
1789 nReturned++;
1790 } else {
1791 this->printMsg("Could not return saddle connector " + birth.to_string()
1792 + " -> " + death.to_string(),
1794 simplifyS2[o] = false;
1795 }
1796 }
1797
1798 this->printMsg("Returned " + std::to_string(nReturned) + " saddle connectors",
1799 1.0, tm.getElapsedTime(), this->threadNumber_);
1800
1801 return 0;
1802}
#define TTK_PSORT(NTHREADS,...)
Parallel sort macro.
Definition: OpenMP.h:46
AbstractTriangulation is an interface class that defines an interface for efficient traversal methods...
int threadNumber_
Definition: BaseClass.h:95
virtual int setThreadNumber(const int threadNumber)
Definition: BaseClass.h:80
Minimalist debugging class.
Definition: Debug.h:88
int debugLevel_
Definition: Debug.h:379
virtual int setDebugLevel(const int &debugLevel)
Definition: Debug.cpp:147
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
TTK DiscreteMorseSandwich processing package.
TTK processing package for the computation of Morse-Smale complexes.
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 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)
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)
int getAscendingSeparatrices1(const std::vector< SimplexId > &saddles, std::vector< Separatrix > &separatrices, const triangulationType &triangulation) const
double getElapsedTime()
Definition: Timer.h:15
TTK discreteGradient processing package.
void setInputOffsets(const SimplexId *const data)
int setManifoldSize(const std::array< std::vector< SimplexId >, 4 > &criticalCellsByDim, const SimplexId *const ascendingManifold, const SimplexId *const descendingManifold, std::vector< SimplexId > &manifoldSize) const
void setInputScalarField(const void *const data, const size_t mTime)
int getCriticalPoints(std::array< std::vector< SimplexId >, 4 > &criticalCellsByDim, const triangulationType &triangulation) const
void preconditionTriangulation(AbstractTriangulation *const data)
int buildGradient(const triangulationType &triangulation, bool bypassCache=false)
int setCriticalPoints(const std::array< std::vector< SimplexId >, 4 > &criticalCellsByDim, std::vector< std::array< float, 3 > > &points, std::vector< char > &cellDimensions, std::vector< SimplexId > &cellIds, std::vector< char > &isOnBoundary, std::vector< SimplexId > &PLVertexIdentifiers, const triangulationType &triangulation) const
The Topology ToolKit.
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::@1 cl
std::vector< ttk::SimplexId > connectivity_
std::vector< ttk::SimplexId > sourceIds_
struct ttk::MorseSmaleComplex::Output1Separatrices::@0 pt
2-Separatrices point and cell data arrays
struct ttk::MorseSmaleComplex::Output2Separatrices::@2 pt
std::vector< ttk::SimplexId > sourceIds_
struct ttk::MorseSmaleComplex::Output2Separatrices::@3 cl
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.
std::vector< dcg::Cell > geometry_
Auto-cleaning re-usable graph propagations data structure.
Definition: VisitedMask.h:27
printMsg(debug::output::GREEN+"                           "+debug::output::ENDCOLOR+debug::output::GREEN+"▒"+debug::output::ENDCOLOR+debug::output::GREEN+"▒▒▒▒▒▒▒▒▒▒▒▒▒░"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)