TTK
Loading...
Searching...
No Matches
ttkMergeTreeBase.h
Go to the documentation of this file.
1
70
71#pragma once
72
73// ttk code includes
75
76class ttkMergeTreeBase : virtual public ttk::Debug {
77
78public:
80 ~ttkMergeTreeBase() override = default;
81
82 int getSkeletonNodes(vtkUnstructuredGrid *outputSkeletonNodes);
83
85 const int cc,
86 vtkPoints *points,
87 vtkUnstructuredGrid *skeletonArcs,
88 ttk::ftm::ArcData &arcData);
89
91 const int cc,
92 vtkPoints *points,
93 vtkUnstructuredGrid *skeletonArcs,
94 ttk::ftm::ArcData &arcData);
95
97 const int cc,
98 vtkPoints *points,
99 vtkUnstructuredGrid *skeletonArcs,
100 ttk::ftm::ArcData &arcData);
101
102 int getSkeletonArcs(vtkUnstructuredGrid *outputSkeletonArcs);
103
104 int getSegmentation(vtkDataSet *outputSegmentation);
105
107 return params_.treeType;
108 }
109
110 void identify(vtkDataSet *ds) const;
111
112#ifdef TTK_ENABLE_FTM_TREE_STATS_TIME
113 void printCSVStats();
114 void printCSVTree(const ttk::ftm::FTMTree_MT *const tree) const;
115#endif
116
117protected:
120 int nbCC_;
121 std::vector<vtkSmartPointer<vtkDataSet>> connected_components_;
122 std::vector<ttk::Triangulation *> triangulation_;
123 std::vector<ttk::ftm::LocalFTM> ftmTree_;
124 std::vector<vtkDataArray *> inputScalars_;
125 std::vector<std::vector<ttk::SimplexId>> offsets_;
126
127 // TODO: future optimization: flag critical vertices based on if they appear
128 // in the merge tree, count the number of critical points and allocate the
129 // arrays accordingly
130 template <class triangulationType>
131 int getMergeTree(vtkUnstructuredGrid *outputSkeletonArcs,
132 std::vector<ttk::ExTreeM::Branch> &mergeTree,
133 vtkDataArray *inputScalars,
134 const triangulationType *triangulation) {
135 vtkNew<vtkUnstructuredGrid> skeletonArcs{};
136 ttk::SimplexId pointIds[2];
137 ttk::SimplexId pointOrders[2];
138 vtkNew<vtkPoints> points{};
139 vtkNew<vtkLongLongArray> data{};
140 data->SetNumberOfComponents(1);
141 data->SetName("Order");
142 vtkNew<vtkIdTypeArray> gIdArray{};
143 gIdArray->SetNumberOfComponents(1);
144 gIdArray->SetName("GlobalPointIds");
145 vtkNew<vtkIdTypeArray> downId{};
146 downId->SetNumberOfComponents(1);
147 downId->SetName("downNodeId");
148 vtkNew<vtkIdTypeArray> upId{};
149 upId->SetNumberOfComponents(1);
150 upId->SetName("upNodeId");
151 float point[3];
152 vtkSmartPointer<vtkDataArray> const scalarArray
153 = vtkSmartPointer<vtkDataArray>::Take(inputScalars->NewInstance());
154 scalarArray->SetNumberOfComponents(1);
155 scalarArray->SetName("Scalar");
156 std::map<ttk::SimplexId, ttk::SimplexId> addedPoints;
157 ttk::SimplexId currentId = 0;
158 for(auto const &b : mergeTree) {
159 auto &vertices = b.vertices;
160 for(size_t p = 0; p < vertices.size() - 1; p++) {
161 pointIds[0] = vertices[p].second;
162 pointIds[1] = vertices[p + 1].second;
163 pointOrders[0] = vertices[p].first;
164 pointOrders[1] = vertices[p + 1].first;
165 // add each point only once to the vtkPoints
166 // addedPoints.insert(x).second inserts x and is true if x was not in
167 // addedPoints beforehand
168 if(addedPoints.insert({pointIds[0], currentId}).second) {
169 triangulation->getVertexPoint(
170 pointIds[0], point[0], point[1], point[2]);
171 points->InsertNextPoint(point);
172 data->InsertNextTuple1(pointOrders[0]);
173 gIdArray->InsertNextTuple1(pointIds[0]);
174 scalarArray->InsertNextTuple1(inputScalars->GetTuple1(pointIds[0]));
175 currentId++;
176 }
177 if(addedPoints.insert({pointIds[1], currentId}).second) {
178 triangulation->getVertexPoint(
179 pointIds[1], point[0], point[1], point[2]);
180 points->InsertNextPoint(point);
181 data->InsertNextTuple1(pointOrders[1]);
182 gIdArray->InsertNextTuple1(pointIds[1]);
183 scalarArray->InsertNextTuple1(inputScalars->GetTuple1(pointIds[1]));
184 currentId++;
185 }
186 downId->InsertNextTuple1(pointIds[0]);
187 upId->InsertNextTuple1(pointIds[1]);
188 vtkIdType pointIdsASVTKIDTYPE[2];
189 pointIdsASVTKIDTYPE[0] = addedPoints.at(pointIds[0]);
190 pointIdsASVTKIDTYPE[1] = addedPoints.at(pointIds[1]);
191 skeletonArcs->InsertNextCell(VTK_LINE, 2, pointIdsASVTKIDTYPE);
192 }
193 }
194 skeletonArcs->SetPoints(points);
195 outputSkeletonArcs->ShallowCopy(skeletonArcs);
196 outputSkeletonArcs->GetPointData()->AddArray(data);
197 outputSkeletonArcs->GetPointData()->AddArray(gIdArray);
198 outputSkeletonArcs->GetPointData()->AddArray(scalarArray);
199 outputSkeletonArcs->GetCellData()->AddArray(upId);
200 outputSkeletonArcs->GetCellData()->AddArray(downId);
201
202 return 1;
203 }
204
205 template <class triangulationType>
207 vtkUnstructuredGrid *outputSkeletonNodes,
208 std::map<ttk::SimplexId, int> cpMap,
209 std::vector<std::pair<ttk::SimplexId, ttk::SimplexId>> &persistencePairs,
210 vtkDataArray *inputScalars,
211 const triangulationType *triangulation) {
212 vtkNew<vtkUnstructuredGrid> skeletonNodes{};
213 vtkNew<vtkPoints> points{};
214 vtkNew<vtkCellArray> cells{};
215 vtkNew<vtkLongLongArray> gIdArray{};
216 gIdArray->SetNumberOfComponents(1);
217 gIdArray->SetName("VertexId");
218 vtkSmartPointer<vtkDataArray> const scalarArray
219 = vtkSmartPointer<vtkDataArray>::Take(inputScalars->NewInstance());
220 scalarArray->SetNumberOfComponents(1);
221 scalarArray->SetName("Scalar");
222 vtkNew<vtkIntArray> cpArray{};
223 cpArray->SetNumberOfComponents(1);
224 cpArray->SetName("CriticalType");
225
226 float point[3];
227 long long pointId = 0;
228
229 for(auto const &pair : persistencePairs) {
230 triangulation->getVertexPoint(pair.first, point[0], point[1], point[2]);
231 points->InsertNextPoint(point);
232 gIdArray->InsertNextTuple1(pair.first);
233 scalarArray->InsertNextTuple1(inputScalars->GetTuple1(pair.first));
234 cpArray->InsertNextTuple1(cpMap[pair.first]);
235 triangulation->getVertexPoint(pair.second, point[0], point[1], point[2]);
236 points->InsertNextPoint(point);
237 gIdArray->InsertNextTuple1(pair.second);
238 scalarArray->InsertNextTuple1(inputScalars->GetTuple1(pair.second));
239 cpArray->InsertNextTuple1(cpMap[pair.second]);
240 skeletonNodes->InsertNextCell(VTK_VERTEX, 1, &pointId);
241 pointId++;
242 skeletonNodes->InsertNextCell(VTK_VERTEX, 1, &pointId);
243 pointId++;
244 }
245 skeletonNodes->SetPoints(points);
246 outputSkeletonNodes->ShallowCopy(skeletonNodes);
247 outputSkeletonNodes->GetPointData()->AddArray(gIdArray);
248 outputSkeletonNodes->GetPointData()->AddArray(scalarArray);
249 outputSkeletonNodes->GetPointData()->AddArray(cpArray);
250
251 return 1;
252 }
253};
ttk::ftm::TreeType GetTreeType() const
std::vector< ttk::Triangulation * > triangulation_
int addCompleteSkeletonArc(const ttk::ftm::idSuperArc arcId, const int cc, vtkPoints *points, vtkUnstructuredGrid *skeletonArcs, ttk::ftm::ArcData &arcData)
int addSampledSkeletonArc(const ttk::ftm::idSuperArc arcId, const int cc, vtkPoints *points, vtkUnstructuredGrid *skeletonArcs, ttk::ftm::ArcData &arcData)
std::vector< vtkDataArray * > inputScalars_
int getSegmentation(vtkDataSet *outputSegmentation)
ttk::ftm::Params params_
int getSkeletonNodes(vtkUnstructuredGrid *outputSkeletonNodes)
std::vector< std::vector< ttk::SimplexId > > offsets_
int addDirectSkeletonArc(const ttk::ftm::idSuperArc arcId, const int cc, vtkPoints *points, vtkUnstructuredGrid *skeletonArcs, ttk::ftm::ArcData &arcData)
int getSkeletonArcs(vtkUnstructuredGrid *outputSkeletonArcs)
void identify(vtkDataSet *ds) const
std::vector< vtkSmartPointer< vtkDataSet > > connected_components_
int getMergeTreePoints(vtkUnstructuredGrid *outputSkeletonNodes, std::map< ttk::SimplexId, int > cpMap, std::vector< std::pair< ttk::SimplexId, ttk::SimplexId > > &persistencePairs, vtkDataArray *inputScalars, const triangulationType *triangulation)
std::vector< ttk::ftm::LocalFTM > ftmTree_
int getMergeTree(vtkUnstructuredGrid *outputSkeletonArcs, std::vector< ttk::ExTreeM::Branch > &mergeTree, vtkDataArray *inputScalars, const triangulationType *triangulation)
~ttkMergeTreeBase() override=default
Minimalist debugging class.
Definition Debug.h:88
long unsigned int idSuperArc
SuperArc index in vect_superArcs_.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22