VTK  9.3.0
vtkSurfaceNets3D.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
145 #ifndef vtkSurfaceNets3D_h
146 #define vtkSurfaceNets3D_h
147 
148 #include "vtkConstrainedSmoothingFilter.h" // Perform mesh smoothing
149 #include "vtkContourValues.h" // Needed for direct access to ContourValues
150 #include "vtkFiltersCoreModule.h" // For export macro
151 #include "vtkPolyData.h" // To support data caching
152 #include "vtkPolyDataAlgorithm.h"
153 
154 #include <vector> // For selected seeds
155 
156 VTK_ABI_NAMESPACE_BEGIN
157 
158 class vtkImageData;
159 
160 class VTKFILTERSCORE_EXPORT vtkSurfaceNets3D : public vtkPolyDataAlgorithm
161 {
162 public:
164 
170  void PrintSelf(ostream& os, vtkIndent indent) override;
172 
177  vtkMTimeType GetMTime() override;
178 
180 
193  void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
194  void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
196 
198 
201  double GetValue(int i) { return this->Labels->GetValue(i); }
202  double GetLabel(int i) { return this->Labels->GetValue(i); }
204 
206 
210  double* GetValues() { return this->Labels->GetValues(); }
211  double* GetLabels() { return this->Labels->GetValues(); }
213 
215 
220  void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
221  void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
223 
225 
232  void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
233  void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
235 
237 
240  vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
241  vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
243 
245 
249  void GenerateLabels(int numLabels, double range[2])
250  {
251  this->Labels->GenerateValues(numLabels, range);
252  }
253  void GenerateValues(int numContours, double range[2])
254  {
255  this->Labels->GenerateValues(numContours, range);
256  }
257  void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
258  {
259  this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
260  }
261  void GenerateValues(int numContours, double rangeStart, double rangeEnd)
262  {
263  this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
264  }
266 
268 
278  vtkSetMacro(BackgroundLabel, double);
279  vtkGetMacro(BackgroundLabel, double);
281 
283 
287  vtkSetMacro(ArrayComponent, int);
288  vtkGetMacro(ArrayComponent, int);
290 
294  enum MeshType
295  {
296  MESH_TYPE_DEFAULT = 0,
298  MESH_TYPE_QUADS
299  };
300 
302 
312  vtkSetClampMacro(OutputMeshType, int, MESH_TYPE_DEFAULT, MESH_TYPE_QUADS);
313  vtkGetMacro(OutputMeshType, int);
314  void SetOutputMeshTypeToDefault() { this->SetOutputMeshType(MESH_TYPE_DEFAULT); }
315  void SetOutputMeshTypeToTriangles() { this->SetOutputMeshType(MESH_TYPE_TRIANGLES); }
316  void SetOutputMeshTypeToQuads() { this->SetOutputMeshType(MESH_TYPE_QUADS); }
318 
319  // The following code is used to control the smoothing process. Internally
320  // there is a vtkConstrainedSmoothingFilter that can be directly
321  // manipulated. In addition, methods that delegate to this filter are also
322  // provided.
323 
325 
330  vtkSetMacro(Smoothing, bool);
331  vtkGetMacro(Smoothing, bool);
332  vtkBooleanMacro(Smoothing, bool);
334 
336 
341  void SetNumberOfIterations(int n) { this->Smoother->SetNumberOfIterations(n); }
342  int GetNumberOfIterations() { return this->Smoother->GetNumberOfIterations(); }
343  void SetRelaxationFactor(double f) { this->Smoother->SetRelaxationFactor(f); }
344  double GetRelaxationFactor() { return this->Smoother->GetRelaxationFactor(); }
345  void SetConstraintDistance(double d) { this->Smoother->SetConstraintDistance(d); }
346  double GetConstraintDistance() { return this->Smoother->GetConstraintDistance(); }
347  void SetConstraintBox(double sx, double sy, double sz)
348  {
349  this->Smoother->SetConstraintBox(sx, sy, sz);
350  };
351  void SetConstraintBox(double s[3]) { this->Smoother->SetConstraintBox(s); };
352  double* GetConstraintBox() VTK_SIZEHINT(3) { return this->Smoother->GetConstraintBox(); }
353  void GetConstraintBox(double s[3]) { this->Smoother->GetConstraintBox(s); }
355  {
356  this->Smoother->SetConstraintStrategyToConstraintDistance();
357  }
359  {
360  this->Smoother->SetConstraintStrategyToConstraintBox();
361  }
362  int GetConstraintStrategy() { return this->Smoother->GetConstraintStrategy(); }
364 
366 
378  vtkSetMacro(AutomaticSmoothingConstraints, bool);
379  vtkGetMacro(AutomaticSmoothingConstraints, bool);
380  vtkBooleanMacro(AutomaticSmoothingConstraints, bool);
381  vtkSetClampMacro(ConstraintScale, double, 0, 100);
382  vtkGetMacro(ConstraintScale, double);
384 
386 
392  vtkSetMacro(OptimizedSmoothingStencils, bool);
393  vtkGetMacro(OptimizedSmoothingStencils, bool);
394  vtkBooleanMacro(OptimizedSmoothingStencils, bool);
396 
398 
407 
408  // The following code is used to control what is produced for output.
409 
425  {
426  OUTPUT_STYLE_DEFAULT = 0,
428  OUTPUT_STYLE_SELECTED
429  };
430 
432 
444  vtkSetClampMacro(OutputStyle, int, OUTPUT_STYLE_DEFAULT, OUTPUT_STYLE_SELECTED);
445  vtkGetMacro(OutputStyle, int);
446  void SetOutputStyleToDefault() { this->SetOutputStyle(OUTPUT_STYLE_DEFAULT); }
447  void SetOutputStyleToBoundary() { this->SetOutputStyle(OUTPUT_STYLE_BOUNDARY); }
448  void SetOutputStyleToSelected() { this->SetOutputStyle(OUTPUT_STYLE_SELECTED); }
450 
452 
457  void AddSelectedLabel(double label);
458  void DeleteSelectedLabel(double label);
460  double GetSelectedLabel(vtkIdType ithLabel);
462 
467  {
468  TRIANGULATION_GREEDY = 0,
470  TRIANGULATION_MIN_AREA
471  };
472 
474 
484  vtkSetClampMacro(TriangulationStrategy, int, TRIANGULATION_GREEDY, TRIANGULATION_MIN_AREA);
485  vtkGetMacro(TriangulationStrategy, int);
486  void SetTriangulationStrategyToGreedy() { this->SetTriangulationStrategy(TRIANGULATION_GREEDY); }
488  {
489  this->SetTriangulationStrategy(TRIANGULATION_MIN_EDGE);
490  }
492  {
493  this->SetTriangulationStrategy(TRIANGULATION_MIN_AREA);
494  }
496 
498 
508  vtkSetMacro(DataCaching, bool);
509  vtkGetMacro(DataCaching, bool);
510  vtkBooleanMacro(DataCaching, bool);
512 
513 protected:
515  ~vtkSurfaceNets3D() override = default;
516 
517  // Support visualization pipeline operations.
520 
521  // Support the contouring operation.
527 
528  // Support smoothing.
529  bool Smoothing;
534 
535  // Support data caching of the extracted surface nets. This is used to
536  // avoid repeated surface extraction when only smoothing filter
537  // parameters are modified.
542  bool IsCacheEmpty();
544 
545  // Support output style
547  std::vector<double> SelectedLabels;
549 
550  // Support triangulation strategy
552 
553 private:
554  vtkSurfaceNets3D(const vtkSurfaceNets3D&) = delete;
555  void operator=(const vtkSurfaceNets3D&) = delete;
556 };
557 
558 VTK_ABI_NAMESPACE_END
559 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:185
adjust point positions using constrained smoothing
topologically and geometrically regular array of data
Definition: vtkImageData.h:52
a simple class to control print indentation
Definition: vtkIndent.h:38
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:89
generate smoothed isocontours from segmented 3D image data (i.e., "label maps")
double * GetValues()
Get a pointer to an array of labels.
void GetConstraintBox(double s[3])
Convenience methods that delegate to the internal smoothing filter follow below.
void SetTriangulationStrategyToMinArea()
Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_...
double GetRelaxationFactor()
Convenience methods that delegate to the internal smoothing filter follow below.
void SetLabel(int i, double value)
Set a particular label value at label number i.
double GetSelectedLabel(vtkIdType ithLabel)
When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled r...
void SetOutputMeshTypeToTriangles()
Control the type of output mesh.
double GetLabel(int i)
Get the ith label value.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
double * GetLabels()
Get a pointer to an array of labels.
static vtkSurfaceNets3D * New()
Standard methods for instantiation, printing, and obtaining type information.
void InitializeSelectedLabelsList()
When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled r...
double * GetConstraintBox()
Convenience methods that delegate to the internal smoothing filter follow below.
void GetValues(double *contourValues)
Fill a supplied list with label values.
void SetTriangulationStrategyToMinEdge()
Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_...
vtkGetSmartPointerMacro(Smoother, vtkConstrainedSmoothingFilter)
Get the instance of vtkConstrainedSmoothingFilter used to smooth the extracted surface net.
void SetConstraintDistance(double d)
Convenience methods that delegate to the internal smoothing filter follow below.
int GetNumberOfIterations()
Convenience methods that delegate to the internal smoothing filter follow below.
MeshType
This enum is used to control the type of the output polygonal mesh.
void SetOutputMeshTypeToQuads()
Control the type of output mesh.
vtkSmartPointer< vtkContourValues > Labels
void GenerateLabels(int numLabels, double range[2])
Generate numLabels equally spaced labels between the specified range.
void SetConstraintStrategyToConstraintDistance()
Convenience methods that delegate to the internal smoothing filter follow below.
vtkSmartPointer< vtkConstrainedSmoothingFilter > Smoother
std::vector< double > SelectedLabels
vtkIdType GetNumberOfSelectedLabels()
When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled r...
~vtkSurfaceNets3D() override=default
void SetConstraintBox(double sx, double sy, double sz)
Convenience methods that delegate to the internal smoothing filter follow below.
void SetNumberOfIterations(int n)
Convenience methods that delegate to the internal smoothing filter follow below.
vtkTimeStamp SelectedLabelsTime
void SetOutputStyleToDefault()
Specify the form (i.e., the style) of the output.
void SetConstraintBox(double s[3])
Convenience methods that delegate to the internal smoothing filter follow below.
TriangulationType
This enum is used to control how quadrilaterals are triangulated.
void SetRelaxationFactor(double f)
Convenience methods that delegate to the internal smoothing filter follow below.
void GetLabels(double *contourValues)
Fill a supplied list with label values.
void GenerateValues(int numContours, double range[2])
Generate numLabels equally spaced labels between the specified range.
void GenerateValues(int numContours, double rangeStart, double rangeEnd)
Generate numLabels equally spaced labels between the specified range.
void SetOutputStyleToBoundary()
Specify the form (i.e., the style) of the output.
void CacheData(vtkPolyData *pd, vtkCellArray *ca)
void SetOutputMeshTypeToDefault()
Control the type of output mesh.
vtkIdType GetNumberOfLabels()
Get the number of labels in the list of label values.
int GetConstraintStrategy()
Convenience methods that delegate to the internal smoothing filter follow below.
double GetConstraintDistance()
Convenience methods that delegate to the internal smoothing filter follow below.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkTimeStamp SmoothingTime
vtkIdType GetNumberOfContours()
Get the number of labels in the list of label values.
vtkMTimeType GetMTime() override
The modified time is also a function of the label values and the smoothing filter.
void AddSelectedLabel(double label)
When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled r...
void SetNumberOfLabels(int number)
Set the number of labels to place into the list.
void SetOutputStyleToSelected()
Specify the form (i.e., the style) of the output.
void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
Generate numLabels equally spaced labels between the specified range.
void SetTriangulationStrategyToGreedy()
Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_...
void SetValue(int i, double value)
Set a particular label value at label number i.
void SetConstraintStrategyToConstraintBox()
Convenience methods that delegate to the internal smoothing filter follow below.
vtkSmartPointer< vtkPolyData > GeometryCache
bool AutomaticSmoothingConstraints
OutputType
This enum is used to control the production of the filter output.
void DeleteSelectedLabel(double label)
When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled r...
void SetNumberOfContours(int number)
Set the number of labels to place into the list.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, printing, and obtaining type information.
double GetValue(int i)
Get the ith label value.
vtkSmartPointer< vtkCellArray > StencilsCache
record modification and/or execution time
Definition: vtkTimeStamp.h:34
@ info
Definition: vtkX3D.h:376
@ value
Definition: vtkX3D.h:220
@ port
Definition: vtkX3D.h:447
@ range
Definition: vtkX3D.h:238
int vtkIdType
Definition: vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270
#define VTK_SIZEHINT(...)