VTK  9.6.20260605
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
181
182#ifndef vtkSurfaceNets3D_h
183#define vtkSurfaceNets3D_h
184
185#include "vtkConstrainedSmoothingFilter.h" // Perform mesh smoothing
186#include "vtkContourValues.h" // Needed for direct access to ContourValues
187#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_7_0
188#include "vtkFiltersCoreModule.h" // For export macro
189#include "vtkPolyData.h" // To support data caching
190#include "vtkPolyDataAlgorithm.h"
191#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
192
193#include <vector> // For selected seeds
194
195VTK_ABI_NAMESPACE_BEGIN
196
197class vtkImageData;
198
199class VTKFILTERSCORE_EXPORT VTK_MARSHALAUTO vtkSurfaceNets3D : public vtkPolyDataAlgorithm
200{
201public:
203
209 void PrintSelf(ostream& os, vtkIndent indent) override;
211
217
219
232 void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
233 void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
235
237
240 double GetValue(int i) { return this->Labels->GetValue(i); }
241 double GetLabel(int i) { return this->Labels->GetValue(i); }
243
245
249 double* GetValues() { return this->Labels->GetValues(); }
250 double* GetLabels() { return this->Labels->GetValues(); }
252
254
259 void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
260 void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
262
264
271 void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
272 void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
274
276
279 vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
280 vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
282
284
288 void GenerateLabels(int numLabels, double range[2])
289 {
290 this->Labels->GenerateValues(numLabels, range);
291 }
292 void GenerateValues(int numContours, double range[2])
293 {
294 this->Labels->GenerateValues(numContours, range);
295 }
296 void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
297 {
298 this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
299 }
300 void GenerateValues(int numContours, double rangeStart, double rangeEnd)
301 {
302 this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
303 }
304
305
307
317 vtkSetMacro(BackgroundLabel, double);
318 vtkGetMacro(BackgroundLabel, double);
320
322
326 vtkSetMacro(ArrayComponent, int);
327 vtkGetMacro(ArrayComponent, int);
329
339
341
352 vtkGetMacro(OutputMeshType, int);
357
358 // The following code is used to control the smoothing process. Internally
359 // there is a vtkConstrainedSmoothingFilter that can be directly
360 // manipulated. In addition, methods that delegate to this filter are also
361 // provided.
362
364
369 vtkSetMacro(Smoothing, bool);
370 vtkGetMacro(Smoothing, bool);
371 vtkBooleanMacro(Smoothing, bool);
373
375
380 void SetNumberOfIterations(int n) { this->Smoother->SetNumberOfIterations(n); }
381 int GetNumberOfIterations() { return this->Smoother->GetNumberOfIterations(); }
382 void SetRelaxationFactor(double f) { this->Smoother->SetRelaxationFactor(f); }
383 double GetRelaxationFactor() { return this->Smoother->GetRelaxationFactor(); }
384 void SetConstraintDistance(double d) { this->Smoother->SetConstraintDistance(d); }
385 double GetConstraintDistance() { return this->Smoother->GetConstraintDistance(); }
386 void SetConstraintBox(double sx, double sy, double sz)
387 {
388 this->Smoother->SetConstraintBox(sx, sy, sz);
389 }
390 void SetConstraintBox(double s[3]) { this->Smoother->SetConstraintBox(s); }
391 double* GetConstraintBox() VTK_SIZEHINT(3) { return this->Smoother->GetConstraintBox(); }
392 void GetConstraintBox(double s[3]) { this->Smoother->GetConstraintBox(s); }
394 {
395 this->Smoother->SetConstraintStrategyToConstraintDistance();
396 }
398 {
399 this->Smoother->SetConstraintStrategyToConstraintBox();
400 }
401 int GetConstraintStrategy() { return this->Smoother->GetConstraintStrategy(); }
403
405
419 vtkBooleanMacro(AutomaticSmoothingConstraints, bool);
420 vtkSetClampMacro(ConstraintScale, double, 0, 100);
421 vtkGetMacro(ConstraintScale, double);
423
425
431 VTK_DEPRECATED_IN_9_7_0("No longer used")
432 virtual void SetOptimizedSmoothingStencils(bool) {}
433 VTK_DEPRECATED_IN_9_7_0("No longer used")
434 virtual bool GetOptimizedSmoothingStencils() { return true; }
435 VTK_DEPRECATED_IN_9_7_0("No longer used")
437 VTK_DEPRECATED_IN_9_7_0("No longer used")
440
442
449 vtkGetSmartPointerMacro(Smoother, vtkConstrainedSmoothingFilter);
451
452 // The following code is used to control what is produced for output.
453
474
476
489 vtkGetMacro(OutputStyle, int);
494
496
501 void AddSelectedLabel(double label);
502 void DeleteSelectedLabel(double label);
504 double GetSelectedLabel(vtkIdType ithLabel);
506
516
518
529 vtkGetMacro(TriangulationStrategy, int);
539
540
542
552 vtkSetMacro(DataCaching, bool);
553 vtkGetMacro(DataCaching, bool);
554 vtkBooleanMacro(DataCaching, bool);
556
557protected:
559 ~vtkSurfaceNets3D() override = default;
560
561 // Support visualization pipeline operations.
563 int FillInputPortInformation(int port, vtkInformation* info) override;
564
565 // Support the contouring operation.
571
572 // Support smoothing.
577
578 // Support data caching of the extracted surface nets. This is used to
579 // avoid repeated surface extraction when only smoothing filter
580 // parameters are modified.
587
588 // Support output style
590 std::vector<double> SelectedLabels;
592
593 // Support triangulation strategy
595
596private:
597 vtkSurfaceNets3D(const vtkSurfaceNets3D&) = delete;
598 void operator=(const vtkSurfaceNets3D&) = delete;
599};
600
601VTK_ABI_NAMESPACE_END
602#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.
concrete dataset represents vertices, lines, polygons, and triangle strips
Hold a reference to a vtkObjectBase instance.
virtual void SetOutputStyle(int)
Specify the form (i.e., the style) of the output.
double * GetConstraintBox()
Convenience methods that delegate to the internal smoothing filter follow below.
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_...
virtual bool GetOptimizedSmoothingStencils()
Indicate whether to use optimized smoothing stencils.
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.
virtual void OptimizedSmoothingStencilsOn()
Indicate whether to use optimized smoothing stencils.
double * GetLabels()
Get a pointer to an array of labels.
void InitializeSelectedLabelsList()
When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled r...
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_...
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
virtual void OptimizedSmoothingStencilsOff()
Indicate whether to use optimized smoothing stencils.
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.
virtual void SetTriangulationStrategy(int)
Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_...
static vtkSurfaceNets3D * New()
Standard methods for instantiation, printing, and obtaining type information.
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.
virtual void SetOutputMeshType(int)
Control the type of output mesh.
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.
virtual void SetOptimizedSmoothingStencils(bool)
Indicate whether to use optimized smoothing stencils.
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
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...
double * GetValues()
Get a pointer to an array of labels.
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
#define VTK_DEPRECATED_IN_9_7_0(reason)
int vtkIdType
Definition vtkType.h:363
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:318
#define VTK_SIZEHINT(...)
#define VTK_MARSHALAUTO