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