TTK
Loading...
Searching...
No Matches
ttkMergeTreeTemporalReductionDecoding.cpp
Go to the documentation of this file.
2
3#include <ttkMergeTreeUtils.h>
5
6#include <vtkInformation.h>
7
8#include <vtkDataArray.h>
9#include <vtkDataSet.h>
10#include <vtkObjectFactory.h>
11#include <vtkPointData.h>
12#include <vtkTable.h>
13
14#include <ttkUtils.h>
15
16using namespace ttk;
17using namespace ttk::ftm;
18
19// A VTK macro that enables the instantiation of this class via ::New()
20// You do not have to modify this
22
36 this->SetNumberOfInputPorts(2);
37 this->SetNumberOfOutputPorts(2);
38}
39
41 = default;
42
51 int port, vtkInformation *info) {
52 if(port == 0) {
53 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkMultiBlockDataSet");
54 } else if(port == 1) {
55 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTable");
56 } else
57 return 0;
58 return 1;
59}
60
77 int port, vtkInformation *info) {
78 if(port == 0 or port == 1) {
79 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
80 return 1;
81 }
82 return 0;
83}
84
99 vtkInformation *ttkNotUsed(request),
100 vtkInformationVector **inputVector,
101 vtkInformationVector *outputVector) {
102 // ------------------------------------------------------------------------------------
103 // --- Get input object from input vector
104 // ------------------------------------------------------------------------------------
105 printMsg("Get input object from input vector", debug::Priority::VERBOSE);
106 auto blocks = vtkMultiBlockDataSet::GetData(inputVector[0], 0);
107 auto table = vtkTable::GetData(inputVector[1], 0);
108
109 // ------------------------------------------------------------------------------------
110 // --- Load blocks
111 // ------------------------------------------------------------------------------------
112 printMsg("Load Blocks", debug::Priority::VERBOSE);
113 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> inputTrees;
114 loadBlocks(inputTrees, blocks);
115 printMsg("Load Blocks done.", debug::Priority::VERBOSE);
116 auto assignmentSolverArray
117 = blocks->GetFieldData()->GetArray("AssignmentSolver");
118 if(assignmentSolverArray)
119 assignmentSolverID_ = assignmentSolverArray->GetTuple1(0);
120
121 // If we have already computed once but the input has changed
122 if(treesNodes.size() != 0 and inputTrees[0]->GetBlock(0) != treesNodes[0])
123 resetDataVisualization();
124
125 std::vector<std::tuple<double, int, int, int, int>> coefs;
126 int const totalInputs = inputTrees.size() + table->GetNumberOfRows();
127 std::vector<bool> interpolatedTrees(totalInputs, false);
128 for(int i = 0; i < table->GetNumberOfRows(); ++i) {
129 double const alpha
130 = vtkDataArray::SafeDownCast(table->GetColumnByName("Alpha"))
131 ->GetTuple1(i);
132 int const index1R
133 = vtkDataArray::SafeDownCast(table->GetColumnByName("Index1_R"))
134 ->GetTuple1(i);
135 int const index2R
136 = vtkDataArray::SafeDownCast(table->GetColumnByName("Index2_R"))
137 ->GetTuple1(i);
138 int const index1
139 = vtkDataArray::SafeDownCast(table->GetColumnByName("Index1"))
140 ->GetTuple1(i);
141 int const index2
142 = vtkDataArray::SafeDownCast(table->GetColumnByName("Index2"))
143 ->GetTuple1(i);
144
145 coefs.emplace_back(alpha, index1R, index2R, index1, index2);
146 for(int j = index1 + 1; j < index2; ++j)
147 interpolatedTrees[j] = true;
148 }
149
150 return run<float>(outputVector, inputTrees, coefs, interpolatedTrees);
151}
152
153template <class dataType>
155 vtkInformationVector *outputVector,
156 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
157 std::vector<std::tuple<double, int, int, int, int>> &coefs,
158 std::vector<bool> &interpolatedTrees) {
159 if(not isDataVisualizationFilled())
160 runCompute<dataType>(inputTrees, coefs);
161 runOutput<dataType>(outputVector, inputTrees, coefs, interpolatedTrees);
162 return 1;
163}
164
165template <class dataType>
167 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
168 std::vector<std::tuple<double, int, int, int, int>> &coefs) {
169 // ------------------------------------------------------------------------------------
170 // --- Construct trees
171 // ------------------------------------------------------------------------------------
172 printMsg("Construct trees", debug::Priority::VERBOSE);
173
174 const int numInputs = inputTrees.size();
175
176 setDataVisualization(numInputs);
177
178 std::vector<MergeTree<dataType>> intermediateMTrees(numInputs);
179 constructTrees<dataType>(
180 inputTrees, intermediateMTrees, treesNodes, treesArcs, treesSegmentation);
181
182 // ------------------------------------------------------------------------------------
183 // --- Call base
184 // ------------------------------------------------------------------------------------
185 printMsg("Call base", debug::Priority::VERBOSE);
186
187 std::vector<MergeTree<dataType>> allMT_T;
188
189 execute<dataType>(intermediateMTrees, coefs, allMT_T, allMatching);
190 treesNodeCorrMesh = getTreesNodeCorr();
191
192 for(size_t i = 0; i < allMT_T.size(); ++i) {
194 mergeTreeTemplateToDouble<dataType>(allMT_T[i], tree);
195 intermediateSTrees.emplace_back(tree);
196 }
197
198 return 1;
199}
200
201template <class dataType>
203 vtkInformationVector *outputVector,
204 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
205 std::vector<std::tuple<double, int, int, int, int>> &coefs,
206 std::vector<bool> &interpolatedTrees) {
207 bool const OutputSegmentation = false; // TODO
208 // ------------------------------------------------------------------------------------
209 // --- Create output
210 // ------------------------------------------------------------------------------------
211 auto output_sequence = vtkMultiBlockDataSet::GetData(outputVector, 0);
212 auto output_matchings = vtkMultiBlockDataSet::GetData(outputVector, 1);
213
214 // ------------------------------------------
215 // --- Trees
216 // -----------------------------------------
217 std::vector<vtkUnstructuredGrid *> treesNodesT;
218 std::vector<vtkUnstructuredGrid *> treesArcsT;
219 std::vector<vtkDataSet *> treesSegmentationT;
220 std::vector<std::vector<int>> treesNodeCorrMeshT;
221 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> inputTreesT;
222 int index = 0;
223 size_t cpt = 0;
224 while(cpt < coefs.size()) {
225 while(cpt < coefs.size() and std::get<2>(coefs[cpt]) <= index) {
226 treesNodesT.emplace_back(nullptr);
227 treesArcsT.emplace_back(nullptr);
228 treesSegmentationT.emplace_back(nullptr);
229 treesNodeCorrMeshT.emplace_back();
230 inputTreesT.emplace_back(nullptr);
231 ++cpt;
232 }
233 treesNodesT.emplace_back(treesNodes[index]);
234 treesArcsT.emplace_back(treesArcs[index]);
235 treesSegmentationT.emplace_back(treesSegmentation[index]);
236 treesNodeCorrMeshT.emplace_back(treesNodeCorrMesh[index]);
237 inputTreesT.emplace_back(inputTrees[index]);
238 ++index;
239 }
240 treesNodes.swap(treesNodesT);
241 treesArcs.swap(treesArcsT);
242 treesSegmentation.swap(treesSegmentationT);
243 treesNodeCorrMesh.swap(treesNodeCorrMeshT);
244 inputTrees.swap(inputTreesT);
245
246 std::vector<MergeTree<dataType>> intermediateMTrees;
247 mergeTreesDoubleToTemplate<dataType>(intermediateSTrees, intermediateMTrees);
248 std::vector<FTMTree_MT *> trees;
249 mergeTreeToFTMTree<dataType>(intermediateMTrees, trees);
250
251 output_sequence->SetNumberOfBlocks((OutputSegmentation ? 3 : 2));
252 vtkSmartPointer<vtkMultiBlockDataSet> const vtkBlockNodes
254 vtkBlockNodes->SetNumberOfBlocks(intermediateMTrees.size());
255 output_sequence->SetBlock(0, vtkBlockNodes);
258 vtkBlockArcs->SetNumberOfBlocks(intermediateMTrees.size());
259 output_sequence->SetBlock(1, vtkBlockArcs);
260 if(OutputSegmentation) {
263 vtkBlockSegs->SetNumberOfBlocks(intermediateMTrees.size());
264 output_sequence->SetBlock(2, vtkBlockSegs);
265 }
266
267 double prevXMax = 0;
268 std::vector<std::vector<SimplexId>> nodeCorr(intermediateMTrees.size());
269 int cptInterpolatedTree = 0;
270 for(size_t i = 0; i < intermediateMTrees.size(); ++i) {
271 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNode1
273 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputArc1
275 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputSegmentation1
277
279 visuMaker.setPlanarLayout(PlanarLayout);
281 BranchDecompositionPlanarLayout);
282 visuMaker.setRescaleTreesIndividually(RescaleTreesIndividually);
283 visuMaker.setOutputSegmentation(OutputSegmentation);
284 visuMaker.setDimensionSpacing(DimensionSpacing);
285 visuMaker.setDimensionToShift(DimensionToShift);
286 visuMaker.setImportantPairs(ImportantPairs);
287 visuMaker.setMaximumImportantPairs(MaximumImportantPairs);
288 visuMaker.setMinimumImportantPairs(MinimumImportantPairs);
289 visuMaker.setImportantPairsSpacing(ImportantPairsSpacing);
290 visuMaker.setNonImportantPairsSpacing(NonImportantPairsSpacing);
291 visuMaker.setNonImportantPairsProximity(NonImportantPairsProximity);
292 visuMaker.setExcludeImportantPairsHigher(ExcludeImportantPairsHigher);
293 visuMaker.setExcludeImportantPairsLower(ExcludeImportantPairsLower);
294 // visuMaker.setShiftMode(3); // Double Line
295 visuMaker.setShiftMode(2); // Line
296 visuMaker.setVtkOutputNode(vtkOutputNode1);
297 visuMaker.setVtkOutputArc(vtkOutputArc1);
298 visuMaker.setVtkOutputSegmentation(vtkOutputSegmentation1);
299 visuMaker.setTreesNodes(treesNodes);
300 visuMaker.copyPointData(treesNodes[i], treesNodeCorrMesh[i]);
301 visuMaker.setTreesNodeCorrMesh(treesNodeCorrMesh);
302 visuMaker.setTreesSegmentation(treesSegmentation);
303 visuMaker.setInterpolatedTrees(interpolatedTrees);
304 visuMaker.setPrintTreeId(i);
305 visuMaker.setPrintClusterId(0);
306 visuMaker.setPrevXMaxOffset(prevXMax);
307 visuMaker.setDebugLevel(this->debugLevel_);
308 visuMaker.makeTreesOutput<dataType>(trees);
309 prevXMax = visuMaker.getPrevXMax();
310 nodeCorr[i] = visuMaker.getNodeCorr()[i];
311
312 // Field data
313 if(treesNodes[i])
314 vtkOutputNode1->GetFieldData()->ShallowCopy(
315 treesNodes[i]->GetFieldData());
316 if(treesArcs[i])
317 vtkOutputArc1->GetFieldData()->ShallowCopy(treesArcs[i]->GetFieldData());
318 if(treesSegmentation[i] and OutputSegmentation)
319 vtkOutputSegmentation1->GetFieldData()->ShallowCopy(
320 treesSegmentation[i]->GetFieldData());
321
322 // Field data
323 if(interpolatedTrees[i]) {
324 // Distance previous key frame
325 vtkNew<vtkDoubleArray> vtkDistancePreviousKey{};
326 vtkDistancePreviousKey->SetName("DistancePreviousKeyFrame");
327 vtkDistancePreviousKey->SetNumberOfTuples(1);
328 vtkDistancePreviousKey->SetTuple1(
329 0, distancesToKeyFrames_[cptInterpolatedTree * 2]);
330 vtkOutputNode1->GetFieldData()->AddArray(vtkDistancePreviousKey);
331 // Distance next key frame
332 vtkNew<vtkDoubleArray> vtkDistanceNextKey{};
333 vtkDistanceNextKey->SetName("DistanceNextKeyFrame");
334 vtkDistanceNextKey->SetNumberOfTuples(1);
335 vtkDistanceNextKey->SetTuple1(
336 0, distancesToKeyFrames_[cptInterpolatedTree * 2 + 1]);
337 vtkOutputNode1->GetFieldData()->AddArray(vtkDistanceNextKey);
338 // Index previous key frame
339 vtkNew<vtkIntArray> vtkIndexPreviousKey{};
340 vtkIndexPreviousKey->SetName("IndexPreviousKeyFrame");
341 vtkIndexPreviousKey->SetNumberOfTuples(1);
342 vtkIndexPreviousKey->SetTuple1(
343 0, std::get<3>(coefs[cptInterpolatedTree]));
344 vtkOutputNode1->GetFieldData()->AddArray(vtkIndexPreviousKey);
345 // Index previous key frame
346 vtkNew<vtkIntArray> vtkIndexNextKey{};
347 vtkIndexNextKey->SetName("IndexNextKeyFrame");
348 vtkIndexNextKey->SetNumberOfTuples(1);
349 vtkIndexNextKey->SetTuple1(0, std::get<4>(coefs[cptInterpolatedTree]));
350 vtkOutputNode1->GetFieldData()->AddArray(vtkIndexNextKey);
351 // Increment interpolated tree counter
352 ++cptInterpolatedTree;
353 }
354
355 // Construct multiblock
356 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(0))
357 ->SetBlock(i, vtkOutputNode1);
358 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(1))
359 ->SetBlock(i, vtkOutputArc1);
360 if(OutputSegmentation)
361 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(2))
362 ->SetBlock(i, vtkOutputSegmentation1);
363 }
364
365 // -----------------------------------------
366 // --- Matching
367 // -----------------------------------------
368 output_matchings->SetNumberOfBlocks(intermediateMTrees.size() - 1);
369 for(size_t i = 0; i < intermediateMTrees.size() - 1; ++i) {
370 // Declare vtk objects
371 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputMatching
373 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNode1
374 = vtkUnstructuredGrid::SafeDownCast(
375 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(0))
376 ->GetBlock(i));
377 vtkSmartPointer<vtkUnstructuredGrid> const vtkOutputNode2
378 = vtkUnstructuredGrid::SafeDownCast(
379 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(0))
380 ->GetBlock(i + 1));
381 std::vector<std::vector<SimplexId>> nodeCorrTemp;
382 nodeCorrTemp.emplace_back(nodeCorr[i]);
383 nodeCorrTemp.emplace_back(nodeCorr[i + 1]);
384
385 // Fill vtk objects
386 ttkMergeTreeVisualization visuMakerMatching;
387 visuMakerMatching.setVtkOutputMatching(vtkOutputMatching);
388 visuMakerMatching.setOutputMatching(allMatching[i]);
389 visuMakerMatching.setVtkOutputNode1(vtkOutputNode2);
390 visuMakerMatching.setVtkOutputNode2(vtkOutputNode1);
391 visuMakerMatching.setNodeCorr1(nodeCorrTemp);
392 visuMakerMatching.setDebugLevel(this->debugLevel_);
393
394 visuMakerMatching.makeMatchingOutput<dataType>(trees[i], trees[i + 1]);
395
396 // Field data
397 vtkNew<vtkDoubleArray> vtkDistance{};
398 vtkDistance->SetName("Distance");
399 vtkDistance->SetNumberOfTuples(1);
400 vtkDistance->SetTuple1(0, finalDistances_[i]);
401 vtkOutputMatching->GetFieldData()->AddArray(vtkDistance);
402
403 // Construct multiblock
404 output_matchings->SetBlock(i, vtkOutputMatching);
405 }
406
407 // -----------------------------------------
408 treesNodes.swap(treesNodesT);
409 treesArcs.swap(treesArcsT);
410 treesSegmentation.swap(treesSegmentationT);
411 treesNodeCorrMesh.swap(treesNodeCorrMeshT);
412 inputTrees.swap(inputTreesT);
413
414 return 1;
415}
#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::MergeTreeTemporalReductionDecoding module.
int runOutput(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< std::tuple< double, int, int, int, int > > &coefs, std::vector< bool > &interpolatedTrees)
int runCompute(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< std::tuple< double, int, int, int, int > > &coefs)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int run(vtkInformationVector *outputVector, std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< std::tuple< double, int, int, int, int > > &coefs, std::vector< bool > &interpolatedTrees)
void setTreesSegmentation(std::vector< vtkDataSet * > &segmentation)
void setInterpolatedTrees(std::vector< bool > &isInterpolatedTrees)
void copyPointData(vtkUnstructuredGrid *treeNodes, std::vector< int > &nodeCorrT)
void setVtkOutputNode1(vtkUnstructuredGrid *vtkNode1)
void setNodeCorr1(std::vector< std::vector< SimplexId > > &nodeCorrT)
void setVtkOutputSegmentation(vtkDataSet *vtkSegmentation)
void setTreesNodeCorrMesh(std::vector< std::vector< int > > &nodeCorrMesh)
void setVtkOutputMatching(vtkUnstructuredGrid *vtkMatching)
std::vector< std::vector< SimplexId > > getNodeCorr()
void makeTreesOutput(FTMTree_MT *tree1)
void setVtkOutputNode2(vtkUnstructuredGrid *vtkNode2)
void makeMatchingOutput(FTMTree_MT *tree1, FTMTree_MT *tree2)
void setVtkOutputArc(vtkUnstructuredGrid *vtkArc)
void setVtkOutputNode(vtkUnstructuredGrid *vtkNode)
void setOutputMatching(std::vector< std::tuple< idNode, idNode, double > > &matching)
void setTreesNodes(std::vector< vtkUnstructuredGrid * > &nodes)
int debugLevel_
Definition Debug.h:379
virtual int setDebugLevel(const int &debugLevel)
Definition Debug.cpp:147
std::vector< std::vector< int > > getTreesNodeCorr()
void setExcludeImportantPairsLower(std::string &d)
void setExcludeImportantPairsHigher(std::string &d)
void loadBlocks(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, vtkMultiBlockDataSet *blocks)
The Topology ToolKit.
vtkStandardNewMacro(ttkMergeTreeTemporalReductionDecoding)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)