VTK  9.2.6
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 class vtkImageData;
132 
133 class VTKFILTERSCORE_EXPORT vtkSurfaceNets2D : public vtkPolyDataAlgorithm
134 {
135 public:
137 
140  static vtkSurfaceNets2D* New();
142  void PrintSelf(ostream& os, vtkIndent indent) override;
144 
149  vtkMTimeType GetMTime() override;
150 
152 
162  void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
163  void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
165 
167 
170  double GetValue(int i) { return this->Labels->GetValue(i); }
171  double GetLabel(int i) { return this->Labels->GetValue(i); }
173 
175 
179  double* GetValues() { return this->Labels->GetValues(); }
180  double* GetLabels() { return this->Labels->GetValues(); }
182 
184 
189  void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
190  void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
192 
194 
201  void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
202  void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
204 
206 
209  vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
210  vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
212 
214 
218  void GenerateLabels(int numLabels, double range[2])
219  {
220  this->Labels->GenerateValues(numLabels, range);
221  }
222  void GenerateValues(int numContours, double range[2])
223  {
224  this->Labels->GenerateValues(numContours, range);
225  }
226  void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
227  {
228  this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
229  }
230  void GenerateValues(int numContours, double rangeStart, double rangeEnd)
231  {
232  this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
233  }
235 
237 
245  vtkSetMacro(ComputeScalars, bool);
246  vtkGetMacro(ComputeScalars, bool);
247  vtkBooleanMacro(ComputeScalars, bool);
249 
251 
261  vtkSetMacro(BackgroundLabel, double);
262  vtkGetMacro(BackgroundLabel, double);
264 
266 
270  vtkSetMacro(ArrayComponent, int);
271  vtkGetMacro(ArrayComponent, int);
273 
275 
280  vtkSetMacro(Smoothing, bool);
281  vtkGetMacro(Smoothing, bool);
282  vtkBooleanMacro(Smoothing, bool);
284 
286 
293  vtkGetSmartPointerMacro(Smoother, vtkConstrainedSmoothingFilter);
295 
297 
307  vtkSetMacro(DataCaching, bool);
308  vtkGetMacro(DataCaching, bool);
309  vtkBooleanMacro(DataCaching, bool);
311 
312 protected:
314  ~vtkSurfaceNets2D() override = default;
315 
317  int FillInputPortInformation(int port, vtkInformation* info) override;
318 
323 
324  bool Smoothing;
326 
327  // Support data caching of the extracted surface nets. This is used to
328  // avoid repeated surface extraction when only smoothing filter
329  // parameters are modified.
334  bool IsCacheEmpty();
335  void CacheData(vtkPolyData* pd, vtkCellArray* ca);
336 
337 private:
338  vtkSurfaceNets2D(const vtkSurfaceNets2D&) = delete;
339  void operator=(const vtkSurfaceNets2D&) = delete;
340 };
341 
342 #endif
adjust point positions using constrained smoothing
Store vtkAlgorithm input/output information.
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:287
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:35
vtkSmartPointer< vtkConstrainedSmoothingFilter > Smoother
generate smoothed constours from segmented 2D image data (i.e., "label maps")
int vtkIdType
Definition: vtkType.h:332
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:90
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:39
topologically and geometrically regular array of data
Definition: vtkImageData.h:53
virtual vtkMTimeType GetMTime()
Return this object's modified time.
object to represent cell connectivity
Definition: vtkCellArray.h:186
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.