TTK
Loading...
Searching...
No Matches
ttkMergeTreeUtils.h
Go to the documentation of this file.
1
7
8#pragma once
9
10#include <FTMTree.h>
11#include <FTMTreeUtils.h>
12
13#include <ttkUtils.h>
14
15#include <vtkCellData.h>
16#include <vtkDoubleArray.h>
17#include <vtkIntArray.h>
18#include <vtkMultiBlockDataSet.h>
19#include <vtkPointData.h>
20#include <vtkUnstructuredGrid.h>
21
22namespace ttk {
23 namespace ftm {
24
25 template <class dataType>
26 MergeTree<dataType> makeTree(vtkUnstructuredGrid *treeNodes,
27 vtkUnstructuredGrid *treeArcs) {
28 auto treeNodeIdArray = treeNodes->GetPointData()->GetArray("TreeNodeId");
29
30 // Init Scalars
31 auto scalars = std::make_shared<Scalars>();
32 vtkSmartPointer<vtkDataArray> const nodesScalar
33 = treeNodes->GetPointData()->GetArray("Scalar"); // 1: Scalar
34 scalars->size = nodesScalar->GetNumberOfTuples();
35 if(treeNodeIdArray)
36 scalars->size = std::max(
37 (ttk::SimplexId)treeNodeIdArray->GetRange()[1] + 1, scalars->size);
38 auto scalarsValues
39 = std::make_shared<std::vector<dataType>>(scalars->size);
40 for(int i = 0; i < nodesScalar->GetNumberOfTuples(); ++i) {
41 int const index = (treeNodeIdArray ? treeNodeIdArray->GetTuple1(i) : i);
42 (*scalarsValues)[index] = nodesScalar->GetTuple1(i);
43 }
44 scalars->values = (void *)(scalarsValues->data());
45
46 // Init Tree
47 auto params = std::make_shared<Params>();
48 params->treeType = Join_Split;
49 MergeTree<dataType> mergeTree(scalars, scalarsValues, params);
50
51 // Add Nodes
52 for(vtkIdType i = 0; i < scalars->size; ++i) {
53 mergeTree.tree.makeNode(i);
54 }
55
56 // Add Arcs
58 = treeArcs->GetCellData()->GetArray("upNodeId"); // 1: upNodeId
60 = treeArcs->GetCellData()->GetArray("downNodeId"); // 2: downNodeId
61 vtkIdType const arcsNumTuples = arcsUp->GetNumberOfTuples();
62 vtkSmartPointer<vtkDataArray> const dummyArcArray
63 = treeArcs->GetCellData()->GetArray("isDummyArc");
64 std::set<std::tuple<double, double>> added_arcs; // Avoid duplicates
65 for(vtkIdType i = 0; i < arcsNumTuples; ++i) {
66 if(dummyArcArray != nullptr and dummyArcArray->GetTuple1(i) == 1)
67 continue;
68 double downId = arcsDown->GetTuple1(i);
69 double upId = arcsUp->GetTuple1(i);
70 if(treeNodeIdArray) {
71 downId = treeNodeIdArray->GetTuple1(downId);
72 upId = treeNodeIdArray->GetTuple1(upId);
73 }
74 auto it = added_arcs.find(std::make_tuple(downId, upId));
75 if(it == added_arcs.end()) { // arc not added yet
76 mergeTree.tree.makeSuperArc(downId, upId); // (down, Up)
77 added_arcs.insert(std::make_tuple(downId, upId));
78 }
79 }
80
81 // Manage inconsistent arcs
83
84 // Remove self link
85 removeSelfLink(&(mergeTree.tree));
86
87 return mergeTree;
88 }
89
90 // Returns a branch decomposition tree
91 template <class dataType>
93 makeBDTreeFromPDGrid(vtkUnstructuredGrid *persistenceDiagram,
94 bool useSadMaxPairs = true) {
95 auto birthArray
96 = persistenceDiagram->GetCellData()->GetArray(PersistenceBirthName);
97 auto persArray
98 = persistenceDiagram->GetCellData()->GetArray(PersistenceName);
99 auto pairTypeArray
100 = persistenceDiagram->GetCellData()->GetArray(PersistencePairTypeName);
101 auto criticalTypeArray = persistenceDiagram->GetPointData()->GetArray(
103
104 auto treeNodeIdArray
105 = persistenceDiagram->GetPointData()->GetArray("TreeNodeId");
106 auto treeNodeId2Array
107 = persistenceDiagram->GetPointData()->GetArray("TreeNodeIdOrigin");
108 bool const gotNodeArrays = (treeNodeIdArray and treeNodeId2Array);
109
110 auto noPairs = birthArray->GetNumberOfTuples();
111 int noNodes = noPairs * 2;
112 if(gotNodeArrays) {
113 for(vtkIdType i = 0; i < treeNodeIdArray->GetNumberOfTuples(); ++i) {
114 int const val = std::max(treeNodeIdArray->GetTuple1(i),
115 treeNodeId2Array->GetTuple1(i))
116 + 1;
117 noNodes = std::max(noNodes, val);
118 }
119 }
120 std::vector<dataType> scalarsVector(noNodes);
121
122 // Init Tree
123 MergeTree<dataType> mergeTree
124 = ttk::ftm::createEmptyMergeTree<dataType>(scalarsVector.size());
125
126 // Init critical type enum values
127 auto locMin = static_cast<int>(CriticalType::Local_minimum);
128 auto locMax = static_cast<int>(CriticalType::Local_maximum);
129
130 // Get Min-Max pair index
131 int minMaxPairIndex = -1;
132 for(vtkIdType i = 0; i < noPairs; ++i) {
133 vtkIdType npts;
134 vtkIdType const *pts;
135 persistenceDiagram->GetCellPoints(i, npts, pts);
136 auto ct1 = criticalTypeArray->GetTuple1(pts[0]);
137 auto ct2 = criticalTypeArray->GetTuple1(pts[1]);
138 if((ct1 == locMin and ct2 == locMax)
139 or (ct1 == locMax and ct2 == locMin))
140 minMaxPairIndex = i;
141 }
142
143 // Init Tree Structure
144 for(vtkIdType i = 0; i < noNodes; ++i)
145 mergeTree.tree.makeNode(i);
146 for(vtkIdType i = 0; i < noPairs; ++i) {
147 auto pairType = pairTypeArray->GetTuple1(i);
148 vtkIdType npts;
149 vtkIdType const *pts;
150 persistenceDiagram->GetCellPoints(i, npts, pts);
151 auto ct1 = criticalTypeArray->GetTuple1(pts[0]);
152 auto ct2 = criticalTypeArray->GetTuple1(pts[1]);
153 if((pairType == -1
154 or (useSadMaxPairs and ct1 != locMax and ct2 != locMax)
155 or (not useSadMaxPairs and ct1 != locMin and ct2 != locMin))
156 and i != minMaxPairIndex)
157 continue;
158 int const index1
159 = (gotNodeArrays ? treeNodeId2Array->GetTuple1(pts[0]) : i * 2);
160 int const index2
161 = (gotNodeArrays ? treeNodeIdArray->GetTuple1(pts[0]) : i * 2 + 1);
162 mergeTree.tree.getNode(index1)->setOrigin(index2);
163 mergeTree.tree.getNode(index2)->setOrigin(index1);
164 scalarsVector[index1] = birthArray->GetTuple1(i);
165 scalarsVector[index2]
166 = birthArray->GetTuple1(i) + persArray->GetTuple1(i);
167
168 if(i != minMaxPairIndex) {
169 auto up = (gotNodeArrays ? treeNodeIdArray->GetTuple1(minMaxPairIndex)
170 : minMaxPairIndex * 2 + 1);
171 // mergeTree.tree.makeSuperArc(index1, up);
172 mergeTree.tree.makeSuperArc(index2, up);
173 }
174 }
175
176 // Init scalars
177 ttk::ftm::setTreeScalars<dataType>(mergeTree, scalarsVector);
178
179 return mergeTree;
180 }
181
182 inline void
184 vtkMultiBlockDataSet *blocks) {
185 if(blocks != nullptr) {
186 if(blocks->GetBlock(0)->IsA("vtkMultiBlockDataSet"))
187 inputTrees.resize(
188 vtkMultiBlockDataSet::SafeDownCast(blocks->GetBlock(0))
189 ->GetNumberOfBlocks());
190 else if(blocks->GetBlock(0)->IsA("vtkUnstructuredGrid"))
191 inputTrees.resize(blocks->GetNumberOfBlocks());
192 for(size_t i = 0; i < inputTrees.size(); ++i) {
193 if(blocks->GetBlock(0)->IsA("vtkMultiBlockDataSet")) {
196 vtkBlock->SetNumberOfBlocks(blocks->GetNumberOfBlocks());
197 for(unsigned int j = 0; j < blocks->GetNumberOfBlocks(); ++j)
198 vtkBlock->SetBlock(
199 j, vtkMultiBlockDataSet::SafeDownCast(blocks->GetBlock(j))
200 ->GetBlock(i));
201 inputTrees[i] = vtkBlock;
202 } else if(blocks->GetBlock(0)->IsA("vtkUnstructuredGrid")) {
205 vtkBlock->SetNumberOfBlocks(1);
206 vtkBlock->SetBlock(
207 0, vtkUnstructuredGrid::SafeDownCast(blocks->GetBlock(i)));
208 inputTrees[i] = vtkBlock;
209 }
210 }
211 }
212 }
213
214 template <class dataType>
216 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
217 std::vector<MergeTree<dataType>> &intermediateTrees,
218 std::vector<vtkUnstructuredGrid *> &treesNodes,
219 std::vector<vtkUnstructuredGrid *> &treesArcs,
220 std::vector<vtkDataSet *> &treesSegmentation,
221 std::vector<bool> useSadMaxPairs) {
222 bool isPersistenceDiagram = false;
223 const int numInputs = inputTrees.size();
224 intermediateTrees.resize(numInputs);
225 treesNodes.resize(numInputs);
226 treesArcs.resize(numInputs);
227 treesSegmentation.resize(numInputs);
228 for(int i = 0; i < numInputs; i++) {
229 if(inputTrees[i]->GetNumberOfBlocks() >= 2) {
230 treesNodes[i]
231 = vtkUnstructuredGrid::SafeDownCast(inputTrees[i]->GetBlock(0));
232 treesArcs[i]
233 = vtkUnstructuredGrid::SafeDownCast(inputTrees[i]->GetBlock(1));
234 if(inputTrees[i]->GetNumberOfBlocks() > 2)
235 treesSegmentation[i]
236 = vtkDataSet::SafeDownCast(inputTrees[i]->GetBlock(2));
237 intermediateTrees[i]
238 = makeTree<dataType>(treesNodes[i], treesArcs[i]);
239 } else {
240 treesNodes[i]
241 = vtkUnstructuredGrid::SafeDownCast(inputTrees[i]->GetBlock(0));
242 vtkUnstructuredGrid *persistenceDiagram
243 = vtkUnstructuredGrid::SafeDownCast(inputTrees[i]->GetBlock(0));
244 intermediateTrees[i] = makeBDTreeFromPDGrid<dataType>(
245 persistenceDiagram, useSadMaxPairs[i]);
246 isPersistenceDiagram = true;
247 }
248 }
249 return isPersistenceDiagram;
250 }
251
252 template <class dataType>
254 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
255 std::vector<MergeTree<dataType>> &intermediateTrees,
256 std::vector<vtkUnstructuredGrid *> &treesNodes,
257 std::vector<vtkUnstructuredGrid *> &treesArcs,
258 std::vector<vtkDataSet *> &treesSegmentation,
259 bool useSadMaxPairs = true) {
260 std::vector<bool> const useSadMaxPairsVec(
261 inputTrees.size(), useSadMaxPairs);
262 return constructTrees(inputTrees, intermediateTrees, treesNodes,
263 treesArcs, treesSegmentation, useSadMaxPairsVec);
264 }
265
266 template <class dataType>
268 std::vector<vtkSmartPointer<vtkMultiBlockDataSet>> &inputTrees,
269 std::vector<MergeTree<dataType>> &intermediateTrees,
270 bool useSadMaxPairs = true) {
271 std::vector<vtkUnstructuredGrid *> treesNodes;
272 std::vector<vtkUnstructuredGrid *> treesArcs;
273 std::vector<vtkDataSet *> treesSegmentation;
274 return constructTrees(inputTrees, intermediateTrees, treesNodes,
275 treesArcs, treesSegmentation, useSadMaxPairs);
276 }
277 } // namespace ftm
278} // namespace ttk
TTK container representing a node of the MergeTree.
idNode makeNode(SimplexId vertexId, SimplexId linked=nullVertex)
idSuperArc makeSuperArc(idNode downNodeId, idNode upNodeId)
Node * getNode(idNode nodeId)
Definition FTMTree_MT.h:393
void setOrigin(SimplexId linked)
Definition FTMNode.h:72
bool constructTrees(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, std::vector< MergeTree< dataType > > &intermediateTrees, std::vector< vtkUnstructuredGrid * > &treesNodes, std::vector< vtkUnstructuredGrid * > &treesArcs, std::vector< vtkDataSet * > &treesSegmentation, std::vector< bool > useSadMaxPairs)
MergeTree< dataType > makeBDTreeFromPDGrid(vtkUnstructuredGrid *persistenceDiagram, bool useSadMaxPairs=true)
void removeSelfLink(FTMTree_MT *tree)
void loadBlocks(std::vector< vtkSmartPointer< vtkMultiBlockDataSet > > &inputTrees, vtkMultiBlockDataSet *blocks)
void manageInconsistentArcsMultiParent(FTMTree_MT *tree)
MergeTree< dataType > makeTree(vtkUnstructuredGrid *treeNodes, vtkUnstructuredGrid *treeArcs)
The Topology ToolKit.
const char PersistenceName[]
Definition DataTypes.h:72
const char PersistencePairTypeName[]
Definition DataTypes.h:73
const char PersistenceCriticalTypeName[]
Definition DataTypes.h:67
const char PersistenceBirthName[]
Definition DataTypes.h:68
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
ftm::FTMTree_MT tree
Definition FTMTree_MT.h:903