VTK  9.3.0
vtkSurfaceNets2D.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
114#ifndef vtkSurfaceNets2D_h
115#define vtkSurfaceNets2D_h
116
117#include "vtkConstrainedSmoothingFilter.h" // Perform mesh smoothing
118#include "vtkContourValues.h" // Needed for direct access to ContourValues
119#include "vtkFiltersCoreModule.h" // For export macro
120#include "vtkPolyData.h" // To support data caching
121#include "vtkPolyDataAlgorithm.h"
122
123VTK_ABI_NAMESPACE_BEGIN
124
125class vtkImageData;
126
127class VTKFILTERSCORE_EXPORT vtkSurfaceNets2D : public vtkPolyDataAlgorithm
128{
129public:
131
136 void PrintSelf(ostream& os, vtkIndent indent) override;
138
144
146
156 void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
157 void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
159
161
164 double GetValue(int i) { return this->Labels->GetValue(i); }
165 double GetLabel(int i) { return this->Labels->GetValue(i); }
167
169
173 double* GetValues() { return this->Labels->GetValues(); }
174 double* GetLabels() { return this->Labels->GetValues(); }
176
178
183 void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
184 void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
186
188
195 void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
196 void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
198
200
203 vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
204 vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
206
208
212 void GenerateLabels(int numLabels, double range[2])
213 {
214 this->Labels->GenerateValues(numLabels, range);
215 }
216 void GenerateValues(int numContours, double range[2])
217 {
218 this->Labels->GenerateValues(numContours, range);
219 }
220 void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
221 {
222 this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
223 }
224 void GenerateValues(int numContours, double rangeStart, double rangeEnd)
225 {
226 this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
227 }
229
231
239 vtkSetMacro(ComputeScalars, bool);
240 vtkGetMacro(ComputeScalars, bool);
241 vtkBooleanMacro(ComputeScalars, bool);
243
245
255 vtkSetMacro(BackgroundLabel, double);
256 vtkGetMacro(BackgroundLabel, double);
258
260
264 vtkSetMacro(ArrayComponent, int);
265 vtkGetMacro(ArrayComponent, int);
267
269
274 vtkSetMacro(Smoothing, bool);
275 vtkGetMacro(Smoothing, bool);
276 vtkBooleanMacro(Smoothing, bool);
278
280
289
291
301 vtkSetMacro(DataCaching, bool);
302 vtkGetMacro(DataCaching, bool);
303 vtkBooleanMacro(DataCaching, bool);
305
306protected:
308 ~vtkSurfaceNets2D() override = default;
309
311 int FillInputPortInformation(int port, vtkInformation* info) override;
312
317
320
321 // Support data caching of the extracted surface nets. This is used to
322 // avoid repeated surface extraction when only smoothing filter
323 // parameters are modified.
330
331private:
332 vtkSurfaceNets2D(const vtkSurfaceNets2D&) = delete;
333 void operator=(const vtkSurfaceNets2D&) = delete;
334};
335
336VTK_ABI_NAMESPACE_END
337#endif
object to represent cell connectivity
adjust point positions using constrained smoothing
topologically and geometrically regular array of data
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:80
Hold a reference to a vtkObjectBase instance.
generate smoothed constours from segmented 2D image data (i.e., "label maps")
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, printing, and type information.
void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
Generate numLabels equally spaced labels between the specified range.
void SetNumberOfLabels(int number)
Set the number of labels to place into the list.
vtkMTimeType GetMTime() override
The modified time is also a function of the label values and the smoothing filter.
vtkSmartPointer< vtkConstrainedSmoothingFilter > Smoother
vtkSmartPointer< vtkContourValues > Labels
vtkSmartPointer< vtkPolyData > GeometryCache
void GetValues(double *contourValues)
Fill a supplied list with label values.
vtkTimeStamp SmoothingTime
void GenerateValues(int numContours, double range[2])
Generate numLabels equally spaced labels between the specified range.
void SetLabel(int i, double value)
Set a particular label value at label number i.
vtkGetSmartPointerMacro(Smoother, vtkConstrainedSmoothingFilter)
Get the instance of vtkConstrainedSmoothingFilter used to smooth the extracted surface net.
double * GetLabels()
Get a pointer to an array of labels.
~vtkSurfaceNets2D() override=default
double GetLabel(int i)
Get the ith label value.
double GetValue(int i)
Get the ith label value.
void GenerateLabels(int numLabels, double range[2])
Generate numLabels equally spaced labels between the specified range.
void SetNumberOfContours(int number)
Set the number of labels to place into the list.
void GetLabels(double *contourValues)
Fill a supplied list with label values.
double * GetValues()
Get a pointer to an array of labels.
static vtkSurfaceNets2D * New()
Standard methods for instantiation, printing, and type information.
void CacheData(vtkPolyData *pd, vtkCellArray *ca)
void GenerateValues(int numContours, double rangeStart, double rangeEnd)
Generate numLabels equally spaced labels between the specified range.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkIdType GetNumberOfLabels()
Get the number of labels in the list of label values.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkIdType GetNumberOfContours()
Get the number of labels in the list of label values.
vtkSmartPointer< vtkCellArray > StencilsCache
void SetValue(int i, double value)
Set a particular label value at label number i.
record modification and/or execution time
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270