VTK  9.4.20241217
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
119#ifndef vtkSurfaceNets2D_h
120#define vtkSurfaceNets2D_h
121
122#include "vtkConstrainedSmoothingFilter.h" // Perform mesh smoothing
123#include "vtkContourValues.h" // Needed for direct access to ContourValues
124#include "vtkFiltersCoreModule.h" // For export macro
125#include "vtkPolyData.h" // To support data caching
126#include "vtkPolyDataAlgorithm.h"
127
128VTK_ABI_NAMESPACE_BEGIN
129
130class vtkImageData;
131
132class VTKFILTERSCORE_EXPORT vtkSurfaceNets2D : public vtkPolyDataAlgorithm
133{
134public:
136
141 void PrintSelf(ostream& os, vtkIndent indent) override;
143
149
151
161 void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
162 void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
164
166
169 double GetValue(int i) { return this->Labels->GetValue(i); }
170 double GetLabel(int i) { return this->Labels->GetValue(i); }
172
174
178 double* GetValues() { return this->Labels->GetValues(); }
179 double* GetLabels() { return this->Labels->GetValues(); }
181
183
188 void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
189 void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
191
193
200 void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
201 void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
203
205
208 vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
209 vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
211
213
217 void GenerateLabels(int numLabels, double range[2])
218 {
219 this->Labels->GenerateValues(numLabels, range);
220 }
221 void GenerateValues(int numContours, double range[2])
222 {
223 this->Labels->GenerateValues(numContours, range);
224 }
225 void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
226 {
227 this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
228 }
229 void GenerateValues(int numContours, double rangeStart, double rangeEnd)
230 {
231 this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
232 }
234
236
244 vtkSetMacro(ComputeScalars, bool);
245 vtkGetMacro(ComputeScalars, bool);
246 vtkBooleanMacro(ComputeScalars, bool);
248
250
260 vtkSetMacro(BackgroundLabel, double);
261 vtkGetMacro(BackgroundLabel, double);
263
265
269 vtkSetMacro(ArrayComponent, int);
270 vtkGetMacro(ArrayComponent, int);
272
274
279 vtkSetMacro(Smoothing, bool);
280 vtkGetMacro(Smoothing, bool);
281 vtkBooleanMacro(Smoothing, bool);
283
285
294
296
306 vtkSetMacro(DataCaching, bool);
307 vtkGetMacro(DataCaching, bool);
308 vtkBooleanMacro(DataCaching, bool);
310
311protected:
313 ~vtkSurfaceNets2D() override = default;
314
316 int FillInputPortInformation(int port, vtkInformation* info) override;
317
322
325
326 // Support data caching of the extracted surface nets. This is used to
327 // avoid repeated surface extraction when only smoothing filter
328 // parameters are modified.
335
336private:
337 vtkSurfaceNets2D(const vtkSurfaceNets2D&) = delete;
338 void operator=(const vtkSurfaceNets2D&) = delete;
339};
340
341VTK_ABI_NAMESPACE_END
342#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:108
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
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