TTK
Loading...
Searching...
No Matches
ttkMergeTreeVisualization.h
Go to the documentation of this file.
1
7
8#pragma once
9
10#include <FTMTree.h>
12
13#include <ttkAlgorithm.h>
14
15// VTK Includes
16#include <vtkAppendFilter.h>
17#include <vtkCellData.h>
18#include <vtkDoubleArray.h>
19#include <vtkFloatArray.h>
20#include <vtkImageData.h>
21#include <vtkIntArray.h>
22#include <vtkPointData.h>
23#include <vtkPoints.h>
24#include <vtkStringArray.h>
25#include <vtkUnstructuredGrid.h>
26
28private:
29 // redefine types
30 using idNode = ttk::ftm::idNode;
31 using SimplexId = ttk::SimplexId;
33
34 // Visualization parameters
35 bool PlanarLayout = false;
36 double DimensionSpacing = 1.;
37 int DimensionToShift = 0;
38 bool OutputSegmentation = false;
39 int MaximumImportantPairs = 0;
40 int MinimumImportantPairs = 0;
41 bool outputTreeNodeIndex = false;
42 bool isPersistenceDiagram = false;
43 bool convertedToDiagram = false;
44 bool isPDSadMax = true;
45
46 // Shift mode
47 // -1: None ; 0: Star ; 1: Star Barycenter ; 2: Line ; 3: Double Line
48 int ShiftMode = 0;
49
50 // Offset
51 int iSampleOffset = 0;
52 int noSampleOffset = 0;
53 double prevXMaxOffset = 0;
54
55 // Print only one tree
56 int printTreeId = -1; // -1 for all
57 int printClusterId = -1; // -1 for all
58
59 // Barycenter position according alpha
60 bool BarycenterPositionAlpha = false;
61 double Alpha = 0.5;
62
63 // Used for critical type and for point coordinates
64 std::vector<vtkUnstructuredGrid *> treesNodes;
65 std::vector<std::vector<int>>
66 treesNodeCorrMesh; // used to access treesNodes given input trees
67
68 // Segmentation
69 std::vector<vtkDataSet *> treesSegmentation;
70
71 // Clustering output
72 std::vector<int> clusteringAssignment;
73 std::vector<std::vector<std::vector<std::tuple<idNode, idNode, double>>>>
74 outputMatchingBarycenter;
75
76 // Barycenter output
77 std::vector<std::vector<float>> allBaryPercentMatch;
78
79 // Temporal Subsampling Output
80 std::vector<bool> interpolatedTrees;
81
82 // Output
83 vtkUnstructuredGrid *vtkOutputNode{};
84 vtkUnstructuredGrid *vtkOutputArc{};
85 vtkDataSet *vtkOutputSegmentation{};
86
87 // Matching output
88 vtkUnstructuredGrid *vtkOutputNode1{}, *vtkOutputNode2{}; // input data
89 std::vector<std::vector<SimplexId>> nodeCorr1, nodeCorr2;
90 vtkUnstructuredGrid *vtkOutputMatching{}; // output
91
92 // Custom array
93 std::vector<std::tuple<std::string, std::vector<double>>> customArrays;
94 std::vector<std::tuple<std::string, std::vector<int>>> customIntArrays;
95 std::vector<std::tuple<std::string, std::vector<std::string>>>
96 customStringArrays;
97
98 // Filled by the algorithm
99 std::vector<std::vector<SimplexId>> nodeCorr;
100 std::vector<double> clusterShift;
101 double prevXMax = 0;
102
103public:
105 ;
106 ~ttkMergeTreeVisualization() override = default;
107 ;
108
109 // ==========================================================================
110 // Getter / Setter
111 // ==========================================================================
112 // Visualization parameters
113 void setPlanarLayout(bool b) {
114 PlanarLayout = b;
115 }
116 void setDimensionSpacing(double d) {
117 DimensionSpacing = d;
118 }
120 DimensionToShift = i;
121 }
123 OutputSegmentation = b;
124 }
125 void setShiftMode(int mode) {
126 ShiftMode = mode;
127 }
128 void setMaximumImportantPairs(int maxPairs) {
129 MaximumImportantPairs = maxPairs;
130 }
131 void setMinimumImportantPairs(int minPairs) {
132 MinimumImportantPairs = minPairs;
133 }
134
135 void setOutputTreeNodeId(int doOutput) {
136 outputTreeNodeIndex = doOutput;
137 }
138
139 void setIsPersistenceDiagram(bool isPD) {
140 isPersistenceDiagram = isPD;
141 }
142 void setConvertedToDiagram(bool converted) {
143 convertedToDiagram = converted;
144 }
145 void setIsPDSadMax(bool isSadMax) {
146 isPDSadMax = isSadMax;
147 }
148
149 // Offset
150 void setISampleOffset(int offset) {
151 iSampleOffset = offset;
152 }
153 void setNoSampleOffset(int offset) {
154 noSampleOffset = offset;
155 }
156 void setPrevXMaxOffset(double offset) {
157 prevXMaxOffset = offset;
158 }
159
160 // Print only one tree
161 void setPrintTreeId(int id) {
162 printTreeId = id;
163 }
164 void setPrintClusterId(int id) {
165 printClusterId = id;
166 }
167
168 // Barycenter position according alpha
170 BarycenterPositionAlpha = pos;
171 }
172 void setAlpha(double alpha) {
173 Alpha = alpha;
174 }
175
176 // Used for critical type and if not planar layout for point coordinates
177 void setTreesNodes(std::vector<vtkUnstructuredGrid *> &nodes) {
178 treesNodes = nodes;
179 }
180 void setTreesNodeCorrMesh(std::vector<std::vector<int>> &nodeCorrMesh) {
181 treesNodeCorrMesh = nodeCorrMesh;
182 }
183 void setTreesNodes(vtkUnstructuredGrid *nodes) {
184 treesNodes.clear();
185 treesNodes.emplace_back(nodes);
186 }
187 void setTreesNodeCorrMesh(std::vector<int> &nodeCorrMesh) {
188 treesNodeCorrMesh.clear();
189 treesNodeCorrMesh.emplace_back(nodeCorrMesh);
190 }
191
192 // Segmentation
193 void setTreesSegmentation(std::vector<vtkDataSet *> &segmentation) {
194 treesSegmentation = segmentation;
195 }
196 void setTreesSegmentation(vtkDataSet *segmentation) {
197 treesSegmentation.clear();
198 treesSegmentation.emplace_back(segmentation);
199 }
200
201 // Clustering output
202 void setClusteringAssignment(std::vector<int> &asgn) {
203 clusteringAssignment = asgn;
204 }
206 std::vector<std::vector<std::vector<std::tuple<idNode, idNode, double>>>>
207 &matching) {
208 outputMatchingBarycenter = matching;
209 }
210
211 // Barycenter output
212 std::vector<std::vector<float>> getAllBaryPercentMatch() {
213 return allBaryPercentMatch;
214 }
215 void
216 setAllBaryPercentMatch(std::vector<std::vector<float>> &baryPercentMatch) {
217 allBaryPercentMatch = baryPercentMatch;
218 }
219
220 // Temporal Subsampling Output
221 void setInterpolatedTrees(std::vector<bool> &isInterpolatedTrees) {
222 interpolatedTrees = isInterpolatedTrees;
223 }
224
225 // Output
226 void setVtkOutputNode(vtkUnstructuredGrid *vtkNode) {
227 vtkOutputNode = vtkNode;
228 }
229 void setVtkOutputArc(vtkUnstructuredGrid *vtkArc) {
230 vtkOutputArc = vtkArc;
231 }
232 void setVtkOutputSegmentation(vtkDataSet *vtkSegmentation) {
233 vtkOutputSegmentation = vtkSegmentation;
234 }
235
236 // Matching output
237 void setVtkOutputNode1(vtkUnstructuredGrid *vtkNode1) {
238 vtkOutputNode1 = vtkNode1;
239 }
240 void setVtkOutputNode2(vtkUnstructuredGrid *vtkNode2) {
241 vtkOutputNode2 = vtkNode2;
242 }
243 void setNodeCorr1(std::vector<std::vector<SimplexId>> &nodeCorrT) {
244 nodeCorr1 = nodeCorrT;
245 }
246 void setNodeCorr2(std::vector<std::vector<SimplexId>> &nodeCorrT) {
247 nodeCorr2 = nodeCorrT;
248 }
249 void setVtkOutputMatching(vtkUnstructuredGrid *vtkMatching) {
250 vtkOutputMatching = vtkMatching;
251 }
253 std::vector<std::tuple<idNode, idNode, double>> &matching) {
254 outputMatchingBarycenter.resize(1);
255 outputMatchingBarycenter[0].resize(1);
256 outputMatchingBarycenter[0][0] = matching;
257 }
258
259 // Custom array
260 void addCustomArray(std::string &name, std::vector<double> &vec) {
261 customArrays.emplace_back(name, vec);
262 }
263 void addCustomIntArray(std::string &name, std::vector<int> &vec) {
264 customIntArrays.emplace_back(name, vec);
265 }
266 void addCustomStringArray(std::string &name, std::vector<std::string> &vec) {
267 customStringArrays.emplace_back(name, vec);
268 }
270 customArrays.clear();
271 }
273 customIntArrays.clear();
274 }
276 customStringArrays.clear();
277 }
283
284 template <class dataType>
286 std::vector<std::tuple<std::string, std::vector<dataType>>> &cArrays,
287 std::vector<std::vector<dataType>> &cArraysValues,
288 vtkUnstructuredGrid *vtkOutput,
289 int type,
290 int output) {
291 for(unsigned int i = 0; i < cArrays.size(); ++i) {
292 vtkNew<vtkDoubleArray> customDoubleArrayVtk;
293 vtkNew<vtkIntArray> customIntArrayVtk;
294 vtkNew<vtkStringArray> customStringArrayVtk;
295 vtkAbstractArray *customArrayVtk;
296 if(type == 0)
297 customArrayVtk = customDoubleArrayVtk;
298 else if(type == 1)
299 customArrayVtk = customIntArrayVtk;
300 else
301 customArrayVtk = customStringArrayVtk;
302 customArrayVtk->SetName(std::get<0>(cArrays[i]).c_str());
303 customArrayVtk->SetNumberOfTuples(cArraysValues[i].size());
304 for(unsigned int j = 0; j < cArraysValues[i].size(); ++j) {
305 // Add value depending on type (vtkAbstractArray can not be used here)
306 if(type == 0) {
307 double const doubleValue
308 = (*(std::vector<double> *)&(cArraysValues[i]))[j];
309 customDoubleArrayVtk->SetValue(j, doubleValue);
310 } else if(type == 1) {
311 int const intValue = (*(std::vector<int> *)&(cArraysValues[i]))[j];
312 customIntArrayVtk->SetValue(j, intValue);
313 } else {
314 std::string const stringValue
315 = (*(std::vector<std::string> *)&(cArraysValues[i]))[j];
316 customStringArrayVtk->SetValue(j, stringValue);
317 }
318 }
319 if(output == 0)
320 vtkOutput->GetPointData()->AddArray(customArrayVtk);
321 else
322 vtkOutput->GetCellData()->AddArray(customArrayVtk);
323 }
324 }
325
326 void getTreeNodeIdRev(vtkDataArray *treeNodeIdArray,
327 std::vector<int> &treeNodeIdRev) {
328 double valueRange[2];
329 treeNodeIdArray->GetRange(valueRange);
330 int const maxValue = valueRange[1];
331 treeNodeIdRev.clear();
332 treeNodeIdRev.resize(maxValue + 1);
333 for(int i = 0; i < treeNodeIdArray->GetNumberOfValues(); ++i)
334 treeNodeIdRev[treeNodeIdArray->GetTuple1(i)] = i;
335 }
336
337 void copyPointData(vtkUnstructuredGrid *treeNodes,
338 std::vector<int> &nodeCorrT) {
339 if(!treeNodes)
340 return;
341
342 auto treeNodeIdArray = treeNodes->GetPointData()->GetArray("TreeNodeId");
343 std::vector<int> treeNodeIdRev;
344 if(treeNodeIdArray)
345 getTreeNodeIdRev(treeNodeIdArray, treeNodeIdRev);
346
347 for(int i = 0; i < treeNodes->GetPointData()->GetNumberOfArrays(); ++i) {
348 auto dataArray
349 = vtkDataArray::SafeDownCast(treeNodes->GetPointData()->GetArray(i));
350 auto stringArray
351 = vtkStringArray::SafeDownCast(treeNodes->GetPointData()->GetArray(i));
352 vtkAbstractArray *array;
353 if(dataArray)
354 array = dataArray;
355 else if(stringArray)
356 array = stringArray;
357 else
358 continue;
359 auto vecSize = (nodeCorrT.size() == 0 ? array->GetNumberOfValues()
360 : nodeCorrT.size());
361 std::vector<double> vec(vecSize);
362 std::vector<std::string> vecString(vecSize);
363 for(unsigned int j = 0; j < vec.size(); ++j) {
364 int toGet = (nodeCorrT.size() == 0 ? j : nodeCorrT[j]);
365 if(treeNodeIdArray)
366 toGet = (nodeCorrT.size() == 0 ? treeNodeIdRev[j]
367 : treeNodeIdRev[nodeCorrT[j]]);
368 auto value = array->GetVariantValue(toGet);
369 if(dataArray)
370 vec[j] = value.ToDouble();
371 else
372 vecString[j] = value.ToString();
373 }
374 std::string name{array->GetName()};
375 if(dataArray)
376 addCustomArray(name, vec);
377 else
378 addCustomStringArray(name, vecString);
379 }
380 }
381 void copyPointData(vtkUnstructuredGrid *treeNodes) {
382 std::vector<int> nodeCorrT;
383 copyPointData(treeNodes, nodeCorrT);
384 }
385
386 // Filled by the algorithm
387 std::vector<std::vector<SimplexId>> getNodeCorr() {
388 return nodeCorr;
389 }
390 std::vector<double> getClusterShift() {
391 return clusterShift;
392 }
393 double getPrevXMax() {
394 return prevXMax;
395 }
396
397 // ==========================================================================
398 // Matching Visualization
399 // ==========================================================================
400 template <class dataType>
402 std::vector<FTMTree_MT *> trees{tree1, tree2};
403 std::vector<FTMTree_MT *> barycenters;
404
405 makeMatchingOutput<dataType>(trees, barycenters);
406 }
407
408 template <class dataType>
409 void makeMatchingOutput(std::vector<FTMTree_MT *> &trees,
410 std::vector<FTMTree_MT *> &barycenters) {
411 int numInputs = trees.size();
412 int NumberOfBarycenters = barycenters.size();
413 bool const clusteringOutput = (NumberOfBarycenters != 0);
414 NumberOfBarycenters
415 = std::max(NumberOfBarycenters, 1); // to always enter the outer loop
416 if(not clusteringOutput)
417 numInputs = 1;
418
419 vtkNew<vtkUnstructuredGrid> vtkMatching{};
420 vtkNew<vtkPoints> pointsM{};
421
422 // Fields
423 vtkNew<vtkIntArray> matchingID{};
424 matchingID->SetName("MatchingID");
425 vtkNew<vtkIntArray> matchingType{};
426 matchingType->SetName("MatchingType");
427 vtkNew<vtkDoubleArray> matchPers{};
428 matchPers->SetName("MeanMatchedPersistence");
429 vtkNew<vtkDoubleArray> costArray{};
430 costArray->SetName("Cost");
431 vtkNew<vtkIntArray> tree1NodeIdField{};
432 tree1NodeIdField->SetName("tree1NodeId");
433 vtkNew<vtkIntArray> tree2NodeIdField{};
434 tree2NodeIdField->SetName("tree2NodeId");
435
436 vtkNew<vtkFloatArray> matchingPercentMatch{};
437 matchingPercentMatch->SetName("MatchingPercentMatch");
438
439 // Iterate through clusters and trees
440 printMsg(
441 "// Iterate through clusters and trees", ttk::debug::Priority::VERBOSE);
442 int count = 0;
443 for(int c = 0; c < NumberOfBarycenters; ++c) {
444 for(int i = 0; i < numInputs; ++i) {
445 if((printTreeId == -1 and printClusterId != -1 and c != printClusterId)
446 or (printTreeId != -1 and printClusterId == -1 and i != printTreeId)
447 or (printTreeId != -1 and printClusterId != -1
448 and (c != printClusterId or i != printTreeId)))
449 continue;
450 for(std::tuple<idNode, idNode, double> match :
451 outputMatchingBarycenter[c][i]) {
452 vtkIdType pointIds[2];
453 idNode tree1NodeId = std::get<0>(match);
454 idNode tree2NodeId = std::get<1>(match);
455 double const cost = std::get<2>(match);
456 FTMTree_MT *tree1;
457 FTMTree_MT *tree2 = trees[i];
458 if(not clusteringOutput) {
459 tree1 = trees[0];
460 tree2 = trees[1];
461 } else
462 tree1 = barycenters[c];
463
464 // Get first point
465 printMsg("// Get first point", ttk::debug::Priority::VERBOSE);
466 SimplexId const pointToGet1 = clusteringOutput
467 ? nodeCorr2[c][tree1NodeId]
468 : nodeCorr1[0][tree1NodeId];
469 double *point1 = vtkOutputNode2->GetPoints()->GetPoint(pointToGet1);
470 const SimplexId nextPointId1 = pointsM->InsertNextPoint(point1);
471 pointIds[0] = nextPointId1;
472
473 // Get second point
474 printMsg("// Get second point", ttk::debug::Priority::VERBOSE);
475 SimplexId const pointToGet2 = clusteringOutput
476 ? nodeCorr1[i][tree2NodeId]
477 : nodeCorr1[1][tree2NodeId];
478 double *point2 = vtkOutputNode1->GetPoints()->GetPoint(pointToGet2);
479 const SimplexId nextPointId2 = pointsM->InsertNextPoint(point2);
480 pointIds[1] = nextPointId2;
481
482 // Add cell
484 vtkMatching->InsertNextCell(VTK_LINE, 2, pointIds);
485
486 // Add arc matching percentage
487 printMsg(
488 "// Add arc matching percentage", ttk::debug::Priority::VERBOSE);
489 if(allBaryPercentMatch.size() != 0)
490 matchingPercentMatch->InsertNextTuple1(
491 allBaryPercentMatch[c][tree1NodeId]);
492
493 // Add tree1 and tree2 node ids
494 printMsg(
495 "// Add tree1 and tree2 node ids", ttk::debug::Priority::VERBOSE);
496 tree1NodeIdField->InsertNextTuple1(pointToGet1);
497 tree2NodeIdField->InsertNextTuple1(pointToGet2);
498
499 // Add matching ID
500 matchingID->InsertNextTuple1(count);
501
502 // Add matching type
503 printMsg("// Add matching type", ttk::debug::Priority::VERBOSE);
504 int thisType = 0;
505 int const tree1NodeDown
506 = tree1->getNode(tree1NodeId)->getNumberOfDownSuperArcs();
507 int const tree1NodeUp
508 = tree1->getNode(tree1NodeId)->getNumberOfUpSuperArcs();
509 int const tree2NodeDown
510 = tree2->getNode(tree2NodeId)->getNumberOfDownSuperArcs();
511 int const tree2NodeUp
512 = tree2->getNode(tree2NodeId)->getNumberOfUpSuperArcs();
513 if(tree1NodeDown != 0 and tree1NodeUp != 0 and tree2NodeDown != 0
514 and tree2NodeUp != 0)
515 thisType = 1; // Saddle to Saddle
516 if(tree1NodeDown == 0 and tree1NodeUp != 0 and tree2NodeDown == 0
517 and tree2NodeUp != 0)
518 thisType = 2; // Leaf to leaf
519 if(tree1NodeDown != 0 and tree1NodeUp == 0 and tree2NodeDown != 0
520 and tree2NodeUp == 0)
521 thisType = 3; // Root to root
522 matchingType->InsertNextTuple1(thisType);
523
524 // Add mean matched persistence
525 printMsg(
526 "// Add mean matched persistence", ttk::debug::Priority::VERBOSE);
527 double tree1Pers = tree1->getNodePersistence<dataType>(tree1NodeId);
528 double tree2Pers = tree2->getNodePersistence<dataType>(tree2NodeId);
529 double const meanPersistence = (tree1Pers + tree2Pers) / 2;
530 matchPers->InsertNextTuple1(meanPersistence);
531
532 // Add cost
533 costArray->InsertNextTuple1(cost);
534
535 count++;
536 }
537 }
538 }
539 vtkMatching->SetPoints(pointsM);
540 vtkMatching->GetCellData()->AddArray(matchingType);
541 vtkMatching->GetCellData()->AddArray(matchPers);
542 vtkMatching->GetCellData()->AddArray(matchingID);
543 vtkMatching->GetCellData()->AddArray(costArray);
544 vtkMatching->GetCellData()->AddArray(tree1NodeIdField);
545 vtkMatching->GetCellData()->AddArray(tree2NodeIdField);
546 if(allBaryPercentMatch.size() != 0)
547 vtkMatching->GetCellData()->AddArray(matchingPercentMatch);
548 vtkOutputMatching->ShallowCopy(vtkMatching);
549 }
550
551 // ==========================================================================
552 // Trees Visualization
553 // ==========================================================================
554 template <class dataType>
556 std::vector<FTMTree_MT *> trees{tree1};
557
558 makeTreesOutput<dataType>(trees);
559 }
560
561 template <class dataType>
562 void makeTreesOutput(FTMTree_MT *tree1, FTMTree_MT *tree2) {
563 std::vector<FTMTree_MT *> trees{tree1, tree2};
564
565 makeTreesOutput<dataType>(trees);
566 }
567
568 template <class dataType>
569 void makeTreesOutput(std::vector<FTMTree_MT *> &trees) {
570 std::vector<FTMTree_MT *> barycenters;
571 clusteringAssignment.clear();
572 clusteringAssignment.resize(trees.size(), 0);
573
574 makeTreesOutput<dataType>(trees, barycenters);
575 }
576
577 template <class dataType>
578 void makeTreesOutput(std::vector<FTMTree_MT *> &trees,
579 std::vector<FTMTree_MT *> &barycenters) {
580 int numInputs = trees.size();
581 int const numInputsOri = numInputs;
582 int NumberOfBarycenters = barycenters.size();
583 bool const clusteringOutput = (NumberOfBarycenters != 0);
584 NumberOfBarycenters
585 = std::max(NumberOfBarycenters, 1); // to always enter the outer loop
586 bool const embeddedDiagram = not PlanarLayout and isPersistenceDiagram;
587
588 // TreeNodeIdRev
589 for(int i = 0; i < numInputs; ++i) {
590 if(i < (int)treesNodes.size() and treesNodes[i]) {
591 auto treeNodeIdArray
592 = treesNodes[i]->GetPointData()->GetArray("TreeNodeId");
593 if(treeNodeIdArray) {
594 std::vector<int> treeNodeIdRev;
595 getTreeNodeIdRev(treeNodeIdArray, treeNodeIdRev);
596 for(unsigned int j = 0; j < treesNodeCorrMesh[i].size(); ++j)
597 treesNodeCorrMesh[i][j] = treeNodeIdRev[treesNodeCorrMesh[i][j]];
598 }
599 }
600 }
601
602 // Bounds
603 printMsg("Bounds and branching", ttk::debug::Priority::VERBOSE);
604 std::vector<std::tuple<double, double, double, double, double, double>>
605 allBounds(numInputs);
606 for(int i = 0; i < numInputs; ++i) {
607 if(OutputSegmentation and treesSegmentation[i]) {
608 double *tBounds = treesSegmentation[i]->GetBounds();
609 allBounds[i] = std::make_tuple(tBounds[0], tBounds[1], tBounds[2],
610 tBounds[3], tBounds[4], tBounds[5]);
611 } else if(treesNodes.size() != 0 and treesNodes[i] != nullptr) {
612 if(not isPersistenceDiagram)
613 allBounds[i]
614 = getRealBounds(treesNodes[i], trees[i], treesNodeCorrMesh[i]);
615 else {
616 double bounds[6];
617 treesNodes[i]->GetBounds(bounds);
618 allBounds[i] = std::make_tuple(
619 bounds[0], bounds[1], bounds[2], bounds[3], bounds[4], bounds[5]);
620 }
621 } else {
622 allBounds[i] = allBounds[0];
623 }
624 }
625
626 std::vector<std::tuple<double, double, double, double, double, double>>
627 allBaryBounds(barycenters.size());
628 std::vector<std::vector<idNode>> allBaryBranching(barycenters.size());
629 std::vector<std::vector<int>> allBaryBranchingID(barycenters.size());
630 for(size_t c = 0; c < barycenters.size(); ++c) {
631 allBaryBounds[c] = getMaximalBounds(allBounds, clusteringAssignment, c);
632 if(not isPersistenceDiagram)
633 barycenters[c]->getTreeBranching(
634 allBaryBranching[c], allBaryBranchingID[c]);
635 }
636 if(not clusteringOutput)
637 allBaryBounds.emplace_back(
638 getMaximalBounds(allBounds, clusteringAssignment, 0));
639
640 // ----------------------------------------------------------------------
641 // Make Trees Output
642 // ----------------------------------------------------------------------
643 printMsg("--- Make Trees Output", ttk::debug::Priority::VERBOSE);
644 std::vector<FTMTree_MT *> const treesOri(trees);
645 if(ShiftMode == 1) { // Star Barycenter
646 trees.clear();
647 clusteringAssignment.clear();
648 for(unsigned int j = 0; j < barycenters.size(); ++j) {
649 trees.emplace_back(barycenters[j]);
650 clusteringAssignment.emplace_back(j);
651 }
652 numInputs = trees.size();
653 }
654 // - Declare VTK arrays
655 vtkNew<vtkUnstructuredGrid> vtkArcs{};
656 vtkNew<vtkPoints> points{};
657
658 // Node fields
659 vtkNew<vtkIntArray> criticalType{};
660 criticalType->SetName(ttk::PersistenceCriticalTypeName);
661 vtkNew<vtkDoubleArray> persistenceNode{};
662 persistenceNode->SetName(ttk::PersistenceName);
663 vtkNew<vtkIntArray> clusterIDNode{};
664 clusterIDNode->SetName("ClusterID");
665 vtkNew<vtkIntArray> isDummyNode{};
666 isDummyNode->SetName("isDummyNode");
667 vtkNew<vtkIntArray> branchNodeID{};
668 branchNodeID->SetName("BranchNodeID");
669 vtkNew<vtkFloatArray> scalar{};
670 scalar->SetName("Scalar");
671 vtkNew<vtkIntArray> isImportantPairsNode{};
672 isImportantPairsNode->SetName("isImportantPair");
673 vtkNew<vtkIntArray> nodeID{};
674 nodeID->SetName("NodeId"); // Simplex Id
675 vtkNew<vtkIntArray> trueNodeID{};
676 trueNodeID->SetName("TrueNodeId");
677 vtkNew<vtkIntArray> vertexID{};
678 vertexID->SetName(
679 (isPersistenceDiagram ? ttk::VertexScalarFieldName : "VertexId"));
680
681 vtkNew<vtkIntArray> treeIDNode{};
682 treeIDNode->SetName("TreeID");
683 vtkNew<vtkIntArray> branchBaryNodeID{};
684 branchBaryNodeID->SetName("BranchBaryNodeID");
685 vtkNew<vtkIntArray> isInterpolatedTreeNode{};
686 isInterpolatedTreeNode->SetName("isInterpolatedTree");
687
688 vtkNew<vtkFloatArray> percentMatch{};
689 percentMatch->SetName("PercentMatchNode");
690 vtkNew<vtkFloatArray> persistenceBaryNode{};
691 persistenceBaryNode->SetName("PersistenceBarycenter");
692 vtkNew<vtkIntArray> persistenceBaryOrderNode{};
693 persistenceBaryOrderNode->SetName("PersistenceBarycenterOrder");
694
695 vtkNew<vtkDoubleArray> pairBirthNode{};
696 pairBirthNode->SetName(ttk::PersistenceBirthName);
697
698 vtkNew<vtkFloatArray> treeNodeId{};
699 treeNodeId->SetName("TreeNodeId");
700
701 vtkNew<vtkFloatArray> treeNodeIdOrigin{};
702 treeNodeIdOrigin->SetName("TreeNodeIdOrigin");
703 vtkNew<vtkDoubleArray> coordinates{};
704 coordinates->SetName(ttk::PersistenceCoordinatesName);
705 coordinates->SetNumberOfComponents(3);
706
707 vtkNew<vtkIntArray> isMultiPersPairNode{};
708 isMultiPersPairNode->SetName("isMultiPersPairNode");
709
710 std::vector<std::vector<double>> customArraysValues(customArrays.size());
711 std::vector<std::vector<int>> customIntArraysValues(customIntArrays.size());
712 std::vector<std::vector<std::string>> customStringArraysValues(
713 customStringArrays.size());
714
715 // Arc fields
716 vtkNew<vtkDoubleArray> persistenceArc{};
717 persistenceArc->SetName(ttk::PersistenceName);
718 vtkNew<vtkIntArray> clusterIDArc{};
719 clusterIDArc->SetName("ClusterID");
720 vtkNew<vtkIntArray> isImportantPairsArc{};
721 isImportantPairsArc->SetName("isImportantPair");
722 vtkNew<vtkIntArray> isDummyArc{};
723 isDummyArc->SetName("isDummyArc");
724 vtkNew<vtkIntArray> branchID{};
725 branchID->SetName("BranchID");
726 vtkNew<vtkIntArray> upNodeId{};
727 upNodeId->SetName("upNodeId");
728 vtkNew<vtkIntArray> downNodeId{};
729 downNodeId->SetName("downNodeId");
730
731 vtkNew<vtkIntArray> treeIDArc{};
732 treeIDArc->SetName((isPersistenceDiagram ? "DiagramID" : "TreeID"));
733 vtkNew<vtkIntArray> branchBaryID{};
734 branchBaryID->SetName("BranchBaryNodeID");
735 vtkNew<vtkIntArray> isInterpolatedTreeArc{};
736 isInterpolatedTreeArc->SetName("isInterpolatedTree");
737
738 vtkNew<vtkFloatArray> percentMatchArc{};
739 percentMatchArc->SetName("PercentMatchArc");
740 vtkNew<vtkFloatArray> persistenceBaryArc{};
741 persistenceBaryArc->SetName("PersistenceBarycenter");
742 vtkNew<vtkIntArray> persistenceBaryOrderArc{};
743 persistenceBaryOrderArc->SetName("PersistenceBarycenterOrder");
744
745 vtkNew<vtkIntArray> pairIdentifier{};
746 pairIdentifier->SetName(ttk::PersistencePairIdentifierName);
747 vtkNew<vtkIntArray> pairType{};
748 pairType->SetName(ttk::PersistencePairTypeName);
749 vtkNew<vtkIntArray> pairIsFinite{};
750 pairIsFinite->SetName(ttk::PersistenceIsFinite);
751 vtkNew<vtkDoubleArray> pairBirth{};
752 pairBirth->SetName(ttk::PersistenceBirthName);
753
754 vtkNew<vtkIntArray> isMultiPersPairArc{};
755 isMultiPersPairArc->SetName("isMultiPersPairArc");
756
757 std::vector<std::vector<double>> customCellArraysValues(
758 customArrays.size());
759 std::vector<std::vector<int>> customCellIntArraysValues(
760 customIntArrays.size());
761 std::vector<std::vector<std::string>> customCellStringArraysValues(
762 customStringArrays.size());
763
764 // Segmentation
765 vtkNew<vtkAppendFilter> appendFilter{};
766
767 // Internal data
768 int cellCount = 0;
769 int pointCount = 0;
770 bool foundOneInterpolatedTree = false;
771 nodeCorr.clear();
772 nodeCorr.resize(numInputs);
773 clusterShift.clear();
774 clusterShift.resize(NumberOfBarycenters, 0);
775 allBaryPercentMatch.clear();
776 allBaryPercentMatch.resize(NumberOfBarycenters);
777
778 // --------------------------------------------------------
779 // Iterate through all clusters
780 // --------------------------------------------------------
781 printMsg("Iterate through all clusters", ttk::debug::Priority::VERBOSE);
782 double const importantPairsOriginal = importantPairs_;
783 for(int c = 0; c < NumberOfBarycenters; ++c) {
784
785 // Get persistence order
786 std::vector<int> baryPersistenceOrder;
787 if(clusteringOutput and ShiftMode != 1) {
788 baryPersistenceOrder.resize(barycenters[c]->getNumberOfNodes(), -1);
789 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, dataType>>
790 pairsBary;
791 barycenters[c]->getPersistencePairsFromTree<dataType>(pairsBary, false);
792 for(unsigned int j = 0; j < pairsBary.size(); ++j) {
793 int const index = pairsBary.size() - 1 - j;
794 baryPersistenceOrder[std::get<0>(pairsBary[j])] = index;
795 baryPersistenceOrder[std::get<1>(pairsBary[j])] = index;
796 }
797 }
798
799 // Get radius
800 printMsg("// Get radius", ttk::debug::Priority::VERBOSE);
801 double delta_max = 1.0;
802 int noSample = 0 + noSampleOffset;
803 for(int i = 0; i < numInputsOri; ++i) {
804 delta_max = std::max(
805 (std::get<3>(allBounds[i]) - std::get<2>(allBounds[i])), delta_max);
806 delta_max = std::max(
807 (std::get<1>(allBounds[i]) - std::get<0>(allBounds[i])), delta_max);
808 if(clusteringAssignment[i] != c)
809 continue;
810 noSample += 1;
811 }
812 double const radius = delta_max * 2 * DimensionSpacing;
813 int iSample = 0 + iSampleOffset - 1;
814
815 if(c < NumberOfBarycenters - 1)
816 clusterShift[c + 1] = radius * 4 + clusterShift[c];
817
818 // Line/Double line attributes
819 prevXMax = 0 + prevXMaxOffset;
820 std::vector<double> allPrevXMax;
821 double prevYMax = std::numeric_limits<double>::lowest();
822
823 // ------------------------------------------
824 // Iterate through all trees of this cluster
825 // ------------------------------------------
826 printMsg("Iterate through all trees of this cluster",
828 for(int i = 0; i < numInputs; ++i) {
829 if(clusteringAssignment[i] != c)
830 continue;
831
832 iSample += 1;
833
834 if((printTreeId == -1 and printClusterId != -1 and c != printClusterId)
835 or (printTreeId != -1 and printClusterId == -1 and i != printTreeId)
836 or (printTreeId != -1 and printClusterId != -1
837 and (c != printClusterId or i != printTreeId)))
838 continue;
839
840 // Manage important pairs threshold
841 importantPairs_ = importantPairsOriginal;
842 if(MaximumImportantPairs > 0 or MinimumImportantPairs > 0) {
843 std::vector<std::tuple<idNode, idNode, dataType>> pairs;
844 trees[i]->getPersistencePairsFromTree(pairs, false);
845 if(MaximumImportantPairs > 0) {
846 int firstIndex = pairs.size() - MaximumImportantPairs;
847 firstIndex
848 = std::max(std::min(firstIndex, int(pairs.size()) - 1), 0);
849 double tempThreshold = 0.999 * std::get<2>(pairs[firstIndex])
850 / std::get<2>(pairs[pairs.size() - 1]);
851 tempThreshold *= 100;
852 importantPairs_ = std::max(importantPairs_, tempThreshold);
853 }
854 if(MinimumImportantPairs > 0) {
855 int firstIndex = pairs.size() - MinimumImportantPairs;
856 firstIndex
857 = std::max(std::min(firstIndex, int(pairs.size()) - 1), 0);
858 double tempThreshold = 0.999 * std::get<2>(pairs[firstIndex])
859 / std::get<2>(pairs[pairs.size() - 1]);
860 tempThreshold *= 100;
861 importantPairs_ = std::min(importantPairs_, tempThreshold);
862 }
863 }
864
865 // Get is interpolated tree (temporal subsampling)
866 bool isInterpolatedTree = false;
867 if(interpolatedTrees.size() != 0)
868 isInterpolatedTree = interpolatedTrees[i];
869 foundOneInterpolatedTree |= isInterpolatedTree;
870
871 // Get branching
872 printMsg("// Get branching", ttk::debug::Priority::VERBOSE);
873 std::vector<idNode> treeBranching;
874 std::vector<int> treeBranchingID;
875 if(not isPersistenceDiagram)
876 trees[i]->getTreeBranching(treeBranching, treeBranchingID);
877
878 // Get shift
880 double const angle = 360.0 / noSample * iSample;
881 double const pi = M_PI;
882 double diff_x = 0, diff_y = 0;
883 double const alphaShift
884 = BarycenterPositionAlpha ? (-radius + 2 * radius * Alpha) * -1 : 0;
885 switch(ShiftMode) {
886 case -1:
887 diff_x = 0.0;
888 diff_y = 0.0;
889 break;
890 case 0: // Star
891 diff_x
892 = -1 * radius * std::cos(-1 * angle * pi / 180) + clusterShift[c];
893 diff_y = -1 * radius * std::sin(-1 * angle * pi / 180);
894 break;
895 case 1: // Star Barycenter
896 diff_x = clusterShift[c] + alphaShift;
897 diff_y = 0;
898 break;
899 case 2: // Line
900 diff_x = prevXMax + radius;
901 break;
902 case 3: // Double Line
903 diff_x = prevXMax + radius;
904 if(i >= numInputs / 2) {
905 diff_y = -(prevYMax + radius / 2);
906 diff_x = allPrevXMax[i - int(numInputs / 2)] + radius;
907 } else
908 allPrevXMax.emplace_back(prevXMax);
909 break;
910 default:
911 break;
912 }
913
914 // Get dimension shift
915 printMsg("// Get dimension shift", ttk::debug::Priority::VERBOSE);
916 double diff_z = PlanarLayout ? 0 : -std::get<4>(allBounds[i]);
917 // TODO DimensionToShift for Planar Layout
918 if(not PlanarLayout)
919 if(DimensionToShift != 0) { // is not X
920 if(DimensionToShift == 2) // is Z
921 diff_z = diff_x;
922 else if(DimensionToShift == 1) // is Y
923 diff_y = diff_x;
924 diff_x = -std::get<0>(allBounds[i]);
925 }
926
927 // Planar layout
928 printMsg("// Planar Layout", ttk::debug::Priority::VERBOSE);
929 std::vector<float> layout;
930 if(PlanarLayout) {
931 double refPersistence;
932 if(clusteringOutput)
933 refPersistence = barycenters[0]->getNodePersistence<dataType>(
934 barycenters[0]->getRoot());
935 else
936 refPersistence
937 = trees[0]->getNodePersistence<dataType>(trees[0]->getRoot());
938 if(not isPersistenceDiagram)
939 treePlanarLayout<dataType>(
940 trees[i], allBaryBounds[c], refPersistence, layout);
941 else {
942 persistenceDiagramPlanarLayout<dataType>(trees[i], layout);
943 }
944 }
945
946 // Internal arrays
947 printMsg("// Internal arrays", ttk::debug::Priority::VERBOSE);
948 int cptNode = 0;
949 nodeCorr[i].resize(trees[i]->getNumberOfNodes());
950 std::vector<SimplexId> treeSimplexId(trees[i]->getNumberOfNodes());
951 std::vector<SimplexId> treeDummySimplexId(trees[i]->getNumberOfNodes());
952 std::vector<SimplexId> layoutCorr(trees[i]->getNumberOfNodes());
953 std::vector<idNode> treeMatching(trees[i]->getNumberOfNodes(), -1);
954 if(clusteringOutput and ShiftMode != 1)
955 for(auto match : outputMatchingBarycenter[c][i])
956 treeMatching[std::get<1>(match)] = std::get<0>(match);
957 // _ m[i][j] contains the node in treesOri[j] matched to the node i in
958 // the barycenter
959 std::vector<std::vector<idNode>> baryMatching(
960 trees[i]->getNumberOfNodes(),
961 std::vector<idNode>(
962 numInputsOri, std::numeric_limits<idNode>::max()));
963 if(ShiftMode == 1) {
964 for(size_t j = 0; j < outputMatchingBarycenter[c].size(); ++j)
965 for(auto match : outputMatchingBarycenter[c][j])
966 baryMatching[std::get<0>(match)][j] = std::get<1>(match);
967 allBaryPercentMatch[c].resize(trees[i]->getNumberOfNodes(), 100.0);
968 }
969 double minBirth = std::numeric_limits<double>::max(),
970 maxBirth = std::numeric_limits<double>::lowest();
971 SimplexId minBirthNode = 0, maxBirthNode = 0;
972
973 // ----------------------------
974 // Tree traversal
975 // ----------------------------
976 printMsg("// Tree traversal", ttk::debug::Priority::VERBOSE);
977 std::queue<idNode> queue;
978 queue.emplace(trees[i]->getRoot());
979 while(!queue.empty()) {
980 idNode const node = queue.front();
981 queue.pop();
982 idNode const nodeOrigin = trees[i]->getNode(node)->getOrigin();
983
984 // Push children to the queue
985 printMsg(
986 "// Push children to the queue", ttk::debug::Priority::VERBOSE);
987 std::vector<idNode> children;
988 trees[i]->getChildren(node, children);
989 for(auto child : children)
990 queue.emplace(child);
991
992 // --------------
993 // Insert point
994 // --------------
995 auto getPoint
996 = [&](vtkUnstructuredGrid *vtu, int pointID, double(&point)[3]) {
997 if(not vtu)
998 return;
999 if(not isPersistenceDiagram or convertedToDiagram) {
1000 double *pointTemp = vtu->GetPoints()->GetPoint(pointID);
1001 for(int k = 0; k < 3; ++k)
1002 point[k] += pointTemp[k];
1003 } else {
1004 for(int k = 0; k < 3; ++k) {
1005 auto array = vtu->GetPointData()->GetArray(
1007 if(array)
1008 point[k] += array->GetComponent(pointID, k);
1009 }
1010 }
1011 };
1012
1013 printMsg("// Get and insert point", ttk::debug::Priority::VERBOSE);
1014 int nodeMesh = -1;
1015 int nodeMeshTreeIndex = -1;
1016 double noMatched = 0.0;
1017 double point[3] = {0, 0, 0};
1018 if(ShiftMode == 1) { // Star barycenter
1019 for(int j = 0; j < numInputsOri; ++j) {
1020 if(baryMatching[node][j] != std::numeric_limits<idNode>::max()) {
1021 nodeMesh = treesNodeCorrMesh[j][baryMatching[node][j]];
1022 if(not PlanarLayout)
1023 getPoint(treesNodes[j], nodeMesh, point);
1024 noMatched += 1;
1025 nodeMeshTreeIndex = j;
1026 }
1027 }
1028 for(int k = 0; k < 3; ++k)
1029 point[k] /= noMatched;
1030 } else if(not isInterpolatedTree and treesNodes.size() != 0
1031 and treesNodes[i] != nullptr) {
1032 nodeMesh = treesNodeCorrMesh[i][node];
1033 if(not PlanarLayout)
1034 getPoint(treesNodes[i], nodeMesh, point);
1035 }
1036 if(PlanarLayout) {
1037 layoutCorr[node] = cptNode;
1038 point[0] = layout[cptNode];
1039 point[1] = layout[cptNode + 1];
1040 point[2] = 0;
1041 cptNode += 2;
1042 }
1043 point[0] += diff_x;
1044 point[1] += diff_y;
1045 point[2] += diff_z;
1046
1047 // Bary percentage matching
1048 if(ShiftMode == 1) { // Star Barycenter
1049 float const percentMatchT = noMatched * 100 / numInputs;
1050 allBaryPercentMatch[c][node] = percentMatchT;
1051 }
1052
1053 // Get x Max and y Min for next iteration if needed (double line mode)
1054 prevXMax = std::max(prevXMax, point[0]);
1055 if(ShiftMode == 3) { // Double line
1056 if(i < numInputs / 2)
1057 prevYMax = std::max(prevYMax, point[1]);
1058 if(i == int(numInputs / 2) - 1)
1059 prevXMax = 0;
1060 }
1061
1062 // TODO too many dummy nodes are created
1063 bool dummyNode
1064 = PlanarLayout and not branchDecompositionPlanarLayout_
1065 and (!trees[i]->isRoot(node) or isPersistenceDiagram);
1066 dummyNode = dummyNode or embeddedDiagram;
1067 if(dummyNode) {
1068 double pointToAdd[3] = {0, 0, 0};
1069 if(embeddedDiagram) {
1070 if(i < (int)treesNodes.size()
1071 and i < (int)treesNodeCorrMesh.size()
1072 and nodeOrigin < treesNodeCorrMesh[i].size())
1073 getPoint(
1074 treesNodes[i], treesNodeCorrMesh[i][nodeOrigin], pointToAdd);
1075 } else {
1076 if(not isPersistenceDiagram) {
1077 // will be modified when processing son
1078 std::copy(
1079 std::begin(point), std::end(point), std::begin(pointToAdd));
1080 } else {
1081 double pdPoint[3] = {
1082 point[0],
1083 point[1]
1084 - (layout[layoutCorr[node] + 1] - layout[layoutCorr[node]]),
1085 0};
1086 std::copy(std::begin(pdPoint), std::end(pdPoint),
1087 std::begin(pointToAdd));
1088 }
1089 }
1090 treeDummySimplexId[node] = points->InsertNextPoint(pointToAdd);
1091 if(not embeddedDiagram) {
1092 if(isPersistenceDiagram) {
1093 if(layout[layoutCorr[node]] < minBirth) {
1094 minBirth = layout[layoutCorr[node]];
1095 minBirthNode = treeDummySimplexId[node];
1096 }
1097 if(layout[layoutCorr[node]] > maxBirth) {
1098 maxBirth = layout[layoutCorr[node]];
1099 maxBirthNode = treeDummySimplexId[node];
1100 }
1101 }
1102 }
1103 }
1104 SimplexId const nextPointId = points->InsertNextPoint(point);
1105 treeSimplexId[node] = nextPointId;
1106 nodeCorr[i][node] = nextPointId;
1107 if(dummyNode)
1108 nodeCorr[i][node] = treeDummySimplexId[node];
1109 if(isPersistenceDiagram)
1110 nodeCorr[i][node] = nextPointId;
1111
1112 idNode const nodeBranching
1113 = ((PlanarLayout and branchDecompositionPlanarLayout_)
1114 or isPersistenceDiagram
1115 ? node
1116 : treeBranching[node]);
1117
1118 // --------------
1119 // Insert cell connecting parent
1120 // --------------
1121 printMsg(
1122 "// Add cell connecting parent", ttk::debug::Priority::VERBOSE);
1123 if(!trees[i]->isRoot(node) or isPersistenceDiagram) {
1124 vtkIdType pointIds[2];
1125 pointIds[0] = treeSimplexId[node];
1126
1127 idNode const nodeParent = trees[i]->getParentSafe(node);
1128 // TODO too many dummy cells are created
1129 bool const dummyCell
1130 = PlanarLayout and not branchDecompositionPlanarLayout_
1131 and (node < treeBranching.size()
1132 and treeBranching[node] == nodeParent)
1133 and !trees[i]->isRoot(nodeParent) and not isPersistenceDiagram;
1134 if(isPersistenceDiagram) {
1135 pointIds[1] = treeDummySimplexId[node];
1136 } else if(PlanarLayout and branchDecompositionPlanarLayout_) {
1137 pointIds[1] = treeSimplexId[treeBranching[node]];
1138 } else if(dummyCell) {
1139 double dummyPoint[3]
1140 = {point[0], layout[layoutCorr[nodeParent] + 1] + diff_y,
1141 0. + diff_z};
1142 SimplexId const dummyPointId = treeDummySimplexId[nodeParent];
1143 points->SetPoint(dummyPointId, dummyPoint);
1144 vtkIdType dummyPointIds[2];
1145 dummyPointIds[0] = dummyPointId;
1146 dummyPointIds[1] = treeSimplexId[nodeParent];
1147 vtkArcs->InsertNextCell(VTK_LINE, 2, dummyPointIds);
1148 pointIds[1] = dummyPointId;
1149 } else
1150 pointIds[1] = treeSimplexId[nodeParent];
1151
1152 vtkArcs->InsertNextCell(VTK_LINE, 2, pointIds);
1153
1154 // --------------
1155 // Arc field
1156 // --------------
1157 int const toAdd = (dummyCell ? 2 : 1);
1158 for(int toAddT = 0; toAddT < toAdd; ++toAddT) {
1159 // Add arc matching percentage
1160 if(ShiftMode == 1) { // Star Barycenter
1161 auto nodeToGet
1162 = (!isPersistenceDiagram ? allBaryBranching[c][node] : node);
1163 percentMatchArc->InsertNextTuple1(
1164 allBaryPercentMatch[c][nodeToGet]);
1165 }
1166
1167 // Add branch bary ID
1168 printMsg(
1169 "// Push arc bary branch id", ttk::debug::Priority::VERBOSE);
1170 if(clusteringOutput and ShiftMode != 1) {
1171 int tBranchID = -1;
1172 auto nodeToGet = node;
1173 if(treeMatching[nodeToGet] < allBaryBranchingID[c].size())
1174 tBranchID = allBaryBranchingID[c][treeMatching[nodeToGet]];
1175 branchBaryID->InsertNextTuple1(tBranchID);
1176 }
1177
1178 // Add branch ID
1179 if(not isPersistenceDiagram) {
1180 int const tBranchID = treeBranchingID[node];
1181 branchID->InsertNextTuple1(tBranchID);
1182 }
1183
1184 // Add up and down nodeId
1185 if(not isPersistenceDiagram) {
1186 upNodeId->InsertNextTuple1(treeSimplexId[nodeParent]);
1187 downNodeId->InsertNextTuple1(treeSimplexId[node]);
1188 }
1189
1190 // Add arc persistence
1191 printMsg(
1192 "// Push arc persistence", ttk::debug::Priority::VERBOSE);
1193 idNode const nodeToGetPers
1194 = (isPersistenceDiagram ? node : nodeBranching);
1195 double const persToAdd
1196 = trees[i]->getNodePersistence<dataType>(nodeToGetPers);
1197 persistenceArc->InsertNextTuple1(persToAdd);
1198
1199 // Add birth
1200 auto birthDeath
1201 = trees[i]->getBirthDeath<dataType>(nodeToGetPers);
1202 pairBirth->InsertNextTuple1(std::get<0>(birthDeath));
1203
1204 // Add arc persistence barycenter and order
1205 if(clusteringOutput and ShiftMode != 1) {
1206 idNode const nodeToGet = nodeBranching;
1207 if(treeMatching[nodeToGet] < allBaryBranchingID[c].size()) {
1208 persistenceBaryArc->InsertTuple1(
1209 cellCount, barycenters[c]->getNodePersistence<dataType>(
1210 treeMatching[nodeToGet]));
1211 persistenceBaryOrderArc->InsertTuple1(
1212 cellCount, baryPersistenceOrder[treeMatching[nodeToGet]]);
1213 }
1214 }
1215
1216 // Add arc cluster ID
1217 clusterIDArc->InsertNextTuple1(clusteringAssignment[i]);
1218
1219 // Add arc tree ID
1220 treeIDArc->InsertNextTuple1(i + iSampleOffset);
1221
1222 // Add isImportantPair
1223 bool isImportant = false;
1224 idNode const nodeToGetImportance = nodeBranching;
1225 isImportant = trees[i]->isImportantPair<dataType>(
1226 nodeToGetImportance, importantPairs_,
1229 isImportantPairsArc->InsertNextTuple1(isImportant);
1230
1231 // Add isDummyArc
1232 bool const isDummy = toAdd == 2 and toAddT == 0;
1233 isDummyArc->InsertNextTuple1(isDummy);
1234
1235 // Add isInterpolatedTree
1236 isInterpolatedTreeArc->InsertNextTuple1(isInterpolatedTree);
1237
1238 // Add pairIdentifier
1239 pairIdentifier->InsertNextTuple1(treeSimplexId[node]);
1240
1241 // Add isMinMaxPair
1242 bool const isMinMaxPair
1243 = (trees[i]->isRoot(node) and not trees[i]->isLeaf(node))
1244 or (trees[i]->isRoot(nodeOrigin)
1245 and not trees[i]->isLeaf(nodeOrigin));
1246 pairIsFinite->InsertNextTuple1(!isMinMaxPair);
1247
1248 // Add pairType TODO
1249 pairType->InsertNextTuple1(0);
1250
1251 // Add isMultiPersPairArc
1252 bool const isMultiPersPair
1253 = (trees[i]->isMultiPersPair(nodeBranching)
1254 or trees[i]->isMultiPersPair(
1255 trees[i]->getNode(nodeBranching)->getOrigin()));
1256 isMultiPersPairArc->InsertNextTuple1(isMultiPersPair);
1257
1258 // Add custom point arrays to cells
1259 for(unsigned int ca = 0; ca < customArrays.size(); ++ca)
1260 customCellArraysValues[ca].push_back(
1261 std::get<1>(customArrays[ca])[nodeBranching]);
1262 for(unsigned int ca = 0; ca < customIntArrays.size(); ++ca)
1263 customCellIntArraysValues[ca].push_back(
1264 std::get<1>(customIntArrays[ca])[nodeBranching]);
1265 for(unsigned int ca = 0; ca < customStringArrays.size(); ++ca)
1266 customCellStringArraysValues[ca].push_back(
1267 std::get<1>(customStringArrays[ca])[nodeBranching]);
1268
1269 cellCount++;
1270 }
1271 }
1272
1273 // --------------
1274 // Node field
1275 // --------------
1276 int const toAdd = (dummyNode ? 2 : 1);
1277 for(int toAddT = 0; toAddT < toAdd; ++toAddT) {
1278 // Add node id
1279 nodeID->InsertNextTuple1(treeSimplexId[node]);
1280
1281 // Add trueNodeId
1282 trueNodeID->InsertNextTuple1(node);
1283
1284 // Add VertexId
1285 int nodeVertexId = -1;
1286 if(i < int(treesNodes.size()) and treesNodes[i]) {
1287 auto vertexIdArray = treesNodes[i]->GetPointData()->GetArray(
1288 (isPersistenceDiagram ? ttk::VertexScalarFieldName
1289 : "VertexId"));
1290 if(vertexIdArray and nodeMesh != -1)
1291 nodeVertexId = vertexIdArray->GetTuple1(nodeMesh);
1292 }
1293 vertexID->InsertNextTuple1(nodeVertexId);
1294
1295 // Add node scalar
1296 scalar->InsertNextTuple1(trees[i]->getValue<dataType>(node));
1297
1298 // Add criticalType
1299 printMsg("// Add criticalType", ttk::debug::Priority::VERBOSE);
1300 int criticalTypeT = -1;
1301 if(not isPersistenceDiagram) {
1302 if(not isInterpolatedTree) {
1303 if(ShiftMode == 1) {
1304 if(nodeMeshTreeIndex != -1) {
1305 auto array
1306 = treesNodes[nodeMeshTreeIndex]->GetPointData()->GetArray(
1307 "CriticalType");
1308 if(array)
1309 criticalTypeT = array->GetTuple1(nodeMesh);
1310 }
1311 } else if(treesNodes.size() != 0 and treesNodes[i] != nullptr) {
1312 auto array
1313 = treesNodes[i]->GetPointData()->GetArray("CriticalType");
1314 if(array)
1315 criticalTypeT = array->GetTuple1(nodeMesh);
1316 }
1317 } else {
1318 // TODO critical type for interpolated trees
1319 }
1320 } else {
1321 auto locMin = static_cast<int>(ttk::CriticalType::Local_minimum);
1322 auto saddle1 = static_cast<int>(ttk::CriticalType::Saddle1);
1323 auto locMax = static_cast<int>(ttk::CriticalType::Local_maximum);
1324 auto saddle2 = static_cast<int>(ttk::CriticalType::Saddle2);
1325 auto nodeIsRoot = trees[i]->isRoot(node);
1326 criticalTypeT
1327 = (toAddT == 1
1328 ? (isPDSadMax or nodeIsRoot ? locMax : saddle1)
1329 : (not isPDSadMax or nodeIsRoot ? locMin : saddle2));
1330 if(embeddedDiagram) {
1331 bool const nodeSup = trees[i]->getValue<dataType>(node)
1332 > trees[i]->getValue<dataType>(nodeOrigin);
1333 criticalTypeT
1334 = ((nodeSup and toAddT == 1) or (not nodeSup and toAddT == 0)
1335 ? (isPDSadMax or nodeIsRoot ? locMax : saddle1)
1336 : (not isPDSadMax or nodeIsRoot ? locMin : saddle2));
1337 }
1338 }
1339 criticalType->InsertNextTuple1(criticalTypeT);
1340
1341 // Add node matching percentage
1342 if(ShiftMode == 1) { // Star Barycenter
1343 percentMatch->InsertNextTuple1(allBaryPercentMatch[c][node]);
1344 }
1345
1346 // Add node branch bary id
1347 printMsg(
1348 "// Add node bary branch id", ttk::debug::Priority::VERBOSE);
1349 if(clusteringOutput and ShiftMode != 1) {
1350 int tBranchID = -1;
1351 if(treeMatching[node] < allBaryBranchingID[c].size()) {
1352 tBranchID = allBaryBranchingID[c][treeMatching[node]];
1353 if(!trees[i]->isLeaf(node)
1354 && treeMatching[nodeOrigin] < allBaryBranchingID[c].size())
1355 tBranchID = allBaryBranchingID[c][treeMatching[nodeOrigin]];
1356 }
1357 branchBaryNodeID->InsertNextTuple1(tBranchID);
1358 }
1359
1360 // Add node branch id
1361 if(not isPersistenceDiagram) {
1362 int tBranchID = treeBranchingID[node];
1363 if(not trees[i]->isLeaf(node))
1364 tBranchID = treeBranchingID[nodeOrigin];
1365 branchNodeID->InsertNextTuple1(tBranchID);
1366 }
1367
1368 // Add node persistence
1369 printMsg("// Push node persistence", ttk::debug::Priority::VERBOSE);
1370 persistenceNode->InsertNextTuple1(
1371 trees[i]->getNodePersistence<dataType>(node));
1372
1373 // Add birth
1374 auto birthDeath = trees[i]->getBirthDeath<dataType>(node);
1375 pairBirthNode->InsertNextTuple1(std::get<0>(birthDeath));
1376
1377 // Add node persistence barycenter
1378 if(clusteringOutput and ShiftMode != 1) {
1379 if(treeMatching[node] < allBaryBranchingID[c].size()) {
1380 persistenceBaryNode->InsertTuple1(
1381 pointCount, barycenters[c]->getNodePersistence<dataType>(
1382 treeMatching[node]));
1383 persistenceBaryOrderNode->InsertTuple1(
1384 pointCount, baryPersistenceOrder[treeMatching[node]]);
1385 }
1386 }
1387
1388 // Add node clusterID
1389 clusterIDNode->InsertNextTuple1(clusteringAssignment[i]);
1390
1391 // Add node tree ID
1392 treeIDNode->InsertNextTuple1(i + iSampleOffset);
1393
1394 // Add isDummyNode
1395 bool const isDummy
1396 = toAdd == 2 and toAddT == 1 and !trees[i]->isRoot(node);
1397 isDummyNode->InsertNextTuple1(isDummy);
1398
1399 // Add isInterpolatedTree
1400 isInterpolatedTreeNode->InsertNextTuple1(isInterpolatedTree);
1401
1402 // Add isImportantPair
1403 bool isImportant = false;
1404 isImportant = trees[i]->isImportantPair<dataType>(
1407 isImportantPairsNode->InsertNextTuple1(isImportant);
1408
1409 // Add treeNodeId
1410 treeNodeId->InsertNextTuple1(node);
1411
1412 // Add treeNodeIdOrigin
1413 treeNodeIdOrigin->InsertNextTuple1(nodeOrigin);
1414
1415 // Add coordinates
1416 printMsg("// Add coordinates", ttk::debug::Priority::VERBOSE);
1417 if(isPersistenceDiagram and !treesNodes.empty()
1418 and ShiftMode != 1) {
1419 double coord[3] = {0.0, 0.0, 0.0};
1420 getPoint(treesNodes[i], treesNodeCorrMesh[i][node], coord);
1421 coordinates->InsertNextTuple3(coord[0], coord[1], coord[2]);
1422 }
1423
1424 // Add isMultiPersPairArc
1425 printMsg("// isMultiPersPairArc", ttk::debug::Priority::VERBOSE);
1426 bool const isMultiPersPair
1427 = (trees[i]->isMultiPersPair(node)
1428 or trees[i]->isMultiPersPair(
1429 trees[i]->getNode(node)->getOrigin()));
1430 isMultiPersPairNode->InsertNextTuple1(isMultiPersPair);
1431
1432 // Add custom arrays
1433 for(unsigned int ca = 0; ca < customArrays.size(); ++ca)
1434 customArraysValues[ca].emplace_back(
1435 std::get<1>(customArrays[ca])[node]);
1436 for(unsigned int ca = 0; ca < customIntArrays.size(); ++ca)
1437 customIntArraysValues[ca].emplace_back(
1438 std::get<1>(customIntArrays[ca])[node]);
1439 for(unsigned int ca = 0; ca < customStringArrays.size(); ++ca)
1440 customStringArraysValues[ca].emplace_back(
1441 std::get<1>(customStringArrays[ca])[node]);
1442
1443 pointCount++;
1444 }
1445
1447 } // end tree traversal
1448
1449 // Add diagonal if isPersistenceDiagram
1450 if(isPersistenceDiagram and not embeddedDiagram) {
1451 vtkIdType pointIds[2];
1452 pointIds[0] = minBirthNode;
1453 pointIds[1] = maxBirthNode;
1454 vtkArcs->InsertNextCell(VTK_LINE, 2, pointIds);
1455 cellCount++;
1456
1457 pairIdentifier->InsertNextTuple1(-1);
1458 pairType->InsertNextTuple1(-1);
1459 persistenceArc->InsertNextTuple1(-1);
1460 pairIsFinite->InsertNextTuple1(0);
1461 pairBirth->InsertNextTuple1(0);
1462
1463 for(unsigned int ca = 0; ca < customArrays.size(); ++ca)
1464 customCellArraysValues[ca].push_back(-1);
1465 for(unsigned int ca = 0; ca < customIntArrays.size(); ++ca)
1466 customCellIntArraysValues[ca].push_back(-1);
1467 for(unsigned int ca = 0; ca < customStringArrays.size(); ++ca)
1468 customCellStringArraysValues[ca].emplace_back("");
1469
1470 isMultiPersPairArc->InsertNextTuple1(0);
1471 clusterIDArc->InsertNextTuple1(clusteringAssignment[i]);
1472 treeIDArc->InsertNextTuple1(i + iSampleOffset);
1473 isImportantPairsArc->InsertNextTuple1(0);
1474 branchBaryID->InsertNextTuple1(-1);
1475 percentMatchArc->InsertNextTuple1(100);
1476 }
1477
1478 // --------------
1479 // Manage segmentation
1480 // --------------
1481 // Use TransformFilter (see commit
1482 // 85600763a8907674b8e57d6ad77ca97640725b30) when issue #513 is
1483 // solved.
1484 printMsg("// Shift segmentation", ttk::debug::Priority::VERBOSE);
1485 if(OutputSegmentation and not PlanarLayout and treesSegmentation[i]) {
1486 vtkNew<vtkUnstructuredGrid> iTreesSegmentationCopy{};
1487 if(ShiftMode != -1)
1488 iTreesSegmentationCopy->DeepCopy(treesSegmentation[i]);
1489 else
1490 iTreesSegmentationCopy->ShallowCopy(treesSegmentation[i]);
1491 auto iVkOutputSegmentationTemp
1492 = vtkUnstructuredGrid::SafeDownCast(iTreesSegmentationCopy);
1493 if(!iVkOutputSegmentationTemp
1494 or !iVkOutputSegmentationTemp->GetPoints()) {
1495 printWrn("Convert segmentation to vtkUnstructuredGrid.");
1496 vtkNew<vtkAppendFilter> appendFilter2{};
1497 appendFilter2->AddInputData(treesSegmentation[i]);
1498 appendFilter2->Update();
1499 iVkOutputSegmentationTemp->ShallowCopy(appendFilter2->GetOutput());
1500 }
1501 if(ShiftMode != -1) {
1502 for(int p = 0;
1503 p < iVkOutputSegmentationTemp->GetPoints()->GetNumberOfPoints();
1504 ++p) {
1505 double *point
1506 = iVkOutputSegmentationTemp->GetPoints()->GetPoint(p);
1507 point[0] += diff_x;
1508 point[1] += diff_y;
1509 point[2] += diff_z;
1510 iVkOutputSegmentationTemp->GetPoints()->SetPoint(p, point);
1511 }
1512 }
1513 appendFilter->AddInputData(iVkOutputSegmentationTemp);
1514 }
1515 printMsg("// Shift segmentation DONE", ttk::debug::Priority::VERBOSE);
1516 }
1517 }
1518 for(int i = persistenceBaryNode->GetNumberOfTuples(); i < pointCount; ++i)
1519 persistenceBaryNode->InsertNextTuple1(0);
1520 for(int i = persistenceBaryArc->GetNumberOfTuples(); i < cellCount; ++i)
1521 persistenceBaryArc->InsertNextTuple1(0);
1522 for(int i = persistenceBaryOrderNode->GetNumberOfTuples(); i < pointCount;
1523 ++i)
1524 persistenceBaryOrderNode->InsertNextTuple1(0);
1525 for(int i = persistenceBaryOrderArc->GetNumberOfTuples(); i < cellCount;
1526 ++i)
1527 persistenceBaryOrderArc->InsertNextTuple1(0);
1528
1529 // --- Add VTK arrays to output
1530 printMsg("// Add VTK arrays to output", ttk::debug::Priority::VERBOSE);
1531 // - Manage node output
1532 // Custom arrays
1533 addVtkCustomArrays(customArrays, customArraysValues, vtkOutputNode, 0, 0);
1535 customIntArrays, customIntArraysValues, vtkOutputNode, 1, 0);
1537 customStringArrays, customStringArraysValues, vtkOutputNode, 2, 0);
1538
1539 // Classical arrays
1540 vtkOutputNode->SetPoints(points);
1541 vtkOutputNode->GetPointData()->AddArray(criticalType);
1542 vtkOutputNode->GetPointData()->AddArray(persistenceNode);
1543 vtkOutputNode->GetPointData()->AddArray(clusterIDNode);
1544 vtkOutputNode->GetPointData()->AddArray(treeIDNode);
1545 vtkOutputNode->GetPointData()->AddArray(trueNodeID);
1546 vtkOutputNode->GetPointData()->AddArray(vertexID);
1547 vtkOutputNode->GetPointData()->AddArray(isImportantPairsNode);
1548 vtkOutputNode->GetPointData()->AddArray(isMultiPersPairNode);
1549 if(not isPersistenceDiagram) {
1550 vtkOutputNode->GetPointData()->AddArray(nodeID);
1551 vtkOutputNode->GetPointData()->AddArray(branchNodeID);
1552 vtkOutputNode->GetPointData()->AddArray(isDummyNode);
1553 }
1554 if(not branchDecompositionPlanarLayout_ and not isPersistenceDiagram)
1555 vtkOutputNode->GetPointData()->AddArray(scalar);
1556 if(clusteringOutput and ShiftMode != 1) {
1557 vtkOutputNode->GetPointData()->AddArray(branchBaryNodeID);
1558 vtkOutputNode->GetPointData()->AddArray(persistenceBaryNode);
1559 vtkOutputNode->GetPointData()->AddArray(persistenceBaryOrderNode);
1560 }
1561 if(foundOneInterpolatedTree)
1562 vtkOutputNode->GetPointData()->AddArray(isInterpolatedTreeNode);
1563 if(ShiftMode == 1) // Star Barycenter
1564 vtkOutputNode->GetPointData()->AddArray(percentMatch);
1565 if(outputTreeNodeIndex)
1566 vtkOutputNode->GetPointData()->AddArray(treeNodeId);
1567 if(isPersistenceDiagram) {
1568 vtkOutputNode->GetPointData()->AddArray(treeNodeIdOrigin);
1569 if(!treesNodes.empty() and ShiftMode != 1)
1570 vtkOutputNode->GetPointData()->AddArray(coordinates);
1571 }
1572 vtkOutputNode->GetPointData()->AddArray(pairBirthNode);
1573
1574 // - Manage arc output
1575 // Custom arrays
1576 addVtkCustomArrays(customArrays, customCellArraysValues, vtkArcs, 0, 1);
1578 customIntArrays, customCellIntArraysValues, vtkArcs, 1, 1);
1580 customStringArrays, customCellStringArraysValues, vtkArcs, 2, 1);
1581
1582 // Classical arrays
1583 vtkArcs->SetPoints(points);
1584 vtkArcs->GetCellData()->AddArray(persistenceArc);
1585 vtkArcs->GetCellData()->AddArray(clusterIDArc);
1586 vtkArcs->GetCellData()->AddArray(treeIDArc);
1587 vtkArcs->GetCellData()->AddArray(isImportantPairsArc);
1588 vtkArcs->GetCellData()->AddArray(isMultiPersPairArc);
1589 if(not isPersistenceDiagram) {
1590 vtkArcs->GetCellData()->AddArray(isDummyArc);
1591 vtkArcs->GetCellData()->AddArray(branchID);
1592 vtkArcs->GetCellData()->AddArray(upNodeId);
1593 vtkArcs->GetCellData()->AddArray(downNodeId);
1594 }
1595 if(clusteringOutput and ShiftMode != 1) {
1596 vtkArcs->GetCellData()->AddArray(branchBaryID);
1597 vtkArcs->GetCellData()->AddArray(persistenceBaryArc);
1598 vtkArcs->GetCellData()->AddArray(persistenceBaryOrderArc);
1599 }
1600 if(foundOneInterpolatedTree)
1601 vtkArcs->GetCellData()->AddArray(isInterpolatedTreeArc);
1602 if(ShiftMode == 1) // Star Barycenter
1603 vtkArcs->GetCellData()->AddArray(percentMatchArc);
1604 if(isPersistenceDiagram) {
1605 vtkArcs->GetCellData()->AddArray(pairIdentifier);
1606 vtkArcs->GetCellData()->AddArray(pairType);
1607 vtkArcs->GetCellData()->AddArray(pairIsFinite);
1608 }
1609 vtkArcs->GetCellData()->AddArray(pairBirth);
1610 if(vtkOutputArc == vtkOutputNode)
1611 vtkArcs->GetPointData()->ShallowCopy(vtkOutputNode->GetPointData());
1612 if(not branchDecompositionPlanarLayout_ and not isPersistenceDiagram)
1613 vtkArcs->GetPointData()->AddArray(scalar);
1614 vtkOutputArc->ShallowCopy(vtkArcs);
1615
1616 // - Manage segmentation output
1617 if(OutputSegmentation and not PlanarLayout
1618 and appendFilter->GetNumberOfInputConnections(0) != 0) {
1619 appendFilter->SetMergePoints(false);
1620 appendFilter->Update();
1621 vtkOutputSegmentation->ShallowCopy(appendFilter->GetOutput());
1622 }
1623
1624 //
1625 if(ShiftMode == 1) // Star Barycenter
1626 trees = treesOri;
1627 }
1628
1629 // ==========================================================================
1630 // Bounds Utils
1631 // ==========================================================================
1632 std::tuple<double, double, double, double, double, double>
1633 getRealBounds(vtkUnstructuredGrid *treeNodes,
1634 FTMTree_MT *tree,
1635 std::vector<int> &nodeCorrT) {
1636 double x_min = std::numeric_limits<double>::max();
1637 double y_min = std::numeric_limits<double>::max();
1638 double z_min = std::numeric_limits<double>::max();
1639 double x_max = std::numeric_limits<double>::lowest();
1640 double y_max = std::numeric_limits<double>::lowest();
1641 double z_max = std::numeric_limits<double>::lowest();
1642 std::queue<idNode> queue;
1643 queue.emplace(tree->getRoot());
1644 while(!queue.empty()) {
1645 idNode const node = queue.front();
1646 queue.pop();
1647 double *point = treeNodes->GetPoints()->GetPoint(nodeCorrT[node]);
1648 x_min = std::min(x_min, point[0]);
1649 x_max = std::max(x_max, point[0]);
1650 y_min = std::min(y_min, point[1]);
1651 y_max = std::max(y_max, point[1]);
1652 z_min = std::min(z_min, point[2]);
1653 z_max = std::max(z_max, point[2]);
1654 std::vector<idNode> children;
1655 tree->getChildren(node, children);
1656 for(auto child : children)
1657 queue.emplace(child);
1658 }
1659 return std::make_tuple(x_min, x_max, y_min, y_max, z_min, z_max);
1660 }
1661
1662 std::tuple<double, double, double, double, double, double>
1663 getRealBounds(vtkUnstructuredGrid *treeNodes, FTMTree_MT *tree) {
1664 std::vector<int> nodeCorrT(tree->getNumberOfNodes());
1665 for(size_t i = 0; i < nodeCorrT.size(); ++i)
1666 nodeCorrT[i] = i;
1667
1668 return getRealBounds(treeNodes, tree, nodeCorrT);
1669 }
1670};
#define M_PI
Definition Os.h:50
void makeTreesOutput(FTMTree_MT *tree1, FTMTree_MT *tree2)
void setTreesSegmentation(std::vector< vtkDataSet * > &segmentation)
void setInterpolatedTrees(std::vector< bool > &isInterpolatedTrees)
std::vector< double > getClusterShift()
void copyPointData(vtkUnstructuredGrid *treeNodes, std::vector< int > &nodeCorrT)
void setVtkOutputNode1(vtkUnstructuredGrid *vtkNode1)
void setNodeCorr1(std::vector< std::vector< SimplexId > > &nodeCorrT)
void setOutputMatchingBarycenter(std::vector< std::vector< std::vector< std::tuple< idNode, idNode, double > > > > &matching)
void setVtkOutputSegmentation(vtkDataSet *vtkSegmentation)
void setTreesNodeCorrMesh(std::vector< std::vector< int > > &nodeCorrMesh)
ttkMergeTreeVisualization()=default
void addCustomIntArray(std::string &name, std::vector< int > &vec)
void setVtkOutputMatching(vtkUnstructuredGrid *vtkMatching)
void setClusteringAssignment(std::vector< int > &asgn)
void addCustomStringArray(std::string &name, std::vector< std::string > &vec)
void makeTreesOutput(std::vector< FTMTree_MT * > &trees, std::vector< FTMTree_MT * > &barycenters)
std::vector< std::vector< float > > getAllBaryPercentMatch()
std::vector< std::vector< SimplexId > > getNodeCorr()
std::tuple< double, double, double, double, double, double > getRealBounds(vtkUnstructuredGrid *treeNodes, FTMTree_MT *tree)
std::tuple< double, double, double, double, double, double > getRealBounds(vtkUnstructuredGrid *treeNodes, FTMTree_MT *tree, std::vector< int > &nodeCorrT)
void setTreesNodes(vtkUnstructuredGrid *nodes)
void copyPointData(vtkUnstructuredGrid *treeNodes)
void setTreesSegmentation(vtkDataSet *segmentation)
void makeTreesOutput(FTMTree_MT *tree1)
~ttkMergeTreeVisualization() override=default
void setConvertedToDiagram(bool converted)
void makeMatchingOutput(std::vector< FTMTree_MT * > &trees, std::vector< FTMTree_MT * > &barycenters)
void addVtkCustomArrays(std::vector< std::tuple< std::string, std::vector< dataType > > > &cArrays, std::vector< std::vector< dataType > > &cArraysValues, vtkUnstructuredGrid *vtkOutput, int type, int output)
void getTreeNodeIdRev(vtkDataArray *treeNodeIdArray, std::vector< int > &treeNodeIdRev)
void setVtkOutputNode2(vtkUnstructuredGrid *vtkNode2)
void setTreesNodeCorrMesh(std::vector< int > &nodeCorrMesh)
void makeTreesOutput(std::vector< FTMTree_MT * > &trees)
void makeMatchingOutput(FTMTree_MT *tree1, FTMTree_MT *tree2)
void setVtkOutputArc(vtkUnstructuredGrid *vtkArc)
void setNodeCorr2(std::vector< std::vector< SimplexId > > &nodeCorrT)
void setVtkOutputNode(vtkUnstructuredGrid *vtkNode)
void setOutputMatching(std::vector< std::tuple< idNode, idNode, double > > &matching)
void addCustomArray(std::string &name, std::vector< double > &vec)
void setTreesNodes(std::vector< vtkUnstructuredGrid * > &nodes)
void setAllBaryPercentMatch(std::vector< std::vector< float > > &baryPercentMatch)
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
std::vector< double > excludeImportantPairsLowerValues_
std::vector< double > excludeImportantPairsHigherValues_
std::tuple< double, double, double, double, double, double > getMaximalBounds(std::vector< std::tuple< double, double, double, double, double, double > > &allBounds, std::vector< int > &clusteringAssignmentT, int clusterID)
dataType getNodePersistence(idNode nodeId)
idNode getNumberOfNodes() const
Definition FTMTree_MT.h:389
void getChildren(idNode nodeId, std::vector< idNode > &res)
Node * getNode(idNode nodeId)
Definition FTMTree_MT.h:393
idSuperArc getNumberOfDownSuperArcs() const
Definition FTMNode.h:82
idSuperArc getNumberOfUpSuperArcs() const
Definition FTMNode.h:86
unsigned int idNode
Node index in vect_nodes_.
const char PersistenceName[]
Definition DataTypes.h:72
const char PersistenceCoordinatesName[]
Definition DataTypes.h:70
const char PersistencePairTypeName[]
Definition DataTypes.h:73
const char PersistenceCriticalTypeName[]
Definition DataTypes.h:67
const char PersistenceIsFinite[]
Definition DataTypes.h:74
const char VertexScalarFieldName[]
default name for vertex scalar field
Definition DataTypes.h:35
const char PersistencePairIdentifierName[]
Definition DataTypes.h:71
const char PersistenceBirthName[]
Definition DataTypes.h:68
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)