TTK
Loading...
Searching...
No Matches
ttkMergeTreeAutoencoder.cpp
Go to the documentation of this file.
4#include <ttkMergeTreeUtils.h>
6
7#include <vtkInformation.h>
8
9#include <vtkDataArray.h>
10#include <vtkDataSet.h>
11#include <vtkFloatArray.h>
12#include <vtkObjectFactory.h>
13#include <vtkPointData.h>
14#include <vtkSmartPointer.h>
15#include <vtkTable.h>
16#include <vtkUnsignedIntArray.h>
17
18#include <ttkMacros.h>
19#include <ttkUtils.h>
20
21// A VTK macro that enables the instantiation of this class via ::New()
22// You do not have to modify this
24
38 this->SetNumberOfInputPorts(3);
39 this->SetNumberOfOutputPorts(4);
40}
41
50 vtkInformation *info) {
51 if(port == 0) {
52 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
53 } else if(port == 1) {
54 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
55 info->Set(vtkAlgorithm::INPUT_IS_OPTIONAL(), 1);
56 } else if(port == 2) {
57 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
58 info->Set(vtkAlgorithm::INPUT_IS_OPTIONAL(), 1);
59 } else
60 return 0;
61 return 1;
62}
63
80 vtkInformation *info) {
81 if(port == 0 or port == 1 or port == 2 or port == 3) {
82 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
83 } else {
84 return 0;
85 }
86 return 1;
87}
88
103 vtkInformationVector **inputVector,
104 vtkInformationVector *outputVector) {
105#ifndef TTK_ENABLE_TORCH
106 TTK_FORCE_USE(inputVector);
107 TTK_FORCE_USE(outputVector);
108 printErr("This filter requires Torch.");
109 return 0;
110#else
111 // ------------------------------------------------------------------------------------
112 // --- Get input object from input vector
113 // ------------------------------------------------------------------------------------
114 auto blocks = vtkMultiBlockDataSet::GetData(inputVector[0], 0);
115 auto blocks2 = vtkMultiBlockDataSet::GetData(inputVector[1], 0);
116 auto table = vtkTable::GetData(inputVector[2], 0);
117
118 // ------------------------------------------------------------------------------------
119 // --- Load blocks
120 // ------------------------------------------------------------------------------------
121 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> inputTrees, inputTrees2;
122 ttk::ftm::loadBlocks(inputTrees, blocks);
123 ttk::ftm::loadBlocks(inputTrees2, blocks2);
124
125 // Load table
126 clusterAsgn_.clear();
127 vtkAbstractArray *clusterAsgn;
128 if(table) {
129 clusterAsgn = this->GetInputArrayToProcess(0, inputVector);
130 if(clusterAsgn) {
131 clusterAsgn_.resize(clusterAsgn->GetNumberOfValues());
132 for(unsigned int i = 0; i < clusterAsgn_.size(); ++i)
133 clusterAsgn_[i] = clusterAsgn->GetVariantValue(i).ToInt();
134 }
135 }
136 if((not table or not clusterAsgn) and clusteringLossWeight_ != 0) {
137 printErr(
138 "You must provide a table column in info input to use clustering loss");
139 return 0;
140 }
141 if(clusteringLossWeight_ != 0) {
142 std::stringstream ss;
143 for(auto &e : clusterAsgn_)
144 ss << e << " ";
145 printMsg(ss.str());
146 }
147
148 // ------------------------------------------------------------------------------------
149 // If we have already computed once but the input has changed
150 if((treesNodes.size() != 0 and inputTrees[0]->GetBlock(0) != treesNodes[0])
151 or (treesNodes2.size() != inputTrees2.size()))
152 resetDataVisualization();
153
154 // Parameters
156 if(not normalizedWasserstein_) {
157 oldEpsilonTree1 = epsilonTree1_;
158 epsilonTree1_ = 100;
159 } else
160 epsilonTree1_ = oldEpsilonTree1;
162 printMsg("Computation with normalized Wasserstein.");
163 else
164 printMsg("Computation without normalized Wasserstein.");
165
166 return run(outputVector, inputTrees, inputTrees2);
167#endif
168}
169
170#ifdef TTK_ENABLE_TORCH
172 vtkInformationVector *outputVector,
173 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
174 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2) {
175 runCompute(outputVector, inputTrees, inputTrees2);
176 runOutput(outputVector, inputTrees, inputTrees2);
177 return 1;
178}
179
181 vtkInformationVector *ttkNotUsed(outputVector),
182 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
183 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2) {
184 // ------------------------------------------------------------------------------------
185 // --- Construct trees
186 // ------------------------------------------------------------------------------------
187 std::vector<ttk::ftm::MergeTree<float>> intermediateMTrees,
188 intermediateMTrees2;
189
190 bool useSadMaxPairs = (mixtureCoefficient_ == 0);
191 isPersistenceDiagram_ = ttk::ftm::constructTrees<float>(
192 inputTrees, intermediateMTrees, treesNodes, treesArcs, treesSegmentation,
193 useSadMaxPairs);
194 // If merge trees are provided in input and normalization is not asked
198 or (mixtureCoefficient_ != 0 and mixtureCoefficient_ != 1)) {
199 auto &inputTrees2ToUse
200 = (not isPersistenceDiagram_ ? inputTrees2 : inputTrees);
201 ttk::ftm::constructTrees<float>(inputTrees2ToUse, intermediateMTrees2,
202 treesNodes2, treesArcs2, treesSegmentation2,
203 !useSadMaxPairs);
204 }
206
207 const int numInputs = intermediateMTrees.size();
208 const int numInputs2 = intermediateMTrees2.size();
209 setDataVisualization(numInputs, numInputs2);
210
211 // ------------------------------------------------------------------------------------
212 // --- Call base
213 // ------------------------------------------------------------------------------------
214 execute(intermediateMTrees, intermediateMTrees2);
215
216 ttk::ftm::mergeTreesTemplateToDouble<float>(
217 intermediateMTrees, intermediateDTrees);
218
219 return 1;
220}
221
222// TODO manage double input
224 vtkInformationVector *outputVector,
225 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
226 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &ttkNotUsed(inputTrees2)) {
227 if(not createOutput_)
228 return 1;
229 // ------------------------------------------------------------------------------------
230 // --- Create output
231 // ------------------------------------------------------------------------------------
232 auto output_origins = vtkMultiBlockDataSet::GetData(outputVector, 0);
233 auto output_vectors = vtkMultiBlockDataSet::GetData(outputVector, 1);
234 auto output_coef = vtkMultiBlockDataSet::GetData(outputVector, 2);
235 auto output_data = vtkMultiBlockDataSet::GetData(outputVector, 3);
236
237 // ------------------------------------------
238 // --- Tracking information
239 // ------------------------------------------
240 auto originsMatchingSize = originsMatchings_.size();
241 std::vector<std::vector<ttk::ftm::idNode>> originsMatchingVectorT(
242 originsMatchingSize),
243 invOriginsMatchingVectorT = originsMatchingVectorT;
244 for(unsigned int l = 0; l < originsMatchingVectorT.size(); ++l) {
245 auto &tree1 = (l == 0 ? origins_[0] : originsPrime_[l - 1]);
246 auto &tree2 = (l == 0 ? originsPrime_[0] : originsPrime_[l]);
247 getMatchingVector(tree1.mTree, tree2.mTree, originsMatchings_[l],
248 originsMatchingVectorT[l]);
249 getInverseMatchingVector(tree1.mTree, tree2.mTree, originsMatchings_[l],
250 invOriginsMatchingVectorT[l]);
251 }
252 std::vector<std::vector<ttk::ftm::idNode>> originsMatchingVector;
253 std::vector<std::vector<double>> originsPersPercent, originsPersDiff;
254 std::vector<double> originPersPercent, originPersDiff;
255 std::vector<int> originPersistenceOrder;
256 ttk::wae::computeTrackingInformation(
257 origins_, originsPrime_, originsMatchingVectorT, invOriginsMatchingVectorT,
258 isPersistenceDiagram_, originsMatchingVector, originsPersPercent,
259 originsPersDiff, originPersPercent, originPersDiff, originPersistenceOrder);
260
261 std::vector<std::vector<std::vector<ttk::ftm::idNode>>>
262 invDataMatchingVectorT(dataMatchings_.size());
263 for(unsigned int l = 0; l < invDataMatchingVectorT.size(); ++l) {
264 invDataMatchingVectorT[l].resize(dataMatchings_[l].size());
265 for(unsigned int i = 0; i < invDataMatchingVectorT[l].size(); ++i)
266 getInverseMatchingVector(origins_[l].mTree, recs_[i][l].mTree,
267 dataMatchings_[l][i],
268 invDataMatchingVectorT[l][i]);
269 }
270 std::vector<std::vector<ttk::ftm::idNode>> invReconstMatchingVectorT(
271 reconstMatchings_.size());
272 for(unsigned int i = 0; i < invReconstMatchingVectorT.size(); ++i) {
273 auto l = recs_[i].size() - 1;
274 getInverseMatchingVector(recs_[i][0].mTree, recs_[i][l].mTree,
276 invReconstMatchingVectorT[i]);
277 }
278
279 // ------------------------------------------
280 // --- Data
281 // ------------------------------------------
282 output_data->SetNumberOfBlocks(1);
285 data->SetNumberOfBlocks(1);
288 dataSeg->SetNumberOfBlocks(recs_.size());
289 bool outputSegmentation = !treesSegmentation.empty() and treesSegmentation[0];
290 for(unsigned int l = 0; l < 1; ++l) {
293 out_layer_i->SetNumberOfBlocks(recs_.size());
294 std::vector<ttk::ftm::MergeTree<float> *> trees(recs_.size());
295 for(unsigned int i = 0; i < recs_.size(); ++i)
296 trees[i] = &(recs_[i][l].mTree);
297
298 // Custom arrays
299 std::vector<std::vector<std::tuple<std::string, std::vector<int>>>>
300 customIntArrays(recs_.size());
301 std::vector<std::vector<std::tuple<std::string, std::vector<double>>>>
302 customDoubleArrays(recs_.size());
303 unsigned int lShift = 0;
304 ttk::wae::computeCustomArrays(
305 recs_, persCorrelationMatrix_, invDataMatchingVectorT,
306 invReconstMatchingVectorT, originsMatchingVector, originsMatchingVectorT,
307 originsPersPercent, originsPersDiff, originPersistenceOrder, l, lShift,
308 customIntArrays, customDoubleArrays);
309
310 // Create output
311 ttk::wae::makeManyOutput(trees, treesNodes, treesNodeCorr_, out_layer_i,
312 customIntArrays, customDoubleArrays,
315 if(outputSegmentation and l == 0) {
316 ttk::wae::makeManyOutput(
317 trees, treesNodes, treesNodeCorr_, treesSegmentation, dataSeg,
318 customIntArrays, customDoubleArrays, mixtureCoefficient_,
320 }
321 data->SetBlock(l, out_layer_i);
322 std::stringstream ss;
323 ss << (l == 0 ? "Input" : "Layer") << l;
324 data->GetMetaData(l)->Set(vtkCompositeDataSet::NAME(), ss.str());
325 }
326 output_data->SetBlock(0, data);
327 unsigned int num = 0;
328 output_data->GetMetaData(num)->Set(
329 vtkCompositeDataSet::NAME(), "layersTrees");
330 if(outputSegmentation)
331 output_data->SetBlock(1, dataSeg);
332 vtkNew<vtkFloatArray> lossArray{};
333 lossArray->SetName("Loss");
334 lossArray->InsertNextTuple1(bestLoss_);
335 output_data->GetFieldData()->AddArray(lossArray);
336
337 // ------------------------------------------
338 // --- Origins
339 // ------------------------------------------
340 output_origins->SetNumberOfBlocks(2);
341 // Origins
346 origins->SetNumberOfBlocks(noLayers_);
347 originsP->SetNumberOfBlocks(noLayers_);
348 std::vector<ttk::ftm::MergeTree<float> *> trees(noLayers_);
349 std::vector<std::vector<std::tuple<std::string, std::vector<int>>>>
350 customIntArrays(noLayers_);
351 std::vector<std::vector<std::tuple<std::string, std::vector<double>>>>
352 customDoubleArrays(noLayers_);
353 for(unsigned int l = 0; l < noLayers_; ++l) {
354 trees[l] = &(origins_[l].mTree);
355 if(l == 0) {
356 std::string name2{"OriginPersPercent"};
357 customDoubleArrays[l].emplace_back(
358 std::make_tuple(name2, originPersPercent));
359 std::string name3{"OriginPersDiff"};
360 customDoubleArrays[l].emplace_back(
361 std::make_tuple(name3, originPersDiff));
362 std::string nameOrder{"OriginPersOrder"};
363 customIntArrays[l].emplace_back(
364 std::make_tuple(nameOrder, originPersistenceOrder));
365 }
366 }
367 ttk::wae::makeManyOutput(trees, origins, customIntArrays, customDoubleArrays,
370
371 customIntArrays.clear();
372 customIntArrays.resize(noLayers_);
373 customDoubleArrays.clear();
374 customDoubleArrays.resize(noLayers_);
375 for(unsigned int l = 0; l < noLayers_; ++l) {
376 trees[l] = &(originsPrime_[l].mTree);
377 if(l < originsMatchingVector.size()) {
378 std::vector<int> customArrayMatching,
379 originPersOrder(trees[l]->tree.getNumberOfNodes(), -1);
380 for(unsigned int i = 0; i < originsMatchingVector[l].size(); ++i) {
381 customArrayMatching.emplace_back(originsMatchingVector[l][i]);
382 if(originsMatchingVector[l][i] < originPersistenceOrder.size())
383 originPersOrder[i]
384 = originPersistenceOrder[originsMatchingVector[l][i]];
385 }
386 std::string name{"OriginTrueNodeId"};
387 customIntArrays[l].emplace_back(
388 std::make_tuple(name, customArrayMatching));
389 std::string nameOrder{"OriginPersOrder"};
390 customIntArrays[l].emplace_back(
391 std::make_tuple(nameOrder, originPersOrder));
392 std::string name2{"OriginPersPercent"};
393 customDoubleArrays[l].emplace_back(
394 std::make_tuple(name2, originsPersPercent[l]));
395 std::string name3{"OriginPersDiff"};
396 customDoubleArrays[l].emplace_back(
397 std::make_tuple(name3, originsPersDiff[l]));
398 }
399 }
400 ttk::wae::makeManyOutput(trees, originsP, customIntArrays, customDoubleArrays,
403 output_origins->SetBlock(0, origins);
404 output_origins->SetBlock(1, originsP);
405 // for(unsigned int l = 0; l < 2; ++l) {
406 for(unsigned int l = 0; l < noLayers_; ++l) {
407 if(l >= 2)
408 break;
409 std::stringstream ss;
410 ss << (l == 0 ? "InputOrigin" : "LayerOrigin") << l;
411 auto originsMetaData = origins->GetMetaData(l);
412 if(originsMetaData)
413 originsMetaData->Set(vtkCompositeDataSet::NAME(), ss.str());
414 ss.str("");
415 ss << (l == 0 ? "InputOriginPrime" : "LayerOriginPrime") << l;
416 auto originsPMetaData = originsP->GetMetaData(l);
417 if(originsPMetaData)
418 originsPMetaData->Set(vtkCompositeDataSet::NAME(), ss.str());
419 }
420 num = 0;
421 output_origins->GetMetaData(num)->Set(
422 vtkCompositeDataSet::NAME(), "layersOrigins");
423 num = 1;
424 output_origins->GetMetaData(num)->Set(
425 vtkCompositeDataSet::NAME(), "layersOriginsPrime");
426
427 // ------------------------------------------
428 // --- Coefficients
429 // ------------------------------------------
430 output_coef->SetNumberOfBlocks(allAlphas_[0].size());
431 for(unsigned int l = 0; l < allAlphas_[0].size(); ++l) {
433 vtkNew<vtkIntArray> treeIDArray{};
434 treeIDArray->SetName("TreeID");
435 treeIDArray->SetNumberOfTuples(inputTrees.size());
436 for(unsigned int i = 0; i < inputTrees.size(); ++i)
437 treeIDArray->SetTuple1(i, i);
438 coef_table->AddColumn(treeIDArray);
439 auto noVec = allAlphas_[0][l].sizes()[0];
440 for(unsigned int v = 0; v < noVec; ++v) {
441 // Alphas
442 vtkNew<vtkFloatArray> tArray{};
443 std::string name = ttk::axa::getTableCoefficientName(noVec, v);
444 tArray->SetName(name.c_str());
445 tArray->SetNumberOfTuples(allAlphas_.size());
446 // Act Alphas
447 vtkNew<vtkFloatArray> actArray{};
448 std::string actName = "Act" + name;
449 actArray->SetName(actName.c_str());
450 actArray->SetNumberOfTuples(allAlphas_.size());
451 // Scaled Alphas
452 vtkNew<vtkFloatArray> tArrayNorm{};
453 std::string nameNorm = ttk::axa::getTableCoefficientNormName(noVec, v);
454 tArrayNorm->SetName(nameNorm.c_str());
455 tArrayNorm->SetNumberOfTuples(allAlphas_.size());
456 // Act Scaled Alphas
457 vtkNew<vtkFloatArray> actArrayNorm{};
458 std::string actNameNorm = "Act" + nameNorm;
459 actArrayNorm->SetName(actNameNorm.c_str());
460 actArrayNorm->SetNumberOfTuples(allAlphas_.size());
461 // Fill Arrays
462 for(unsigned int i = 0; i < allAlphas_.size(); ++i) {
463 tArray->SetTuple1(i, allAlphas_[i][l][v].item<float>());
464 actArray->SetTuple1(i, allActAlphas_[i][l][v].item<float>());
465 tArrayNorm->SetTuple1(i, allScaledAlphas_[i][l][v].item<float>());
466 actArrayNorm->SetTuple1(i, allActScaledAlphas_[i][l][v].item<float>());
467 }
468 coef_table->AddColumn(tArray);
469 coef_table->AddColumn(actArray);
470 coef_table->AddColumn(tArrayNorm);
471 coef_table->AddColumn(actArrayNorm);
472 }
473 if(!clusterAsgn_.empty()) {
474 vtkNew<vtkIntArray> clusterArray{};
475 clusterArray->SetName("ClusterAssignment");
476 clusterArray->SetNumberOfTuples(inputTrees.size());
477 for(unsigned int i = 0; i < clusterAsgn_.size(); ++i)
478 clusterArray->SetTuple1(i, clusterAsgn_[i]);
479 coef_table->AddColumn(clusterArray);
480 }
481 if(l == 0) {
482 vtkNew<vtkIntArray> treesNoNodesArray{};
483 treesNoNodesArray->SetNumberOfTuples(recs_.size());
484 treesNoNodesArray->SetName("treeNoNodes");
485 for(unsigned int i = 0; i < recs_.size(); ++i)
486 treesNoNodesArray->SetTuple1(
487 i, recs_[i][0].mTree.tree.getNumberOfNodes());
488 coef_table->AddColumn(treesNoNodesArray);
489 }
490 output_coef->SetBlock(l, coef_table);
491 std::stringstream ss;
492 ss << "Coef" << l;
493 output_coef->GetMetaData(l)->Set(vtkCompositeDataSet::NAME(), ss.str());
494 }
495
496 // Copy Field Data
497 // - aggregate input field data
498 for(unsigned int b = 0; b < inputTrees[0]->GetNumberOfBlocks(); ++b) {
499 vtkNew<vtkFieldData> fd{};
500 fd->CopyStructure(inputTrees[0]->GetBlock(b)->GetFieldData());
501 fd->SetNumberOfTuples(inputTrees.size());
502 for(size_t i = 0; i < inputTrees.size(); ++i) {
503 fd->SetTuple(i, 0, inputTrees[i]->GetBlock(b)->GetFieldData());
504 }
505
506 // - copy input field data to output row data
507 for(int i = 0; i < fd->GetNumberOfArrays(); ++i) {
508 auto array = fd->GetAbstractArray(i);
509 array->SetName(array->GetName());
510 vtkTable::SafeDownCast(output_coef->GetBlock(0))->AddColumn(array);
511 }
512 }
513
514 // Field Data Input Parameters
515 std::vector<std::string> paramNames;
516 getParamNames(paramNames);
517 for(auto paramName : paramNames) {
518 vtkNew<vtkDoubleArray> array{};
519 array->SetName(paramName.c_str());
520 array->InsertNextTuple1(getParamValueFromName(paramName));
521 output_coef->GetFieldData()->AddArray(array);
522 }
523 vtkNew<vtkDoubleArray> arrayActivate{};
524 arrayActivate->SetName("activate");
525 arrayActivate->InsertNextTuple1(activate_);
526 output_coef->GetFieldData()->AddArray(arrayActivate);
527 vtkNew<vtkDoubleArray> arrayActivateFunction{};
528 arrayActivateFunction->SetName("activationFunction");
529 arrayActivateFunction->InsertNextTuple1(activationFunction_);
530 output_coef->GetFieldData()->AddArray(arrayActivateFunction);
531
532 // ------------------------------------------
533 // --- Axes Vectors
534 // ------------------------------------------
535 std::vector<std::vector<std::vector<ttk::ftm::idNode>>> dataMatchingVectorT(
536 dataMatchings_.size());
537 for(unsigned int l = 0; l < dataMatchingVectorT.size(); ++l) {
538 dataMatchingVectorT[l].resize(dataMatchings_[l].size());
539 for(unsigned int i = 0; i < dataMatchingVectorT[l].size(); ++i) {
540 auto &origin = (l == 0 ? origins_[0] : originsPrime_[l - 1]);
541 getMatchingVector(origin.mTree, recs_[i][l].mTree, dataMatchings_[l][i],
542 dataMatchingVectorT[l][i]);
543 }
544 }
545 output_vectors->SetNumberOfBlocks(2);
548 vectors->SetNumberOfBlocks(vSTensor_.size());
551 vectorsPrime->SetNumberOfBlocks(vSTensor_.size());
552 for(unsigned int l = 0; l < vSTensor_.size(); ++l) {
554 vtkSmartPointer<vtkTable> vectorsPrimeTable
556 for(unsigned int v = 0; v < vSTensor_[l].sizes()[1]; ++v) {
557 // Vs
558 vtkNew<vtkFloatArray> vectorArray{};
559 std::string name
560 = ttk::axa::getTableVectorName(vSTensor_[l].sizes()[1], v, 0, 0, false);
561 vectorArray->SetName(name.c_str());
562 vectorArray->SetNumberOfTuples(vSTensor_[l].sizes()[0]);
563 for(unsigned int i = 0; i < vSTensor_[l].sizes()[0]; ++i)
564 vectorArray->SetTuple1(i, vSTensor_[l][i][v].item<float>());
565 vectorsTable->AddColumn(vectorArray);
566 // Vs Prime
567 vtkNew<vtkFloatArray> vectorPrimeArray{};
568 std::string name2
569 = ttk::axa::getTableVectorName(vSTensor_[l].sizes()[1], v, 0, 0, false);
570 vectorPrimeArray->SetName(name2.c_str());
571 vectorPrimeArray->SetNumberOfTuples(vSPrimeTensor_[l].sizes()[0]);
572 for(unsigned int i = 0; i < vSPrimeTensor_[l].sizes()[0]; ++i)
573 vectorPrimeArray->SetTuple1(i, vSPrimeTensor_[l][i][v].item<float>());
574 vectorsPrimeTable->AddColumn(vectorPrimeArray);
575 }
576 // Rev node corr
577 vtkNew<vtkUnsignedIntArray> revNodeCorrArray{};
578 revNodeCorrArray->SetName("revNodeCorr");
579 revNodeCorrArray->SetNumberOfTuples(vSTensor_[l].sizes()[0]);
580 std::vector<unsigned int> revNodeCorr;
581 getReverseTorchNodeCorr(origins_[l], revNodeCorr);
582 for(unsigned int i = 0; i < vSTensor_[l].sizes()[0]; ++i)
583 revNodeCorrArray->SetTuple1(i, revNodeCorr[i]);
584 vectorsTable->AddColumn(revNodeCorrArray);
585 // Rev node corr prime
586 vtkNew<vtkUnsignedIntArray> revNodeCorrPrimeArray{};
587 revNodeCorrPrimeArray->SetNumberOfTuples(vSPrimeTensor_[l].sizes()[0]);
588 revNodeCorrPrimeArray->SetName("revNodeCorr");
589 std::vector<unsigned int> revNodeCorrPrime;
590 getReverseTorchNodeCorr(originsPrime_[l], revNodeCorrPrime);
591 for(unsigned int i = 0; i < vSPrimeTensor_[l].sizes()[0]; ++i)
592 revNodeCorrPrimeArray->SetTuple1(i, revNodeCorrPrime[i]);
593 vectorsPrimeTable->AddColumn(revNodeCorrPrimeArray);
594 // Origins Matchings
595 auto addOriginMatchingArray
596 = [&](vtkSmartPointer<vtkTable> &vectorsTableT,
597 std::vector<ttk::ftm::idNode> &originMatchingVector) {
598 vtkNew<vtkIntArray> matchingArray{};
599 matchingArray->SetNumberOfTuples(originMatchingVector.size());
600 matchingArray->SetName("nextOriginMatching");
601 for(unsigned int i = 0; i < originMatchingVector.size(); ++i)
602 matchingArray->SetTuple1(i, (int)originMatchingVector[i]);
603 vectorsTableT->AddColumn(matchingArray);
604 };
605 if(l == 0)
606 addOriginMatchingArray(vectorsTable, originsMatchingVectorT[l]);
607 if(l < originsMatchingVectorT.size() - 1)
608 addOriginMatchingArray(vectorsPrimeTable, originsMatchingVectorT[l + 1]);
609 // Data Matchings
610 auto addDataMatchingArray
611 = [&](vtkSmartPointer<vtkTable> &vectorsTableT,
612 std::vector<std::vector<ttk::ftm::idNode>> &dataMatchingVector) {
613 for(unsigned int i = 0; i < dataMatchingVector.size(); ++i) {
614 vtkNew<vtkIntArray> matchingArray{};
615 matchingArray->SetNumberOfTuples(dataMatchingVector[i].size());
616 std::stringstream ss;
617 ss << "matching"
618 << ttk::axa::getTableTreeName(dataMatchingVector.size(), i);
619 matchingArray->SetName(ss.str().c_str());
620 for(unsigned int j = 0; j < dataMatchingVector[i].size(); ++j)
621 matchingArray->SetTuple1(j, (int)dataMatchingVector[i][j]);
622 vectorsTableT->AddColumn(matchingArray);
623 }
624 };
625 if(l == 0)
626 addDataMatchingArray(vectorsTable, dataMatchingVectorT[l]);
627 if(l < dataMatchingVectorT.size() - 1)
628 addDataMatchingArray(vectorsPrimeTable, dataMatchingVectorT[l + 1]);
629 // Reconst Matchings
630 if(l == vSTensor_.size() - 1) {
631 for(unsigned int i = 0; i < invReconstMatchingVectorT.size(); ++i) {
632 vtkNew<vtkIntArray> matchingArray{};
633 matchingArray->SetNumberOfTuples(invReconstMatchingVectorT[i].size());
634 std::stringstream ss;
635 ss << "reconstMatching"
636 << ttk::axa::getTableTreeName(invReconstMatchingVectorT.size(), i);
637 matchingArray->SetName(ss.str().c_str());
638 for(unsigned int j = 0; j < invReconstMatchingVectorT[i].size(); ++j)
639 matchingArray->SetTuple1(j, (int)invReconstMatchingVectorT[i][j]);
640 vectorsPrimeTable->AddColumn(matchingArray);
641 }
642 }
643
644 // Add new block
645 vectors->SetBlock(l, vectorsTable);
646 std::stringstream ss;
647 ss << "Vectors" << l;
648 vectors->GetMetaData(l)->Set(vtkCompositeDataSet::NAME(), ss.str());
649 vectorsPrime->SetBlock(l, vectorsPrimeTable);
650 ss.str("");
651 ss << "VectorsPrime" << l;
652 vectorsPrime->GetMetaData(l)->Set(vtkCompositeDataSet::NAME(), ss.str());
653 }
654 output_vectors->SetBlock(0, vectors);
655 output_vectors->SetBlock(1, vectorsPrime);
656 num = 0;
657 output_vectors->GetMetaData(num)->Set(vtkCompositeDataSet::NAME(), "Vectors");
658 num = 1;
659 output_vectors->GetMetaData(num)->Set(
660 vtkCompositeDataSet::NAME(), "VectorsPrime");
661
662 // return success
663 return 1;
664}
665#endif
#define TTK_FORCE_USE(x)
Force the compiler to use the function/method parameter.
Definition BaseClass.h:57
#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::MergeTreeAutoencoder module.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int run(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees2)
int FillInputPortInformation(int port, vtkInformation *info) override
int runCompute(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 debugLevel_
Definition Debug.h:379
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
void execute(std::vector< ftm::MergeTree< float > > &trees, std::vector< ftm::MergeTree< float > > &trees2)
std::vector< unsigned int > clusterAsgn_
std::vector< std::vector< std::tuple< ftm::idNode, ftm::idNode, double > > > originsMatchings_
std::vector< std::vector< double > > persCorrelationMatrix_
std::vector< std::vector< std::tuple< ftm::idNode, ftm::idNode, double > > > reconstMatchings_
std::vector< std::vector< std::vector< std::tuple< ftm::idNode, ftm::idNode, double > > > > dataMatchings_
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)
void getMatchingVector(ftm::MergeTree< dataType > &barycenter, ftm::MergeTree< dataType > &tree, std::vector< std::tuple< ftm::idNode, ftm::idNode, double > > &matchings, std::vector< ftm::idNode > &matchingVector)
void getParamNames(std::vector< std::string > &paramNames)
std::vector< std::vector< int > > treesNodeCorr_
double getParamValueFromName(std::string &paramName)
std::string getTableCoefficientNormName(int noAxes, int axeNum)
std::string getTableVectorName(int noAxes, int axeNum, int vId, int vComp, bool isSecondInput)
std::string getTableTreeName(int noTrees, int treeNum)
std::string getTableCoefficientName(int noAxes, int axeNum)
void loadBlocks(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, vtkMultiBlockDataSet *blocks)
vtkStandardNewMacro(ttkMergeTreeAutoencoder)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)