TTK
Loading...
Searching...
No Matches
ttkFTRGraph.cpp
Go to the documentation of this file.
1#include "ttkFTRGraph.h"
2#include <ttkMacros.h>
3#include <ttkUtils.h>
4
5// only used on the cpp
6#include <vtkConnectivityFilter.h>
7#include <vtkDataObject.h>
8#include <vtkInformation.h>
9
10using namespace ttk::ftr;
11
13
15 this->setDebugMsgPrefix("FTRGraph");
16 this->SetNumberOfInputPorts(1);
17 this->SetNumberOfOutputPorts(3);
18}
19
20int ttkFTRGraph::FillInputPortInformation(int port, vtkInformation *info) {
21 if(port == 0) {
22 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
23 return 1;
24 }
25 return 0;
26}
27
28int ttkFTRGraph::FillOutputPortInformation(int port, vtkInformation *info) {
29 if(port == 0 || port == 1) {
30 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
31 return 1;
32 } else if(port == 2) {
34 return 1;
35 }
36 return 0;
37}
38
40 const idSuperArc arcId,
41 vtkPoints *points,
42 vtkUnstructuredGrid *skeletonArcs,
43 ttk::ftr::ArcData &arcData) {
44 const SuperArc &arc = graph.getArc(arcId);
45 float pointCoord[3];
46 const idNode upNodeId = arc.getUpNodeId();
47 const idNode downNodeId = arc.getDownNodeId();
48
49#ifndef TTK_ENABLE_KAMIKAZE
50 if(upNodeId == nullNode || downNodeId == nullNode) {
51 this->printErr("NULL NODES IN SKELETON " + graph.printArc(arcId));
52 return 1;
53 }
54#endif
55
56 const idVertex upVertId = graph.getNode(upNodeId).getVertexIdentifier();
57 const idVertex downVertId = graph.getNode(downNodeId).getVertexIdentifier();
58
59 vtkIdType pointIds[2];
60 if(arcData.points.count(downVertId)) {
61 pointIds[0] = arcData.points[downVertId];
62 } else {
63 triangulation_->getVertexPoint(
64 downVertId, pointCoord[0], pointCoord[1], pointCoord[2]);
65 pointIds[0] = points->InsertNextPoint(pointCoord);
66 arcData.points.emplace(downVertId, pointIds[0]);
67 arcData.setPointInfo(graph, arcId, pointIds[0]);
68 }
69
70 for(const idVertex regV : arc.segmentation()) {
71
72 // for regular mask to be set by the boundary vertices correctly
73 if(regV == downVertId || regV == upVertId)
74 continue;
75
76 if(arcData.points.count(regV)) {
77 pointIds[1] = arcData.points[regV];
78 } else {
79 triangulation_->getVertexPoint(
80 regV, pointCoord[0], pointCoord[1], pointCoord[2]);
81 pointIds[1] = points->InsertNextPoint(pointCoord);
82 arcData.points.emplace(regV, pointIds[1]);
83 arcData.setPointInfo(graph, arcId, pointIds[1], true);
84 }
85 if(pointIds[0] != pointIds[1]) {
86 vtkIdType nextCellId
87 = skeletonArcs->InsertNextCell(VTK_LINE, 2, pointIds);
88 arcData.setArcInfo(graph, arcId, nextCellId);
89 }
90 pointIds[0] = pointIds[1];
91 }
92
93 if(arcData.points.count(upVertId)) {
94 pointIds[1] = arcData.points[upVertId];
95 } else {
96 triangulation_->getVertexPoint(
97 upVertId, pointCoord[0], pointCoord[1], pointCoord[2]);
98 pointIds[1] = points->InsertNextPoint(pointCoord);
99 arcData.points.emplace(upVertId, pointIds[1]);
100 arcData.setPointInfo(graph, arcId, pointIds[1]);
101 }
102 if(pointIds[0] != pointIds[1]) {
103 vtkIdType nextCellId = skeletonArcs->InsertNextCell(VTK_LINE, 2, pointIds);
104 arcData.setArcInfo(graph, arcId, nextCellId);
105 }
106
107 return 0;
108}
109
111 const idSuperArc arcId,
112 vtkPoints *points,
113 vtkUnstructuredGrid *skeletonArcs,
114 ttk::ftr::ArcData &arcData) {
115 float pointCoord[3];
116 const idNode upNodeId = graph.getArc(arcId).getUpNodeId();
117 const idNode downNodeId = graph.getArc(arcId).getDownNodeId();
118
119#ifndef TTK_ENABLE_KAMIKAZE
120 if(upNodeId == nullNode || downNodeId == nullNode) {
121 this->printErr("NULL NODES IN SKELETON " + graph.printArc(arcId));
122 return 1;
123 }
124#endif
125
126 const idVertex upVertId = graph.getNode(upNodeId).getVertexIdentifier();
127 const idVertex downVertId = graph.getNode(downNodeId).getVertexIdentifier();
128
129 // Mesh skeleton
130 vtkIdType pointIds[2];
131 if(arcData.points.count(downVertId)) {
132 pointIds[0] = arcData.points[downVertId];
133 } else {
134 triangulation_->getVertexPoint(
135 downVertId, pointCoord[0], pointCoord[1], pointCoord[2]);
136 pointIds[0] = points->InsertNextPoint(pointCoord);
137 arcData.points.emplace(downVertId, pointIds[0]);
138 arcData.setPointInfo(graph, arcId, pointIds[0]);
139 }
140 if(arcData.points.count(upVertId)) {
141 pointIds[1] = arcData.points[upVertId];
142 } else {
143 triangulation_->getVertexPoint(
144 upVertId, pointCoord[0], pointCoord[1], pointCoord[2]);
145 pointIds[1] = points->InsertNextPoint(pointCoord);
146 arcData.points.emplace(upVertId, pointIds[1]);
147 arcData.setPointInfo(graph, arcId, pointIds[1]);
148 }
149 vtkIdType nextCellId = skeletonArcs->InsertNextCell(VTK_LINE, 2, pointIds);
150
151 // Scalar arrays
152 arcData.setArcInfo(graph, arcId, nextCellId);
153 return 0;
154}
155
157 const idSuperArc arcId,
158 vtkPoints *points,
159 vtkUnstructuredGrid *skeletonArcs,
160 ttk::ftr::ArcData &arcData) {
161 const SuperArc &arc = graph.getArc(arcId);
162 float pointCoord[3];
163 const idNode upNodeId = arc.getUpNodeId();
164 const idNode downNodeId = arc.getDownNodeId();
165
166 const double scalarMin
167 = inputScalars_->GetTuple1(graph.getNode(downNodeId).getVertexIdentifier());
168 const double scalarMax
169 = inputScalars_->GetTuple1(graph.getNode(upNodeId).getVertexIdentifier());
170 const double delta = (scalarMax - scalarMin) / (params_.samplingLvl + 1);
171
172#ifndef TTK_ENABLE_KAMIKAZE
173 if(upNodeId == nullNode || downNodeId == nullNode) {
174 this->printErr("NULL NODES IN SKELETON " + graph.printArc(arcId));
175 return 1;
176 }
177#endif
178
179 const idVertex upVertId = graph.getNode(upNodeId).getVertexIdentifier();
180 const idVertex downVertId = graph.getNode(downNodeId).getVertexIdentifier();
181
182 vtkIdType pointIds[2];
183 if(arcData.points.count(downVertId)) {
184 pointIds[0] = arcData.points[downVertId];
185 } else {
186 triangulation_->getVertexPoint(
187 downVertId, pointCoord[0], pointCoord[1], pointCoord[2]);
188 pointIds[0] = points->InsertNextPoint(pointCoord);
189 arcData.points.emplace(downVertId, pointIds[0]);
190 arcData.setPointInfo(graph, arcId, pointIds[0]);
191 }
192
193 float sum[3]{};
194 int chunk = 0;
195 double scalarLimit = scalarMin + delta;
196
197 for(const idVertex regV : arc.segmentation()) {
198
199 // for regular mask to be set by the boundary vertices correctly
200 if(regV == downVertId || regV == upVertId)
201 continue;
202
203 triangulation_->getVertexPoint(
204 regV, pointCoord[0], pointCoord[1], pointCoord[2]);
205 const double scalar = inputScalars_->GetTuple1(regV);
206 if(scalar < scalarLimit) {
207 sum[0] += pointCoord[0];
208 sum[1] += pointCoord[1];
209 sum[2] += pointCoord[2];
210 ++chunk;
211 } else {
212 if(chunk) {
213 sum[0] /= chunk;
214 sum[1] /= chunk;
215 sum[2] /= chunk;
216
217 // do not use memorized points as even if a point already exist it is
218 // from another vertex and should not be used
219 pointIds[1] = points->InsertNextPoint(sum);
220 arcData.points.emplace(regV, pointIds[1]);
221 arcData.setPointInfo(graph, arcId, pointIds[1], true);
222
223 if(pointIds[0] != pointIds[1]) {
224 vtkIdType nextCellId
225 = skeletonArcs->InsertNextCell(VTK_LINE, 2, pointIds);
226 arcData.setArcInfo(graph, arcId, nextCellId);
227 }
228
229 pointIds[0] = pointIds[1];
230 }
231
232 // reset for next iteration
233 scalarLimit += delta;
234 sum[0] = 0;
235 sum[1] = 0;
236 sum[2] = 0;
237 chunk = 0;
238 }
239 }
240
241 if(arcData.points.count(upVertId)) {
242 pointIds[1] = arcData.points[upVertId];
243 } else {
244 triangulation_->getVertexPoint(
245 upVertId, pointCoord[0], pointCoord[1], pointCoord[2]);
246 pointIds[1] = points->InsertNextPoint(pointCoord);
247 arcData.points.emplace(upVertId, pointIds[1]);
248 arcData.setPointInfo(graph, arcId, pointIds[1]);
249 }
250 if(pointIds[0] != pointIds[1]) {
251 vtkIdType nextCellId = skeletonArcs->InsertNextCell(VTK_LINE, 2, pointIds);
252 arcData.setArcInfo(graph, arcId, nextCellId);
253 }
254
255 return 0;
256}
257
258template <typename VTK_TT, typename TTK_TT>
261 static_cast<TTK_TT *>(triangulation_->getData()));
262
263 // common parameters
264 ftrGraph_.setParams(params_);
265 ftrGraph_.setDebugLevel(this->debugLevel_);
266 // reeb graph parameters
267 ftrGraph_.setScalars(ttkUtils::GetVoidPointer(inputScalars_));
268 // TODO: SimplexId -> template to int + long long?
269 ftrGraph_.setVertexSoSoffsets(
270 static_cast<ttk::SimplexId *>(ttkUtils::GetVoidPointer(offsets_)));
271 // build
272 this->printMsg("Starting computation on field: "
273 + std::string{inputScalars_->GetName()});
274 ftrGraph_.build();
275 // get output
276 graph = std::move(ftrGraph_.extractOutputGraph());
277
278 return 0;
279}
280
281int ttkFTRGraph::RequestData(vtkInformation *ttkNotUsed(request),
282 vtkInformationVector **inputVector,
283 vtkInformationVector *outputVector) {
284
285 // Input
286 mesh_ = vtkDataSet::GetData(inputVector[0]);
287
288 // Output
289 auto outputSkeletonNodes = vtkUnstructuredGrid::GetData(outputVector, 0);
290 auto outputSkeletonArcs = vtkUnstructuredGrid::GetData(outputVector, 1);
291 auto outputSegmentation = vtkDataSet::GetData(outputVector, 2);
292
293 // Skeleton
294 Graph graph;
295
296 triangulation_ = ttkAlgorithm::GetTriangulation(mesh_);
297#ifndef TTK_ENABLE_KAMIKAZE
298 if(!triangulation_) {
299 this->printErr("ttkTriangulation::getTriangulation() is null.");
300 return -1;
301 }
302 if(triangulation_->isEmpty()) {
303 this->printErr(
304 "ttkTriangulation on connected component allocation problem.");
305 return -1;
306 }
307#endif
308
309 // gives parameters to the params_ structure
310 params_.debugLevel = debugLevel_;
311 params_.threadNumber = threadNumber_;
312
313 // Scalar field related
314 inputScalars_ = this->GetInputArrayToProcess(0, inputVector);
315 if(inputScalars_ == nullptr) {
316 this->printErr("input scalar field pointer is null.");
317 return -3;
318 }
319
320 offsets_ = this->GetOrderArray(mesh_, 0, 1, ForceInputOffsetScalarField);
321
322 // compute graph
323 ttkVtkTemplateMacro(inputScalars_->GetDataType(), triangulation_->getType(),
324 (dispatch<VTK_TT, TTK_TT>(graph)));
325
326 UpdateProgress(0.50);
327
328 // Construct output
329 if(getSkeletonNodes(graph, outputSkeletonNodes)) {
330#ifndef TTK_ENABLE_KAMIKAZE
331 this->printErr("wrong properties on skeleton nodes.");
332 return -7;
333#endif
334 }
335
336 if(getSkeletonArcs(graph, outputSkeletonArcs)) {
337#ifndef TTK_ENABLE_KAMIKAZE
338 this->printErr("wrong properties on skeleton arcs.");
339 return -8;
340#endif
341 }
342
343 if(GetWithSegmentation() && getSegmentation(graph, outputSegmentation)) {
344#ifndef TTK_ENABLE_KAMIKAZE
345 this->printErr("wrong properties on segmentation.");
346 return -9;
347#endif
348 }
349
350 UpdateProgress(1);
351
352 return 1;
353}
354
356 vtkDataSet *outputSegmentation) {
357 outputSegmentation->ShallowCopy(mesh_);
358
359 const idVertex numberOfVertices = mesh_->GetNumberOfPoints();
360 ttk::ftr::VertData vertData(numberOfVertices);
361
362 // TODO parallel
363 for(idVertex v = 0; v < numberOfVertices; ++v) {
364 vertData.setVertexInfo(graph, v);
365 }
366
367 vertData.addArrays(outputSegmentation, params_);
368
369 return 0;
370}
371
373 vtkUnstructuredGrid *outputSkeletonArcs) {
374 const idSuperArc nbArcs = graph.getNumberOfArcs();
375
376 idSuperArc nbFinArc = 0;
377 switch(params_.samplingLvl) {
378 case -1:
379 // loops may create more arcs
380 nbFinArc = triangulation_->getNumberOfVertices() * 1.5;
381 break;
382 case 0:
383 nbFinArc = graph.getNumberOfVisibleArcs();
384 break;
385 default:
386 nbFinArc = graph.getNumberOfVisibleArcs() * (params_.samplingLvl + 1);
387 break;
388 }
389
390 ttk::ftr::ArcData arcData(nbFinArc);
391 vtkNew<vtkUnstructuredGrid> arcs{};
392 vtkNew<vtkPoints> points{};
393
394 for(idSuperArc arcId = 0; arcId < nbArcs; ++arcId) {
395 if(!graph.getArc(arcId).isVisible())
396 continue;
397
398 switch(params_.samplingLvl) {
399 case -1:
400 addCompleteSkeletonArc(graph, arcId, points, arcs, arcData);
401 break;
402 case 0:
403 addDirectSkeletonArc(graph, arcId, points, arcs, arcData);
404 break;
405 default:
406 addSampledSkeletonArc(graph, arcId, points, arcs, arcData);
407 break;
408 }
409 }
410
411 arcs->SetPoints(points);
412 arcData.addArrays(arcs, params_);
413 outputSkeletonArcs->ShallowCopy(arcs);
414
415 return 0;
416}
417
419 vtkUnstructuredGrid *outputSkeletonNodes) {
420 const idNode nbNodes = graph.getNumberOfNodes();
421
422 ttk::ftr::NodeData nodeData(nbNodes);
423 vtkNew<vtkUnstructuredGrid> nodes{};
424 vtkNew<vtkPoints> points{};
425
426 for(idNode nodeId = 0; nodeId < nbNodes; nodeId++) {
427 const ttk::ftr::idVertex vertId
428 = graph.getNode(nodeId).getVertexIdentifier();
429 float point[3];
430 triangulation_->getVertexPoint(vertId, point[0], point[1], point[2]);
431 points->InsertNextPoint(point);
432
433 double scalar = inputScalars_->GetTuple1(vertId);
434 nodeData.addNode(graph, nodeId, scalar);
435 }
436
437 ttkUtils::CellVertexFromPoints(nodes, points);
438 nodeData.addArrays(nodes->GetPointData(), params_);
439
440 outputSkeletonNodes->ShallowCopy(nodes);
441
442 return 0;
443}
444
445// protected
446
447void ttkFTRGraph::identify(vtkDataSet *ds) const {
448 vtkNew<vtkIntArray> identifiers{};
449 const vtkIdType nbPoints = ds->GetNumberOfPoints();
450 identifiers->SetName("VertexIdentifier");
451 identifiers->SetNumberOfComponents(1);
452 identifiers->SetNumberOfTuples(nbPoints);
453
454 for(int i = 0; i < nbPoints; i++) {
455 identifiers->SetTuple1(i, i);
456 }
457
458 ds->GetPointData()->AddArray(identifiers);
459}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition: BaseClass.h:47
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
TTK VTK-filter for the computation of Reeb Graphs.
Definition: ttkFTRGraph.h:63
int FillInputPortInformation(int port, vtkInformation *info) override
Definition: ttkFTRGraph.cpp:20
int addDirectSkeletonArc(const ttk::ftr::Graph &graph, const ttk::ftr::idSuperArc arcId, vtkPoints *points, vtkUnstructuredGrid *skeletonArcs, ttk::ftr::ArcData &arcData)
int FillOutputPortInformation(int port, vtkInformation *info) override
Definition: ttkFTRGraph.cpp:28
int getSkeletonNodes(const ttk::ftr::Graph &graph, vtkUnstructuredGrid *outputSkeletonNodes)
bool GetWithSegmentation() const
Definition: ttkFTRGraph.h:90
void identify(vtkDataSet *ds) const
int getSegmentation(const ttk::ftr::Graph &graph, vtkDataSet *outputSegmentation)
int getSkeletonArcs(const ttk::ftr::Graph &graph, vtkUnstructuredGrid *outputSkeletonArcs)
int addSampledSkeletonArc(const ttk::ftr::Graph &graph, const ttk::ftr::idSuperArc arcId, vtkPoints *points, vtkUnstructuredGrid *skeletonArcs, ttk::ftr::ArcData &arcData)
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int dispatch(ttk::ftr::Graph &graph)
int addCompleteSkeletonArc(const ttk::ftr::Graph &graph, const ttk::ftr::idSuperArc arcId, vtkPoints *points, vtkUnstructuredGrid *skeletonArcs, ttk::ftr::ArcData &arcData)
Definition: ttkFTRGraph.cpp:39
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition: ttkUtils.cpp:225
static int CellVertexFromPoints(vtkDataSet *const dataSet, vtkPoints *const points)
Definition: ttkUtils.cpp:327
int threadNumber_
Definition: BaseClass.h:95
int debugLevel_
Definition: Debug.h:379
void setDebugMsgPrefix(const std::string &prefix)
Definition: Debug.h:364
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
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition: Debug.h:149
AbstractTriangulation * getData()
Triangulation::Type getType() const
int getVertexPoint(const SimplexId &vertexId, float &x, float &y, float &z) const override
bool isEmpty() const override
SimplexId getNumberOfVertices() const override
TTK FTRGraph processing package.
Definition: FTRGraph.h:77
void setVertexSoSoffsets(SimplexId *sos)
Definition: FTRGraph.h:198
void setScalars(const void *scalars)
Scalar field used to compute the Reeb Graph.
Definition: FTRGraph.h:188
int setDebugLevel(const int &lvl) override
Control the verbosity of the base code.
Definition: FTRGraph.h:177
Graph && extractOutputGraph()
Definition: FTRGraph.h:160
void setParams(const Params &p)
Definition: FTRGraph.h:182
TTK FTRGraph graph skeleton.
Definition: Graph.h:38
idNode getNumberOfNodes() const
Definition: Graph.h:84
idSuperArc getNumberOfArcs() const
Definition: Graph.h:88
std::string printArc(const idSuperArc arcId) const
Definition: Graph.cpp:54
const SuperArc & getArc(const idSuperArc id) const
Definition: Graph.h:121
const Node & getNode(const idNode id) const
Definition: Graph.h:113
idSuperArc getNumberOfVisibleArcs() const
Definition: Graph.h:92
idVertex getVertexIdentifier() const
Definition: FTRNode.h:34
TTK FTRGraph graph arc.
Definition: FTRSuperArc.h:30
idNode getUpNodeId() const
Definition: FTRSuperArc.h:48
idNode getDownNodeId() const
Definition: FTRSuperArc.h:57
const decltype(segmentation_) & segmentation() const
Definition: FTRSuperArc.h:147
bool isVisible() const
Definition: FTRSuperArc.h:89
long unsigned int idSuperArc
SuperArc index in vect_superArcs_.
Definition: FTRDataTypes.h:25
SimplexId idVertex
Vertex index in scalars_.
Definition: FTRDataTypes.h:29
unsigned int idNode
Node index in vect_nodes_.
Definition: FTRDataTypes.h:27
int SimplexId
Identifier type for simplices of any dimension.
Definition: DataTypes.h:22
std::map< ttk::ftr::idVertex, vtkIdType > points
void setPointInfo(const ttk::ftr::Graph &ttkNotUsed(graph), const ttk::ftr::idSuperArc ttkNotUsed(a), const vtkIdType skeletonVert, bool r=false)
void addArrays(vtkUnstructuredGrid *arcs, ttk::ftr::Params ttkNotUsed(params))
void setArcInfo(const ttk::ftr::Graph &graph, const ttk::ftr::idSuperArc a, const vtkIdType skeletonCell)
void addNode(const ttk::ftr::Graph &graph, const ttk::ftr::idNode n, const double scalar)
void addArrays(vtkPointData *pointData, ttk::ftr::Params ttkNotUsed(params))
void setVertexInfo(const ttk::ftr::Graph &graph, const ttk::ftr::idVertex v)
void addArrays(vtkDataSet *segmentation, ttk::ftr::Params ttkNotUsed(params))
vtkStandardNewMacro(ttkFTRGraph)
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition: ttkMacros.h:69