TTK
Loading...
Searching...
No Matches
ttkMergeTreePrincipalGeodesicsDecoding.cpp
Go to the documentation of this file.
2#include <ttkMergeTreeUtils.h>
3
4#include <vtkInformation.h>
5
6#include <vtkDataArray.h>
7#include <vtkDataSet.h>
8#include <vtkObjectFactory.h>
9#include <vtkPointData.h>
10#include <vtkSmartPointer.h>
11#include <vtkVariantArray.h>
12
13#include <ttkUtils.h>
14
15// A VTK macro that enables the instantiation of this class via ::New()
16// You do not have to modify this
18
33 this->SetNumberOfInputPorts(5);
34 this->SetNumberOfOutputPorts(1);
35}
36
45 int port, vtkInformation *info) {
46 if(port == 0)
47 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
48 else if(port == 1 or port == 2 or port == 3) {
49 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
50 if(port == 3) // correlation table
51 info->Set(vtkAlgorithm::INPUT_IS_OPTIONAL(), 1);
52 } else if(port == 4) {
53 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
54 info->Set(vtkAlgorithm::INPUT_IS_OPTIONAL(), 1);
55 } else
56 return 0;
57 return 1;
58}
59
76 int port, vtkInformation *info) {
77 if(port == 0) {
78 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
79 } else
80 return 0;
81 return 1;
82}
83
98 vtkInformation *ttkNotUsed(request),
99 vtkInformationVector **inputVector,
100 vtkInformationVector *outputVector) {
101 // ------------------------------------------------------------------------------------
102 // --- Get input object from input vector
103 // ------------------------------------------------------------------------------------
104 auto blockBary = vtkMultiBlockDataSet::GetData(inputVector[0], 0);
105 auto tableCoefficients = vtkTable::GetData(inputVector[1], 0);
106 auto tableVectors = vtkTable::GetData(inputVector[2], 0);
107 auto tableCorrelation = vtkTable::GetData(inputVector[3], 0);
108 auto blockInputTrees = vtkMultiBlockDataSet::GetData(inputVector[4], 0);
109
110 // ------------------------------------------------------------------------------------
111 // --- Load blocks
112 // ------------------------------------------------------------------------------------
113 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> inputBary, inputTrees;
114 ttk::ftm::loadBlocks(inputBary, blockBary);
115 ttk::ftm::loadBlocks(inputTrees, blockInputTrees);
116
117 // If we have already computed once but the input has changed
118 if((baryTreeNodes.size() != 0
119 and inputBary[0]->GetBlock(0) != baryTreeNodes[0])
120 or (allTs_.size() != 0
121 and int(allTs_[0].size()) != tableCoefficients->GetNumberOfRows())
122 or tableCoefficientsMTime != tableCoefficients->GetMTime()
123 or tableVectorsMTime != tableVectors->GetMTime()
124 or (tableCorrelation
125 and tableCorrelationMTime != tableCorrelation->GetMTime())
126 or (blockInputTrees
127 and blockInputTreesMTime != blockInputTrees->GetMTime())) {
128 resetDataVisualization();
129 }
130 tableCoefficientsMTime = tableCoefficients->GetMTime();
131 tableVectorsMTime = tableVectors->GetMTime();
132 if(tableCorrelation)
133 tableCorrelationMTime = tableCorrelation->GetMTime();
134 if(blockInputTrees)
135 blockInputTreesMTime = blockInputTrees->GetMTime();
136
137 // Parameters
138 printMsg("Load parameters from field data.");
139 std::vector<std::string> paramNames;
140 getParamNames(paramNames);
141 for(auto paramName : paramNames) {
142 auto array = tableCoefficients->GetFieldData()->GetArray(paramName.c_str());
143 if(array) {
144 double const value = array->GetTuple1(0);
145 setParamValueFromName(paramName, value);
146 } else
147 printMsg(" - " + paramName + " was not found in the field data.");
148 printMsg(" - " + paramName + " = "
150 }
152 printMsg("Computation with normalized Wasserstein.");
153 else
154 printMsg("Computation without normalized Wasserstein.");
155
156 // ------------------------------------------------------------------------------------
157 // --- Load tables
158 // ------------------------------------------------------------------------------------
159 auto tableVNoCols = tableVectors->GetNumberOfColumns() / inputBary.size();
160 unsigned int const numberOfGeodesics = tableVNoCols / 4;
161 auto numberOfInputs = tableCoefficients->GetNumberOfRows();
162
163 auto baryNoNodes = tableVectors->GetNumberOfRows();
164 auto baryNoNodes2 = baryNoNodes;
165 if(inputBary.size() > 1) {
166 for(unsigned int i = 0; i < 2; ++i) {
167 auto &baryNoNodesT = (i == 0 ? baryNoNodes : baryNoNodes2);
168 while(std::isnan(
169 tableVectors
170 ->GetColumnByName(
171 getTableVectorName(numberOfGeodesics, 0, 0, 0, (i == 1)).c_str())
172 ->GetVariantValue(baryNoNodesT - 1)
173 .ToDouble()))
174 --baryNoNodesT;
175 }
176 }
177
178 // -----------------
179 // Coefficients
180 // -----------------
181 // Pointers can't be used because we need to access the rows (VariantArray)
182 std::vector<std::string> nonFieldDataNames{"TreeID"};
183 allTs_.resize(numberOfGeodesics, std::vector<double>(numberOfInputs, 0.0));
184 for(unsigned int i = 0; i < numberOfGeodesics; ++i)
185 for(unsigned int j = 0; j < numberOfInputs; ++j) {
186 std::string const name = getTableCoefficientName(numberOfGeodesics, i);
187 allTs_[i][j] = tableCoefficients->GetColumnByName(name.c_str())
188 ->GetVariantValue(j)
189 .ToDouble();
190 nonFieldDataNames.push_back(name);
191 std::string const normName
192 = getTableCoefficientNormName(numberOfGeodesics, i);
193 nonFieldDataNames.push_back(normName);
194 }
196
197 // - aggregate input field data
198 vtkNew<vtkFieldData> fd{};
199 for(int i = 0; i < tableCoefficients->GetNumberOfColumns(); ++i) {
200 std::string const name{tableCoefficients->GetColumnName(i)};
201 if(std::count(nonFieldDataNames.begin(), nonFieldDataNames.end(), name)
202 == 0) {
203 auto array = tableCoefficients->GetColumn(i);
204 array->SetName(name.c_str());
205 fd->AddArray(array);
206 }
207 }
208 inputFieldData = vtkFieldData::New();
209 inputFieldData->ShallowCopy(fd);
210
211 // -----------------
212 // Vectors
213 // -----------------
214 pVS_.resize(numberOfGeodesics, std::vector<double *>(2));
215 pV2s_ = pVS_;
216 if(inputBary.size() > 1) {
217 pTrees2Vs_.resize(numberOfGeodesics, std::vector<double *>(2));
219 }
220 for(unsigned int h = 0; h < inputBary.size(); ++h) {
221 auto &pVS = (h == 0 ? pVS_ : pTrees2Vs_);
222 auto &pV2s = (h == 0 ? pV2s_ : pTrees2V2s_);
223 bool const secondInput = (h == 1);
224 for(unsigned int i = 0; i < numberOfGeodesics; ++i) {
225 for(unsigned int k = 0; k < 2; ++k) {
226 std::string const name1
227 = getTableVectorName(numberOfGeodesics, i, 0, k, secondInput);
228 std::string const name2
229 = getTableVectorName(numberOfGeodesics, i, 1, k, secondInput);
230 pVS[i][k] = ttkUtils::GetPointer<double>(vtkDataArray::SafeDownCast(
231 tableVectors->GetColumnByName(name1.c_str())));
232 pV2s[i][k] = ttkUtils::GetPointer<double>(vtkDataArray::SafeDownCast(
233 tableVectors->GetColumnByName(name2.c_str())));
234 }
235 }
236 }
237 vSize_ = baryNoNodes;
238 vSize2_ = baryNoNodes2;
239
240 // -----------------
241 // Correlation
242 // -----------------
243 if(tableCorrelation) {
244 pBranchesCorrelationMatrix_.resize(numberOfGeodesics);
246
247 auto tableCorrNoRows = tableCorrelation->GetNumberOfRows();
248 auto baryNodeIdArray = tableCorrelation->GetColumnByName("BaryNodeId");
249
250 for(unsigned int j = 0; j < numberOfGeodesics; ++j) {
251 std::string const name = getTableCorrelationName(numberOfGeodesics, j);
252 std::string const name2
253 = getTableCorrelationPersName(numberOfGeodesics, j);
255 = ttkUtils::GetPointer<double>(vtkDataArray::SafeDownCast(
256 tableCorrelation->GetColumnByName(name.c_str())));
258 = ttkUtils::GetPointer<double>(vtkDataArray::SafeDownCast(
259 tableCorrelation->GetColumnByName(name2.c_str())));
260 }
261
262 // Pointers can't be used because processing on these data will be done
263 // later
264 baryMatchings_.resize(numberOfInputs);
265 for(unsigned int i = 0; i < numberOfInputs; ++i) {
266 std::string const name = getTableCorrelationTreeName(numberOfInputs, i);
267 auto array = tableCorrelation->GetColumnByName(name.c_str());
268 if(array) {
269 baryMatchings_[i].resize(tableCorrNoRows);
270 for(unsigned int j = 0; j < tableCorrNoRows; ++j) {
271 auto baryNodeId = baryNodeIdArray->GetVariantValue(j).ToInt();
272 baryMatchings_[i][j] = std::make_tuple(
273 baryNodeId, array->GetVariantValue(j).ToDouble(), 0.0);
274 }
275 }
276 }
277 } else {
280 baryMatchings_.clear();
281 }
282
283 // ------------------------------------------------------------------------------------
284
285 return run<float>(outputVector, inputBary, inputTrees);
286}
287
288template <class dataType>
290 vtkInformationVector *outputVector,
291 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputBary,
292 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees) {
293 if(not isDataVisualizationFilled())
294 runCompute<dataType>(outputVector, inputBary, inputTrees);
295 runOutput<dataType>(outputVector, inputBary, inputTrees);
296 return 1;
297}
298
299template <class dataType>
301 vtkInformationVector *ttkNotUsed(outputVector),
302 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputBary,
303 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees) {
304 useDoubleInput_ = (inputBary.size() > 1);
305 processFirstInput = (not useDoubleInput_ or not ProcessSecondInput);
306 // ------------------------------------------------------------------------------------
307 // --- Construct trees
308 // ------------------------------------------------------------------------------------
309 const int numTrees = inputTrees.size();
310
311 setDataVisualization(numTrees);
312
313 std::vector<ttk::ftm::MergeTree<dataType>> baryDTree, inputDTrees;
314
315 std::vector<bool> useSadMaxPairsVec{false, true};
316 if(not useDoubleInput_ and mixtureCoefficient_ == 0)
317 useSadMaxPairsVec.erase(useSadMaxPairsVec.begin()); // {true}
318 ttk::ftm::constructTrees<dataType>(inputBary, baryDTree, baryTreeNodes,
319 baryTreeArcs, baryTreeSegmentation,
320 useSadMaxPairsVec);
321
322 if(OutputInputTrees
323 or (ReconstructInputTrees
325 bool const useSadMaxPairs
326 = (useDoubleInput_ and not processFirstInput) or mixtureCoefficient_ == 0;
327 bool const isInputPD = ttk::ftm::constructTrees<dataType>(
328 inputTrees, inputDTrees, inputTreesNodes, inputTreesArcs,
329 inputTreesSegmentation, useSadMaxPairs);
330 if(not isInputPD and isPersistenceDiagram_)
331 mtsFlattening(inputDTrees);
332 }
333
334 //------------------------------------------------------------------------------------
335 // --- Call base
336 //------------------------------------------------------------------------------------
337 std::string const preFix{(processFirstInput ? "" : "Second Input ")};
338 auto &baryToUse = (processFirstInput ? baryDTree[0] : baryDTree[1]);
339
340 // Geodesics distances
341 if(geodesicsDistances_.size() == 0) {
343 printMsg("Compute Geodesics Distance...");
344 computeGeodesicsDistance<dataType>(baryDTree);
345 }
346
347 // Input Trees
348 if(OutputInputTrees) {
350 printMsg("Process Input Trees...");
351 processInputTrees<dataType>(inputDTrees);
352 }
353
354 // Reconstructed Trees
355 if(ReconstructInputTrees) {
357 printMsg("Reconstruct " + preFix + "Input Trees...");
358 std::vector<ttk::ftm::MergeTree<dataType>> reconstructedTreesT;
359 reconstruction<dataType>(baryToUse, inputDTrees, reconstructedTreesT,
360 reconstructionErrors, recInputMatchings,
361 recBaryMatchings, !processFirstInput);
362 ttk::ftm::mergeTreesTemplateToDouble<dataType>(
363 reconstructedTreesT, reconstructedTrees);
364 } else
365 reconstructedTrees.clear();
366
367 // Geodesics Trees
368 if(ConstructGeodesicsTrees) {
370 printMsg("Construct Geodesics Trees...");
371 std::vector<std::vector<ttk::ftm::MergeTree<dataType>>> geodesicsTreesT;
372 constructGeodesicsTrees<dataType>(
373 baryToUse, geodesicsTreesT, !processFirstInput);
374 ttk::ftm::mergeTreesTemplateToDouble<dataType>(
375 geodesicsTreesT, geodesicsTrees);
376 } else
377 geodesicsTrees.clear();
378
379 // Ellipses Trees
380 if(ConstructEllipses) {
382 printMsg("Construct Ellipses Trees...");
383 std::vector<ttk::ftm::MergeTree<dataType>> geodesicsEllipsesT;
384 constructGeodesicsEllipses<dataType>(
385 baryToUse, geodesicsEllipsesT, !processFirstInput);
386 ttk::ftm::mergeTreesTemplateToDouble<dataType>(
387 geodesicsEllipsesT, geodesicsEllipses);
388 } else
389 geodesicsEllipses.clear();
390
391 // Rectangle Trees
392 if(ConstructRectangle) {
394 printMsg("Construct Rectangle Trees...");
395 std::vector<ttk::ftm::MergeTree<dataType>> geodesicsRectangleT;
396 constructGeodesicsRectangle<dataType>(
397 baryToUse, geodesicsRectangleT, RectangleMultiplier, !processFirstInput);
398 ttk::ftm::mergeTreesTemplateToDouble<dataType>(
399 geodesicsRectangleT, geodesicsRectangle);
400 } else
401 geodesicsRectangle.clear();
402
403 // Surface Trees
404 if(ConstructSurface) {
406 printMsg("Construct " + preFix + "Surface Trees...");
407 std::vector<ttk::ftm::MergeTree<dataType>> geodesicsSurfaceT;
408 constructGeodesicsSurface<dataType>(
409 baryToUse, geodesicsSurfaceT, !processFirstInput);
410 ttk::ftm::mergeTreesTemplateToDouble<dataType>(
411 geodesicsSurfaceT, geodesicsSurface);
412 } else
413 geodesicsSurface.clear();
414
415 ttk::ftm::mergeTreesTemplateToDouble<dataType>(baryDTree, baryMTree);
416 ttk::ftm::mergeTreesTemplateToDouble<dataType>(inputDTrees, inputMTrees);
417
418 return 1;
419}
420
421template <class dataType>
423 vtkInformationVector *outputVector,
424 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &ttkNotUsed(inputBary),
425 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees) {
427 printMsg("Make output...");
428 // ------------------------------------------------------------------------------------
429 // --- Create output
430 // ------------------------------------------------------------------------------------
431 auto output = vtkMultiBlockDataSet::GetData(outputVector, 0);
432
433 unsigned int const noInput = (OutputInputTrees ? inputTrees.size() : 0);
434 unsigned int const noReconst = reconstructedTrees.size();
435 unsigned int const noBary = (OutputBarycenter ? 1 : 0);
436 unsigned int const noGeod
437 = geodesicsTrees.size()
438 * (geodesicsTrees.size() == 0 ? 0 : geodesicsTrees[0].size());
439 unsigned int const noEllipses = geodesicsEllipses.size();
440 unsigned int const noRectangle = geodesicsRectangle.size();
441 unsigned int const noSurface = geodesicsSurface.size();
442
443 printMsg("noInput = " + std::to_string(noInput));
444 printMsg("noReconst = " + std::to_string(noReconst));
445 printMsg("noBary = " + std::to_string(noBary));
446 printMsg("noGeod = " + std::to_string(noGeod));
447 printMsg("noEllipses = " + std::to_string(noEllipses));
448 printMsg("noRectangle = " + std::to_string(noRectangle));
449 printMsg("noSurface = " + std::to_string(noSurface));
450
451 std::vector<unsigned int> allNo{
452 noInput, noReconst, noBary, noGeod, noEllipses, noRectangle, noSurface};
453 std::vector<unsigned int> steps(allNo.size());
454 for(unsigned int i = 0; i < allNo.size(); ++i)
455 steps[i] = (i == 0 ? 0 : steps[i - 1]) + allNo[i];
456
457 unsigned int const noTreesTotal = steps[steps.size() - 1];
458
459 unsigned int const noBlocks
460 = 2 - isPersistenceDiagram_ + OutputInputTreesSegmentation;
461 output->SetNumberOfBlocks(noBlocks);
462 vtkSmartPointer<vtkMultiBlockDataSet> const vtkBlockNodes
464 vtkBlockNodes->SetNumberOfBlocks(noTreesTotal);
465 output->SetBlock(0, vtkBlockNodes);
466 if(not isPersistenceDiagram_) {
469 vtkBlockArcs->SetNumberOfBlocks(noTreesTotal);
470 output->SetBlock(1, vtkBlockArcs);
471 }
472 if(OutputInputTreesSegmentation) {
475 vtkBlockSegs->SetNumberOfBlocks(noTreesTotal);
476 output->SetBlock(2 - isPersistenceDiagram_, vtkBlockSegs);
477 }
478
479 int const geodesicsTreesOffset = std::max((int)geodesicsTrees.size() - 1, 0);
480
481 // ------------------------------------------
482 // --- Trees
483 // ------------------------------------------
484 std::vector<std::vector<ttk::ftm::idNode>> matchingMatrix;
485 if(!baryMatchings_.empty())
486 getMatchingMatrix<double>(
487 baryMTree[0], inputMTrees, baryMatchings_, matchingMatrix);
488 // TODO compute matching to barycenter if correlation matrix is not provided
490 and (inputMTrees.empty() or baryMatchings_.empty()))
491 printWrn("Please provide input trees and correlation matrix to transfer "
492 "input trees information.");
493 // TODO fix if an interpolation is empty
494#ifdef TTK_ENABLE_OPENMP
495#pragma omp parallel for schedule(dynamic) num_threads(this->threadNumber_)
496#endif
497 for(unsigned int i = 0; i < noTreesTotal; ++i) {
498 int index = 0;
499 int treeType = 0;
501 std::vector<double> ts;
502 if(i < steps[0]) {
503 // Input Trees
504 index = i;
505 ttk::ftm::mergeTreeDoubleToTemplate<dataType>(inputMTrees[index], mt);
506 treeType = -1;
507 ts = allTreesTs_[index];
508 } else if(i < steps[1]) {
509 // Reconstructed trees
510 index = i - steps[0];
511 ttk::ftm::mergeTreeDoubleToTemplate<dataType>(
512 reconstructedTrees[index], mt);
513 treeType = 0;
514 ts = allTreesTs_[index];
515 } else if(i < steps[2]) {
516 // Barycenter
517 ttk::ftm::mergeTreeDoubleToTemplate<dataType>(baryMTree[0], mt);
518 treeType = 1;
519 std::array<double, 2> middle;
520 getGeodesicsMiddle<dataType>(mt, pVS_, pV2s_, vSize_, middle);
521 ts.resize(2);
522 ts[0] = middle[0];
523 ts[1] = middle[1];
524 } else if(i < steps[3]) {
525 // Geodesics Trees
526 index = i - steps[2];
527 int const i0 = index / geodesicsTrees[0].size();
528 int const i1 = index % geodesicsTrees[0].size();
529 ttk::ftm::mergeTreeDoubleToTemplate<dataType>(geodesicsTrees[i0][i1], mt);
530 treeType = 2 + i0;
531 ts = tGeodesics_[i0][i1];
532 } else if(i < steps[4]) {
533 // Ellipses Trees
534 index = i - steps[3];
535 ttk::ftm::mergeTreeDoubleToTemplate<dataType>(
536 geodesicsEllipses[index], mt);
537 treeType = 3 + geodesicsTreesOffset;
538 ts = tEllipses_[index];
539 } else if(i < steps[5]) {
540 // Rectangle Trees
541 index = i - steps[4];
542 ttk::ftm::mergeTreeDoubleToTemplate<dataType>(
543 geodesicsRectangle[index], mt);
544 treeType = 4 + geodesicsTreesOffset;
545 ts = tRectangle_[index];
546 } else if(i < steps[6]) {
547 // Surface Trees
548 index = i - steps[5];
549 ttk::ftm::mergeTreeDoubleToTemplate<dataType>(
550 geodesicsSurface[index], mt);
551 treeType = 5 + geodesicsTreesOffset;
552 ts = tSurface_[index];
553 }
554
555 bool const isInputTree = (i < steps[0]);
556 bool const isInputTreeAndGotMesh
557 = (isInputTree and inputTrees[i]->GetNumberOfBlocks() >= 3);
558 bool const isReconstructedTree = (i >= steps[0] and i < steps[1]);
559 bool const isInputOrReconstructedTree
560 = (isInputTree or isReconstructedTree);
561
562 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNode
564 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputArc
566 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputSegmentation
568
570
571 if(isInputTree) {
572 for(unsigned int g = 0; g < pBranchesCorrelationMatrix_.size(); ++g) {
573 std::vector<double> corr(mt.tree.getNumberOfNodes(), 0.0),
574 corrPers = corr;
575 for(unsigned int n = 0; n < vSize_; ++n) {
576 if((int)matchingMatrix[n][i] != -1) {
577 corr[matchingMatrix[n][i]] = pBranchesCorrelationMatrix_[g][n];
578 corrPers[matchingMatrix[n][i]] = pPersCorrelationMatrix_[g][n];
579 }
580 }
581 std::string name
583 visuMaker.addCustomArray(name, corr);
584 std::string name2
586 visuMaker.addCustomArray(name2, corrPers);
587 }
588 if(!baryMatchings_.empty()) {
589 std::vector<double> baryNodePers(mt.tree.getNumberOfNodes(), -1);
590 std::vector<int> baryNodeID(mt.tree.getNumberOfNodes(), -1);
591 for(unsigned int n = 0; n < vSize_; ++n) {
592 if((int)matchingMatrix[n][i] != -1) {
593 baryNodeID[matchingMatrix[n][i]] = n;
594 baryNodePers[matchingMatrix[n][i]]
595 = baryMTree[0].tree.template getNodePersistence<double>(n);
596 }
597 }
598 std::string nameBaryNodePers{"BaryNodePers"};
599 visuMaker.addCustomArray(nameBaryNodePers, baryNodePers);
600 std::string nameBaryNodeId{"BaryNodeId"};
601 visuMaker.addCustomIntArray(nameBaryNodeId, baryNodeID);
602 }
603 // visuMaker.setShiftMode(0); // Star
604 visuMaker.setShiftMode(-1); // None
605 visuMaker.setISampleOffset(i);
606 visuMaker.setNoSampleOffset(noInput - 1);
607 visuMaker.setTreesNodes(inputTreesNodes[index]);
608 visuMaker.setTreesNodeCorrMesh(treesNodeCorr_[index]);
609 visuMaker.setDimensionSpacing(2);
611 if(isInputTreeAndGotMesh and OutputInputTreesSegmentation) {
612 visuMaker.setOutputSegmentation(true);
613 visuMaker.setVtkOutputSegmentation(vtkOutputSegmentation);
614 visuMaker.setTreesSegmentation(inputTreesSegmentation[index]);
615 }
616 } else
617 visuMaker.setShiftMode(2); // Line
619 visuMaker.setPlanarLayout(true);
620
621 // TODO transfer other information?
622 if(isReconstructedTree) {
625 ttk::ftm::mergeTreeDoubleToTemplate<dataType>(baryMTree[0], baryMT);
626 std::vector<ttk::ftm::idNode> matchingVector;
628 mt, baryMT, recBaryMatchings[index], matchingVector);
629 std::vector<int> baryNodeID(mt.tree.getNumberOfNodes(), -1);
630 for(unsigned int n = 0; n < vSize_; ++n) {
631 auto match = matchingVector[n];
632 if((int)match != -1)
633 baryNodeID[match] = n;
634 }
635 std::string nameBaryNodeId{"BaryNodeId"};
636 visuMaker.addCustomIntArray(nameBaryNodeId, baryNodeID);
637 }
638 if(transferInputTreesInformation_ and !inputMTrees.empty()
639 and !baryMatchings_.empty()) {
641 ttk::ftm::mergeTreeDoubleToTemplate<dataType>(
642 inputMTrees[index], inputMT);
643 std::vector<ttk::ftm::idNode> matchingVector;
645 mt, inputMT, recInputMatchings[index], matchingVector);
646 std::vector<int> baryNodeID(mt.tree.getNumberOfNodes(), -1);
647 for(unsigned int n = 0; n < vSize_; ++n) {
648 auto match = matchingMatrix[n][index];
649 if((int)match != -1) {
650 match = matchingVector[match];
651 if((int)match != -1)
652 baryNodeID[match] = n;
653 }
654 }
655 std::string nameBaryNodeId{"BaryNodeId"};
656 visuMaker.addCustomIntArray(nameBaryNodeId, baryNodeID);
657 }
658 }
659
660 visuMaker.setVtkOutputNode(vtkOutputNode);
662 visuMaker.setVtkOutputArc(vtkOutputArc);
663 else {
664 visuMaker.setVtkOutputArc(vtkOutputNode);
666 visuMaker.setIsPDSadMax(!processFirstInput);
667 else
668 visuMaker.setIsPDSadMax(mixtureCoefficient_ == 0);
669 }
670 visuMaker.setDebugLevel(this->debugLevel_);
672 visuMaker.makeTreesOutput<dataType>(&(mt.tree));
673
674 // Copy Input Field Data
675 if(noReconst != 0 or noInput != 0) {
676 vtkNew<vtkFieldData> fd{};
677 if(isInputOrReconstructedTree)
678 fd->DeepCopy(inputFieldData);
679 else
680 fd->CopyStructure(inputFieldData);
681 for(int a = 0; a < inputFieldData->GetNumberOfArrays(); ++a) {
682 auto array = fd->GetAbstractArray(a);
683 array->SetNumberOfTuples(1);
684 if(isInputOrReconstructedTree) {
685 vtkNew<vtkIdList> id;
686 id->InsertNextId(i % (noInput == 0 ? noReconst : noInput));
687 inputFieldData->GetAbstractArray(a)->GetTuples(id, array);
688 } else {
689 auto downArray = vtkDataArray::SafeDownCast(array);
690 if(downArray)
691 downArray->SetTuple1(0, std::nan(""));
692 }
693 }
694 vtkOutputNode->GetFieldData()->ShallowCopy(fd);
695 }
696
697 // Field Data
698 vtkNew<vtkIntArray> treeTypeArray{};
699 treeTypeArray->SetName("TreeType");
700 treeTypeArray->InsertNextTuple1(treeType);
701 vtkOutputNode->GetFieldData()->AddArray(treeTypeArray);
702
703 vtkNew<vtkIntArray> treeIDArray{};
704 treeIDArray->SetName("TreeID");
705 treeIDArray->InsertNextTuple1(i);
706 vtkOutputNode->GetFieldData()->AddArray(treeIDArray);
707
708 for(unsigned int j = 0; j < ts.size(); ++j) {
709 vtkNew<vtkDoubleArray> tArray{};
710 std::string const name = getTableCoefficientName(ts.size(), j);
711 tArray->SetName(name.c_str());
712 tArray->InsertNextTuple1(ts[j]);
713 vtkOutputNode->GetFieldData()->AddArray(tArray);
714
715 vtkNew<vtkDoubleArray> scaledTArray{};
716 std::string const name2 = getTableCoefficientNormName(ts.size(), j);
717 scaledTArray->SetName(name2.c_str());
718 scaledTArray->InsertNextTuple1(ts[j] * geodesicsDistances_[j]);
719 vtkOutputNode->GetFieldData()->AddArray(scaledTArray);
720 }
721
722 bool const isSurfaceTree = (treeType == 5 + geodesicsTreesOffset);
723 vtkNew<vtkIntArray> isSurfaceArray{};
724 isSurfaceArray->SetName("isSurface");
725 isSurfaceArray->InsertNextTuple1(isSurfaceTree);
726 vtkOutputNode->GetFieldData()->AddArray(isSurfaceArray);
727 bool const isBoundary = (isSurfaceTree ? surfaceIsBoundary_[index] : false);
728 vtkNew<vtkIntArray> isBoundaryArray{};
729 isBoundaryArray->SetName("isBoundary");
730 isBoundaryArray->InsertNextTuple1(isBoundary);
731 vtkOutputNode->GetFieldData()->AddArray(isBoundaryArray);
732 int const boundaryID = (isSurfaceTree ? surfaceBoundaryID_[index] : -1);
733 vtkNew<vtkIntArray> BoundaryIDArray{};
734 BoundaryIDArray->SetName("BoundaryID");
735 BoundaryIDArray->InsertNextTuple1(boundaryID);
736 vtkOutputNode->GetFieldData()->AddArray(BoundaryIDArray);
737
738 vtkNew<vtkDoubleArray> reconstructionErrorArray{};
739 reconstructionErrorArray->SetName("ReconstructionError");
740 auto err = (isReconstructedTree and reconstructionErrors.size() != 0
741 ? reconstructionErrors[index]
742 : std::nan(""));
743 reconstructionErrorArray->InsertNextTuple1(err);
744 vtkOutputNode->GetFieldData()->AddArray(reconstructionErrorArray);
745
746 // Add new block
747 vtkMultiBlockDataSet::SafeDownCast(output->GetBlock(0))
748 ->SetBlock(i, vtkOutputNode);
750 vtkMultiBlockDataSet::SafeDownCast(output->GetBlock(1))
751 ->SetBlock(i, vtkOutputArc);
752 if(OutputInputTreesSegmentation)
753 vtkMultiBlockDataSet::SafeDownCast(
754 output->GetBlock(2 - isPersistenceDiagram_))
755 ->SetBlock(i, vtkOutputSegmentation);
756 }
757
758 // Field Data Input Parameters
759 std::vector<std::string> paramNames;
760 getParamNames(paramNames);
761 for(auto paramName : paramNames) {
762 vtkNew<vtkDoubleArray> array{};
763 array->SetName(paramName.c_str());
764 array->InsertNextTuple1(getParamValueFromName(paramName));
765 output->GetFieldData()->AddArray(array);
766 }
767
768 return 1;
769}
#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::MergeTreePrincipalGeodesicsDecoding module.
int runOutput(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputBary, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees)
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int runCompute(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputBary, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int run(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputBary, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees)
void setTreesSegmentation(std::vector< vtkDataSet * > &segmentation)
void setVtkOutputSegmentation(vtkDataSet *vtkSegmentation)
void setTreesNodeCorrMesh(std::vector< std::vector< int > > &nodeCorrMesh)
void addCustomIntArray(std::string &name, std::vector< int > &vec)
void makeTreesOutput(FTMTree_MT *tree1)
void setVtkOutputArc(vtkUnstructuredGrid *vtkArc)
void setVtkOutputNode(vtkUnstructuredGrid *vtkNode)
void addCustomArray(std::string &name, std::vector< double > &vec)
void setTreesNodes(std::vector< vtkUnstructuredGrid * > &nodes)
int debugLevel_
Definition Debug.h:379
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
virtual int setDebugLevel(const int &debugLevel)
Definition Debug.cpp:147
std::string getTableCorrelationTreeName(int noTrees, int treeNum)
std::string getTableVectorName(int noGeodesics, int geodesicNum, int vId, int vComp, bool isSecondInput=false)
std::string getTableCorrelationPersName(int noGeodesics, int geodesicNum)
std::string getTableCoefficientNormName(int noGeodesics, int geodesicNum)
void getInverseMatchingVector(ftm::MergeTree< dataType > &barycenter, ftm::MergeTree< dataType > &tree, std::vector< std::tuple< ftm::idNode, ftm::idNode, double > > &matchings, std::vector< ftm::idNode > &matchingVector)
std::string getTableCorrelationName(int noGeodesics, int geodesicNum)
std::string getTableCoefficientName(int noGeodesics, int geodesicNum)
void setParamValueFromName(std::string &paramName, double value)
void getParamNames(std::vector< std::string > &paramNames)
std::vector< std::vector< int > > treesNodeCorr_
void mtsFlattening(std::vector< ftm::MergeTree< dataType > > &mts)
double getParamValueFromName(std::string &paramName)
std::vector< std::vector< double > > allTreesTs_
std::vector< std::vector< double > > allTs_
std::vector< std::vector< std::tuple< ftm::idNode, ftm::idNode, double > > > baryMatchings_
std::vector< std::vector< std::vector< double > > > tGeodesics_
idNode getNumberOfNodes() const
Definition FTMTree_MT.h:389
std::string to_string(__int128)
Definition ripserpy.cpp:99
void transposeMatrix(const std::vector< std::vector< T > > &a, std::vector< std::vector< T > > &out)
Definition Geometry.cpp:764
void loadBlocks(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, vtkMultiBlockDataSet *blocks)
ftm::FTMTree_MT tree
Definition FTMTree_MT.h:901
vtkStandardNewMacro(ttkMergeTreePrincipalGeodesicsDecoding)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)