TTK
Loading...
Searching...
No Matches
ttkMergeTreePrincipalGeodesics.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 <vtkDoubleArray.h>
9#include <vtkObjectFactory.h>
10#include <vtkPointData.h>
11#include <vtkSmartPointer.h>
12#include <vtkTable.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
33 this->SetNumberOfInputPorts(2);
34 this->SetNumberOfOutputPorts(4);
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) {
49 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
50 info->Set(vtkAlgorithm::INPUT_IS_OPTIONAL(), 1);
51 } else
52 return 0;
53 return 1;
54}
55
72 int port, vtkInformation *info) {
73 if(port == 0) {
74 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
75 } else if(port == 1 or port == 2 or port == 3) {
76 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkTable");
77 } else {
78 return 0;
79 }
80 return 1;
81}
82
97 vtkInformation *ttkNotUsed(request),
98 vtkInformationVector **inputVector,
99 vtkInformationVector *outputVector) {
100 // ------------------------------------------------------------------------------------
101 // --- Get input object from input vector
102 // ------------------------------------------------------------------------------------
103 auto blocks = vtkMultiBlockDataSet::GetData(inputVector[0], 0);
104 auto blocks2 = vtkMultiBlockDataSet::GetData(inputVector[1], 0);
105
106 // ------------------------------------------------------------------------------------
107 // --- Load blocks
108 // ------------------------------------------------------------------------------------
109 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> inputTrees, inputTrees2;
110 ttk::ftm::loadBlocks(inputTrees, blocks);
111 ttk::ftm::loadBlocks(inputTrees2, blocks2);
112
113 // ------------------------------------------------------------------------------------
114 // If we have already computed once but the input has changed
115 if((treesNodes.size() != 0 and inputTrees[0]->GetBlock(0) != treesNodes[0])
116 or (treesNodes2.size() != inputTrees2.size()))
117 resetDataVisualization();
118
119 // Parameters
121 if(not normalizedWasserstein_) {
122 oldEpsilonTree1 = epsilonTree1_;
123 epsilonTree1_ = 100;
124 } else
125 epsilonTree1_ = oldEpsilonTree1;
127 printMsg("Computation with normalized Wasserstein.");
128 else
129 printMsg("Computation without normalized Wasserstein.");
130
131 return run<float>(outputVector, inputTrees, inputTrees2);
132}
133
134template <class dataType>
136 vtkInformationVector *outputVector,
137 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
138 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2) {
139 if(not isDataVisualizationFilled())
140 runCompute<dataType>(outputVector, inputTrees, inputTrees2);
141 runOutput<dataType>(outputVector, inputTrees, inputTrees2);
142 return 1;
143}
144
145template <class dataType>
147 vtkInformationVector *ttkNotUsed(outputVector),
148 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
149 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2) {
150 // ------------------------------------------------------------------------------------
151 // --- Construct trees
152 // ------------------------------------------------------------------------------------
153 std::vector<ttk::ftm::MergeTree<dataType>> intermediateMTrees,
154 intermediateMTrees2;
155
156 bool const useSadMaxPairs = (mixtureCoefficient_ == 0);
157 isPersistenceDiagram_ = ttk::ftm::constructTrees<dataType>(
158 inputTrees, intermediateMTrees, treesNodes, treesArcs, treesSegmentation,
159 useSadMaxPairs);
160 // If merge trees are provided in input and normalization is not asked
165 = (intermediateMTrees[0].tree.template isJoinTree<dataType>() ? 1 : 0);
166 }
168 or (mixtureCoefficient_ != 0 and mixtureCoefficient_ != 1)) {
169 auto &inputTrees2ToUse
170 = (not isPersistenceDiagram_ ? inputTrees2 : inputTrees);
171 ttk::ftm::constructTrees<dataType>(inputTrees2ToUse, intermediateMTrees2,
172 treesNodes2, treesArcs2,
173 treesSegmentation2, !useSadMaxPairs);
174 }
176
177 const int numInputs = intermediateMTrees.size();
178 const int numInputs2 = intermediateMTrees2.size();
179 setDataVisualization(numInputs, numInputs2);
180
181 // ------------------------------------------------------------------------------------
182 // --- Call base
183 // ------------------------------------------------------------------------------------
184 execute<dataType>(intermediateMTrees, intermediateMTrees2);
185
186 ttk::ftm::mergeTreesTemplateToDouble<dataType>(
187 intermediateMTrees, intermediateDTrees);
188
189 return 1;
190}
191
192template <class dataType>
195 int blockId,
196 vtkMultiBlockDataSet *output_barycenter) {
197 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNodeBary
199 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputArcBary
201
202 ttkMergeTreeVisualization visuMakerBary;
203 visuMakerBary.setShiftMode(2); // Line
204 visuMakerBary.setVtkOutputNode(vtkOutputNodeBary);
206 visuMakerBary.setVtkOutputArc(vtkOutputArcBary);
207 else {
208 visuMakerBary.setVtkOutputArc(vtkOutputNodeBary);
210 visuMakerBary.setIsPDSadMax(blockId);
211 else
212 visuMakerBary.setIsPDSadMax(mixtureCoefficient_ == 0);
213 visuMakerBary.setPlanarLayout(true);
214 }
215 visuMakerBary.setDebugLevel(this->debugLevel_);
216 visuMakerBary.setOutputTreeNodeId(true);
219 visuMakerBary.makeTreesOutput<dataType>(&(barycenter.tree));
220
221 // Construct multiblock
222 if(not isPersistenceDiagram_) {
223 vtkMultiBlockDataSet::SafeDownCast(output_barycenter->GetBlock(0))
224 ->SetBlock(blockId, vtkOutputNodeBary);
225 vtkMultiBlockDataSet::SafeDownCast(output_barycenter->GetBlock(1))
226 ->SetBlock(blockId, vtkOutputArcBary);
227 } else
228 output_barycenter->SetBlock(blockId, vtkOutputNodeBary);
229}
230
231template <class dataType>
233 vtkInformationVector *outputVector,
234 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
235 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees2) {
236 // ------------------------------------------------------------------------------------
237 // --- Create output
238 // ------------------------------------------------------------------------------------
239 auto output_barycenter = vtkMultiBlockDataSet::GetData(outputVector, 0);
240 auto output_coef = vtkTable::GetData(outputVector, 1);
241 auto output_geodVector = vtkTable::GetData(outputVector, 2);
242 auto output_corrMatrix = vtkTable::GetData(outputVector, 3);
243
244 // ------------------------------------------
245 // --- Barycenter
246 // ------------------------------------------
247 auto noBlockBary
248 = 1
250 or (mixtureCoefficient_ != 0 and mixtureCoefficient_ != 1));
251 output_barycenter->SetNumberOfBlocks(noBlockBary);
252 if(not isPersistenceDiagram_) {
253 vtkSmartPointer<vtkMultiBlockDataSet> const vtkBlockNodes
257 int const noBlock = (inputTrees2.size() == 0 ? 1 : 2);
258 vtkBlockNodes->SetNumberOfBlocks(noBlock);
259 vtkBlockArcs->SetNumberOfBlocks(noBlock);
260 output_barycenter->SetBlock(0, vtkBlockNodes);
261 output_barycenter->SetBlock(1, vtkBlockArcs);
262 }
263
265 mergeTreeDoubleToTemplate(barycenter_, barycenterT);
266 makeBarycenterOutput<dataType>(barycenterT, 0, output_barycenter);
267
268 if(inputTrees2.size() != 0
270 and mixtureCoefficient_ != 1)) {
272 mergeTreeDoubleToTemplate(barycenterInput2_, barycenter2T);
273 makeBarycenterOutput<dataType>(barycenter2T, 1, output_barycenter);
274 }
275
276 // ------------------------------------------
277 // --- Coefficients
278 // ------------------------------------------
279 std::vector<std::vector<std::vector<double>> *> allTs{&allTs_, &allScaledTs_};
280
281 vtkNew<vtkIntArray> treeIDArray{};
282 treeIDArray->SetName("TreeID");
283 treeIDArray->SetNumberOfTuples(inputTrees.size());
284 for(unsigned int i = 0; i < inputTrees.size(); ++i) {
285 treeIDArray->SetTuple1(i, i);
286 }
287 output_coef->AddColumn(treeIDArray);
288
289 for(unsigned int t = 0; t < 2; ++t) {
290 for(unsigned int i = 0; i < (*allTs[t]).size(); ++i) {
291 vtkNew<vtkDoubleArray> tArray{};
292 auto size = (*allTs[t]).size();
293 std::string const name = (t == 0 ? getTableCoefficientName(size, i)
294 : getTableCoefficientNormName(size, i));
295 tArray->SetName(name.c_str());
296 tArray->SetNumberOfTuples(inputTrees.size());
297 for(unsigned int j = 0; j < (*allTs[t])[i].size(); ++j)
298 tArray->SetTuple1(j, (*allTs[t])[i][j]);
299 output_coef->AddColumn(tArray);
300 }
301 }
302
303 // Copy Field Data
304 // - aggregate input field data
305 for(unsigned int b = 0; b < inputTrees[0]->GetNumberOfBlocks(); ++b) {
306 vtkNew<vtkFieldData> fd{};
307 fd->CopyStructure(inputTrees[0]->GetBlock(b)->GetFieldData());
308 fd->SetNumberOfTuples(inputTrees.size());
309 for(size_t i = 0; i < inputTrees.size(); ++i) {
310 fd->SetTuple(i, 0, inputTrees[i]->GetBlock(b)->GetFieldData());
311 }
312
313 // - copy input field data to output row data
314 for(int i = 0; i < fd->GetNumberOfArrays(); ++i) {
315 auto array = fd->GetAbstractArray(i);
316 array->SetName(array->GetName());
317 output_coef->AddColumn(array);
318 }
319 }
320
321 // Field Data Input Parameters
322 std::vector<std::string> paramNames;
323 getParamNames(paramNames);
324 for(auto paramName : paramNames) {
325 vtkNew<vtkDoubleArray> array{};
326 array->SetName(paramName.c_str());
327 array->InsertNextTuple1(getParamValueFromName(paramName));
328 output_coef->GetFieldData()->AddArray(array);
329 }
330
331 // ------------------------------------------
332 // --- Geodesics Vectors
333 // ------------------------------------------
334 std::vector<std::vector<std::vector<std::vector<double>>> *> vectors{
335 &vS_, &v2s_};
336 unsigned int maxSize = vS_[0].size();
337 if(inputTrees2.size() != 0
339 and mixtureCoefficient_ != 1)) {
340 vectors = {&vS_, &v2s_, &trees2Vs_, &trees2V2s_};
341 maxSize = std::max(vS_[0].size(), trees2Vs_[0].size());
342 }
343 for(unsigned int v = 0; v < vectors.size(); ++v) {
344 auto vector = vectors[v];
345 for(unsigned int i = 0; i < (*vector).size(); ++i)
346 for(unsigned int k = 0; k < 2; ++k) {
347 vtkNew<vtkDoubleArray> vectorArray{};
348 bool const isSecondInput = (v >= 2);
349 std::string const name
350 = getTableVectorName((*vector).size(), i, v % 2, k, isSecondInput);
351 vectorArray->SetName(name.c_str());
352 vectorArray->SetNumberOfTuples(maxSize);
353 for(unsigned int j = 0; j < (*vector)[i].size(); ++j)
354 vectorArray->SetTuple1(j, (*vector)[i][j][k]);
355 for(unsigned int j = (*vector)[i].size(); j < maxSize; ++j)
356 vectorArray->SetTuple1(j, std::nan(""));
357 output_geodVector->AddColumn(vectorArray);
358 }
359 }
360
361 // ------------------------------------------
362 // --- Correlation Matrix
363 // ------------------------------------------
364 // TODO manage second input
365 // Correlation coefficients
366 unsigned int const noCols = branchesCorrelationMatrix_[0].size();
367 unsigned int const noRows = branchesCorrelationMatrix_.size();
368 for(unsigned int j = 0; j < noCols; ++j) {
369 vtkNew<vtkDoubleArray> corrArray{}, persArray{};
370 std::string const name = getTableCorrelationName(noCols, j);
371 corrArray->SetName(name.c_str());
372 corrArray->SetNumberOfTuples(noRows);
373 std::string const name2 = getTableCorrelationPersName(noCols, j);
374 persArray->SetName(name2.c_str());
375 persArray->SetNumberOfTuples(noRows);
376 for(unsigned int i = 0; i < noRows; ++i) {
377 corrArray->SetTuple1(i, branchesCorrelationMatrix_[i][j]);
378 persArray->SetTuple1(i, persCorrelationMatrix_[i][j]);
379 }
380 output_corrMatrix->AddColumn(corrArray);
381 output_corrMatrix->AddColumn(persArray);
382 }
383
384 // Tree matching
385 std::vector<std::vector<ttk::ftm::idNode>> matchingMatrix;
387 barycenter_, intermediateDTrees, baryMatchings_, matchingMatrix);
389 for(unsigned int j = 0; j < inputTrees.size(); ++j)
391 ->getOrigin()][j]
392 = intermediateDTrees[j]
393 .tree.getNode(intermediateDTrees[j].tree.getRoot())
394 ->getOrigin();
395 for(unsigned int j = 0; j < inputTrees.size(); ++j) {
396 vtkNew<vtkIntArray> corrArray{};
397 std::string const name = getTableCorrelationTreeName(inputTrees.size(), j);
398 corrArray->SetName(name.c_str());
399 corrArray->SetNumberOfTuples(noRows);
400 for(unsigned int i = 0; i < noRows; ++i)
401 corrArray->SetTuple1(i, (int)matchingMatrix[i][j]);
402 output_corrMatrix->AddColumn(corrArray);
403 }
404
405 // Barycenter node id
406 vtkNew<vtkIntArray> baryNodeId{};
407 baryNodeId->SetName("BaryNodeId");
408 baryNodeId->SetNumberOfTuples(noRows);
409 for(unsigned int i = 0; i < noRows; ++i)
410 baryNodeId->SetTuple1(i, i);
411 output_corrMatrix->AddColumn(baryNodeId);
412
413 // Barycenter node isLeaf
414 vtkNew<vtkIntArray> baryNodeIsLeaf{};
415 baryNodeIsLeaf->SetName("BaryNodeIsLeaf");
416 baryNodeIsLeaf->SetNumberOfTuples(noRows);
417 for(unsigned int i = 0; i < noRows; ++i)
418 baryNodeIsLeaf->SetTuple1(i, barycenter_.tree.isLeaf(i));
419 output_corrMatrix->AddColumn(baryNodeIsLeaf);
420
421 // Barycenter node persistence
422 vtkNew<vtkDoubleArray> baryNodePers{};
423 baryNodePers->SetName("BaryNodePers");
424 baryNodePers->SetNumberOfTuples(noRows);
425 for(unsigned int i = 0; i < noRows; ++i) {
426 auto nodePers = barycenter_.tree.template getNodePersistence<double>(i);
427 baryNodePers->SetTuple1(i, nodePers);
428 }
429 output_corrMatrix->AddColumn(baryNodePers);
430
431 // Trees persistence order
432 std::vector<std::vector<int>> treesOrder(intermediateDTrees.size());
433 for(unsigned int i = 0; i < intermediateDTrees.size(); ++i) {
434 treesOrder[i].resize(intermediateDTrees[i].tree.getNumberOfNodes());
435 std::vector<std::tuple<ttk::ftm::idNode, ttk::ftm::idNode, double>> pairs;
436 intermediateDTrees[i].tree.template getPersistencePairsFromTree<double>(
437 pairs, false);
438 for(unsigned int j = 0; j < pairs.size(); ++j) {
439 int const index = pairs.size() - 1 - j;
440 treesOrder[i][std::get<0>(pairs[j])] = index;
441 treesOrder[i][std::get<1>(pairs[j])] = index;
442 }
443 }
444
445 vtkNew<vtkDoubleArray> treesPersOrder{};
446 treesPersOrder->SetName("TreesPersOrder");
447 treesPersOrder->SetNumberOfTuples(noRows);
448 for(unsigned int i = 0; i < noRows; ++i) {
449 int minIndex = std::numeric_limits<int>::max();
450 for(unsigned int j = 0; j < intermediateDTrees.size(); ++j) {
451 if(int(matchingMatrix[i][j]) != -1)
452 minIndex = std::min(minIndex, treesOrder[j][matchingMatrix[i][j]]);
453 }
454 treesPersOrder->SetTuple1(i, minIndex);
455 }
456 output_corrMatrix->AddColumn(treesPersOrder);
457
458 return 1;
459}
#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::MergeTreePrincipalGeodesics module.
void makeBarycenterOutput(ttk::ftm::MergeTree< dataType > &barycenter, int blockId, vtkMultiBlockDataSet *output_barycenter)
int runCompute(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees2)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int runOutput(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees2)
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
void makeTreesOutput(FTMTree_MT *tree1)
void setConvertedToDiagram(bool converted)
void setVtkOutputArc(vtkUnstructuredGrid *vtkArc)
void setVtkOutputNode(vtkUnstructuredGrid *vtkNode)
int debugLevel_
Definition Debug.h:379
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)
std::string getTableCorrelationName(int noGeodesics, int geodesicNum)
std::string getTableCoefficientName(int noGeodesics, int geodesicNum)
void getMatchingMatrix(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)
void getParamNames(std::vector< std::string > &paramNames)
double getParamValueFromName(std::string &paramName)
std::vector< std::vector< std::vector< double > > > v2s_
std::vector< std::vector< std::vector< double > > > trees2V2s_
std::vector< std::vector< std::vector< double > > > vS_
std::vector< std::vector< double > > allTs_
std::vector< std::vector< double > > branchesCorrelationMatrix_
std::vector< std::vector< double > > persCorrelationMatrix_
std::vector< std::vector< std::tuple< ftm::idNode, ftm::idNode, double > > > baryMatchings_
std::vector< std::vector< std::vector< double > > > trees2Vs_
std::vector< std::vector< double > > allScaledTs_
bool isLeaf(idNode nodeId)
Node * getNode(idNode nodeId)
Definition FTMTree_MT.h:393
SimplexId getOrigin() const
Definition FTMNode.h:64
void loadBlocks(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, vtkMultiBlockDataSet *blocks)
ftm::FTMTree_MT tree
Definition FTMTree_MT.h:901
vtkStandardNewMacro(ttkMergeTreePrincipalGeodesics)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)