TTK
Loading...
Searching...
No Matches
FTMTree_MT.h
Go to the documentation of this file.
1
14
15#pragma once
16
17#include <functional>
18#include <map>
19#include <numeric>
20#include <queue>
21#include <set>
22#include <vector>
23
24#include <Geometry.h>
25#include <Triangulation.h>
26#include <Wrapper.h>
27
28#include "FTMAtomicUF.h"
29#include "FTMAtomicVector.h"
30#include "FTMDataTypes.h"
31#include "FTMNode.h"
32#include "FTMStructures.h"
33#include "FTMSuperArc.h"
34
35static ttk::Timer _launchGlobalTime;
36
37namespace ttk {
38 namespace ftm {
39 using UF = AtomicUF *;
40
41 // Tree data ( 1 per tree )
42 struct TreeData {
44
45 // components : tree / nodes / extrema
46 std::shared_ptr<FTMAtomicVector<SuperArc>> superArcs;
47 std::shared_ptr<FTMAtomicVector<Node>> nodes;
48 std::shared_ptr<FTMAtomicVector<idNode>> roots;
49 std::vector<idNode> leaves;
50
51 // vertex 2 node / superarc
52 std::vector<idCorresp> vert2tree;
53 std::vector<SimplexId> visitOrder;
54 std::vector<std::list<std::vector<SimplexId>>> trunkSegments;
55
56 // Track information
57 std::vector<AtomicUF> storage;
58 std::vector<UF> ufs;
59 std::vector<UF> propagation;
60 std::shared_ptr<FTMAtomicVector<CurrentState>> states;
61 // valences
62 std::vector<valence> valences;
63 // opened nodes
64 std::vector<char> openedNodes;
65
66 // current nb of tasks
68
69 // Segmentation, stay empty for Contour tree as
70 // they are created by Merge Tree
72
73#ifdef TTK_ENABLE_FTM_TREE_STATS_TIME
74 std::vector<ActiveTask> activeTasksStats;
75#endif
76
77#ifdef TTK_ENABLE_OMP_PRIORITY
78 // Is this MT to be computed with greater task priority than others
79 bool prior = false;
80#endif
81 };
82
83 class FTMTree_MT : virtual public Debug {
84 protected:
85 // global
86 std::shared_ptr<Params> params_;
87 std::shared_ptr<Scalars> scalars_;
88
89 // local
92
93 public:
94 // -----------
95 // CONSTRUCT
96 // -----------
97
98 // Tree with global data and partition number
99 FTMTree_MT(const std::shared_ptr<Params> &params,
100 const std::shared_ptr<Scalars> &scalars,
101 TreeType type);
102
103 ~FTMTree_MT() override;
104
105 void clear();
106
107 // --------------------
108 // Init
109 // --------------------
110
111 inline void setParamsScalars(const std::shared_ptr<Params> &params,
112 const std::shared_ptr<Scalars> &scalars) {
113 this->scalars_ = scalars;
114 this->params_ = params;
115 this->mt_data_.treeType = params->treeType;
116 }
117
118 template <class triangulationType>
119 void initNbScalars(const triangulationType *triangulation) {
120 scalars_->size = triangulation->getNumberOfVertices();
121 }
122
123 void initComp() {
124 if(isST()) {
125 comp_.vertLower
126 = [this](const SimplexId a, const SimplexId b) -> bool {
127 return this->scalars_->isHigher(a, b);
128 };
129 comp_.vertHigher
130 = [this](const SimplexId a, const SimplexId b) -> bool {
131 return this->scalars_->isLower(a, b);
132 };
133 } else {
134 comp_.vertLower
135 = [this](const SimplexId a, const SimplexId b) -> bool {
136 return this->scalars_->isLower(a, b);
137 };
138 comp_.vertHigher
139 = [this](const SimplexId a, const SimplexId b) -> bool {
140 return this->scalars_->isHigher(a, b);
141 };
142 }
143 }
144
145 bool compLower(const SimplexId a, const SimplexId b) {
146 return comp_.vertLower(a, b);
147 }
148
150 template <typename scalarType>
151 void sortInput();
152
154 void makeAlloc() {
156
157 // Stats alloc
158
160 mt_data_.nodes->reserve(scalars_->size / 2);
161
163 mt_data_.roots->reserve(10);
164
166 mt_data_.leaves.reserve(scalars_->size / 3);
167
168 // Known size
169
171 mt_data_.vert2tree.resize(scalars_->size);
172
174
176 mt_data_.visitOrder.resize(scalars_->size);
177
179 mt_data_.ufs.resize(scalars_->size);
180
181 createVector<UF>(mt_data_.propagation);
182 mt_data_.propagation.resize(scalars_->size);
183
185 mt_data_.valences.resize(scalars_->size);
186
187 createVector<char>(mt_data_.openedNodes);
188 mt_data_.openedNodes.resize(scalars_->size);
189
190 mt_data_.segments_.clear();
191 }
192
193 void makeInit() {
194 initVector<idCorresp>(mt_data_.vert2tree, nullCorresp);
195 initVector<SimplexId>(mt_data_.visitOrder, nullVertex);
196 initVector<UF>(mt_data_.ufs, nullptr);
197 initVector<UF>(mt_data_.propagation, nullptr);
198 initVector<valence>(mt_data_.valences, 0);
199 initVector<char>(mt_data_.openedNodes, 0);
200 }
201
202 void initVectStates(const SimplexId nbLeaves) {
203 if(!mt_data_.states) {
204 mt_data_.states = std::make_shared<FTMAtomicVector<CurrentState>>(
205 nbLeaves, comp_.vertHigher);
206 }
207 mt_data_.states->clear();
208 mt_data_.states->reserve(nbLeaves);
209 }
210
211 // -------------------
212 // Process
213 // -------------------
214
216 template <class triangulationType>
217 void build(const triangulationType *mesh, const bool ct);
218
219 // extrema
220
221 template <class triangulationType>
222 int leafSearch(const triangulationType *mesh);
223
224 // skeleton
225
226 template <class triangulationType>
227 void leafGrowth(const triangulationType *mesh);
228
229 template <class triangulationType>
230 void arcGrowth(const triangulationType *mesh,
231 const SimplexId startVert,
232 const SimplexId orig);
233
234 template <class triangulationType>
235 std::tuple<bool, bool> propagate(const triangulationType *mesh,
236 CurrentState &currentState,
237 UF curUF);
238
239 template <class triangulationType>
240 void closeAndMergeOnSaddle(const triangulationType *mesh,
241 SimplexId saddleVert);
242
243 template <class triangulationType>
244 void closeOnBackBone(const triangulationType *mesh, SimplexId saddleVert);
245
246 void closeArcsUF(idNode closeNode, UF uf);
247
248 template <class triangulationType>
249 SimplexId trunk(const triangulationType *mesh, const bool ct);
250
251 virtual SimplexId
252 trunkSegmentation(const std::vector<SimplexId> &pendingNodesVerts,
253 const SimplexId begin,
254 const SimplexId stop);
255
256 // fill treedata_.trunkSegments
258 trunkCTSegmentation(const std::vector<SimplexId> &pendingNodesVerts,
259 const SimplexId begin,
260 const SimplexId stop);
261
262 // segmentation
263
266 void buildSegmentation();
267
268 // Create the segmentation of all arcs by operating the pending operations
270
271 void normalizeIds();
272
273 // -------------
274 // ACCESSOR
275 // ------------
276
277 // Tree info for wrapper
278
279#ifdef TTK_ENABLE_FTM_TREE_STATS_TIME
280 const ActiveTask &getActiveTasks(const idSuperArc taskId) const {
281 return (*mt_data_.activeTasksStats)[taskId];
282 }
283#endif
284
285 inline SimplexId getArcSize(const idSuperArc arcId) {
286 return getSuperArc(arcId)->size();
287 }
288
289 inline bool isJT() const {
290 return mt_data_.treeType == TreeType::Join;
291 }
292
293 inline bool isST() const {
294 return mt_data_.treeType == TreeType::Split;
295 }
296
297 // global
298 // called for the tree used by the wrapper (only).
299 // On this implementation, the warpper communicate with ContourForest
300 // A child class of this one.
301
303 const bool preproc = true) {
304 if(tri && preproc) {
305 // propagate through vertices (build)
307 }
308 }
309
310 inline void setScalars(void *local_scalars) {
311 scalars_->values = local_scalars;
312 }
313
314 inline void setTreeType(const int local_treeType) {
315 params_->treeType = static_cast<TreeType>(local_treeType);
316 }
317
318 inline void setSegmentation(const bool segm) {
319 params_->segm = segm;
320 }
321
322 inline void setNormalizeIds(const bool normalize) {
323 params_->normalize = normalize;
324 }
325
326#ifdef TTK_ENABLE_OMP_PRIORITY
327 inline void setPrior(void) {
328 mt_data_.prior = true;
329 }
330
331 inline bool isPrior(void) const {
332 return mt_data_.prior;
333 }
334#endif
335
336 // scalar
337
338 template <typename scalarType>
339 inline const scalarType &getValue(SimplexId nodeId) const {
340 return (((scalarType *)scalars_->values))[nodeId];
341 }
342
343 template <typename scalarType>
344 inline void setVertexScalars(const scalarType *vals) {
345 scalars_->values = static_cast<void *>(const_cast<scalarType *>(vals));
346 }
347
348 // offset
357 inline void setVertexSoSoffsets(const SimplexId *const sos) {
358 scalars_->offsets = sos;
359 }
360
361 // arcs
362
364 return mt_data_.superArcs->size();
365 }
366
368#ifndef TTK_ENABLE_KAMIKAZE
369 if(i >= mt_data_.superArcs->size()) {
370 std::cout << "[Merge Tree] get superArc on bad id :" << i;
371 std::cout << " / " << mt_data_.superArcs->size() << std::endl;
372 }
373#endif
374 return &((*mt_data_.superArcs)[i]);
375 }
376
377 inline const SuperArc *getSuperArc(idSuperArc i) const {
378#ifndef TTK_ENABLE_KAMIKAZE
379 if(i >= mt_data_.superArcs->size()) {
380 std::cout << "[Merge Tree] get superArc on bad id :" << i;
381 std::cout << " / " << mt_data_.superArcs->size() << std::endl;
382 }
383#endif
384 return &((*mt_data_.superArcs)[i]);
385 }
386
387 // nodes
388
389 inline idNode getNumberOfNodes() const {
390 return mt_data_.nodes->size();
391 }
392
393 inline Node *getNode(idNode nodeId) const {
394 return &((*mt_data_.nodes)[nodeId]);
395 }
396
397 inline void setValence(const SimplexId v, const SimplexId val) {
398 mt_data_.valences[v] = val;
399 }
400
401 // leaves / root
402
403 inline idNode getNumberOfLeaves() const {
404 return mt_data_.leaves.size();
405 }
406
407 inline const std::vector<idNode> &getLeaves() const {
408 // break encapsulation...
409 return mt_data_.leaves;
410 }
411
412 inline idNode getLeave(const idNode id) const {
413#ifndef TTK_ENABLE_KAMIKAZE
414 if(id > mt_data_.leaves.size()) {
415 this->printErr("getLeaves out of bounds: " + std::to_string(id));
416 return mt_data_.leaves[0];
417 }
418#endif
419 return mt_data_.leaves[id];
420 }
421
422 inline const std::vector<idNode> &getRoots() const {
423 // break encapsulation...
424 return (*mt_data_.roots);
425 }
426
427 // vertices
428
430 return scalars_->size;
431 }
432
433 // vert2tree
434
435 inline void setVert2Tree(decltype(mt_data_.vert2tree) const &vect2tree) {
436 mt_data_.vert2tree = vect2tree;
437 }
438
439 // --------------------
440 // VERT 2 TREE Special functions
441 // --------------------
442
443 // test vertex correpondance
444
445 inline bool isCorrespondingArc(const SimplexId val) const {
446 return !isCorrespondingNull(val) && mt_data_.vert2tree[val] >= 0;
447 }
448
449 inline bool isCorrespondingNode(const SimplexId val) const {
450 return mt_data_.vert2tree[val] < 0;
451 }
452
453 inline bool isCorrespondingNull(const SimplexId val) const {
454 return mt_data_.vert2tree[val] == nullCorresp;
455 }
456
457 // Get vertex info
458
459 inline idNode getCorrespondingNodeId(const SimplexId val) const {
460#ifndef TTK_ENABLE_KAMIKAZE
461 if(!isCorrespondingNode(val)) {
462 this->printErr("getCorrespondingNode, Vertex: " + std::to_string(val)
463 + " is not a node: "
464 + std::to_string(mt_data_.vert2tree[val]));
465 }
466#endif
467 return corr2idNode(val);
468 }
469
471#ifndef TTK_ENABLE_KAMIKAZE
472 if(!isCorrespondingArc(val)) {
473 this->printErr(
474 "getCorrespondingSuperArcId, Vertex: " + std::to_string(val)
475 + " is not on an arc: " + std::to_string(mt_data_.vert2tree[val]));
476 }
477#endif
478 return mt_data_.vert2tree[val];
479 }
480
481 // Get corresponding element
482
483 inline SuperArc *vertex2SuperArc(const SimplexId vert) {
484 return &((*mt_data_.superArcs)[getCorrespondingSuperArcId(vert)]);
485 }
486
487 inline Node *vertex2Node(const SimplexId vert) {
488 return &((*mt_data_.nodes)[getCorrespondingNodeId(vert)]);
489 }
490
491 // Update vertex info
492
493 inline void updateCorrespondingArc(const SimplexId vert,
494 const idSuperArc arc) {
495 mt_data_.vert2tree[vert] = arc;
496 }
497
498 inline void updateCorrespondingNode(const SimplexId vert,
499 const idNode node) {
500 mt_data_.vert2tree[vert] = idNode2corr(node);
501 }
502
503 inline idCorresp idNode2corr(const idNode id) const {
504 // transform idNode to special value for the array : -idNode -1
505 return -static_cast<idCorresp>(id + 1);
506 }
507
508 inline idNode corr2idNode(const idCorresp &corr) const {
509 return static_cast<idNode>(-(mt_data_.vert2tree[corr] + 1));
510 }
511
512 // --------------------------------
513 // Arcs and node manipulations
514 // --------------------------------
515 // SuperArcs
516
517 idSuperArc openSuperArc(idNode downNodeId);
518
519 idSuperArc makeSuperArc(idNode downNodeId, idNode upNodeId);
520
521 void closeSuperArc(idSuperArc superArcId, idNode upNodeId);
522
523 // Nodes
524
525 std::vector<idNode> sortedNodes(const bool parallel = false);
526
527 void sortLeaves(const bool parallel = false);
528
535 void sortNodes();
536
544 void sortArcs();
545
546 idNode makeNode(SimplexId vertexId, SimplexId linked = nullVertex);
547
548 idNode makeNode(const Node *const n, SimplexId linked = nullVertex);
549
550 idSuperArc insertNode(Node *node, const bool segm = true);
551
552 // get node starting / ending this arc
553 // orientation depends on Join/Split tree
554 Node *getDownNode(const SuperArc *a);
555 Node *getUpNode(const SuperArc *a);
556 idNode getDownNodeId(const SuperArc *a);
557 idNode getUpNodeId(const SuperArc *a);
558
559 // get node above / below this arc
560 // in term of scalar value
561 Node *getLowerNode(const SuperArc *a);
562 Node *getUpperNode(const SuperArc *a);
565
567 return getSuperArc(getNode(n)->getUpSuperArcId(0))->getUpNodeId();
568 }
569
570 void delNode(idNode node);
571
572 // ---------------------------
573 // Operators : clone/ move & print
574 // ---------------------------
575
576 std::shared_ptr<FTMTree_MT> clone() const;
577
578 void move(FTMTree_MT &mt);
579
580 // Print
581 std::string printArc(idSuperArc a);
582
583 std::string printNode(idNode n);
584
585 void printTree2();
586
587 void printParams() const;
588
589 int printTime(Timer &t,
590 const std::string &s,
591 const int debugLevel = 2) const;
592
593 // ----------------------------------------
594 // Utils functions
595 // Mathieu Pont (mathieu.pont@lip6.fr)
596 // 2021
597 // ----------------------------------------
598
599 // --------------------
600 // Is
601 // --------------------
602 bool isNodeOriginDefined(idNode nodeId) const;
603
604 bool isRoot(idNode nodeId) const;
605
606 bool isLeaf(idNode nodeId) const;
607
608 bool isNodeAlone(idNode nodeId) const;
609
610 bool isFullMerge() const;
611
612 bool isBranchOrigin(idNode nodeId) const;
613
614 template <class dataType>
615 bool isJoinTree() const;
616
617 template <class dataType>
618 bool isImportantPair(idNode nodeId,
619 double threshold,
620 std::vector<double> &excludeLower,
621 std::vector<double> &excludeHigher) const;
622
623 template <class dataType>
624 bool isImportantPair(idNode nodeId, double threshold) const;
625
626 bool isNodeMerged(idNode nodeId) const;
627
628 bool isNodeIdInconsistent(idNode nodeId) const;
629
631
632 // Do not normalize node is if root or son of a merged root
633 bool notNeedToNormalize(idNode nodeId) const;
634
635 bool isMultiPersPair(idNode nodeId) const;
636
637 template <class dataType>
638 bool isParentInconsistent(ftm::idNode nodeId) const;
639
640 template <class dataType>
642
643 // --------------------
644 // Get
645 // --------------------
646 idNode getRoot() const;
647
648 idNode getParentSafe(idNode nodeId) const;
649
650 void getChildren(idNode nodeId, std::vector<idNode> &res) const;
651
652 void getLeavesFromTree(std::vector<idNode> &res) const;
653
654 int getNumberOfLeavesFromTree() const;
655
656 int getNumberOfNodeAlone() const;
657
658 int getRealNumberOfNodes() const;
659
660 template <class dataType>
662
664 idNode node,
665 std::tuple<std::vector<idNode>, std::vector<idNode>> &res) const;
666
667 void
668 getTreeBranching(std::vector<idNode> &branching,
669 std::vector<int> &branchingID,
670 std::vector<std::vector<idNode>> &nodeBranching) const;
671
672 void getTreeBranching(std::vector<idNode> &branching,
673 std::vector<int> &branchingID) const;
674
675 void getAllRoots(std::vector<idNode> &res) const;
676
677 int getNumberOfRoot() const;
678
679 int getNumberOfChildren(idNode nodeId) const;
680
681 int getTreeDepth() const;
682
683 int getNodeLevel(idNode nodeId) const;
684
685 void getAllNodeLevel(std::vector<int> &res) const;
686
687 void getLevelToNode(std::vector<std::vector<idNode>> &res) const;
688
689 void getBranchSubtree(std::vector<idNode> &branching,
690 idNode branchRoot,
691 std::vector<idNode> &res) const;
692
693 template <class dataType>
694 idNode getLowestNode(idNode nodeStart) const;
695
696 // --------------------
697 // Persistence
698 // --------------------
699 template <class dataType>
700 std::tuple<dataType, dataType> getBirthDeathFromIds(idNode nodeId1,
701 idNode nodeId2) const;
702
703 template <class dataType>
704 std::tuple<dataType, dataType>
705 getBirthDeathNodeFromIds(idNode nodeId1, idNode nodeId2) const;
706
707 template <class dataType>
708 std::tuple<dataType, dataType> getBirthDeath(idNode nodeId) const;
709
710 template <class dataType>
711 std::tuple<ftm::idNode, ftm::idNode>
712 getBirthDeathNode(idNode nodeId) const;
713
714 template <class dataType>
715 std::tuple<dataType, dataType> getMergedRootBirthDeath() const;
716
717 template <class dataType>
718 std::tuple<ftm::idNode, ftm::idNode> getMergedRootBirthDeathNode() const;
719
720 template <class dataType>
721 dataType getBirth(idNode nodeId) const;
722
723 template <class dataType>
724 dataType getNodePersistence(idNode nodeId) const;
725
726 template <class dataType>
727 dataType getMaximumPersistence() const;
728
729 template <class dataType>
731
732 template <class dataType>
733 dataType getSecondMaximumPersistence() const;
734
735 template <class dataType>
737 std::vector<std::tuple<ftm::idNode, ftm::idNode, dataType>> &pairs,
738 bool useBD) const;
739
740 template <class dataType>
741 std::vector<ftm::idNode> getMultiPersOrigins(bool useBD) const;
742
744 std::vector<std::vector<idNode>> &res) const;
745
746 // --------------------
747 // Set
748 // --------------------
749 void setParent(idNode nodeId, idNode newParentNodeId);
750
751 // --------------------
752 // Delete
753 // --------------------
754 // Delete node by keeping subtree
755 void deleteNode(idNode nodeId);
756
757 void deleteIthUpArc(idNode nodeId, int arcIth);
758
759 // Delete arc of the node to its parent
760 void deleteParent(idNode nodeId);
761
762 // Delete node without keeping subtree
763 void deleteSubtree(idNode nodeId);
764
765 // --------------------
766 // Create/Delete/Modify Tree
767 // --------------------
768 void copyMergeTreeStructure(const FTMTree_MT *tree);
769
770 // --------------------
771 // Utils
772 // --------------------
773 void printNodeSS(idNode node, std::stringstream &ss) const;
774
775 template <class dataType>
776 std::stringstream printNode2(idNode nodeId, bool doPrint = true) const;
777
778 template <class dataType>
779 std::stringstream printMergedRoot(bool doPrint = true) const;
780
781 std::stringstream printSubTree(idNode subRoot) const;
782
783 std::stringstream printTree(bool doPrint = true) const;
784
785 std::stringstream printTreeStats(bool doPrint = true) const;
786
787 template <class dataType>
788 std::stringstream printTreeScalars(bool printNodeAlone = true,
789 bool doPrint = true) const;
790
791 template <class dataType>
792 std::stringstream printPairsFromTree(bool useBD = false,
793 bool printPairs = true,
794 bool doPrint = true) const;
795
796 std::stringstream printMultiPersOriginsVectorFromTree(bool doPrint
797 = true) const;
798
799 template <class dataType>
800 std::stringstream printMultiPersPairsFromTree(bool useBD = false,
801 bool printPairs = true,
802 bool doPrint = true) const;
803
804 // ----------------------------------------
805 // End of utils functions
806 // ----------------------------------------
807
808 protected:
809 // -----
810 // Tools
811 // -----
812
813 idNode getVertInRange(const std::vector<SimplexId> &range,
814 const SimplexId v,
815 const idNode last = 0) const;
816
817 std::tuple<SimplexId, SimplexId>
818 getBoundsFromVerts(const std::vector<SimplexId> &nodes) const;
819
823
824 inline SimplexId getChunkSize(const SimplexId nbVerts = -1,
825 const SimplexId nbtasks = 100) const {
826 const SimplexId s = (nbVerts == -1) ? scalars_->size : nbVerts;
827#ifndef NDEBUG
828 // Debug mode
829 static const SimplexId minWorks = 1;
830#else
831 // Release mode
832 static const SimplexId minWorks = 10000;
833#endif
834 return std::max(minWorks, 1 + (s / (nbtasks * threadNumber_)));
835 }
836
837 inline SimplexId getChunkCount(const SimplexId nbVerts = -1,
838 const SimplexId nbTasks = 100) const {
839 const SimplexId s = (nbVerts == -1) ? scalars_->size : nbVerts;
840 return 1 + (s / getChunkSize(s, nbTasks));
841 }
842
843 void sortUpArcs(const idNode nid) {
844 auto comp = [&](const idSuperArc a, const idSuperArc b) -> bool {
845 return comp_.vertLower(getUpperNode(getSuperArc(a))->getVertexId(),
846 getUpperNode(getSuperArc(b))->getVertexId());
847 };
848
849 getNode(nid)->sortUpArcs(comp);
850 }
851
852 void sortDownArcs(const idNode nid) {
853 auto comp = [&](const idSuperArc a, const idSuperArc b) -> bool {
854 return comp_.vertHigher(getUpperNode(getSuperArc(a))->getVertexId(),
855 getUpperNode(getSuperArc(b))->getVertexId());
856 };
857
858 getNode(nid)->sortDownArcs(comp);
859 }
860
861 // ------------------
862 // Comparisons
863 // -----------------
864 // Compare using the scalar array : only for sort step
865
866 inline bool isLower(SimplexId a, SimplexId b) const {
867 return scalars_->offsets[a] < scalars_->offsets[b];
868 }
869
870 inline bool isHigher(SimplexId a, SimplexId b) const {
871 return scalars_->offsets[a] > scalars_->offsets[b];
872 }
873
874 template <typename type>
875 void createVector(std::vector<type> &vec) {
876 vec.clear();
877 }
878
879 template <typename type>
880 void createAtomicVector(std::shared_ptr<FTMAtomicVector<type>> &ptr) {
881 if(!ptr)
882 ptr = std::make_shared<FTMAtomicVector<type>>();
883 ptr->clear();
884 }
885
886 template <typename type>
887 void initVector(std::vector<type> &vect, const type val) {
888 auto s = vect.size();
889#ifdef TTK_ENABLE_OPENMP
890#pragma omp parallel for num_threads(threadNumber_) schedule(static)
891#endif
892 for(typename std::vector<type>::size_type i = 0; i < s; i++) {
893 vect[i] = val;
894 }
895 }
896 }; // end of FTMTree_MT class
897
898 std::ostream &operator<<(std::ostream &o, Node const &n);
899 std::ostream &operator<<(std::ostream &o, SuperArc const &a);
900
901 template <typename dataType>
902 struct MergeTree {
903 std::shared_ptr<ftm::Scalars> scalars;
904 std::shared_ptr<std::vector<dataType>> scalarsValues;
905 std::shared_ptr<ftm::Params> params;
907
908 std::shared_ptr<ftm::Scalars> emptyScalars() {
909 auto scalarsT = std::make_shared<ftm::Scalars>();
910 scalarsT->size = 0;
911 scalarsT->values = nullptr;
912 return scalarsT;
913 }
914
915 std::shared_ptr<ftm::Params> emptyParams() {
916 auto paramsT = std::make_shared<ftm::Params>();
917 paramsT->treeType = ftm::Join_Split;
918 return paramsT;
919 }
920
923
924 template <typename T, typename U>
925 MergeTree(const T scalarsT, U paramsT)
926 : scalars(scalarsT), params(paramsT),
927 tree(paramsT, scalarsT, params->treeType) {
928 tree.makeAlloc();
929 scalarsValues = std::make_shared<std::vector<dataType>>();
930 for(unsigned int i = 0; i < tree.getNumberOfNodes(); ++i)
931 scalarsValues->push_back(tree.getValue<dataType>(i));
932 scalars->values = (void *)(scalarsValues->data());
933 }
934
935 MergeTree(const std::shared_ptr<ftm::Scalars> &scalarsT,
936 const std::shared_ptr<std::vector<dataType>> &scalarValuesT,
937 std::shared_ptr<ftm::Params> &paramsT)
938 : scalars(scalarsT), scalarsValues(scalarValuesT), params(paramsT),
939 tree(paramsT, scalarsT, params->treeType) {
940 tree.makeAlloc();
941 scalars->values = (void *)(scalarsValues->data());
942 }
943
944 void copy(const MergeTree<dataType> &mt) {
945 // Copy scalars
946 scalars = std::make_shared<ftm::Scalars>();
947 scalars->size = mt.scalars->size;
949 scalars->values = (void *)(scalarsValues->data());
950
951 // Copy params
952 params = std::make_shared<ftm::Params>();
953 params->treeType = mt.params->treeType;
954
955 // Copy tree
956 tree.clear();
957 tree.setParamsScalars(params, scalars);
958 tree.makeAlloc();
959 tree.copyMergeTreeStructure(const_cast<FTMTree_MT *>(&(mt.tree)));
960 }
961
964 params(mt.params), tree(params, scalars, params->treeType) {
965 copy(mt);
966 }
967
969 if(&mt != this) {
970 copy(mt);
971 }
972 return *this;
973 }
974 };
975
976 } // namespace ftm
977} // namespace ttk
978
980#include <FTMTree_MT_Template.h>
AbstractTriangulation is an interface class that defines an interface for efficient traversal methods...
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
std::stringstream printMultiPersPairsFromTree(bool useBD=false, bool printPairs=true, bool doPrint=true) const
std::stringstream printMultiPersOriginsVectorFromTree(bool doPrint=true) const
std::tuple< dataType, dataType > getMergedRootBirthDeath() const
void getMultiPersOriginsVectorFromTree(std::vector< std::vector< idNode > > &res) const
Node * getUpperNode(const SuperArc *a)
bool isBranchOrigin(idNode nodeId) const
idCorresp idNode2corr(const idNode id) const
Definition FTMTree_MT.h:503
int getNumberOfLeavesFromTree() const
idNode corr2idNode(const idCorresp &corr) const
Definition FTMTree_MT.h:508
bool isCorrespondingNode(const SimplexId val) const
Definition FTMTree_MT.h:449
int getNumberOfNodeAlone() const
Node * getNode(idNode nodeId) const
Definition FTMTree_MT.h:393
std::shared_ptr< Params > params_
Definition FTMTree_MT.h:86
void closeOnBackBone(const triangulationType *mesh, SimplexId saddleVert)
std::stringstream printTree(bool doPrint=true) const
idNode getLeave(const idNode id) const
Definition FTMTree_MT.h:412
const scalarType & getValue(SimplexId nodeId) const
Definition FTMTree_MT.h:339
FTMTree_MT(const std::shared_ptr< Params > &params, const std::shared_ptr< Scalars > &scalars, TreeType type)
std::string printNode(idNode n)
void getChildren(idNode nodeId, std::vector< idNode > &res) const
SuperArc * getSuperArc(idSuperArc i)
Definition FTMTree_MT.h:367
Node * getDownNode(const SuperArc *a)
void setScalars(void *local_scalars)
Definition FTMTree_MT.h:310
std::tuple< ftm::idNode, ftm::idNode > getBirthDeathNode(idNode nodeId) const
void printParams() const
std::tuple< dataType, dataType > getBirthDeathNodeFromIds(idNode nodeId1, idNode nodeId2) const
void createAtomicVector(std::shared_ptr< FTMAtomicVector< type > > &ptr)
Definition FTMTree_MT.h:880
void getTreeBranching(std::vector< idNode > &branching, std::vector< int > &branchingID, std::vector< std::vector< idNode > > &nodeBranching) const
std::stringstream printNode2(idNode nodeId, bool doPrint=true) const
idNode getNumberOfNodes() const
Definition FTMTree_MT.h:389
bool isMultiPersPair(idNode nodeId) const
void setParent(idNode nodeId, idNode newParentNodeId)
idNode getUpNodeId(const SuperArc *a)
void initVector(std::vector< type > &vect, const type val)
Definition FTMTree_MT.h:887
SimplexId getChunkSize(const SimplexId nbVerts=-1, const SimplexId nbtasks=100) const
Definition FTMTree_MT.h:824
bool isParentInconsistent(ftm::idNode nodeId) const
bool isCorrespondingNull(const SimplexId val) const
Definition FTMTree_MT.h:453
void sortInput()
if sortedVertices_ is null, define and fill it
void delNode(idNode node)
SuperArc * vertex2SuperArc(const SimplexId vert)
Definition FTMTree_MT.h:483
std::string printArc(idSuperArc a)
Node * vertex2Node(const SimplexId vert)
Definition FTMTree_MT.h:487
void copyMergeTreeStructure(const FTMTree_MT *tree)
idNode getRoot() const
void move(FTMTree_MT &mt)
idSuperArc upArcFromVert(const SimplexId v)
Definition FTMTree_MT.h:820
std::tuple< SimplexId, SimplexId > getBoundsFromVerts(const std::vector< SimplexId > &nodes) const
int leafSearch(const triangulationType *mesh)
Node * getLowerNode(const SuperArc *a)
idNode getDownNodeId(const SuperArc *a)
void initNbScalars(const triangulationType *triangulation)
Definition FTMTree_MT.h:119
idSuperArc insertNode(Node *node, const bool segm=true)
bool isCorrespondingArc(const SimplexId val) const
Definition FTMTree_MT.h:445
std::stringstream printTreeStats(bool doPrint=true) const
int getNodeLevel(idNode nodeId) const
bool verifyBranchDecompositionInconsistency() const
std::stringstream printTreeScalars(bool printNodeAlone=true, bool doPrint=true) const
void setVertexSoSoffsets(const SimplexId *const sos)
Definition FTMTree_MT.h:357
idNode getLowestNode(idNode nodeStart) const
void closeArcsUF(idNode closeNode, UF uf)
ftm::idNode getSecondMaximumPersistenceNode() const
idNode getParentSafe(idNode nodeId) const
void closeAndMergeOnSaddle(const triangulationType *mesh, SimplexId saddleVert)
idNode getLowerNodeId(const SuperArc *a)
void makeAlloc()
clear local data for new computation
Definition FTMTree_MT.h:154
SimplexId getNumberOfVertices() const
Definition FTMTree_MT.h:429
std::vector< ftm::idNode > getMultiPersOrigins(bool useBD) const
void setNormalizeIds(const bool normalize)
Definition FTMTree_MT.h:322
void sortUpArcs(const idNode nid)
Definition FTMTree_MT.h:843
void getLevelToNode(std::vector< std::vector< idNode > > &res) const
void getAllNodeLevel(std::vector< int > &res) const
virtual SimplexId trunkSegmentation(const std::vector< SimplexId > &pendingNodesVerts, const SimplexId begin, const SimplexId stop)
std::stringstream printSubTree(idNode subRoot) const
bool isLower(SimplexId a, SimplexId b) const
Definition FTMTree_MT.h:866
dataType getMaximumPersistence() const
int getNumberOfChildren(idNode nodeId) const
void sortNodes()
Sort tree nodes according to vertex order.
std::tuple< bool, bool > propagate(const triangulationType *mesh, CurrentState &currentState, UF curUF)
idNode getCorrespondingNodeId(const SimplexId val) const
Definition FTMTree_MT.h:459
idSuperArc getNumberOfSuperArcs() const
Definition FTMTree_MT.h:363
idNode getUpperNodeId(const SuperArc *a)
void deleteIthUpArc(idNode nodeId, int arcIth)
const std::vector< idNode > & getRoots() const
Definition FTMTree_MT.h:422
dataType getSecondMaximumPersistence() const
bool isST() const
Definition FTMTree_MT.h:293
int getRealNumberOfNodes() const
void sortLeaves(const bool parallel=false)
idNode getNumberOfLeaves() const
Definition FTMTree_MT.h:403
std::tuple< ftm::idNode, ftm::idNode > getMergedRootBirthDeathNode() const
void getLeavesFromTree(std::vector< idNode > &res) const
dataType getNodePersistence(idNode nodeId) const
void build(const triangulationType *mesh, const bool ct)
Compute the merge.
idNode getVertInRange(const std::vector< SimplexId > &range, const SimplexId v, const idNode last=0) const
std::stringstream printPairsFromTree(bool useBD=false, bool printPairs=true, bool doPrint=true) const
void setTreeType(const int local_treeType)
Definition FTMTree_MT.h:314
void getBranchOriginsFromThisBranch(idNode node, std::tuple< std::vector< idNode >, std::vector< idNode > > &res) const
bool isNodeOriginDefined(idNode nodeId) const
void sortDownArcs(const idNode nid)
Definition FTMTree_MT.h:852
bool isLeaf(idNode nodeId) const
std::tuple< dataType, dataType > getBirthDeathFromIds(idNode nodeId1, idNode nodeId2) const
void deleteNode(idNode nodeId)
int printTime(Timer &t, const std::string &s, const int debugLevel=2) const
void getAllRoots(std::vector< idNode > &res) const
const std::vector< idNode > & getLeaves() const
Definition FTMTree_MT.h:407
void setVertexScalars(const scalarType *vals)
Definition FTMTree_MT.h:344
bool isRoot(idNode nodeId) const
std::tuple< dataType, dataType > getBirthDeath(idNode nodeId) const
std::stringstream printMergedRoot(bool doPrint=true) const
SimplexId getChunkCount(const SimplexId nbVerts=-1, const SimplexId nbTasks=100) const
Definition FTMTree_MT.h:837
void createVector(std::vector< type > &vec)
Definition FTMTree_MT.h:875
bool notNeedToNormalize(idNode nodeId) const
bool isNodeIdInconsistent(idNode nodeId) const
SimplexId trunkCTSegmentation(const std::vector< SimplexId > &pendingNodesVerts, const SimplexId begin, const SimplexId stop)
idSuperArc openSuperArc(idNode downNodeId)
int getNumberOfRoot() const
std::shared_ptr< FTMTree_MT > clone() const
SimplexId trunk(const triangulationType *mesh, const bool ct)
idNode getParent(const idNode n)
Definition FTMTree_MT.h:566
idNode makeNode(SimplexId vertexId, SimplexId linked=nullVertex)
bool isNodeAlone(idNode nodeId) const
void buildSegmentation()
use vert2tree to compute the segmentation of the fresh built merge tree.
void getPersistencePairsFromTree(std::vector< std::tuple< ftm::idNode, ftm::idNode, dataType > > &pairs, bool useBD) const
const SuperArc * getSuperArc(idSuperArc i) const
Definition FTMTree_MT.h:377
void deleteParent(idNode nodeId)
idSuperArc getCorrespondingSuperArcId(const SimplexId val) const
Definition FTMTree_MT.h:470
void leafGrowth(const triangulationType *mesh)
bool isImportantPair(idNode nodeId, double threshold, std::vector< double > &excludeLower, std::vector< double > &excludeHigher) const
bool isNodeMerged(idNode nodeId) const
void setValence(const SimplexId v, const SimplexId val)
Definition FTMTree_MT.h:397
idSuperArc makeSuperArc(idNode downNodeId, idNode upNodeId)
void updateCorrespondingNode(const SimplexId vert, const idNode node)
Definition FTMTree_MT.h:498
void printNodeSS(idNode node, std::stringstream &ss) const
Node * getUpNode(const SuperArc *a)
bool isJT() const
Definition FTMTree_MT.h:289
void preconditionTriangulation(AbstractTriangulation *tri, const bool preproc=true)
Definition FTMTree_MT.h:302
void getBranchSubtree(std::vector< idNode > &branching, idNode branchRoot, std::vector< idNode > &res) const
void closeSuperArc(idSuperArc superArcId, idNode upNodeId)
bool isFullMerge() const
SimplexId getArcSize(const idSuperArc arcId)
Definition FTMTree_MT.h:285
bool isHigher(SimplexId a, SimplexId b) const
Definition FTMTree_MT.h:870
std::shared_ptr< Scalars > scalars_
Definition FTMTree_MT.h:87
void sortArcs()
Sort tree arcs.
void deleteSubtree(idNode nodeId)
bool isThereOnlyOnePersistencePair() const
void setSegmentation(const bool segm)
Definition FTMTree_MT.h:318
void initVectStates(const SimplexId nbLeaves)
Definition FTMTree_MT.h:202
dataType getBirth(idNode nodeId) const
void updateCorrespondingArc(const SimplexId vert, const idSuperArc arc)
Definition FTMTree_MT.h:493
void arcGrowth(const triangulationType *mesh, const SimplexId startVert, const SimplexId orig)
std::vector< idNode > sortedNodes(const bool parallel=false)
void setParamsScalars(const std::shared_ptr< Params > &params, const std::shared_ptr< Scalars > &scalars)
Definition FTMTree_MT.h:111
void setVert2Tree(decltype(mt_data_.vert2tree) const &vect2tree)
Definition FTMTree_MT.h:435
bool compLower(const SimplexId a, const SimplexId b)
Definition FTMTree_MT.h:145
void sortDownArcs(const std::function< bool(const idSuperArc, const idSuperArc)> &comp)
Definition FTMNode.h:190
idSuperArc getUpSuperArcId(idSuperArc neighborId) const
Definition FTMNode.h:105
void sortUpArcs(const std::function< bool(const idSuperArc, const idSuperArc)> &comp)
Definition FTMNode.h:185
idNode getUpNodeId() const
Definition FTMSuperArc.h:68
size_t size() const
long unsigned int idSuperArc
SuperArc index in vect_superArcs_.
AtomicUF * UF
Definition FTMTree_MT.h:39
std::ostream & operator<<(std::ostream &o, Node const &n)
long long int idCorresp
type used to recover Node/Arc in vert2tree SIGNED ONLY
unsigned int idNode
Node index in vect_nodes_.
TTK base package defining the standard types.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
coefficient_t normalize(const coefficient_t n, const coefficient_t modulus)
Definition ripser.cpp:171
T begin(std::pair< T, T > &p)
Definition ripser.cpp:499
std::shared_ptr< ftm::Params > params
Definition FTMTree_MT.h:905
MergeTree(const std::shared_ptr< ftm::Scalars > &scalarsT, const std::shared_ptr< std::vector< dataType > > &scalarValuesT, std::shared_ptr< ftm::Params > &paramsT)
Definition FTMTree_MT.h:935
void copy(const MergeTree< dataType > &mt)
Definition FTMTree_MT.h:944
MergeTree(const T scalarsT, U paramsT)
Definition FTMTree_MT.h:925
MergeTree(const MergeTree< dataType > &mt)
Definition FTMTree_MT.h:962
std::shared_ptr< std::vector< dataType > > scalarsValues
Definition FTMTree_MT.h:904
std::shared_ptr< ftm::Scalars > emptyScalars()
Definition FTMTree_MT.h:908
MergeTree< dataType > & operator=(const MergeTree< dataType > &mt)
Definition FTMTree_MT.h:968
ftm::FTMTree_MT tree
Definition FTMTree_MT.h:906
std::shared_ptr< ftm::Params > emptyParams()
Definition FTMTree_MT.h:915
std::shared_ptr< ftm::Scalars > scalars
Definition FTMTree_MT.h:903
std::shared_ptr< FTMAtomicVector< SuperArc > > superArcs
Definition FTMTree_MT.h:46
std::vector< AtomicUF > storage
Definition FTMTree_MT.h:57
std::shared_ptr< FTMAtomicVector< CurrentState > > states
Definition FTMTree_MT.h:60
std::vector< std::list< std::vector< SimplexId > > > trunkSegments
Definition FTMTree_MT.h:54
std::shared_ptr< FTMAtomicVector< Node > > nodes
Definition FTMTree_MT.h:47
std::vector< SimplexId > visitOrder
Definition FTMTree_MT.h:53
std::vector< UF > propagation
Definition FTMTree_MT.h:59
std::vector< idCorresp > vert2tree
Definition FTMTree_MT.h:52
std::vector< UF > ufs
Definition FTMTree_MT.h:58
std::shared_ptr< FTMAtomicVector< idNode > > roots
Definition FTMTree_MT.h:48
std::vector< valence > valences
Definition FTMTree_MT.h:62
std::vector< idNode > leaves
Definition FTMTree_MT.h:49
std::vector< char > openedNodes
Definition FTMTree_MT.h:64