TTK
Loading...
Searching...
No Matches
ttkMergeTreeClustering.cpp
Go to the documentation of this file.
1#include <FTMStructures.h>
2#include <FTMTree.h>
3#include <FTMTreeUtils.h>
4#include <MergeTreeUtils.h>
7#include <ttkMergeTreeUtils.h>
9
10#include <vtkDataObject.h> // For port information
11#include <vtkObjectFactory.h> // for new macro
12
13#include <vtkCellData.h>
14#include <vtkDataArray.h>
15#include <vtkDataSet.h>
16#include <vtkFloatArray.h>
17#include <vtkInformation.h>
18#include <vtkInformationVector.h>
19#include <vtkMultiBlockDataSet.h>
20#include <vtkPointData.h>
21#include <vtkTable.h>
22
23#include <string>
24
25using namespace ttk;
26using namespace ftm;
27
28// A VTK macro that enables the instantiation of this class via ::New()
29// You do not have to modify this
31
45 this->setDebugMsgPrefix("MergeTreeClustering");
46 this->SetNumberOfInputPorts(2);
47 this->SetNumberOfOutputPorts(3);
48}
49
51
60 vtkInformation *info) {
61 if(port == 0)
62 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
63 else if(port == 1) {
64 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
65 info->Set(vtkAlgorithm::INPUT_IS_OPTIONAL(), 1);
66 } else
67 return 0;
68
69 return 1;
70}
71
88 vtkInformation *info) {
89 if(port == 0 || port == 1 || port == 2) {
90 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
91 } else {
92 return 0;
93 }
94 return 1;
95}
96
111 vtkInformationVector **inputVector,
112 vtkInformationVector *outputVector) {
113 // ------------------------------------------------------------------------------------
114 // --- Get input object from input vector
115 // ------------------------------------------------------------------------------------
116 auto blocks = vtkMultiBlockDataSet::GetData(inputVector[0], 0);
117 auto blocks2 = vtkMultiBlockDataSet::GetData(inputVector[1], 0);
118
119 if(Backend == 0) {
120 baseModule = 0;
121 } else if(Backend == 1) {
122 baseModule = 0;
123 } else if(Backend == 3) {
124 baseModule = 1;
125 } else if(Backend == 4) {
126 baseModule = 2;
127 } else {
128 baseModule = 0;
129 }
130
131 // filter out new backends (not yet supported)
132 if(baseModule != 0) {
133 printErr("Invalid Backend chosen. Path Mapping Distance and Branch Mapping "
134 "Distance not yet supported. Canceling computation.");
135 return 1;
136 }
137
138 // ------------------------------------------------------------------------------------
139 // --- Load blocks
140 // ------------------------------------------------------------------------------------
141 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> inputTrees, inputTrees2;
142 loadBlocks(inputTrees, blocks);
143 loadBlocks(inputTrees2, blocks2);
144
145 // If we have already computed once but the input has changed
146 if(treesNodes.size() != 0 and inputTrees[0]->GetBlock(0) != treesNodes[0])
147 resetDataVisualization();
148
149 return run<float>(outputVector, inputTrees, inputTrees2);
150}
151
152template <class dataType>
154 vtkInformationVector *outputVector,
155 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
156 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2) {
157 if(not isDataVisualizationFilled())
158 runCompute<dataType>(outputVector, inputTrees, inputTrees2);
159 runOutput<dataType>(outputVector, inputTrees, inputTrees2);
160 return 1;
161}
162
163template <class dataType>
165 vtkInformationVector *ttkNotUsed(outputVector),
166 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
167 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2) {
168 // ------------------------------------------------------------------------------------
169 // --- Construct trees
170 // ------------------------------------------------------------------------------------
171 const int numInputs = inputTrees.size();
172 const int numInputs2
173 = (not inputTrees[0]->GetBlock(0)->IsA("vtkUnstructuredGrid")
174 ? inputTrees2.size()
175 : (JoinSplitMixtureCoefficient != 0
176 and JoinSplitMixtureCoefficient != 1
177 ? numInputs
178 : 0));
179
180 setDataVisualization(numInputs, numInputs2);
181
182 std::vector<MergeTree<dataType>> intermediateMTrees(numInputs),
183 intermediateMTrees2(numInputs2), barycenters(NumberOfBarycenters),
184 barycenters2((numInputs2 != 0) * NumberOfBarycenters);
185 std::vector<FTMTree_MT *> intermediateTrees(numInputs),
186 intermediateTrees2(numInputs2);
187
188 bool const useSecondPairsType = (JoinSplitMixtureCoefficient == 0);
189 IsPersistenceDiagram = constructTrees<dataType>(
190 inputTrees, intermediateMTrees, treesNodes, treesArcs, treesSegmentation,
191 useSecondPairsType, DiagramPairTypes);
192 if(not IsPersistenceDiagram
193 or (JoinSplitMixtureCoefficient != 0
194 and JoinSplitMixtureCoefficient != 1)) {
195 auto &inputTrees2ToUse
196 = (not IsPersistenceDiagram ? inputTrees2 : inputTrees);
197 constructTrees<dataType>(inputTrees2ToUse, intermediateMTrees2, treesNodes2,
198 treesArcs2, treesSegmentation2,
199 !useSecondPairsType, DiagramPairTypes);
200 }
201
202 mergeTreeToFTMTree<dataType>(intermediateMTrees, intermediateTrees);
203 mergeTreeToFTMTree<dataType>(intermediateMTrees2, intermediateTrees2);
204
205 // ------------------------------------------------------------------------------------
206 // --- Call base
207 // ------------------------------------------------------------------------------------
208 bool const AddNodes = true;
209 // Classical distance
210 double distance = 0;
211
212 // Verify parameters
213 if(Backend == 0) {
214 BranchDecomposition = true;
215 NormalizedWasserstein = true;
216 KeepSubtree = false;
217 } else if(Backend == 1) {
218 BranchDecomposition = false;
219 NormalizedWasserstein = false;
220 KeepSubtree = true;
221 ComputeBarycenter = false;
222 }
223 if(IsPersistenceDiagram) {
224 BranchDecomposition = true;
225 }
226 if(ComputeBarycenter) {
227 if(not BranchDecomposition)
228 printMsg("BranchDecomposition is set to true since the barycenter "
229 "computation is asked.");
230 BranchDecomposition = true;
231 if(KeepSubtree)
232 printMsg("KeepSubtree is set to false since the barycenter computation "
233 "is asked.");
234 KeepSubtree = false;
235 }
236 if(not BranchDecomposition) {
237 if(NormalizedWasserstein)
238 printMsg("NormalizedWasserstein is set to false since branch "
239 "decomposition is not asked.");
240 NormalizedWasserstein = false;
241 }
242 EpsilonTree2 = EpsilonTree1;
243 Epsilon2Tree2 = Epsilon2Tree1;
244 Epsilon3Tree2 = Epsilon3Tree1;
245 printMsg("BranchDecomposition: " + std::to_string(BranchDecomposition));
246 printMsg("NormalizedWasserstein: " + std::to_string(NormalizedWasserstein));
247 printMsg("KeepSubtree: " + std::to_string(KeepSubtree));
248
249 // Call base
250 if(not ComputeBarycenter) {
251 MergeTreeDistance mergeTreeDistance;
252 mergeTreeDistance.setAssignmentSolver(AssignmentSolver);
253 mergeTreeDistance.setEpsilonTree1(EpsilonTree1);
254 mergeTreeDistance.setEpsilonTree2(EpsilonTree2);
255 mergeTreeDistance.setEpsilon2Tree1(Epsilon2Tree1);
256 mergeTreeDistance.setEpsilon2Tree2(Epsilon2Tree2);
257 mergeTreeDistance.setEpsilon3Tree1(Epsilon3Tree1);
258 mergeTreeDistance.setEpsilon3Tree2(Epsilon3Tree2);
259 mergeTreeDistance.setBranchDecomposition(BranchDecomposition);
260 mergeTreeDistance.setPersistenceThreshold(PersistenceThreshold);
261 mergeTreeDistance.setNormalizedWasserstein(NormalizedWasserstein);
262 mergeTreeDistance.setKeepSubtree(KeepSubtree);
263 mergeTreeDistance.setUseMinMaxPair(UseMinMaxPair);
264 mergeTreeDistance.setCleanTree(true);
265 mergeTreeDistance.setPostprocess(OutputTrees);
266 mergeTreeDistance.setDeleteMultiPersPairs(DeleteMultiPersPairs);
267 mergeTreeDistance.setEpsilon1UseFarthestSaddle(Epsilon1UseFarthestSaddle);
268 mergeTreeDistance.setIsPersistenceDiagram(IsPersistenceDiagram);
269 mergeTreeDistance.setNonMatchingWeight(NonMatchingWeight);
270 mergeTreeDistance.setThreadNumber(this->threadNumber_);
271 mergeTreeDistance.setDebugLevel(this->debugLevel_);
272
273 distance = mergeTreeDistance.execute<dataType>(
274 intermediateMTrees[0], intermediateMTrees[1], outputMatching);
275 trees1NodeCorrMesh = mergeTreeDistance.getTreesNodeCorr();
276 finalDistances = std::vector<double>{distance};
277 } else {
278 if(NumberOfBarycenters == 1) { // and numInputs2==0){
279 MergeTreeBarycenter mergeTreeBarycenter;
280 mergeTreeBarycenter.setAssignmentSolver(AssignmentSolver);
281 mergeTreeBarycenter.setEpsilonTree1(EpsilonTree1);
282 mergeTreeBarycenter.setEpsilonTree2(EpsilonTree2);
283 mergeTreeBarycenter.setEpsilon2Tree1(Epsilon2Tree1);
284 mergeTreeBarycenter.setEpsilon2Tree2(Epsilon2Tree2);
285 mergeTreeBarycenter.setEpsilon3Tree1(Epsilon3Tree1);
286 mergeTreeBarycenter.setEpsilon3Tree2(Epsilon3Tree2);
287 mergeTreeBarycenter.setBranchDecomposition(BranchDecomposition);
288 mergeTreeBarycenter.setPersistenceThreshold(PersistenceThreshold);
289 mergeTreeBarycenter.setNormalizedWasserstein(NormalizedWasserstein);
290 mergeTreeBarycenter.setKeepSubtree(KeepSubtree);
291 mergeTreeBarycenter.setUseMinMaxPair(UseMinMaxPair);
292 mergeTreeBarycenter.setAddNodes(AddNodes);
293 mergeTreeBarycenter.setDeterministic(Deterministic);
294 mergeTreeBarycenter.setBarycenterSizeLimitPercent(
295 BarycenterSizeLimitPercent);
296 mergeTreeBarycenter.setAlpha(Alpha);
297 mergeTreeBarycenter.setPostprocess(OutputTrees);
298 mergeTreeBarycenter.setDeleteMultiPersPairs(DeleteMultiPersPairs);
299 mergeTreeBarycenter.setEpsilon1UseFarthestSaddle(
300 Epsilon1UseFarthestSaddle);
301 mergeTreeBarycenter.setIsPersistenceDiagram(IsPersistenceDiagram);
302 mergeTreeBarycenter.setThreadNumber(this->threadNumber_);
303 mergeTreeBarycenter.setDebugLevel(this->debugLevel_);
304
305 mergeTreeBarycenter.execute<dataType>(
306 intermediateMTrees, outputMatchingBarycenter[0], barycenters[0]);
307 trees1NodeCorrMesh = mergeTreeBarycenter.getTreesNodeCorr();
308 finalDistances = mergeTreeBarycenter.getFinalDistances();
309 } else {
310 MergeTreeClustering<dataType> mergeTreeClustering;
311 mergeTreeClustering.setAssignmentSolver(AssignmentSolver);
312 mergeTreeClustering.setEpsilonTree1(EpsilonTree1);
313 mergeTreeClustering.setEpsilonTree2(EpsilonTree2);
314 mergeTreeClustering.setEpsilon2Tree1(Epsilon2Tree1);
315 mergeTreeClustering.setEpsilon2Tree2(Epsilon2Tree2);
316 mergeTreeClustering.setEpsilon3Tree1(Epsilon3Tree1);
317 mergeTreeClustering.setEpsilon3Tree2(Epsilon3Tree2);
318 mergeTreeClustering.setBranchDecomposition(BranchDecomposition);
319 mergeTreeClustering.setPersistenceThreshold(PersistenceThreshold);
320 mergeTreeClustering.setNormalizedWasserstein(NormalizedWasserstein);
321 mergeTreeClustering.setKeepSubtree(KeepSubtree);
322 mergeTreeClustering.setUseMinMaxPair(UseMinMaxPair);
323 mergeTreeClustering.setAddNodes(AddNodes);
324 mergeTreeClustering.setDeterministic(Deterministic);
325 mergeTreeClustering.setNoCentroids(NumberOfBarycenters);
326 mergeTreeClustering.setBarycenterSizeLimitPercent(
327 BarycenterSizeLimitPercent);
328 mergeTreeClustering.setPostprocess(OutputTrees);
329 mergeTreeClustering.setDeleteMultiPersPairs(DeleteMultiPersPairs);
330 mergeTreeClustering.setEpsilon1UseFarthestSaddle(
331 Epsilon1UseFarthestSaddle);
332 mergeTreeClustering.setIsPersistenceDiagram(IsPersistenceDiagram);
333 mergeTreeClustering.setMixtureCoefficient(JoinSplitMixtureCoefficient);
334 mergeTreeClustering.setThreadNumber(this->threadNumber_);
335 mergeTreeClustering.setDebugLevel(this->debugLevel_);
336
337 mergeTreeClustering.template execute<dataType>(
338 intermediateMTrees, outputMatchingBarycenter, clusteringAssignment,
339 intermediateMTrees2, outputMatchingBarycenter2, barycenters,
340 barycenters2);
341 trees1NodeCorrMesh = mergeTreeClustering.getTreesNodeCorr();
342 trees2NodeCorrMesh = mergeTreeClustering.getTrees2NodeCorr();
343 finalDistances = mergeTreeClustering.getFinalDistances();
344 }
345 }
346
347 mergeTreesTemplateToDouble<dataType>(intermediateMTrees, intermediateSTrees);
348 if(numInputs2 != 0)
350 intermediateMTrees2, intermediateSTrees2);
351 if(ComputeBarycenter) {
352 mergeTreesTemplateToDouble<dataType>(barycenters, barycentersS);
353 if(numInputs2 != 0)
354 mergeTreesTemplateToDouble<dataType>(barycenters2, barycentersS2);
355 }
356
357 return 1;
358}
359
360static void addFieldData(vtkDataSet *in, vtkDataSet *out) {
361 auto inFieldData = in->GetFieldData();
362 auto outFieldData = out->GetFieldData();
363 for(int i = 0; i < inFieldData->GetNumberOfArrays(); ++i) {
364 outFieldData->AddArray(inFieldData->GetAbstractArray(i));
365 }
366}
367
368template <class dataType>
370 ttkMergeTreeVisualization &visuMaker,
372 std::vector<FTMTree_MT *> &intermediateTrees2,
373 std::vector<FTMTree_MT *> &barycentersTree2,
374 std::vector<std::vector<
375 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, double>>>>
376 &outputMatchingBarycenter2,
377 std::vector<std::vector<int>> &trees2NodeCorrMesh,
378 vtkUnstructuredGrid *treeNodes,
379 int i,
380 std::vector<SimplexId> &nodeCorr2) {
381 int const nodeCorrShift = vtkOutputNode->GetNumberOfPoints();
382 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNode2
384 visuMaker.setVtkOutputNode(vtkOutputNode2);
385 visuMaker.setVtkOutputArc(vtkOutputNode2);
386 visuMaker.setOutputMatchingBarycenter(outputMatchingBarycenter2);
387 visuMaker.clearAllCustomArrays();
388 if(treeNodes)
389 visuMaker.copyPointData(treeNodes, trees2NodeCorrMesh[i]);
390 visuMaker.setTreesNodeCorrMesh(trees2NodeCorrMesh);
391 visuMaker.setIsPDSadMax(true);
392 visuMaker.makeTreesOutput<dataType>(intermediateTrees2, barycentersTree2);
393 auto nodeCorrT = visuMaker.getNodeCorr();
394 nodeCorr2 = nodeCorrT[i];
395 for(unsigned int j = 0; j < nodeCorr2.size(); ++j)
396 nodeCorr2[j] += nodeCorrShift;
397
398 vtkNew<vtkAppendFilter> appendFilter{};
399 appendFilter->AddInputData(vtkOutputNode);
400 appendFilter->AddInputData(vtkOutputNode2);
401 appendFilter->SetMergePoints(false);
402 appendFilter->Update();
403 vtkOutputNode->ShallowCopy(appendFilter->GetOutput());
404}
405
406template <class dataType>
408 ttkMergeTreeVisualization &visuMakerMatching,
409 vtkSmartPointer<vtkUnstructuredGrid> &vtkOutputMatching,
410 std::vector<FTMTree_MT *> &intermediateTrees2,
411 std::vector<FTMTree_MT *> &barycentersTree2,
412 std::vector<std::vector<
413 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, double>>>>
414 &outputMatchingBarycenter2,
415 std::vector<std::vector<SimplexId>> &nodeCorr2,
416 std::vector<std::vector<SimplexId>> &nodeCorrBary2,
417 std::vector<std::vector<float>> &allBaryPercentMatch2) {
418 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputMatching2
420 visuMakerMatching.setVtkOutputMatching(vtkOutputMatching2);
421 visuMakerMatching.setOutputMatchingBarycenter(outputMatchingBarycenter2);
422 visuMakerMatching.setAllBaryPercentMatch(allBaryPercentMatch2);
423 visuMakerMatching.setNodeCorr1(nodeCorr2);
424 visuMakerMatching.setNodeCorr2(nodeCorrBary2);
425 visuMakerMatching.makeMatchingOutput<dataType>(
426 intermediateTrees2, barycentersTree2);
427
428 vtkNew<vtkAppendFilter> appendFilter{};
429 appendFilter->AddInputData(vtkOutputMatching);
430 appendFilter->AddInputData(vtkOutputMatching2);
431 appendFilter->SetMergePoints(false);
432 appendFilter->Update();
433 vtkOutputMatching->ShallowCopy(appendFilter->GetOutput());
434}
435
436template <class dataType>
438 vtkInformationVector *outputVector,
439 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
440 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &ttkNotUsed(inputTrees2)) {
441 std::vector<MergeTree<dataType>> intermediateMTrees, intermediateMTrees2;
442 mergeTreesDoubleToTemplate<dataType>(intermediateSTrees, intermediateMTrees);
444 intermediateSTrees2, intermediateMTrees2);
445 std::vector<FTMTree_MT *> intermediateTrees, intermediateTrees2;
446 mergeTreeToFTMTree<dataType>(intermediateMTrees, intermediateTrees);
447 mergeTreeToFTMTree<dataType>(intermediateMTrees2, intermediateTrees2);
448
449 std::vector<MergeTree<dataType>> barycenters, barycenters2;
450 if(ComputeBarycenter) {
451 mergeTreesDoubleToTemplate<dataType>(barycentersS, barycenters);
452 mergeTreesDoubleToTemplate<dataType>(barycentersS2, barycenters2);
453 }
454
455 const int numInputs = inputTrees.size();
456 const int numInputs2 = intermediateMTrees2.size();
457 // ------------------------------------------------------------------------------------
458 // --- Create output
459 // ------------------------------------------------------------------------------------
460 Timer t_makeTreesOutput;
461
462 auto output_clusters = vtkMultiBlockDataSet::GetData(outputVector, 0);
463 auto output_centroids = vtkMultiBlockDataSet::GetData(outputVector, 1);
464 auto output_matchings = vtkMultiBlockDataSet::GetData(outputVector, 2);
465
466 // Declare internal arrays
467 std::vector<std::vector<SimplexId>> nodeCorr(numInputs),
468 nodeCorr2(numInputs2);
469
470 if(not ComputeBarycenter) {
471 // ---------------------------------------------------------------
472 // - Classical output
473 // ---------------------------------------------------------------
474 if(OutputTrees) {
475 FTMTree_MT *tree1 = intermediateTrees[0], *tree2 = intermediateTrees[1];
476 // ------------------------------------------
477 // --- Input trees
478 // ------------------------------------------
479 // Declare vtk objects
480 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNode1
482 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputArc1
484 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputSegmentation1
486
487 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNode2
489 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputArc2
491 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputSegmentation2
493
494 vtkSmartPointer<vtkMultiBlockDataSet> const vtkBlockNodes
500
501 // Fill vtk objects
503 visuMaker.setPlanarLayout(PlanarLayout);
505 BranchDecompositionPlanarLayout);
506 visuMaker.setBranchSpacing(BranchSpacing);
507 visuMaker.setRescaleTreesIndividually(RescaleTreesIndividually);
508 visuMaker.setOutputSegmentation(OutputSegmentation);
509 visuMaker.setDimensionSpacing(DimensionSpacing);
510 visuMaker.setDimensionToShift(DimensionToShift);
511 visuMaker.setImportantPairs(ImportantPairs);
512 visuMaker.setMaximumImportantPairs(MaximumImportantPairs);
513 visuMaker.setMinimumImportantPairs(MinimumImportantPairs);
514 visuMaker.setImportantPairsSpacing(ImportantPairsSpacing);
515 visuMaker.setNonImportantPairsSpacing(NonImportantPairsSpacing);
516 visuMaker.setNonImportantPairsProximity(NonImportantPairsProximity);
517 visuMaker.setExcludeImportantPairsHigher(ExcludeImportantPairsHigher);
518 visuMaker.setExcludeImportantPairsLower(ExcludeImportantPairsLower);
519 visuMaker.setIsPersistenceDiagram(IsPersistenceDiagram);
520
521 nodeCorr.clear();
522 // First tree
523 visuMaker.setNoSampleOffset(1);
524 visuMaker.setTreesNodes(treesNodes[0]);
525 visuMaker.copyPointData(treesNodes[0], trees1NodeCorrMesh[0]);
526 visuMaker.setTreesNodeCorrMesh(trees1NodeCorrMesh[0]);
527 visuMaker.setTreesSegmentation(treesSegmentation[0]);
528 visuMaker.setVtkOutputNode(vtkOutputNode1);
529 if(IsPersistenceDiagram)
530 visuMaker.setVtkOutputArc(vtkOutputNode1);
531 else
532 visuMaker.setVtkOutputArc(vtkOutputArc1);
533 visuMaker.setVtkOutputSegmentation(vtkOutputSegmentation1);
534 visuMaker.setDebugLevel(this->debugLevel_);
535 visuMaker.setIsPersistenceDiagram(IsPersistenceDiagram);
536
537 visuMaker.makeTreesOutput<dataType>(tree1);
538 nodeCorr.push_back(visuMaker.getNodeCorr()[0]);
539
540 // Second tree
541 visuMaker.setISampleOffset(1);
542 visuMaker.setTreesNodes(treesNodes[1]);
543 visuMaker.clearAllCustomArrays();
544 visuMaker.copyPointData(treesNodes[1], trees1NodeCorrMesh[1]);
545 visuMaker.setTreesNodeCorrMesh(trees1NodeCorrMesh[1]);
546 visuMaker.setTreesSegmentation(treesSegmentation[1]);
547 visuMaker.setVtkOutputNode(vtkOutputNode2);
548 if(IsPersistenceDiagram)
549 visuMaker.setVtkOutputArc(vtkOutputNode2);
550 else
551 visuMaker.setVtkOutputArc(vtkOutputArc2);
552 visuMaker.setVtkOutputSegmentation(vtkOutputSegmentation2);
553
554 visuMaker.makeTreesOutput<dataType>(tree2);
555 nodeCorr.push_back(visuMaker.getNodeCorr()[0]);
556
557 // Field data
558 vtkOutputNode1->GetFieldData()->ShallowCopy(
559 treesNodes[0]->GetFieldData());
560 vtkOutputNode2->GetFieldData()->ShallowCopy(
561 treesNodes[1]->GetFieldData());
562 if(not IsPersistenceDiagram) {
563 vtkOutputArc1->GetFieldData()->ShallowCopy(
564 treesArcs[0]->GetFieldData());
565 vtkOutputArc2->GetFieldData()->ShallowCopy(
566 treesArcs[1]->GetFieldData());
567 }
568 if(treesSegmentation[0])
569 addFieldData(treesSegmentation[0], vtkOutputNode1);
570 if(treesSegmentation[1])
571 addFieldData(treesSegmentation[1], vtkOutputNode2);
572 if(OutputSegmentation) {
573 vtkOutputSegmentation1->GetFieldData()->ShallowCopy(
574 treesSegmentation[0]->GetFieldData());
575 vtkOutputSegmentation2->GetFieldData()->ShallowCopy(
576 treesSegmentation[1]->GetFieldData());
577 }
578
579 // Construct multiblock
580 if(IsPersistenceDiagram and not OutputSegmentation) {
581 output_clusters->SetNumberOfBlocks(2);
582 output_clusters->SetBlock(0, vtkOutputNode1);
583 output_clusters->SetBlock(1, vtkOutputNode2);
584 } else {
585 vtkBlockNodes->SetNumberOfBlocks(2);
586 vtkBlockNodes->SetBlock(0, vtkOutputNode1);
587 vtkBlockNodes->SetBlock(1, vtkOutputNode2);
588
589 if(not IsPersistenceDiagram) {
590 vtkBlockArcs->SetNumberOfBlocks(2);
591 vtkBlockArcs->SetBlock(0, vtkOutputArc1);
592 vtkBlockArcs->SetBlock(1, vtkOutputArc2);
593 }
594
595 output_clusters->SetNumberOfBlocks(1 + !IsPersistenceDiagram
596 + OutputSegmentation);
597 output_clusters->SetBlock(0, vtkBlockNodes);
598 if(not IsPersistenceDiagram)
599 output_clusters->SetBlock(1, vtkBlockArcs);
600 if(OutputSegmentation) {
601 vtkBlockSegs->SetNumberOfBlocks(2);
602 vtkBlockSegs->SetBlock(0, vtkOutputSegmentation1);
603 vtkBlockSegs->SetBlock(1, vtkOutputSegmentation2);
604 int const segBlockID = 1 + !IsPersistenceDiagram;
605 output_clusters->SetBlock(segBlockID, vtkBlockSegs);
606 }
607 }
608
609 // ------------------------------------------
610 // --- Matching
611 // ------------------------------------------
612 // Declare vtk objects
613 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputMatching
615
616 // Fill vtk objects
617 ttkMergeTreeVisualization visuMakerMatching;
618 visuMakerMatching.setVtkOutputMatching(vtkOutputMatching);
619 visuMakerMatching.setOutputMatching(outputMatching);
620 visuMakerMatching.setVtkOutputNode2(vtkOutputNode1);
621 visuMakerMatching.setVtkOutputNode1(vtkOutputNode2);
622 visuMakerMatching.setNodeCorr1(nodeCorr);
623 visuMakerMatching.setDebugLevel(this->debugLevel_);
624
625 visuMakerMatching.makeMatchingOutput<dataType>(tree1, tree2);
626
627 // Field data
628 vtkNew<vtkDoubleArray> vtkDistance{};
629 vtkDistance->SetName("Distance");
630 vtkDistance->SetNumberOfTuples(1);
631 vtkDistance->SetTuple1(0, finalDistances[0]);
632 vtkOutputMatching->GetFieldData()->AddArray(vtkDistance);
633
634 // Construct multiblock
635 output_matchings->SetNumberOfBlocks(1);
636 output_matchings->SetBlock(0, vtkOutputMatching);
637 }
638 } else {
639 // ---------------------------------------------------------------
640 // - Barycenter/Clustering output
641 // ---------------------------------------------------------------
642 if(OutputTrees) {
643 // --- Declare internal arrays
644 std::vector<std::vector<SimplexId>> nodeCorrBary(NumberOfBarycenters),
645 nodeCorrBary2((numInputs2 != 0) * NumberOfBarycenters);
646 std::vector<std::vector<float>> allBaryPercentMatch(NumberOfBarycenters),
647 allBaryPercentMatch2((numInputs2 != 0) * NumberOfBarycenters);
648 std::vector<FTMTree_MT *> barycentersTree, barycentersTree2;
649 mergeTreeToFTMTree<dataType>(barycenters, barycentersTree);
650 mergeTreeToFTMTree<dataType>(barycenters2, barycentersTree2);
651
652 // ------------------------------------------
653 // --- Input trees
654 // ------------------------------------------
655 if(IsPersistenceDiagram and not OutputSegmentation) {
656 output_clusters->SetNumberOfBlocks(numInputs);
657 } else {
658 output_clusters->SetNumberOfBlocks((OutputSegmentation ? 3 : 2));
659 vtkSmartPointer<vtkMultiBlockDataSet> const vtkBlockNodes
661 vtkBlockNodes->SetNumberOfBlocks(numInputs);
662 output_clusters->SetBlock(0, vtkBlockNodes);
665 vtkBlockArcs->SetNumberOfBlocks(numInputs);
666 output_clusters->SetBlock(1, vtkBlockArcs);
667 if(OutputSegmentation) {
670 vtkBlockSegs->SetNumberOfBlocks(numInputs);
671 output_clusters->SetBlock(2, vtkBlockSegs);
672 }
673 }
674 for(unsigned int c = 0; c < NumberOfBarycenters; ++c) {
675#ifdef TTK_ENABLE_OPENMP
676#pragma omp parallel for schedule(dynamic) num_threads(this->threadNumber_)
677#endif
678 for(int i = 0; i < numInputs; ++i) {
679 if(clusteringAssignment[i] != (int)c)
680 continue;
681
682 // Declare vtk objects
685 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputArc
687 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputSegmentation
689
690 // Fill vtk objects
692 visuMaker.setPlanarLayout(PlanarLayout);
694 BranchDecompositionPlanarLayout);
695 visuMaker.setBranchSpacing(BranchSpacing);
696 visuMaker.setRescaleTreesIndividually(RescaleTreesIndividually);
697 visuMaker.setOutputSegmentation(OutputSegmentation);
698 visuMaker.setDimensionSpacing(DimensionSpacing);
699 visuMaker.setDimensionToShift(DimensionToShift);
700 visuMaker.setImportantPairs(ImportantPairs);
701 visuMaker.setMaximumImportantPairs(MaximumImportantPairs);
702 visuMaker.setMinimumImportantPairs(MinimumImportantPairs);
703 visuMaker.setImportantPairsSpacing(ImportantPairsSpacing);
704 visuMaker.setNonImportantPairsSpacing(NonImportantPairsSpacing);
705 visuMaker.setNonImportantPairsProximity(NonImportantPairsProximity);
706 visuMaker.setExcludeImportantPairsHigher(ExcludeImportantPairsHigher);
707 visuMaker.setExcludeImportantPairsLower(ExcludeImportantPairsLower);
708 visuMaker.setIsPersistenceDiagram(IsPersistenceDiagram);
709 visuMaker.setTreesNodes(treesNodes);
710 visuMaker.copyPointData(treesNodes[i], trees1NodeCorrMesh[i]);
711 visuMaker.setTreesNodeCorrMesh(trees1NodeCorrMesh);
712 visuMaker.setTreesSegmentation(treesSegmentation);
713 visuMaker.setVtkOutputNode(vtkOutputNode);
714 if(IsPersistenceDiagram)
715 visuMaker.setVtkOutputArc(vtkOutputNode);
716 else
717 visuMaker.setVtkOutputArc(vtkOutputArc);
718 visuMaker.setVtkOutputSegmentation(vtkOutputSegmentation);
719 visuMaker.setClusteringAssignment(clusteringAssignment);
720 visuMaker.setOutputMatchingBarycenter(outputMatchingBarycenter);
721 visuMaker.setPrintTreeId(i);
722 visuMaker.setPrintClusterId(c);
723 visuMaker.setDebugLevel(this->debugLevel_);
724 visuMaker.setIsPersistenceDiagram(IsPersistenceDiagram);
725 visuMaker.setIsPDSadMax(JoinSplitMixtureCoefficient == 0);
726
727 visuMaker.makeTreesOutput<dataType>(
728 intermediateTrees, barycentersTree);
729 auto nodeCorrT = visuMaker.getNodeCorr();
730 nodeCorr[i] = nodeCorrT[i];
731 if(IsPersistenceDiagram and JoinSplitMixtureCoefficient != 0
732 and JoinSplitMixtureCoefficient != 1)
734 visuMaker, vtkOutputNode, intermediateTrees2, barycentersTree2,
735 outputMatchingBarycenter2, trees2NodeCorrMesh, treesNodes[i], i,
736 nodeCorr2[i]);
737
738 // Field data
739 vtkOutputNode->GetFieldData()->ShallowCopy(
740 treesNodes[i]->GetFieldData());
741 if(not IsPersistenceDiagram)
742 vtkOutputArc->GetFieldData()->ShallowCopy(
743 treesArcs[i]->GetFieldData());
744 if(treesSegmentation[i])
745 addFieldData(treesSegmentation[i], vtkOutputNode);
746 if(OutputSegmentation)
747 vtkOutputSegmentation->GetFieldData()->ShallowCopy(
748 treesSegmentation[i]->GetFieldData());
749
750 vtkNew<vtkDoubleArray> vtkClusterAssignment{};
751 vtkClusterAssignment->SetName("ClusterAssignment");
752 vtkClusterAssignment->SetNumberOfTuples(1);
753 vtkClusterAssignment->SetTuple1(0, clusteringAssignment[i]);
754 vtkOutputNode->GetFieldData()->AddArray(vtkClusterAssignment);
755
756 // Construct multiblock
757 if(IsPersistenceDiagram and not OutputSegmentation) {
758 output_clusters->SetBlock(i, vtkOutputNode);
759 } else {
760 vtkMultiBlockDataSet::SafeDownCast(output_clusters->GetBlock(0))
761 ->SetBlock(i, vtkOutputNode);
762 if(not IsPersistenceDiagram)
763 vtkMultiBlockDataSet::SafeDownCast(output_clusters->GetBlock(1))
764 ->SetBlock(i, vtkOutputArc);
765 if(OutputSegmentation) {
766 int const segBlockID = 1 + !IsPersistenceDiagram;
767 vtkMultiBlockDataSet::SafeDownCast(
768 output_clusters->GetBlock(segBlockID))
769 ->SetBlock(i, vtkOutputSegmentation);
770 }
771 }
772 }
773 }
774
775 // ------------------------------------------
776 // --- Barycenter(s)
777 // ------------------------------------------
778 if(IsPersistenceDiagram) {
779 output_centroids->SetNumberOfBlocks(NumberOfBarycenters);
780 } else {
781 output_centroids->SetNumberOfBlocks(2);
782 vtkSmartPointer<vtkMultiBlockDataSet> const vtkBlockNodes2
784 vtkBlockNodes2->SetNumberOfBlocks(NumberOfBarycenters);
785 output_centroids->SetBlock(0, vtkBlockNodes2);
786 vtkSmartPointer<vtkMultiBlockDataSet> const vtkBlockArcs2
788 vtkBlockArcs2->SetNumberOfBlocks(NumberOfBarycenters);
789 output_centroids->SetBlock(1, vtkBlockArcs2);
790 }
791 for(unsigned int c = 0; c < NumberOfBarycenters; ++c) {
792
793 // Declare vtk objects
796 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputArc
798 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputSegmentation
800
801 // Fill vtk objects
802 ttkMergeTreeVisualization visuMakerBary;
803 visuMakerBary.setPlanarLayout(PlanarLayout);
805 BranchDecompositionPlanarLayout);
806 visuMakerBary.setBranchSpacing(BranchSpacing);
807 visuMakerBary.setRescaleTreesIndividually(RescaleTreesIndividually);
808 visuMakerBary.setOutputSegmentation(false);
809 visuMakerBary.setDimensionSpacing(DimensionSpacing);
810 visuMakerBary.setDimensionToShift(DimensionToShift);
811 visuMakerBary.setImportantPairs(ImportantPairs);
812 visuMakerBary.setMaximumImportantPairs(MaximumImportantPairs);
813 visuMakerBary.setMinimumImportantPairs(MinimumImportantPairs);
814 visuMakerBary.setImportantPairsSpacing(ImportantPairsSpacing);
815 visuMakerBary.setNonImportantPairsSpacing(NonImportantPairsSpacing);
816 visuMakerBary.setNonImportantPairsProximity(NonImportantPairsProximity);
817 visuMakerBary.setExcludeImportantPairsHigher(
818 ExcludeImportantPairsHigher);
819 visuMakerBary.setExcludeImportantPairsLower(ExcludeImportantPairsLower);
820 visuMakerBary.setIsPersistenceDiagram(IsPersistenceDiagram);
821 visuMakerBary.setShiftMode(1); // Star Barycenter
822 visuMakerBary.setTreesNodes(treesNodes);
823 visuMakerBary.setTreesNodeCorrMesh(trees1NodeCorrMesh);
824 visuMakerBary.setTreesSegmentation(treesSegmentation);
825 visuMakerBary.setVtkOutputNode(vtkOutputNode);
826 if(IsPersistenceDiagram)
827 visuMakerBary.setVtkOutputArc(vtkOutputNode);
828 else
829 visuMakerBary.setVtkOutputArc(vtkOutputArc);
830 visuMakerBary.setVtkOutputSegmentation(vtkOutputSegmentation);
831 visuMakerBary.setClusteringAssignment(clusteringAssignment);
832 visuMakerBary.setOutputMatchingBarycenter(outputMatchingBarycenter);
833 visuMakerBary.setPrintTreeId(c);
834 visuMakerBary.setPrintClusterId(c);
835 if(numInputs == 2 and NumberOfBarycenters == 1) {
836 visuMakerBary.setBarycenterPositionAlpha(BarycenterPositionAlpha);
837 visuMakerBary.setAlpha(Alpha);
838 }
839 visuMakerBary.setDebugLevel(this->debugLevel_);
840 visuMakerBary.setIsPersistenceDiagram(IsPersistenceDiagram);
841 visuMakerBary.setIsPDSadMax(JoinSplitMixtureCoefficient == 0);
842
843 visuMakerBary.makeTreesOutput<dataType>(
844 intermediateTrees, barycentersTree);
845 auto nodeCorrBaryT = visuMakerBary.getNodeCorr();
846 nodeCorrBary[c] = nodeCorrBaryT[c];
847 auto allBaryPercentMatchT = visuMakerBary.getAllBaryPercentMatch();
848 allBaryPercentMatch[c] = allBaryPercentMatchT[c];
849
850 if(IsPersistenceDiagram and JoinSplitMixtureCoefficient != 0
851 and JoinSplitMixtureCoefficient != 1) {
853 visuMakerBary, vtkOutputNode, intermediateTrees2, barycentersTree2,
854 outputMatchingBarycenter2, trees2NodeCorrMesh, nullptr, c,
855 nodeCorrBary2[c]);
856 allBaryPercentMatchT = visuMakerBary.getAllBaryPercentMatch();
857 allBaryPercentMatch2[c] = allBaryPercentMatchT[c];
858 }
859
860 // Field data
861 vtkNew<vtkDoubleArray> vtkClusterAssignment{};
862 vtkClusterAssignment->SetName("ClusterAssignment");
863 vtkClusterAssignment->SetNumberOfTuples(1);
864 vtkClusterAssignment->SetTuple1(0, c);
865 vtkOutputNode->GetFieldData()->AddArray(vtkClusterAssignment);
866
867 // Construct multiblock
868 if(IsPersistenceDiagram) {
869 output_centroids->SetBlock(c, vtkOutputNode);
870 } else {
871 vtkMultiBlockDataSet::SafeDownCast(output_centroids->GetBlock(0))
872 ->SetBlock(c, vtkOutputNode);
873 vtkMultiBlockDataSet::SafeDownCast(output_centroids->GetBlock(1))
874 ->SetBlock(c, vtkOutputArc);
875 }
876 }
877
878 // ------------------------------------------
879 // --- Matching
880 // ------------------------------------------
881 output_matchings->SetNumberOfBlocks(numInputs);
882 for(unsigned int c = 0; c < NumberOfBarycenters; ++c) {
883 for(int i = 0; i < numInputs; ++i) {
884 if(clusteringAssignment[i] != (int)c)
885 continue;
886
887 // Declare vtk objects
890 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNode1
891 = vtkUnstructuredGrid::SafeDownCast(
892 (IsPersistenceDiagram and not OutputSegmentation
893 ? output_clusters->GetBlock(i)
894 : vtkMultiBlockDataSet::SafeDownCast(
895 output_clusters->GetBlock(0))
896 ->GetBlock(i)));
897 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNode2
898 = vtkUnstructuredGrid::SafeDownCast(
899 (IsPersistenceDiagram ? output_centroids->GetBlock(c)
900 : vtkMultiBlockDataSet::SafeDownCast(
901 output_centroids->GetBlock(0))
902 ->GetBlock(c)));
903
904 // Fill vtk objects
905 ttkMergeTreeVisualization visuMakerMatching;
906 visuMakerMatching.setVtkOutputMatching(vtkOutputMatching);
907 visuMakerMatching.setOutputMatchingBarycenter(
908 outputMatchingBarycenter);
909 visuMakerMatching.setAllBaryPercentMatch(allBaryPercentMatch);
910 visuMakerMatching.setVtkOutputNode1(vtkOutputNode1);
911 visuMakerMatching.setVtkOutputNode2(vtkOutputNode2);
912 visuMakerMatching.setNodeCorr1(nodeCorr);
913 visuMakerMatching.setNodeCorr2(nodeCorrBary);
914 visuMakerMatching.setPrintTreeId(i);
915 visuMakerMatching.setPrintClusterId(c);
916 visuMakerMatching.setDebugLevel(this->debugLevel_);
917
918 visuMakerMatching.makeMatchingOutput<dataType>(
919 intermediateTrees, barycentersTree);
920 if(IsPersistenceDiagram and JoinSplitMixtureCoefficient != 0
921 and JoinSplitMixtureCoefficient != 1)
923 visuMakerMatching, vtkOutputMatching, intermediateTrees2,
924 barycentersTree2, outputMatchingBarycenter2, nodeCorr2,
925 nodeCorrBary2, allBaryPercentMatch2);
926
927 // Field data
928 vtkNew<vtkDoubleArray> vtkDistance{};
929 vtkDistance->SetName("Distance");
930 vtkDistance->SetNumberOfTuples(1);
931 vtkDistance->SetTuple1(0, finalDistances[i]);
932 vtkOutputMatching->GetFieldData()->AddArray(vtkDistance);
933
934 // Construct multiblock
935 output_matchings->SetBlock(i, vtkOutputMatching);
936 }
937 }
938 }
939 }
940
941 printMsg(
942 "Trees output", 1, t_makeTreesOutput.getElapsedTime(), this->threadNumber_);
943
944 return 1;
945}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
TTK VTK-filter that wraps the ttk::MergeTreeClustering module.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
virtual int IsA(const char *type)
int runCompute(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees2)
int run(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees2)
int runOutput(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees2)
int FillInputPortInformation(int port, vtkInformation *info) override
~ttkMergeTreeClustering() override
void setTreesSegmentation(std::vector< vtkDataSet * > &segmentation)
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)
void setVtkOutputMatching(vtkUnstructuredGrid *vtkMatching)
void setClusteringAssignment(std::vector< int > &asgn)
std::vector< std::vector< float > > getAllBaryPercentMatch()
std::vector< std::vector< SimplexId > > getNodeCorr()
void makeTreesOutput(FTMTree_MT *tree1)
void setVtkOutputNode2(vtkUnstructuredGrid *vtkNode2)
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 setTreesNodes(std::vector< vtkUnstructuredGrid * > &nodes)
void setAllBaryPercentMatch(std::vector< std::vector< float > > &baryPercentMatch)
virtual int setThreadNumber(const int threadNumber)
Definition BaseClass.h:80
int debugLevel_
Definition Debug.h:379
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
virtual int setDebugLevel(const int &debugLevel)
Definition Debug.cpp:147
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
void setAddNodes(bool addNodesT)
void setPostprocess(bool postproc)
void setDeterministic(bool deterministicT)
void setBarycenterSizeLimitPercent(double percent)
std::vector< double > getFinalDistances()
void execute(std::vector< ftm::MergeTree< dataType > > &trees, std::vector< double > &alphas, std::vector< std::vector< std::tuple< ftm::idNode, ftm::idNode, double > > > &finalMatchings, ftm::MergeTree< dataType > &baryMergeTree, bool finalAsgnDoubleInput=false, bool finalAsgnFirstInput=true)
void setBranchDecomposition(bool useBD)
void setNormalizedWasserstein(bool normalizedWasserstein)
void setEpsilon3Tree1(double epsilon)
void setEpsilonTree1(double epsilon)
void setAssignmentSolver(int assignmentSolver)
std::vector< std::vector< int > > getTreesNodeCorr()
void setEpsilon2Tree1(double epsilon)
void setEpsilonTree2(double epsilon)
void setPersistenceThreshold(double pt)
void setNonMatchingWeight(double weight)
void setCleanTree(bool clean)
void setEpsilon2Tree2(double epsilon)
void setKeepSubtree(bool keepSubtree)
void setEpsilon1UseFarthestSaddle(bool b)
void setDeleteMultiPersPairs(bool deleteMultiPersPairsT)
void setUseMinMaxPair(bool useMinMaxPair)
void setEpsilon3Tree2(double epsilon)
void setIsPersistenceDiagram(bool isPD)
void setMixtureCoefficient(double coef)
void setNoCentroids(unsigned int noCentroidsT)
std::vector< std::vector< int > > getTrees2NodeCorr()
void setPostprocess(bool postproc)
dataType execute(ftm::MergeTree< dataType > &mTree1, ftm::MergeTree< dataType > &mTree2, std::vector< std::tuple< ftm::idNode, ftm::idNode, double > > &outputMatching)
void setExcludeImportantPairsLower(std::string &d)
void setExcludeImportantPairsHigher(std::string &d)
double getElapsedTime()
Definition Timer.h:15
void mergeTreeToFTMTree(std::vector< MergeTree< dataType > > &trees, std::vector< ftm::FTMTree_MT * > &treesT)
bool constructTrees(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< MergeTree< dataType > > &intermediateTrees, std::vector< vtkUnstructuredGrid * > &treesNodes, std::vector< vtkUnstructuredGrid * > &treesArcs, std::vector< vtkDataSet * > &treesSegmentation, const std::vector< bool > &useSecondPairsTypeVec, int diagramPairTypes=0)
void loadBlocks(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, vtkMultiBlockDataSet *blocks)
void mergeTreesDoubleToTemplate(std::vector< MergeTree< double > > &mts, std::vector< MergeTree< dataType > > &newMts)
void mergeTreesTemplateToDouble(std::vector< MergeTree< dataType > > &mts, std::vector< MergeTree< double > > &newMts)
TTK base package defining the standard types.
vtkStandardNewMacro(ttkMergeTreeClustering)
void makeDoubleInputPersistenceDiagramOutput(ttkMergeTreeVisualization &visuMaker, vtkSmartPointer< vtkUnstructuredGrid > &vtkOutputNode, std::vector< FTMTree_MT * > &intermediateTrees2, std::vector< FTMTree_MT * > &barycentersTree2, std::vector< std::vector< std::vector< std::tuple< ttk::ftm::idNode, ttk::ftm::idNode, double > > > > &outputMatchingBarycenter2, std::vector< std::vector< int > > &trees2NodeCorrMesh, vtkUnstructuredGrid *treeNodes, int i, std::vector< SimplexId > &nodeCorr2)
void makeDoubleInputPersistenceDiagramMatching(ttkMergeTreeVisualization &visuMakerMatching, vtkSmartPointer< vtkUnstructuredGrid > &vtkOutputMatching, std::vector< FTMTree_MT * > &intermediateTrees2, std::vector< FTMTree_MT * > &barycentersTree2, std::vector< std::vector< std::vector< std::tuple< ttk::ftm::idNode, ttk::ftm::idNode, double > > > > &outputMatchingBarycenter2, std::vector< std::vector< SimplexId > > &nodeCorr2, std::vector< std::vector< SimplexId > > &nodeCorrBary2, std::vector< std::vector< float > > &allBaryPercentMatch2)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)