VTK  9.2.20230320
vtkSurfaceNets2D.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkSurfaceNets2D.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
122 #ifndef vtkSurfaceNets2D_h
123 #define vtkSurfaceNets2D_h
124 
125 #include "vtkConstrainedSmoothingFilter.h" // Perform mesh smoothing
126 #include "vtkContourValues.h" // Needed for direct access to ContourValues
127 #include "vtkFiltersCoreModule.h" // For export macro
128 #include "vtkPolyData.h" // To support data caching
129 #include "vtkPolyDataAlgorithm.h"
130 
131 VTK_ABI_NAMESPACE_BEGIN
132 
133 class vtkImageData;
134 
135 class VTKFILTERSCORE_EXPORT vtkSurfaceNets2D : public vtkPolyDataAlgorithm
136 {
137 public:
139 
144  void PrintSelf(ostream& os, vtkIndent indent) override;
146 
151  vtkMTimeType GetMTime() override;
152 
154 
164  void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
165  void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
167 
169 
172  double GetValue(int i) { return this->Labels->GetValue(i); }
173  double GetLabel(int i) { return this->Labels->GetValue(i); }
175 
177 
181  double* GetValues() { return this->Labels->GetValues(); }
182  double* GetLabels() { return this->Labels->GetValues(); }
184 
186 
191  void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
192  void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
194 
196 
203  void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
204  void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
206 
208 
211  vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
212  vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
214 
216 
220  void GenerateLabels(int numLabels, double range[2])
221  {
222  this->Labels->GenerateValues(numLabels, range);
223  }
224  void GenerateValues(int numContours, double range[2])
225  {
226  this->Labels->GenerateValues(numContours, range);
227  }
228  void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
229  {
230  this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
231  }
232  void GenerateValues(int numContours, double rangeStart, double rangeEnd)
233  {
234  this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
235  }
237 
239 
247  vtkSetMacro(ComputeScalars, bool);
248  vtkGetMacro(ComputeScalars, bool);
249  vtkBooleanMacro(ComputeScalars, bool);
251 
253 
263  vtkSetMacro(BackgroundLabel, double);
264  vtkGetMacro(BackgroundLabel, double);
266 
268 
272  vtkSetMacro(ArrayComponent, int);
273  vtkGetMacro(ArrayComponent, int);
275 
277 
282  vtkSetMacro(Smoothing, bool);
283  vtkGetMacro(Smoothing, bool);
284  vtkBooleanMacro(Smoothing, bool);
286 
288 
297 
299 
309  vtkSetMacro(DataCaching, bool);
310  vtkGetMacro(DataCaching, bool);
311  vtkBooleanMacro(DataCaching, bool);
313 
314 protected:
316  ~vtkSurfaceNets2D() override = default;
317 
320 
325 
326  bool Smoothing;
328 
329  // Support data caching of the extracted surface nets. This is used to
330  // avoid repeated surface extraction when only smoothing filter
331  // parameters are modified.
336  bool IsCacheEmpty();
338 
339 private:
340  vtkSurfaceNets2D(const vtkSurfaceNets2D&) = delete;
341  void operator=(const vtkSurfaceNets2D&) = delete;
342 };
343 
344 VTK_ABI_NAMESPACE_END
345 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:297
adjust point positions using constrained smoothing
topologically and geometrically regular array of data
Definition: vtkImageData.h:164
a simple class to control print indentation
Definition: vtkIndent.h:120
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:201
generate smoothed constours from segmented 2D image data (i.e., "label maps")
double * GetValues()
Get a pointer to an array of labels.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, printing, and type information.
static vtkSurfaceNets2D * New()
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
double * GetLabels()
Get a pointer to an array of labels.
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.
~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.
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
Definition: vtkTimeStamp.h:56
@ info
Definition: vtkX3D.h:388
@ value
Definition: vtkX3D.h:232
@ port
Definition: vtkX3D.h:459
@ range
Definition: vtkX3D.h:250
int vtkIdType
Definition: vtkType.h:327
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:282