TTK
Loading...
Searching...
No Matches
DiscreteMorseSandwich.h
Go to the documentation of this file.
1
21
22#pragma once
23
24#include <DiscreteGradient.h>
25
26#include <algorithm>
27#include <numeric>
28
29namespace ttk {
30 class DiscreteMorseSandwich : virtual public Debug {
31 public:
33
43 int type;
44
46 : birth{b}, death{d}, type{t} {
47 }
48 };
49
52 }
53
54 inline void setInputOffsets(const SimplexId *const offsets) {
55 this->dg_.setInputOffsets(offsets);
56 }
57
58 inline void setComputeMinSad(const bool data) {
59 this->ComputeMinSad = data;
60 }
61 inline void setComputeSadSad(const bool data) {
62 this->ComputeSadSad = data;
63 }
64 inline void setComputeSadMax(const bool data) {
65 this->ComputeSadMax = data;
66 }
67
68 template <typename triangulationType>
69 inline int buildGradient(const void *const scalars,
70 const size_t scalarsMTime,
71 const SimplexId *const offsets,
72 const triangulationType &triangulation) {
73 this->dg_.setDebugLevel(this->debugLevel_);
75 this->dg_.setInputOffsets(offsets);
76 this->dg_.setInputScalarField(scalars, scalarsMTime);
77 return this->dg_.buildGradient(triangulation);
78 }
79
92 this->dg_ = std::move(dg);
93 // reset gradient pointer to local storage
94 this->dg_.setLocalGradient();
95 }
97 return std::move(this->dg_);
98 }
99
100 template <typename triangulationType>
101 inline SimplexId
103 const triangulationType &triangulation) {
104 return this->dg_.getCellGreaterVertex(c, triangulation);
105 }
106
107 inline const std::vector<std::vector<SimplexId>> &
109 return this->s2Children_;
110 }
111
126 template <typename triangulationType>
127 int computePersistencePairs(std::vector<PersistencePair> &pairs,
128 const SimplexId *const offsets,
129 const triangulationType &triangulation,
130 const bool ignoreBoundary,
131 const bool compute2SaddlesChildren = false);
132
141 std::vector<SimplexId> boundary;
146 std::array<SimplexId, 2> critVertsIds;
147 };
148
149 protected:
159 template <typename triangulationType>
160 std::vector<std::vector<SimplexId>>
161 getSaddle1ToMinima(const std::vector<SimplexId> &criticalEdges,
162 const triangulationType &triangulation) const;
163
179 template <typename triangulationType,
180 typename GFS,
181 typename GFSN,
182 typename OB>
183 std::vector<std::vector<SimplexId>>
184 getSaddle2ToMaxima(const std::vector<SimplexId> &criticalCells,
185 const GFS &getFaceStar,
186 const GFSN &getFaceStarNumber,
187 const OB &isOnBoundary,
188 const triangulationType &triangulation) const;
189
201 template <typename triangulationType>
202 void getMinSaddlePairs(std::vector<PersistencePair> &pairs,
203 std::vector<bool> &pairedMinima,
204 std::vector<bool> &paired1Saddles,
205 const std::vector<SimplexId> &criticalEdges,
206 const std::vector<SimplexId> &critEdgesOrder,
207 const SimplexId *const offsets,
208 const triangulationType &triangulation) const;
209
221 template <typename triangulationType>
222 void getMaxSaddlePairs(std::vector<PersistencePair> &pairs,
223 std::vector<bool> &pairedMaxima,
224 std::vector<bool> &pairedSaddles,
225 const std::vector<SimplexId> &criticalSaddles,
226 const std::vector<SimplexId> &critSaddlesOrder,
227 const std::vector<SimplexId> &critMaxsOrder,
228 const triangulationType &triangulation) const;
229
243 template <typename triangulationType>
244 void getSaddleSaddlePairs(std::vector<PersistencePair> &pairs,
245 std::vector<bool> &paired1Saddles,
246 std::vector<bool> &paired2Saddles,
247 const bool exportBoundaries,
248 std::vector<GeneratorType> &boundaries,
249 const std::vector<SimplexId> &critical1Saddles,
250 const std::vector<SimplexId> &critical2Saddles,
251 const std::vector<SimplexId> &crit1SaddlesOrder,
252 const triangulationType &triangulation) const;
253
263 template <typename triangulationType>
265 std::array<std::vector<SimplexId>, 4> &criticalCellsByDim,
266 std::array<std::vector<SimplexId>, 4> &critCellsOrder,
267 const SimplexId *const offsets,
268 const triangulationType &triangulation,
269 const bool sortEdges) const;
270
282 void displayStats(
283 const std::vector<PersistencePair> &pairs,
284 const std::array<std::vector<SimplexId>, 4> &criticalCellsByDim,
285 const std::vector<bool> &pairedMinima,
286 const std::vector<bool> &paired1Saddles,
287 const std::vector<bool> &paired2Saddles,
288 const std::vector<bool> &pairedMaxima) const;
289
297 using tripletType = std::array<SimplexId, 3>;
298
311 void tripletsToPersistencePairs(std::vector<PersistencePair> &pairs,
312 std::vector<bool> &pairedExtrema,
313 std::vector<bool> &pairedSaddles,
314 std::vector<SimplexId> &reps,
315 std::vector<tripletType> &triplets,
316 const SimplexId *const saddlesOrder,
317 const SimplexId *const extremaOrder,
318 const SimplexId pairDim) const;
319
338 template <typename triangulationType, typename Container>
341 std::vector<bool> &onBoundary,
342 std::vector<Container> &s2Boundaries,
343 const std::vector<SimplexId> &s2Mapping,
344 const std::vector<SimplexId> &s1Mapping,
345 std::vector<SimplexId> &partners,
346 std::vector<Lock> &s1Locks,
347 std::vector<Lock> &s2Locks,
348 const triangulationType &triangulation) const;
349
355 template <size_t n>
356 struct Simplex {
361 std::array<SimplexId, n> vertsOrder_{};
364 friend bool operator<(const Simplex<n> &lhs, const Simplex<n> &rhs) {
365 return lhs.vertsOrder_ < rhs.vertsOrder_;
366 }
367 };
368
372 struct EdgeSimplex : Simplex<2> {
373 template <typename triangulationType>
374 void fillEdge(const SimplexId id,
375 const SimplexId *const offsets,
376 const triangulationType &triangulation) {
377 this->id_ = id;
378 triangulation.getEdgeVertex(id, 0, this->vertsOrder_[0]);
379 triangulation.getEdgeVertex(id, 1, this->vertsOrder_[1]);
380 this->vertsOrder_[0] = offsets[this->vertsOrder_[0]];
381 this->vertsOrder_[1] = offsets[this->vertsOrder_[1]];
382 // sort vertices in decreasing order
383 std::sort(this->vertsOrder_.rbegin(), this->vertsOrder_.rend());
384 }
385 };
386
391 template <typename triangulationType>
392 void fillTriangle(const SimplexId id,
393 const SimplexId *const offsets,
394 const triangulationType &triangulation) {
395 this->id_ = id;
396 triangulation.getTriangleVertex(id, 0, this->vertsOrder_[0]);
397 triangulation.getTriangleVertex(id, 1, this->vertsOrder_[1]);
398 triangulation.getTriangleVertex(id, 2, this->vertsOrder_[2]);
399 this->vertsOrder_[0] = offsets[this->vertsOrder_[0]];
400 this->vertsOrder_[1] = offsets[this->vertsOrder_[1]];
401 this->vertsOrder_[2] = offsets[this->vertsOrder_[2]];
402 // sort vertices in decreasing order
403 std::sort(this->vertsOrder_.rbegin(), this->vertsOrder_.rend());
404 }
405 };
406
410 struct TetraSimplex : Simplex<4> {
411 template <typename triangulationType>
412 void fillTetra(const SimplexId id,
413 const SimplexId *const offsets,
414 const triangulationType &triangulation) {
415 this->id_ = id;
416 triangulation.getCellVertex(id, 0, this->vertsOrder_[0]);
417 triangulation.getCellVertex(id, 1, this->vertsOrder_[1]);
418 triangulation.getCellVertex(id, 2, this->vertsOrder_[2]);
419 triangulation.getCellVertex(id, 3, this->vertsOrder_[3]);
420 this->vertsOrder_[0] = offsets[this->vertsOrder_[0]];
421 this->vertsOrder_[1] = offsets[this->vertsOrder_[1]];
422 this->vertsOrder_[2] = offsets[this->vertsOrder_[2]];
423 this->vertsOrder_[3] = offsets[this->vertsOrder_[3]];
424 // sort vertices in decreasing order
425 std::sort(this->vertsOrder_.rbegin(), this->vertsOrder_.rend());
426 }
427 };
428
429 template <typename triangulationType>
430 void alloc(const triangulationType &triangulation) {
431 Timer tm{};
432 const auto dim{this->dg_.getDimensionality()};
433 if(dim > 3 || dim < 1) {
434 return;
435 }
436#ifdef TTK_ENABLE_OPENMP
437#pragma omp parallel master num_threads(threadNumber_)
438#endif
439 {
440#ifdef TTK_ENABLE_OPENMP
441#pragma omp task
442#endif // TTK_ENABLE_OPENMP
443 this->firstRepMin_.resize(triangulation.getNumberOfVertices());
444 if(dim > 1) {
445#ifdef TTK_ENABLE_OPENMP
446#pragma omp task
447#endif
448 this->firstRepMax_.resize(triangulation.getNumberOfCells());
449 }
450 if(dim > 2) {
451#ifdef TTK_ENABLE_OPENMP
452#pragma omp task
453#endif
454 this->critEdges_.resize(triangulation.getNumberOfEdges());
455#ifdef TTK_ENABLE_OPENMP
456#pragma omp task
457#endif
458 this->edgeTrianglePartner_.resize(
459 triangulation.getNumberOfEdges(), -1);
460#ifdef TTK_ENABLE_OPENMP
461#pragma omp task
462#endif
463 this->onBoundary_.resize(triangulation.getNumberOfEdges(), false);
464#ifdef TTK_ENABLE_OPENMP
465#pragma omp task
466#endif
467 this->s2Mapping_.resize(triangulation.getNumberOfTriangles(), -1);
468#ifdef TTK_ENABLE_OPENMP
469#pragma omp task
470#endif
471 this->s1Mapping_.resize(triangulation.getNumberOfEdges(), -1);
472 }
473 for(int i = 0; i < dim + 1; ++i) {
474#ifdef TTK_ENABLE_OPENMP
475#pragma omp task
476#endif
477 this->pairedCritCells_[i].resize(
478 this->dg_.getNumberOfCells(i, triangulation), false);
479 }
480 for(int i = 1; i < dim + 1; ++i) {
481#ifdef TTK_ENABLE_OPENMP
482#pragma omp task
483#endif
484 this->critCellsOrder_[i].resize(
485 this->dg_.getNumberOfCells(i, triangulation), -1);
486 }
487 }
488 this->printMsg("Memory allocations", 1.0, tm.getElapsedTime(), 1,
490 }
491
492 void clear() {
493 Timer tm{};
494 this->firstRepMin_ = {};
495 this->firstRepMax_ = {};
496 this->edgeTrianglePartner_ = {};
497 this->s2Mapping_ = {};
498 this->s1Mapping_ = {};
499 this->critEdges_ = {};
500 this->pairedCritCells_ = {};
501 this->onBoundary_ = {};
502 this->critCellsOrder_ = {};
503 this->printMsg("Memory cleanup", 1.0, tm.getElapsedTime(), 1,
505 }
506
508
509 // factor memory allocations outside computation loops
510 mutable std::vector<SimplexId> firstRepMin_{}, firstRepMax_{},
512 mutable std::vector<EdgeSimplex> critEdges_{};
513 mutable std::array<std::vector<bool>, 4> pairedCritCells_{};
514 mutable std::vector<bool> onBoundary_{};
515 mutable std::array<std::vector<SimplexId>, 4> critCellsOrder_{};
516 mutable std::vector<std::vector<SimplexId>> s2Children_{};
517
518 bool ComputeMinSad{true};
519 bool ComputeSadSad{true};
520 bool ComputeSadMax{true};
522 };
523} // namespace ttk
524
525template <typename triangulationType>
526std::vector<std::vector<SimplexId>>
528 const std::vector<SimplexId> &criticalEdges,
529 const triangulationType &triangulation) const {
530
531 Timer tm{};
532
533 std::vector<std::vector<SimplexId>> res(criticalEdges.size());
534
535 // follow vpaths from 1-saddles to minima
536#ifdef TTK_ENABLE_OPENMP
537#pragma omp parallel for num_threads(threadNumber_)
538#endif
539 for(size_t i = 0; i < criticalEdges.size(); ++i) {
540 auto &mins = res[i];
541
542 const auto followVPath = [this, &mins, &triangulation](const SimplexId v) {
543 std::vector<Cell> vpath{};
544 this->dg_.getDescendingPath(Cell{0, v}, vpath, triangulation);
545 const Cell &lastCell = vpath.back();
546 if(lastCell.dim_ == 0 && this->dg_.isCellCritical(lastCell)) {
547 mins.emplace_back(lastCell.id_);
548 }
549 };
550
551 // critical edge vertices
552 SimplexId v0{}, v1{};
553 triangulation.getEdgeVertex(criticalEdges[i], 0, v0);
554 triangulation.getEdgeVertex(criticalEdges[i], 1, v1);
555
556 // follow vpath from each vertex of the critical edge
557 followVPath(v0);
558 followVPath(v1);
559 }
560
561 this->printMsg("Computed the descending 1-separatrices", 1.0,
562 tm.getElapsedTime(), this->threadNumber_, debug::LineMode::NEW,
564
565 return res;
566}
567
568template <typename triangulationType, typename GFS, typename GFSN, typename OB>
569std::vector<std::vector<SimplexId>>
571 const std::vector<SimplexId> &criticalCells,
572 const GFS &getFaceStar,
573 const GFSN &getFaceStarNumber,
574 const OB &isOnBoundary,
575 const triangulationType &triangulation) const {
576
577 Timer tm{};
578
579 const auto dim = this->dg_.getDimensionality();
580 std::vector<std::vector<SimplexId>> res(criticalCells.size());
581
582 // follow vpaths from 2-saddles to maxima
583#ifdef TTK_ENABLE_OPENMP
584#pragma omp parallel for num_threads(threadNumber_)
585#endif
586 for(size_t i = 0; i < criticalCells.size(); ++i) {
587 const auto sid = criticalCells[i];
588 auto &maxs = res[i];
589
590 const auto followVPath
591 = [this, dim, &maxs, &triangulation](const SimplexId v) {
592 std::vector<Cell> vpath{};
593 this->dg_.getAscendingPath(Cell{dim, v}, vpath, triangulation);
594 const Cell &lastCell = vpath.back();
595 if(lastCell.dim_ == dim && this->dg_.isCellCritical(lastCell)) {
596 maxs.emplace_back(lastCell.id_);
597 } else if(lastCell.dim_ == dim - 1) {
598 maxs.emplace_back(-1);
599 }
600 };
601
602 const auto starNumber = getFaceStarNumber(sid);
603
604 for(SimplexId j = 0; j < starNumber; ++j) {
605 SimplexId cellId{};
606 getFaceStar(sid, j, cellId);
607 followVPath(cellId);
608 }
609
610 if(isOnBoundary(sid)) {
611 // critical saddle is on boundary
612 maxs.emplace_back(-1);
613 }
614 }
615
616 this->printMsg("Computed the ascending 1-separatrices", 1.0,
617 tm.getElapsedTime(), this->threadNumber_, debug::LineMode::NEW,
619
620 return res;
621}
622
623template <typename triangulationType>
625 std::vector<PersistencePair> &pairs,
626 std::vector<bool> &pairedMinima,
627 std::vector<bool> &paired1Saddles,
628 const std::vector<SimplexId> &criticalEdges,
629 const std::vector<SimplexId> &critEdgesOrder,
630 const SimplexId *const offsets,
631 const triangulationType &triangulation) const {
632
633 Timer tm{};
634
635 auto saddle1ToMinima = getSaddle1ToMinima(criticalEdges, triangulation);
636
637 Timer tmseq{};
638
639 auto &firstRep{this->firstRepMin_};
640 std::iota(firstRep.begin(), firstRep.end(), 0);
641 std::vector<tripletType> sadMinTriplets{};
642
643 for(size_t i = 0; i < saddle1ToMinima.size(); ++i) {
644 auto &mins = saddle1ToMinima[i];
645 const auto s1 = criticalEdges[i];
646 // remove duplicates
647 TTK_PSORT(this->threadNumber_, mins.begin(), mins.end());
648 const auto last = std::unique(mins.begin(), mins.end());
649 mins.erase(last, mins.end());
650 if(mins.size() != 2) {
651 continue;
652 }
653 sadMinTriplets.emplace_back(tripletType{s1, mins[0], mins[1]});
654 }
655
656 tripletsToPersistencePairs(pairs, pairedMinima, paired1Saddles, firstRep,
657 sadMinTriplets, critEdgesOrder.data(), offsets, 0);
658
659 const auto nMinSadPairs = pairs.size();
660
661 this->printMsg(
662 "Computed " + std::to_string(nMinSadPairs) + " min-saddle pairs", 1.0,
663 tm.getElapsedTime(), this->threadNumber_);
664
665 this->printMsg("min-saddle pairs sequential part", 1.0,
666 tmseq.getElapsedTime(), 1, debug::LineMode::NEW,
668}
669
670template <typename triangulationType>
672 std::vector<PersistencePair> &pairs,
673 std::vector<bool> &pairedMaxima,
674 std::vector<bool> &pairedSaddles,
675 const std::vector<SimplexId> &criticalSaddles,
676 const std::vector<SimplexId> &critSaddlesOrder,
677 const std::vector<SimplexId> &critMaxsOrder,
678 const triangulationType &triangulation) const {
679
680 Timer tm{};
681
682 const auto dim = this->dg_.getDimensionality();
683
684 auto saddle2ToMaxima
685 = dim == 3
686 ? getSaddle2ToMaxima(
687 criticalSaddles,
688 [&triangulation](const SimplexId a, const SimplexId i, SimplexId &r) {
689 return triangulation.getTriangleStar(a, i, r);
690 },
691 [&triangulation](const SimplexId a) {
692 return triangulation.getTriangleStarNumber(a);
693 },
694 [&triangulation](const SimplexId a) {
695 return triangulation.isTriangleOnBoundary(a);
696 },
697 triangulation)
698 : getSaddle2ToMaxima(
699 criticalSaddles,
700 [&triangulation](const SimplexId a, const SimplexId i, SimplexId &r) {
701 return triangulation.getEdgeStar(a, i, r);
702 },
703 [&triangulation](const SimplexId a) {
704 return triangulation.getEdgeStarNumber(a);
705 },
706 [&triangulation](const SimplexId a) {
707 return triangulation.isEdgeOnBoundary(a);
708 },
709 triangulation);
710
711 Timer tmseq{};
712
713 auto &firstRep{this->firstRepMax_};
714 std::iota(firstRep.begin(), firstRep.end(), 0);
715 std::vector<tripletType> sadMaxTriplets{};
716
717 for(size_t i = 0; i < saddle2ToMaxima.size(); ++i) {
718 auto &maxs = saddle2ToMaxima[i];
719 // remove duplicates
720 TTK_PSORT(this->threadNumber_, maxs.begin(), maxs.end(),
721 [](const SimplexId a, const SimplexId b) {
722 // positive values (actual maxima) before negative ones
723 // (boundary component id)
724 if(a * b >= 0) {
725 return a < b;
726 } else {
727 return a > b;
728 }
729 });
730 const auto last = std::unique(maxs.begin(), maxs.end());
731 maxs.erase(last, maxs.end());
732
733 // remove "doughnut" configurations: two ascending separatrices
734 // leading to the same maximum/boundary component
735 if(maxs.size() != 2) {
736 continue;
737 }
738
739 const auto s2 = criticalSaddles[i];
740 if(!pairedSaddles[s2]) {
741 sadMaxTriplets.emplace_back(tripletType{s2, maxs[0], maxs[1]});
742 }
743 }
744
745 const auto nMinSadPairs = pairs.size();
746
747 tripletsToPersistencePairs(pairs, pairedMaxima, pairedSaddles, firstRep,
748 sadMaxTriplets, critSaddlesOrder.data(),
749 critMaxsOrder.data(), dim - 1);
750
751 const auto nSadMaxPairs = pairs.size() - nMinSadPairs;
752
753 this->printMsg(
754 "Computed " + std::to_string(nSadMaxPairs) + " saddle-max pairs", 1.0,
755 tm.getElapsedTime(), this->threadNumber_);
756
757 this->printMsg("saddle-max pairs sequential part", 1.0,
758 tmseq.getElapsedTime(), 1, debug::LineMode::NEW,
760}
761
762template <typename triangulationType, typename Container>
764 const SimplexId s2,
765 std::vector<bool> &onBoundary,
766 std::vector<Container> &s2Boundaries,
767 const std::vector<SimplexId> &s2Mapping,
768 const std::vector<SimplexId> &s1Mapping,
769 std::vector<SimplexId> &partners,
770 std::vector<Lock> &s1Locks,
771 std::vector<Lock> &s2Locks,
772 const triangulationType &triangulation) const {
773
774 auto &boundaryIds{s2Boundaries[s2Mapping[s2]]};
775
776 const auto addBoundary = [&boundaryIds, &onBoundary](const SimplexId e) {
777 // add edge e to boundaryIds/onBoundary modulo 2
778 if(!onBoundary[e]) {
779 boundaryIds.emplace(e);
780 onBoundary[e] = true;
781 } else {
782 const auto it = boundaryIds.find(e);
783 boundaryIds.erase(it);
784 onBoundary[e] = false;
785 }
786 };
787
788 const auto clearOnBoundary = [&boundaryIds, &onBoundary]() {
789 // clear the onBoundary vector (set everything to false)
790 for(const auto e : boundaryIds) {
791 onBoundary[e] = false;
792 }
793 };
794
795 if(!boundaryIds.empty()) {
796 // restore previously computed s2 boundary
797 for(const auto e : boundaryIds) {
798 onBoundary[e] = true;
799 }
800 } else {
801 // init cascade with s2 triangle boundary (3 edges)
802 for(SimplexId i = 0; i < 3; ++i) {
803 SimplexId e{};
804 triangulation.getTriangleEdge(s2, i, e);
805 addBoundary(e);
806 }
807 }
808
809 // lock the 2-saddle to ensure that only one thread can perform the
810 // boundary expansion
811 s2Locks[s2Mapping[s2]].lock();
812
813 while(!boundaryIds.empty()) {
814 // tau: youngest edge on boundary
815 const auto tau{*boundaryIds.begin()};
816 // use the Discrete Gradient to find a triangle paired to tau
817 auto pTau{this->dg_.getPairedCell(Cell{1, tau}, triangulation)};
818 bool critical{false};
819 if(pTau == -1) {
820 // maybe tau is critical and paired to a critical triangle
821 do {
822#ifdef TTK_ENABLE_OPENMP
823#pragma omp atomic read
824#endif // TTK_ENABLE_OPENMP
825 pTau = partners[tau];
826 if(pTau == -1 || s2Boundaries[s2Mapping[pTau]].empty()) {
827 break;
828 }
829 } while(*s2Boundaries[s2Mapping[pTau]].begin() != tau);
830
831 critical = true;
832 }
833 if(pTau == -1) {
834 // tau is critical and not paired
835
836 // compare-and-swap from "Towards Lockfree Persistent Homology"
837 // using locks over 1-saddles instead of atomics (OpenMP compatibility)
838 s1Locks[s1Mapping[tau]].lock();
839 const auto cap = partners[tau];
840 if(partners[tau] == -1) {
841 partners[tau] = s2;
842 }
843 s1Locks[s1Mapping[tau]].unlock();
844
845 // cleanup before exiting
846 clearOnBoundary();
847 s2Locks[s2Mapping[s2]].unlock();
848 if(cap == -1) {
849 return tau;
850 } else {
851 return this->eliminateBoundariesSandwich(
852 s2, onBoundary, s2Boundaries, s2Mapping, s1Mapping, partners, s1Locks,
853 s2Locks, triangulation);
854 }
855
856 } else {
857 // expand boundary
858 if(critical && s2Mapping[pTau] != -1) {
859 if(s2Mapping[pTau] < s2Mapping[s2]) {
860 // pTau is an already-paired 2-saddle
861 // merge pTau boundary into s2 boundary
862
863 // make sure that pTau boundary is not modified by another
864 // thread while we merge the two boundaries...
865 s2Locks[s2Mapping[pTau]].lock();
866 for(const auto e : s2Boundaries[s2Mapping[pTau]]) {
867 addBoundary(e);
868 }
869 s2Locks[s2Mapping[pTau]].unlock();
870 if(this->Compute2SaddlesChildren) {
871 this->s2Children_[s2Mapping[s2]].emplace_back(s2Mapping[pTau]);
872 }
873
874 } else if(s2Mapping[pTau] > s2Mapping[s2]) {
875
876 // compare-and-swap from "Towards Lockfree Persistent
877 // Homology" using locks over 1-saddles
878 s1Locks[s1Mapping[tau]].lock();
879 const auto cap = partners[tau];
880 if(partners[tau] == pTau) {
881 partners[tau] = s2;
882 }
883 s1Locks[s1Mapping[tau]].unlock();
884
885 if(cap == pTau) {
886 // cleanup before exiting
887 clearOnBoundary();
888 s2Locks[s2Mapping[s2]].unlock();
889 return this->eliminateBoundariesSandwich(
890 pTau, onBoundary, s2Boundaries, s2Mapping, s1Mapping, partners,
891 s1Locks, s2Locks, triangulation);
892 }
893 }
894 } else { // pTau is a regular triangle
895 // add pTau triangle boundary (3 edges)
896 for(SimplexId i = 0; i < 3; ++i) {
897 SimplexId e{};
898 triangulation.getTriangleEdge(pTau, i, e);
899 addBoundary(e);
900 }
901 }
902 }
903 }
904
905 // cleanup before exiting
906 clearOnBoundary();
907 s2Locks[s2Mapping[s2]].unlock();
908 return -1;
909}
910
911template <typename triangulationType>
913 std::vector<PersistencePair> &pairs,
914 std::vector<bool> &paired1Saddles,
915 std::vector<bool> &paired2Saddles,
916 const bool exportBoundaries,
917 std::vector<GeneratorType> &boundaries,
918 const std::vector<SimplexId> &critical1Saddles,
919 const std::vector<SimplexId> &critical2Saddles,
920 const std::vector<SimplexId> &crit1SaddlesOrder,
921 const triangulationType &triangulation) const {
922
923 Timer tm2{};
924 const auto nSadExtrPairs = pairs.size();
925
926 // 1- and 2-saddles yet to be paired
927 std::vector<SimplexId> saddles1{}, saddles2{};
928 // filter out already paired 1-saddles (edge id)
929 for(const auto s1 : critical1Saddles) {
930 if(!paired1Saddles[s1]) {
931 saddles1.emplace_back(s1);
932 }
933 }
934 // filter out already paired 2-saddles (triangle id)
935 for(const auto s2 : critical2Saddles) {
936 if(!paired2Saddles[s2]) {
937 saddles2.emplace_back(s2);
938 }
939 }
940
941 if(this->Compute2SaddlesChildren) {
942 this->s2Children_.resize(saddles2.size());
943 }
944
945 // sort every triangulation edges by filtration order
946 const auto &edgesFiltrOrder{crit1SaddlesOrder};
947
948 auto &onBoundary{this->onBoundary_};
949 auto &edgeTrianglePartner{this->edgeTrianglePartner_};
950
951 const auto cmpEdges
952 = [&edgesFiltrOrder](const SimplexId a, const SimplexId b) {
953 return edgesFiltrOrder[a] > edgesFiltrOrder[b];
954 };
955 using Container = std::set<SimplexId, decltype(cmpEdges)>;
956 std::vector<Container> s2Boundaries(saddles2.size(), Container(cmpEdges));
957
958 // unpaired critical triangle id -> index in saddle2 vector
959 auto &s2Mapping{this->s2Mapping_};
960#ifdef TTK_ENABLE_OPENMP
961#pragma omp parallel for num_threads(threadNumber_)
962#endif // TTK_ENABLE_OPENMP
963 for(size_t i = 0; i < saddles2.size(); ++i) {
964 s2Mapping[saddles2[i]] = i;
965 }
966
967 // unpaired critical edge id -> index in saddle1 vector
968 auto &s1Mapping{this->s1Mapping_};
969#ifdef TTK_ENABLE_OPENMP
970#pragma omp parallel for num_threads(threadNumber_)
971#endif // TTK_ENABLE_OPENMP
972 for(size_t i = 0; i < saddles1.size(); ++i) {
973 s1Mapping[saddles1[i]] = i;
974 }
975
976 // one lock per 1-saddle
977 std::vector<Lock> s1Locks(saddles1.size());
978 // one lock per 2-saddle
979 std::vector<Lock> s2Locks(saddles2.size());
980
981 // compute 2-saddles boundaries in parallel
982
983#ifdef TTK_ENABLE_OPENMP4
984#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
985 firstprivate(onBoundary)
986#endif // TTK_ENABLE_OPENMP4
987 for(size_t i = 0; i < saddles2.size(); ++i) {
988 // 2-saddles sorted in increasing order
989 const auto s2 = saddles2[i];
990 this->eliminateBoundariesSandwich(s2, onBoundary, s2Boundaries, s2Mapping,
991 s1Mapping, edgeTrianglePartner, s1Locks,
992 s2Locks, triangulation);
993 }
994
995 Timer tmseq{};
996
997 // extract saddle-saddle pairs from computed boundaries
998 for(size_t i = 0; i < saddles2.size(); ++i) {
999 if(!s2Boundaries[i].empty()) {
1000 const auto s2 = saddles2[i];
1001 const auto s1 = *s2Boundaries[i].begin();
1002 // we found a pair
1003 pairs.emplace_back(s1, s2, 1);
1004 paired1Saddles[s1] = true;
1005 paired2Saddles[s2] = true;
1006 }
1007 }
1008
1009 if(exportBoundaries) {
1010 boundaries.resize(s2Boundaries.size());
1011 for(size_t i = 0; i < boundaries.size(); ++i) {
1012 const auto &boundSet{s2Boundaries[i]};
1013 if(boundSet.empty()) {
1014 continue;
1015 }
1016 boundaries[i] = {
1017 {boundSet.begin(), boundSet.end()},
1018 saddles2[i],
1019 std::array<SimplexId, 2>{
1020 this->dg_.getCellGreaterVertex(Cell{2, saddles2[i]}, triangulation),
1021 this->dg_.getCellGreaterVertex(
1022 Cell{1, *boundSet.begin()}, triangulation),
1023 }};
1024 }
1025 }
1026
1027 const auto nSadSadPairs = pairs.size() - nSadExtrPairs;
1028
1029 this->printMsg(
1030 "Computed " + std::to_string(nSadSadPairs) + " saddle-saddle pairs", 1.0,
1031 tm2.getElapsedTime(), this->threadNumber_);
1032
1033 this->printMsg("saddle-saddle pairs sequential part", 1.0,
1034 tmseq.getElapsedTime(), 1, debug::LineMode::NEW,
1036}
1037
1038template <typename triangulationType>
1040 std::array<std::vector<SimplexId>, 4> &criticalCellsByDim,
1041 std::array<std::vector<SimplexId>, 4> &critCellsOrder,
1042 const SimplexId *const offsets,
1043 const triangulationType &triangulation,
1044 const bool sortEdges) const {
1045
1046 Timer tm{};
1047
1048 this->dg_.getCriticalPoints(criticalCellsByDim, triangulation);
1049
1050 this->printMsg("Extracted critical cells", 1.0, tm.getElapsedTime(),
1051 this->threadNumber_, debug::LineMode::NEW,
1053
1054 // memory allocations
1055 auto &critEdges{this->critEdges_};
1056 if(!sortEdges) {
1057 critEdges.resize(criticalCellsByDim[1].size());
1058 }
1059 std::vector<TriangleSimplex> critTriangles(criticalCellsByDim[2].size());
1060 std::vector<TetraSimplex> critTetras(criticalCellsByDim[3].size());
1061
1062#ifdef TTK_ENABLE_OPENMP
1063#pragma omp parallel num_threads(threadNumber_)
1064#endif // TTK_ENABLE_OPENMP
1065 {
1066 if(sortEdges) {
1067#ifdef TTK_ENABLE_OPENMP
1068#pragma omp for nowait
1069#endif // TTK_ENABLE_OPENMP
1070 for(size_t i = 0; i < critEdges.size(); ++i) {
1071 critEdges[i].fillEdge(i, offsets, triangulation);
1072 }
1073 } else {
1074#ifdef TTK_ENABLE_OPENMP
1075#pragma omp for nowait
1076#endif // TTK_ENABLE_OPENMP
1077 for(size_t i = 0; i < critEdges.size(); ++i) {
1078 critEdges[i].fillEdge(criticalCellsByDim[1][i], offsets, triangulation);
1079 }
1080 }
1081
1082#ifdef TTK_ENABLE_OPENMP
1083#pragma omp for nowait
1084#endif // TTK_ENABLE_OPENMP
1085 for(size_t i = 0; i < critTriangles.size(); ++i) {
1086 critTriangles[i].fillTriangle(
1087 criticalCellsByDim[2][i], offsets, triangulation);
1088 }
1089
1090#ifdef TTK_ENABLE_OPENMP
1091#pragma omp for
1092#endif // TTK_ENABLE_OPENMP
1093 for(size_t i = 0; i < critTetras.size(); ++i) {
1094 critTetras[i].fillTetra(criticalCellsByDim[3][i], offsets, triangulation);
1095 }
1096 }
1097
1098 TTK_PSORT(this->threadNumber_, critEdges.begin(), critEdges.end());
1099 TTK_PSORT(this->threadNumber_, critTriangles.begin(), critTriangles.end());
1100 TTK_PSORT(this->threadNumber_, critTetras.begin(), critTetras.end());
1101
1102#ifdef TTK_ENABLE_OPENMP
1103#pragma omp parallel num_threads(threadNumber_)
1104#endif // TTK_ENABLE_OPENMP
1105 {
1106#ifdef TTK_ENABLE_OPENMP
1107#pragma omp for nowait
1108#endif // TTK_ENABLE_OPENMP
1109 for(size_t i = 0; i < critEdges.size(); ++i) {
1110 critCellsOrder[1][critEdges[i].id_] = i;
1111 }
1112
1113#ifdef TTK_ENABLE_OPENMP
1114#pragma omp for nowait
1115#endif // TTK_ENABLE_OPENMP
1116 for(size_t i = 0; i < critTriangles.size(); ++i) {
1117 criticalCellsByDim[2][i] = critTriangles[i].id_;
1118 critCellsOrder[2][critTriangles[i].id_] = i;
1119 }
1120
1121#ifdef TTK_ENABLE_OPENMP
1122#pragma omp for
1123#endif // TTK_ENABLE_OPENMP
1124 for(size_t i = 0; i < critTetras.size(); ++i) {
1125 criticalCellsByDim[3][i] = critTetras[i].id_;
1126 critCellsOrder[3][critTetras[i].id_] = i;
1127 }
1128 }
1129
1130 if(sortEdges) {
1131 TTK_PSORT(this->threadNumber_, criticalCellsByDim[1].begin(),
1132 criticalCellsByDim[1].end(),
1133 [&critCellsOrder](const SimplexId a, const SimplexId b) {
1134 return critCellsOrder[1][a] < critCellsOrder[1][b];
1135 });
1136 } else {
1137#ifdef TTK_ENABLE_OPENMP
1138#pragma omp parallel for num_threads(threadNumber_)
1139#endif // TTK_ENABLE_OPENMP
1140 for(size_t i = 0; i < critEdges.size(); ++i) {
1141 criticalCellsByDim[1][i] = critEdges[i].id_;
1142 }
1143 }
1144
1145 this->printMsg("Extracted & sorted critical cells", 1.0, tm.getElapsedTime(),
1146 this->threadNumber_, debug::LineMode::NEW,
1148}
1149
1150template <typename triangulationType>
1152 std::vector<PersistencePair> &pairs,
1153 const SimplexId *const offsets,
1154 const triangulationType &triangulation,
1155 const bool ignoreBoundary,
1156 const bool compute2SaddlesChildren) {
1157
1158 // allocate memory
1159 this->alloc(triangulation);
1160
1161 Timer tm{};
1162 pairs.clear();
1163 const auto dim = this->dg_.getDimensionality();
1164 this->Compute2SaddlesChildren = compute2SaddlesChildren;
1165
1166 // get every critical cell sorted them by dimension
1167 std::array<std::vector<SimplexId>, 4> criticalCellsByDim{};
1168 // holds the critical cells order
1169 auto &critCellsOrder{this->critCellsOrder_};
1170
1171 this->extractCriticalCells(
1172 criticalCellsByDim, critCellsOrder, offsets, triangulation, dim == 3);
1173
1174 // if minima are paired
1175 auto &pairedMinima{this->pairedCritCells_[0]};
1176 // if 1-saddles are paired
1177 auto &paired1Saddles{this->pairedCritCells_[1]};
1178 // if 2-saddles are paired
1179 auto &paired2Saddles{this->pairedCritCells_[dim - 1]};
1180 // if maxima are paired
1181 auto &pairedMaxima{this->pairedCritCells_[dim]};
1182
1183 // connected components (global min/max pair)
1184 size_t nConnComp{};
1185
1186 if(this->ComputeMinSad) {
1187 // minima - saddle pairs
1188 this->getMinSaddlePairs(pairs, pairedMinima, paired1Saddles,
1189 criticalCellsByDim[1], critCellsOrder[1], offsets,
1190 triangulation);
1191
1192 // non-paired minima
1193 for(const auto min : criticalCellsByDim[0]) {
1194 if(!pairedMinima[min]) {
1195 pairs.emplace_back(min, -1, 0);
1196 pairedMinima[min] = true;
1197 nConnComp++;
1198 }
1199 }
1200 } else {
1201 // still extract the global pair
1202 const auto globMin{*std::min_element(
1203 criticalCellsByDim[0].begin(), criticalCellsByDim[0].end(),
1204 [offsets](const SimplexId a, const SimplexId b) {
1205 return offsets[a] < offsets[b];
1206 })};
1207 pairs.emplace_back(globMin, -1, 0);
1208 pairedMinima[globMin] = true;
1209 nConnComp++;
1210 }
1211
1212 if(dim > 1 && this->ComputeSadMax) {
1213 // saddle - maxima pairs
1214 this->getMaxSaddlePairs(
1215 pairs, pairedMaxima, paired2Saddles, criticalCellsByDim[dim - 1],
1216 critCellsOrder[dim - 1], critCellsOrder[dim], triangulation);
1217 }
1218
1219 if(ignoreBoundary) {
1220 // post-process saddle-max pairs: remove the one with the global
1221 // maximum (if it exists) to be (more) compatible with FTM
1222 const auto it
1223 = std::find_if(pairs.begin(), pairs.end(), [&](const PersistencePair &p) {
1224 if(p.type < dim - 1) {
1225 return false;
1226 }
1227 const Cell cmax{dim, p.death};
1228 const auto vmax{this->getCellGreaterVertex(cmax, triangulation)};
1229 return offsets[vmax] == triangulation.getNumberOfVertices() - 1;
1230 });
1231
1232 if(it != pairs.end()) {
1233 // remove saddle-max pair with global maximum
1234 paired2Saddles[it->birth] = false;
1235 pairedMaxima[it->death] = false;
1236 pairs.erase(it);
1237 }
1238 }
1239
1240 // saddle - saddle pairs
1241 if(dim == 3 && !criticalCellsByDim[1].empty()
1242 && !criticalCellsByDim[2].empty() && this->ComputeSadSad) {
1243 std::vector<GeneratorType> tmp{};
1244 this->getSaddleSaddlePairs(
1245 pairs, paired1Saddles, paired2Saddles, false, tmp, criticalCellsByDim[1],
1246 criticalCellsByDim[2], critCellsOrder[1], triangulation);
1247 }
1248
1249 if(std::is_same<triangulationType, ttk::ExplicitTriangulation>::value) {
1250 // create infinite pairs from non-paired 1-saddles, 2-saddles and maxima
1251 size_t nHandles{}, nCavities{}, nNonPairedMax{};
1252 if((dim == 2 && !ignoreBoundary && this->ComputeMinSad
1253 && this->ComputeSadMax)
1254 || (dim == 3 && this->ComputeMinSad && this->ComputeSadSad)) {
1255 // non-paired 1-saddles
1256 for(const auto s1 : criticalCellsByDim[1]) {
1257 if(!paired1Saddles[s1]) {
1258 paired1Saddles[s1] = true;
1259 // topological handles
1260 pairs.emplace_back(s1, -1, 1);
1261 nHandles++;
1262 }
1263 }
1264 }
1265 if(dim == 3 && !ignoreBoundary && this->ComputeSadMax
1266 && this->ComputeSadSad) {
1267 // non-paired 2-saddles
1268 for(const auto s2 : criticalCellsByDim[2]) {
1269 if(!paired2Saddles[s2]) {
1270 paired2Saddles[s2] = true;
1271 // cavities
1272 pairs.emplace_back(s2, -1, 2);
1273 nCavities++;
1274 }
1275 }
1276 }
1277 if(dim == 2 && !ignoreBoundary && this->ComputeSadMax) {
1278 // non-paired maxima
1279 for(const auto max : criticalCellsByDim[dim]) {
1280 if(!pairedMaxima[max]) {
1281 pairs.emplace_back(max, -1, 2);
1282 nNonPairedMax++;
1283 }
1284 }
1285 }
1286
1287 int nBoundComp
1288 = (dim == 3 ? nCavities : nHandles) + nConnComp - nNonPairedMax;
1289 nBoundComp = std::max(nBoundComp, 0);
1290
1291 // print Betti numbers
1292 const std::vector<std::vector<std::string>> rows{
1293 {" #Connected components", std::to_string(nConnComp)},
1294 {" #Topological handles", std::to_string(nHandles)},
1295 {" #Cavities", std::to_string(nCavities)},
1296 {" #Boundary components", std::to_string(nBoundComp)},
1297 };
1298
1299 this->printMsg(rows, debug::Priority::DETAIL);
1300 }
1301
1302 this->printMsg(
1303 "Computed " + std::to_string(pairs.size()) + " persistence pairs", 1.0,
1304 tm.getElapsedTime(), this->threadNumber_);
1305
1306 this->displayStats(pairs, criticalCellsByDim, pairedMinima, paired1Saddles,
1307 paired2Saddles, pairedMaxima);
1308
1309 // free memory
1310 this->clear();
1311
1312 return 0;
1313}
#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
Minimalist debugging class.
Definition Debug.h:88
int debugLevel_
Definition Debug.h:379
virtual int setDebugLevel(const int &debugLevel)
Definition Debug.cpp:147
TTK DiscreteMorseSandwich processing package.
void setComputeSadMax(const bool data)
std::array< SimplexId, 3 > tripletType
Triplet type for persistence pairs.
std::vector< SimplexId > s1Mapping_
ttk::dcg::DiscreteGradient && getGradient()
std::vector< SimplexId > firstRepMin_
std::vector< std::vector< SimplexId > > s2Children_
const std::vector< std::vector< SimplexId > > & get2SaddlesChildren() const
std::array< std::vector< SimplexId >, 4 > critCellsOrder_
SimplexId getCellGreaterVertex(const dcg::Cell &c, const triangulationType &triangulation)
void alloc(const triangulationType &triangulation)
std::vector< EdgeSimplex > critEdges_
std::vector< SimplexId > edgeTrianglePartner_
void setGradient(ttk::dcg::DiscreteGradient &&dg)
Ugly hack to avoid a call to buildGradient()
void setComputeMinSad(const bool data)
int buildGradient(const void *const scalars, const size_t scalarsMTime, const SimplexId *const offsets, const triangulationType &triangulation)
void setInputOffsets(const SimplexId *const offsets)
void preconditionTriangulation(AbstractTriangulation *const data)
std::vector< SimplexId > firstRepMax_
void setComputeSadSad(const bool data)
void getMinSaddlePairs(std::vector< PersistencePair > &pairs, std::vector< bool > &pairedMinima, std::vector< bool > &paired1Saddles, const std::vector< SimplexId > &criticalEdges, const std::vector< SimplexId > &critEdgesOrder, const SimplexId *const offsets, const triangulationType &triangulation) const
Compute the pairs of dimension 0.
void extractCriticalCells(std::array< std::vector< SimplexId >, 4 > &criticalCellsByDim, std::array< std::vector< SimplexId >, 4 > &critCellsOrder, const SimplexId *const offsets, const triangulationType &triangulation, const bool sortEdges) const
Extract & sort critical cell from the DiscreteGradient.
void displayStats(const std::vector< PersistencePair > &pairs, const std::array< std::vector< SimplexId >, 4 > &criticalCellsByDim, const std::vector< bool > &pairedMinima, const std::vector< bool > &paired1Saddles, const std::vector< bool > &paired2Saddles, const std::vector< bool > &pairedMaxima) const
Print number of pairs, critical cells per dimension & unpaired cells.
std::vector< SimplexId > s2Mapping_
void getSaddleSaddlePairs(std::vector< PersistencePair > &pairs, std::vector< bool > &paired1Saddles, std::vector< bool > &paired2Saddles, const bool exportBoundaries, std::vector< GeneratorType > &boundaries, const std::vector< SimplexId > &critical1Saddles, const std::vector< SimplexId > &critical2Saddles, const std::vector< SimplexId > &crit1SaddlesOrder, const triangulationType &triangulation) const
Compute the saddle-saddle pairs (in 3D)
int computePersistencePairs(std::vector< PersistencePair > &pairs, const SimplexId *const offsets, const triangulationType &triangulation, const bool ignoreBoundary, const bool compute2SaddlesChildren=false)
Compute the persistence pairs from the discrete gradient.
std::vector< std::vector< SimplexId > > getSaddle2ToMaxima(const std::vector< SimplexId > &criticalCells, const GFS &getFaceStar, const GFSN &getFaceStarNumber, const OB &isOnBoundary, const triangulationType &triangulation) const
Follow the ascending 1-separatrices to compute the saddles -> maxima association.
std::array< std::vector< bool >, 4 > pairedCritCells_
SimplexId eliminateBoundariesSandwich(const SimplexId s2, std::vector< bool > &onBoundary, std::vector< Container > &s2Boundaries, const std::vector< SimplexId > &s2Mapping, const std::vector< SimplexId > &s1Mapping, std::vector< SimplexId > &partners, std::vector< Lock > &s1Locks, std::vector< Lock > &s2Locks, const triangulationType &triangulation) const
Detect 1-saddles paired to a given 2-saddle.
void getMaxSaddlePairs(std::vector< PersistencePair > &pairs, std::vector< bool > &pairedMaxima, std::vector< bool > &pairedSaddles, const std::vector< SimplexId > &criticalSaddles, const std::vector< SimplexId > &critSaddlesOrder, const std::vector< SimplexId > &critMaxsOrder, const triangulationType &triangulation) const
Compute the pairs of dimension dim - 1.
void tripletsToPersistencePairs(std::vector< PersistencePair > &pairs, std::vector< bool > &pairedExtrema, std::vector< bool > &pairedSaddles, std::vector< SimplexId > &reps, std::vector< tripletType > &triplets, const SimplexId *const saddlesOrder, const SimplexId *const extremaOrder, const SimplexId pairDim) const
Compute persistence pairs from triplets.
std::vector< std::vector< SimplexId > > getSaddle1ToMinima(const std::vector< SimplexId > &criticalEdges, const triangulationType &triangulation) const
Follow the descending 1-separatrices to compute the saddles -> minima association.
TTK discreteGradient processing package.
void setInputOffsets(const SimplexId *const data)
void setInputScalarField(const void *const data, const size_t mTime)
SimplexId getNumberOfCells(const int dimension, const triangulationType &triangulation) const
void setLocalGradient()
Use local storage instead of cache.
SimplexId getCellGreaterVertex(const Cell c, const triangulationType &triangulation) const
void preconditionTriangulation(AbstractTriangulation *const data)
int buildGradient(const triangulationType &triangulation, bool bypassCache=false)
int getDescendingPath(const Cell &cell, std::vector< Cell > &vpath, const triangulationType &triangulation) const
std::string to_string(__int128)
Definition ripserpy.cpp:99
The Topology ToolKit.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
T end(std::pair< T, T > &p)
Definition ripserpy.cpp:472
T begin(std::pair< T, T > &p)
Definition ripserpy.cpp:468
void fillEdge(const SimplexId id, const SimplexId *const offsets, const triangulationType &triangulation)
Type for exporting persistent generators.
Persistence pair struct as exported by DiscreteGradient.
PersistencePair(SimplexId b, SimplexId d, int t)
Ad-hoc struct for sorting simplices.
friend bool operator<(const Simplex< n > &lhs, const Simplex< n > &rhs)
void fillTetra(const SimplexId id, const SimplexId *const offsets, const triangulationType &triangulation)
void fillTriangle(const SimplexId id, const SimplexId *const offsets, const triangulationType &triangulation)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)