TTK
Loading...
Searching...
No Matches
ttkMergeTreeTemporalReductionDecoding.cpp
Go to the documentation of this file.
2
3#include <ttkFTMTreeUtils.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 totalInputs = inputTrees.size() + table->GetNumberOfRows();
127 std::vector<bool> interpolatedTrees(totalInputs, false);
128 for(int i = 0; i < table->GetNumberOfRows(); ++i) {
129 double alpha = vtkDataArray::SafeDownCast(table->GetColumnByName("Alpha"))
130 ->GetTuple1(i);
131 int index1R = vtkDataArray::SafeDownCast(table->GetColumnByName("Index1_R"))
132 ->GetTuple1(i);
133 int index2R = vtkDataArray::SafeDownCast(table->GetColumnByName("Index2_R"))
134 ->GetTuple1(i);
135 int index1 = vtkDataArray::SafeDownCast(table->GetColumnByName("Index1"))
136 ->GetTuple1(i);
137 int index2 = vtkDataArray::SafeDownCast(table->GetColumnByName("Index2"))
138 ->GetTuple1(i);
139
140 coefs.emplace_back(alpha, index1R, index2R, index1, index2);
141 for(int j = index1 + 1; j < index2; ++j)
142 interpolatedTrees[j] = true;
143 }
144
145 return run<float>(outputVector, inputTrees, coefs, interpolatedTrees);
146}
147
148template <class dataType>
150 vtkInformationVector *outputVector,
151 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
152 std::vector<std::tuple<double, int, int, int, int>> &coefs,
153 std::vector<bool> &interpolatedTrees) {
154 if(not isDataVisualizationFilled())
155 runCompute<dataType>(inputTrees, coefs);
156 runOutput<dataType>(outputVector, inputTrees, coefs, interpolatedTrees);
157 return 1;
158}
159
160template <class dataType>
162 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
163 std::vector<std::tuple<double, int, int, int, int>> &coefs) {
164 // ------------------------------------------------------------------------------------
165 // --- Construct trees
166 // ------------------------------------------------------------------------------------
167 printMsg("Construct trees", debug::Priority::VERBOSE);
168
169 const int numInputs = inputTrees.size();
170
171 setDataVisualization(numInputs);
172
173 std::vector<MergeTree<dataType>> intermediateMTrees(numInputs);
174 constructTrees<dataType>(
175 inputTrees, intermediateMTrees, treesNodes, treesArcs, treesSegmentation);
176
177 // ------------------------------------------------------------------------------------
178 // --- Call base
179 // ------------------------------------------------------------------------------------
180 printMsg("Call base", debug::Priority::VERBOSE);
181
182 std::vector<MergeTree<dataType>> allMT_T;
183
184 execute<dataType>(intermediateMTrees, coefs, allMT_T, allMatching);
185 treesNodeCorrMesh = getTreesNodeCorr();
186
187 for(size_t i = 0; i < allMT_T.size(); ++i) {
189 mergeTreeTemplateToDouble<dataType>(allMT_T[i], tree);
190 intermediateSTrees.emplace_back(tree);
191 }
192
193 return 1;
194}
195
196template <class dataType>
198 vtkInformationVector *outputVector,
199 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
200 std::vector<std::tuple<double, int, int, int, int>> &coefs,
201 std::vector<bool> &interpolatedTrees) {
202 bool OutputSegmentation = false; // TODO
203 // ------------------------------------------------------------------------------------
204 // --- Create output
205 // ------------------------------------------------------------------------------------
206 auto output_sequence = vtkMultiBlockDataSet::GetData(outputVector, 0);
207 auto output_matchings = vtkMultiBlockDataSet::GetData(outputVector, 1);
208
209 // ------------------------------------------
210 // --- Trees
211 // -----------------------------------------
212 std::vector<vtkUnstructuredGrid *> treesNodesT;
213 std::vector<vtkUnstructuredGrid *> treesArcsT;
214 std::vector<vtkDataSet *> treesSegmentationT;
215 std::vector<std::vector<int>> treesNodeCorrMeshT;
216 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> inputTreesT;
217 int index = 0;
218 size_t cpt = 0;
219 while(cpt < coefs.size()) {
220 while(cpt < coefs.size() and std::get<2>(coefs[cpt]) <= index) {
221 treesNodesT.emplace_back(nullptr);
222 treesArcsT.emplace_back(nullptr);
223 treesSegmentationT.emplace_back(nullptr);
224 treesNodeCorrMeshT.emplace_back();
225 inputTreesT.emplace_back(nullptr);
226 ++cpt;
227 }
228 treesNodesT.emplace_back(treesNodes[index]);
229 treesArcsT.emplace_back(treesArcs[index]);
230 treesSegmentationT.emplace_back(treesSegmentation[index]);
231 treesNodeCorrMeshT.emplace_back(treesNodeCorrMesh[index]);
232 inputTreesT.emplace_back(inputTrees[index]);
233 ++index;
234 }
235 treesNodes.swap(treesNodesT);
236 treesArcs.swap(treesArcsT);
237 treesSegmentation.swap(treesSegmentationT);
238 treesNodeCorrMesh.swap(treesNodeCorrMeshT);
239 inputTrees.swap(inputTreesT);
240
241 std::vector<MergeTree<dataType>> intermediateMTrees;
242 mergeTreesDoubleToTemplate<dataType>(intermediateSTrees, intermediateMTrees);
243 std::vector<FTMTree_MT *> trees;
244 mergeTreeToFTMTree<dataType>(intermediateMTrees, trees);
245
246 output_sequence->SetNumberOfBlocks((OutputSegmentation ? 3 : 2));
249 vtkBlockNodes->SetNumberOfBlocks(intermediateMTrees.size());
250 output_sequence->SetBlock(0, vtkBlockNodes);
253 vtkBlockArcs->SetNumberOfBlocks(intermediateMTrees.size());
254 output_sequence->SetBlock(1, vtkBlockArcs);
255 if(OutputSegmentation) {
258 vtkBlockSegs->SetNumberOfBlocks(intermediateMTrees.size());
259 output_sequence->SetBlock(2, vtkBlockSegs);
260 }
261
262 double prevXMax = 0;
263 std::vector<std::vector<SimplexId>> nodeCorr(intermediateMTrees.size());
264 int cptInterpolatedTree = 0;
265 for(size_t i = 0; i < intermediateMTrees.size(); ++i) {
270 vtkSmartPointer<vtkUnstructuredGrid> vtkOutputSegmentation1
272
274 visuMaker.setPlanarLayout(PlanarLayout);
276 BranchDecompositionPlanarLayout);
277 visuMaker.setRescaleTreesIndividually(RescaleTreesIndividually);
278 visuMaker.setOutputSegmentation(OutputSegmentation);
279 visuMaker.setDimensionSpacing(DimensionSpacing);
280 visuMaker.setDimensionToShift(DimensionToShift);
281 visuMaker.setImportantPairs(ImportantPairs);
282 visuMaker.setMaximumImportantPairs(MaximumImportantPairs);
283 visuMaker.setMinimumImportantPairs(MinimumImportantPairs);
284 visuMaker.setImportantPairsSpacing(ImportantPairsSpacing);
285 visuMaker.setNonImportantPairsSpacing(NonImportantPairsSpacing);
286 visuMaker.setNonImportantPairsProximity(NonImportantPairsProximity);
287 visuMaker.setExcludeImportantPairsHigher(ExcludeImportantPairsHigher);
288 visuMaker.setExcludeImportantPairsLower(ExcludeImportantPairsLower);
289 // visuMaker.setShiftMode(3); // Double Line
290 visuMaker.setShiftMode(2); // Line
291 visuMaker.setVtkOutputNode(vtkOutputNode1);
292 visuMaker.setVtkOutputArc(vtkOutputArc1);
293 visuMaker.setVtkOutputSegmentation(vtkOutputSegmentation1);
294 visuMaker.setTreesNodes(treesNodes);
295 visuMaker.copyPointData(treesNodes[i], treesNodeCorrMesh[i]);
296 visuMaker.setTreesNodeCorrMesh(treesNodeCorrMesh);
297 visuMaker.setTreesSegmentation(treesSegmentation);
298 visuMaker.setInterpolatedTrees(interpolatedTrees);
299 visuMaker.setPrintTreeId(i);
300 visuMaker.setPrintClusterId(0);
301 visuMaker.setPrevXMaxOffset(prevXMax);
302 visuMaker.setDebugLevel(this->debugLevel_);
303 visuMaker.makeTreesOutput<dataType>(trees);
304 prevXMax = visuMaker.getPrevXMax();
305 nodeCorr[i] = visuMaker.getNodeCorr()[i];
306
307 // Field data
308 if(treesNodes[i])
309 vtkOutputNode1->GetFieldData()->ShallowCopy(
310 treesNodes[i]->GetFieldData());
311 if(treesArcs[i])
312 vtkOutputArc1->GetFieldData()->ShallowCopy(treesArcs[i]->GetFieldData());
313 if(treesSegmentation[i] and OutputSegmentation)
314 vtkOutputSegmentation1->GetFieldData()->ShallowCopy(
315 treesSegmentation[i]->GetFieldData());
316
317 // Field data
318 if(interpolatedTrees[i]) {
319 // Distance previous key frame
320 vtkNew<vtkDoubleArray> vtkDistancePreviousKey{};
321 vtkDistancePreviousKey->SetName("DistancePreviousKeyFrame");
322 vtkDistancePreviousKey->SetNumberOfTuples(1);
323 vtkDistancePreviousKey->SetTuple1(
324 0, distancesToKeyFrames_[cptInterpolatedTree * 2]);
325 vtkOutputNode1->GetFieldData()->AddArray(vtkDistancePreviousKey);
326 // Distance next key frame
327 vtkNew<vtkDoubleArray> vtkDistanceNextKey{};
328 vtkDistanceNextKey->SetName("DistanceNextKeyFrame");
329 vtkDistanceNextKey->SetNumberOfTuples(1);
330 vtkDistanceNextKey->SetTuple1(
331 0, distancesToKeyFrames_[cptInterpolatedTree * 2 + 1]);
332 vtkOutputNode1->GetFieldData()->AddArray(vtkDistanceNextKey);
333 // Index previous key frame
334 vtkNew<vtkIntArray> vtkIndexPreviousKey{};
335 vtkIndexPreviousKey->SetName("IndexPreviousKeyFrame");
336 vtkIndexPreviousKey->SetNumberOfTuples(1);
337 vtkIndexPreviousKey->SetTuple1(
338 0, std::get<3>(coefs[cptInterpolatedTree]));
339 vtkOutputNode1->GetFieldData()->AddArray(vtkIndexPreviousKey);
340 // Index previous key frame
341 vtkNew<vtkIntArray> vtkIndexNextKey{};
342 vtkIndexNextKey->SetName("IndexNextKeyFrame");
343 vtkIndexNextKey->SetNumberOfTuples(1);
344 vtkIndexNextKey->SetTuple1(0, std::get<4>(coefs[cptInterpolatedTree]));
345 vtkOutputNode1->GetFieldData()->AddArray(vtkIndexNextKey);
346 // Increment interpolated tree counter
347 ++cptInterpolatedTree;
348 }
349
350 // Construct multiblock
351 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(0))
352 ->SetBlock(i, vtkOutputNode1);
353 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(1))
354 ->SetBlock(i, vtkOutputArc1);
355 if(OutputSegmentation)
356 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(2))
357 ->SetBlock(i, vtkOutputSegmentation1);
358 }
359
360 // -----------------------------------------
361 // --- Matching
362 // -----------------------------------------
363 output_matchings->SetNumberOfBlocks(intermediateMTrees.size() - 1);
364 for(size_t i = 0; i < intermediateMTrees.size() - 1; ++i) {
365 // Declare vtk objects
369 = vtkUnstructuredGrid::SafeDownCast(
370 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(0))
371 ->GetBlock(i));
373 = vtkUnstructuredGrid::SafeDownCast(
374 vtkMultiBlockDataSet::SafeDownCast(output_sequence->GetBlock(0))
375 ->GetBlock(i + 1));
376 std::vector<std::vector<SimplexId>> nodeCorrTemp;
377 nodeCorrTemp.emplace_back(nodeCorr[i]);
378 nodeCorrTemp.emplace_back(nodeCorr[i + 1]);
379
380 // Fill vtk objects
381 ttkMergeTreeVisualization visuMakerMatching;
382 visuMakerMatching.setVtkOutputMatching(vtkOutputMatching);
383 visuMakerMatching.setOutputMatching(allMatching[i]);
384 visuMakerMatching.setVtkOutputNode1(vtkOutputNode2);
385 visuMakerMatching.setVtkOutputNode2(vtkOutputNode1);
386 visuMakerMatching.setNodeCorr1(nodeCorrTemp);
387 visuMakerMatching.setDebugLevel(this->debugLevel_);
388
389 visuMakerMatching.makeMatchingOutput<dataType>(trees[i], trees[i + 1]);
390
391 // Field data
392 vtkNew<vtkDoubleArray> vtkDistance{};
393 vtkDistance->SetName("Distance");
394 vtkDistance->SetNumberOfTuples(1);
395 vtkDistance->SetTuple1(0, finalDistances_[i]);
396 vtkOutputMatching->GetFieldData()->AddArray(vtkDistance);
397
398 // Construct multiblock
399 output_matchings->SetBlock(i, vtkOutputMatching);
400 }
401
402 // -----------------------------------------
403 treesNodes.swap(treesNodesT);
404 treesArcs.swap(treesArcsT);
405 treesSegmentation.swap(treesSegmentationT);
406 treesNodeCorrMesh.swap(treesNodeCorrMeshT);
407 inputTrees.swap(inputTreesT);
408
409 return 1;
410}
#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 setMinimumImportantPairs(int minPairs)
void setMaximumImportantPairs(int maxPairs)
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
int printMsg(const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
Definition: Debug.h:118
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)