TTK
Loading...
Searching...
No Matches
ttkMorseSmaleComplex.h
Go to the documentation of this file.
1
114
115#pragma once
116
117// VTK Module
118#include <ttkMorseSmaleComplexModule.h>
119
120// ttk code includes
121#include <MorseSmaleComplex.h>
122#include <ttkAlgorithm.h>
123
124class vtkPolyData;
125
126class TTKMORSESMALECOMPLEX_EXPORT ttkMorseSmaleComplex
127 : public ttkAlgorithm,
128 protected ttk::MorseSmaleComplex {
129
130public:
132
134
135 vtkSetMacro(ForceInputOffsetScalarField, bool);
136 vtkGetMacro(ForceInputOffsetScalarField, bool);
137
138 vtkSetMacro(ComputeCriticalPoints, bool);
139 vtkGetMacro(ComputeCriticalPoints, bool);
140
143
146
147 vtkSetMacro(ComputeSaddleConnectors, bool);
148 vtkGetMacro(ComputeSaddleConnectors, bool);
149
152
155
158
161
162 vtkSetMacro(ComputeFinalSegmentation, bool);
163 vtkGetMacro(ComputeFinalSegmentation, bool);
164
165 vtkSetMacro(IterationThreshold, int);
166 vtkGetMacro(IterationThreshold, int);
167
168 vtkSetMacro(ReturnSaddleConnectors, bool);
169 vtkGetMacro(ReturnSaddleConnectors, bool);
170
171 vtkSetMacro(DiscreteGradientBackend, int);
172 vtkGetMacro(DiscreteGradientBackend, int);
173
176
177 vtkSetMacro(ThresholdIsAbsolute, bool);
178 vtkGetMacro(ThresholdIsAbsolute, bool);
179
180 vtkSetMacro(ForceLoopFreeGradient, bool);
181 vtkGetMacro(ForceLoopFreeGradient, bool);
182
183 vtkSetMacro(StochasticGradientSeed, unsigned int);
184 vtkGetMacro(StochasticGradientSeed, unsigned int);
185
186protected:
187 template <typename scalarType, typename triangulationType>
188 int dispatch(vtkDataArray *const inputScalars,
189 vtkPolyData *const outputCriticalPoints,
190 vtkPolyData *const outputSeparatrices1,
191 vtkPolyData *const outputSeparatrices2,
192 const SimplexId *const inputOffsets,
193 const triangulationType &triangulation);
194
196
197 int FillInputPortInformation(int port, vtkInformation *info) override;
198 int FillOutputPortInformation(int port, vtkInformation *info) override;
199 int RequestData(vtkInformation *request,
200 vtkInformationVector **inputVector,
201 vtkInformationVector *outputVector) override;
202
203private:
204 bool ForceInputOffsetScalarField{};
205 int IterationThreshold{-1};
206 int DiscreteGradientBackend{0};
207 OutputManifold segmentations_{};
208 unsigned int StochasticGradientSeed{0};
209};
virtual int RequestData(vtkInformation *ttkNotUsed(request), vtkInformationVector **ttkNotUsed(inputVectors), vtkInformationVector *ttkNotUsed(outputVector))
int FillInputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
int FillOutputPortInformation(int ttkNotUsed(port), vtkInformation *ttkNotUsed(info)) override
static ttkMorseSmaleComplex * New()
int dispatch(vtkDataArray *const inputScalars, vtkPolyData *const outputCriticalPoints, vtkPolyData *const outputSeparatrices1, vtkPolyData *const outputSeparatrices2, const SimplexId *const inputOffsets, const triangulationType &triangulation)
TTK processing package for the computation of Morse-Smale complexes.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22