VTK  9.3.1
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 
123 VTK_ABI_NAMESPACE_BEGIN
124 
125 class vtkImageData;
126 
127 class VTKFILTERSCORE_EXPORT vtkSurfaceNets2D : public vtkPolyDataAlgorithm
128 {
129 public:
131 
134  static vtkSurfaceNets2D* New();
136  void PrintSelf(ostream& os, vtkIndent indent) override;
138 
143  vtkMTimeType GetMTime() override;
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 
287  vtkGetSmartPointerMacro(Smoother, vtkConstrainedSmoothingFilter);
289 
291 
301  vtkSetMacro(DataCaching, bool);
302  vtkGetMacro(DataCaching, bool);
303  vtkBooleanMacro(DataCaching, bool);
305 
306 protected:
308  ~vtkSurfaceNets2D() override = default;
309 
311  int FillInputPortInformation(int port, vtkInformation* info) override;
312 
317 
318  bool Smoothing;
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.
328  bool IsCacheEmpty();
329  void CacheData(vtkPolyData* pd, vtkCellArray* ca);
330 
331 private:
332  vtkSurfaceNets2D(const vtkSurfaceNets2D&) = delete;
333  void operator=(const vtkSurfaceNets2D&) = delete;
334 };
335 
336 VTK_ABI_NAMESPACE_END
337 #endif
adjust point positions using constrained smoothing
Store vtkAlgorithm input/output information.
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
vtkSmartPointer< vtkCellArray > StencilsCache
record modification and/or execution time
Definition: vtkTimeStamp.h:24
vtkSmartPointer< vtkConstrainedSmoothingFilter > Smoother
generate smoothed constours from segmented 2D image data (i.e., "label maps")
int vtkIdType
Definition: vtkType.h:315
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:79
void SetValue(int i, double value)
Set a particular label value at label number i.
void GetLabels(double *contourValues)
Fill a supplied list with label values.
void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
Generate numLabels equally spaced labels between the specified range.
vtkIdType GetNumberOfLabels()
Get the number of labels in the list of label values.
void SetNumberOfLabels(int number)
Set the number of labels to place into the list.
static vtkPolyDataAlgorithm * New()
vtkSmartPointer< vtkContourValues > Labels
vtkIdType GetNumberOfContours()
Get the number of labels in the list of label values.
Superclass for algorithms that produce only polydata as output.
double * GetValues()
Get a pointer to an array of labels.
void SetNumberOfContours(int number)
Set the number of labels to place into the list.
double GetValue(int i)
Get the ith label value.
a simple class to control print indentation
Definition: vtkIndent.h:28
topologically and geometrically regular array of data
Definition: vtkImageData.h:42
virtual vtkMTimeType GetMTime()
Return this object's modified time.
object to represent cell connectivity
Definition: vtkCellArray.h:175
double * GetLabels()
Get a pointer to an array of labels.
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.
vtkTimeStamp SmoothingTime
Store zero or more vtkInformation instances.
void SetLabel(int i, double value)
Set a particular label value at label number i.
vtkSmartPointer< vtkPolyData > GeometryCache
double GetLabel(int i)
Get the ith label value.
void GenerateValues(int numContours, double range[2])
Generate numLabels equally spaced labels between the specified range.
void GenerateLabels(int numLabels, double range[2])
Generate numLabels equally spaced labels between the specified range.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void GetValues(double *contourValues)
Fill a supplied list with label values.