TTK
Loading...
Searching...
No Matches
ttkFTRGraph.h
Go to the documentation of this file.
1
48
49#pragma once
50
51// ttk code includes
52#include <FTRGraph.h>
53#include <Graph.h>
54#include <ttkAlgorithm.h>
56
57// VTK includes
58#include <vtkDataArray.h>
59
60// VTK Module
61#include <ttkFTRGraphModule.h>
62
63class TTKFTRGRAPH_EXPORT ttkFTRGraph : public ttkAlgorithm {
64public:
65 static ttkFTRGraph *New();
66 vtkTypeMacro(ttkFTRGraph, ttkAlgorithm);
67
68 vtkSetMacro(ForceInputOffsetScalarField, bool);
69 vtkGetMacro(ForceInputOffsetScalarField, bool);
70
75 void SetSingleSweep(const bool ss) {
76 params_.singleSweep = ss;
77 Modified();
78 }
79 bool GetSingleSweep() const {
80 return params_.singleSweep;
81 }
83
86 void SetWithSegmentation(const bool segm) {
87 params_.segm = segm;
88 Modified();
89 }
90 bool GetWithSegmentation() const {
91 return params_.segm;
92 }
94
98 void SetWithNormalize(const bool norm) {
99 params_.normalize = norm;
100 Modified();
101 }
102 bool GetWithNormalize() const {
103 return params_.normalize;
104 }
106
109 void SetSampling(int lvl) {
110 params_.samplingLvl = lvl;
111 Modified();
112 }
114 return params_.samplingLvl;
115 }
117
118 int getSkeletonNodes(const ttk::ftr::Graph &graph,
119 vtkUnstructuredGrid *outputSkeletonNodes);
120
121 int addDirectSkeletonArc(const ttk::ftr::Graph &graph,
122 const ttk::ftr::idSuperArc arcId,
123 vtkPoints *points,
124 vtkUnstructuredGrid *skeletonArcs,
125 ttk::ftr::ArcData &arcData);
126
127 int addSampledSkeletonArc(const ttk::ftr::Graph &graph,
128 const ttk::ftr::idSuperArc arcId,
129 vtkPoints *points,
130 vtkUnstructuredGrid *skeletonArcs,
131 ttk::ftr::ArcData &arcData);
132
133 int addCompleteSkeletonArc(const ttk::ftr::Graph &graph,
134 const ttk::ftr::idSuperArc arcId,
135 vtkPoints *points,
136 vtkUnstructuredGrid *skeletonArcs,
137 ttk::ftr::ArcData &arcData);
138
139 int getSkeletonArcs(const ttk::ftr::Graph &graph,
140 vtkUnstructuredGrid *outputSkeletonArcs);
141
142 int getSegmentation(const ttk::ftr::Graph &graph,
143 vtkDataSet *outputSegmentation);
144
145 template <typename VTK_TT, typename TTK_TT>
146 int dispatch(ttk::ftr::Graph &graph);
147
148protected:
149 ttkFTRGraph();
150
151 void identify(vtkDataSet *ds) const;
152
153 int FillInputPortInformation(int port, vtkInformation *info) override;
154 int FillOutputPortInformation(int port, vtkInformation *info) override;
155 int RequestData(vtkInformation *request,
156 vtkInformationVector **inputVector,
157 vtkInformationVector *outputVector) override;
158
159private:
160 bool ForceInputOffsetScalarField{};
161 ttk::ftr::Params params_{};
162
163 vtkDataSet *mesh_{};
164 ttk::Triangulation *triangulation_{};
165 vtkDataArray *inputScalars_{};
166 vtkDataArray *offsets_{};
167};
Baseclass of all VTK filters that wrap ttk modules.
Definition: ttkAlgorithm.h:34
virtual int RequestData(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
Definition: ttkAlgorithm.h:390
int FillInputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
Definition: ttkAlgorithm.h:404
int FillOutputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
Definition: ttkAlgorithm.h:419
TTK VTK-filter for the computation of Reeb Graphs.
Definition: ttkFTRGraph.h:63
bool GetSingleSweep() const
Definition: ttkFTRGraph.h:79
int GetSuperArcSamplingLevel() const
Definition: ttkFTRGraph.h:113
static ttkFTRGraph * New()
void SetSampling(int lvl)
control the sampling level of the superarcs
Definition: ttkFTRGraph.h:109
void SetSingleSweep(const bool ss)
control whether the computation should start from min and max (default) or use a single sweep startin...
Definition: ttkFTRGraph.h:75
bool GetWithSegmentation() const
Definition: ttkFTRGraph.h:90
void SetWithNormalize(const bool norm)
if set to true, a post processing pass will be used to enforce consistent node ids between executions
Definition: ttkFTRGraph.h:98
bool GetWithNormalize() const
Definition: ttkFTRGraph.h:102
void SetWithSegmentation(const bool segm)
control if the output should contains the segmentation information
Definition: ttkFTRGraph.h:86
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
Definition: Triangulation.h:48
TTK FTRGraph graph skeleton.
Definition: Graph.h:38
long unsigned int idSuperArc
SuperArc index in vect_superArcs_.
Definition: FTRDataTypes.h:25