VTK  9.6.20260216
vtkGeneralizedSurfaceNets3D.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
128
129#ifndef vtkGeneralizedSurfaceNets3D_h
130#define vtkGeneralizedSurfaceNets3D_h
131
132#include "vtkConstrainedSmoothingFilter.h" // Perform surface smoothing
133#include "vtkContourValues.h" // Manage countour values
134#include "vtkFiltersMeshingModule.h" // For export macro
135#include "vtkPolyDataAlgorithm.h"
136#include "vtkStaticPointLocator.h" // For point locator
137#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
138
139VTK_ABI_NAMESPACE_BEGIN
140
141class VTKFILTERSMESHING_EXPORT VTK_MARSHALAUTO vtkGeneralizedSurfaceNets3D
142 : public vtkPolyDataAlgorithm
143{
144public:
146
151 void PrintSelf(ostream& os, vtkIndent indent) override;
153
154 //---------------The following are methods used to set/get and generate
155 //---------------contour labels. Contour labels are required to be
156 //---------------specified by this filter.
158
170 void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
171 void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
173
175
178 double GetValue(int i) { return this->Labels->GetValue(i); }
179 double GetLabel(int i) { return this->Labels->GetValue(i); }
181
183
187 double* GetValues() { return this->Labels->GetValues(); }
188 double* GetLabels() { return this->Labels->GetValues(); }
190
192
197 void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
198 void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
200
202
209 void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
210 void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
212
214
217 vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
218 vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
220
222
226 void GenerateLabels(int numLabels, double range[2])
227 {
228 this->Labels->GenerateValues(numLabels, range);
229 }
230 void GenerateValues(int numContours, double range[2])
231 {
232 this->Labels->GenerateValues(numContours, range);
233 }
234 void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
235 {
236 this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
237 }
238 void GenerateValues(int numContours, double rangeStart, double rangeEnd)
239 {
240 this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
241 }
242
243
245
252 vtkSetMacro(BackgroundLabel, int);
253 vtkGetMacro(BackgroundLabel, int);
255
256 //---------------Done defining label-related methods.
258
264 vtkBooleanMacro(BoundaryCapping, vtkTypeBool);
266
268
276 vtkBooleanMacro(MergePoints, vtkTypeBool);
278
280
288 vtkSetMacro(Smoothing, vtkTypeBool);
289 vtkGetMacro(Smoothing, vtkTypeBool);
290 vtkBooleanMacro(Smoothing, vtkTypeBool);
292
294
299 void SetNumberOfIterations(int n) { this->Smoother->SetNumberOfIterations(n); }
300 int GetNumberOfIterations() { return this->Smoother->GetNumberOfIterations(); }
301 void SetRelaxationFactor(double f) { this->Smoother->SetRelaxationFactor(f); }
302 double GetRelaxationFactor() { return this->Smoother->GetRelaxationFactor(); }
303 void SetConstraintDistance(double d) { this->Smoother->SetConstraintDistance(d); }
304 double GetConstraintDistance() { return this->Smoother->GetConstraintDistance(); }
306
308
327
329
341 vtkSetVector4Macro(SmoothingConstraints, unsigned char);
342 vtkGetVector4Macro(SmoothingConstraints, unsigned char);
347
349
358 vtkGetSmartPointerMacro(Smoother, vtkConstrainedSmoothingFilter);
360
373
375
386 vtkGetMacro(OutputMeshType, int);
391
393
399 vtkSetClampMacro(Padding, double, 0.001, 0.25);
400 vtkGetMacro(Padding, double);
402
404
411
413
419 vtkSetMacro(Validate, vtkTypeBool);
420 vtkGetMacro(Validate, vtkTypeBool);
421 vtkBooleanMacro(Validate, vtkTypeBool);
423
425
440 vtkSetClampMacro(PointOfInterest, vtkIdType, -1, VTK_ID_MAX);
442 vtkSetObjectMacro(PointsOfInterest, vtkIdTypeArray);
443 vtkGetObjectMacro(PointsOfInterest, vtkIdTypeArray);
447
449
453 vtkSetClampMacro(PruneTolerance, double, 0.0, 0.5);
454 vtkGetMacro(PruneTolerance, double);
456
458
466 vtkSetClampMacro(BatchSize, unsigned int, 1, VTK_INT_MAX);
467 vtkGetMacro(BatchSize, unsigned int);
469
474 int GetNumberOfThreads() { return this->NumberOfThreads; }
475
480 int GetNumberOfPrunes() { return this->NumberOfPrunes; }
481
487
488protected:
490 ~vtkGeneralizedSurfaceNets3D() override = default;
491
492 // Support the contouring operation by defining labels.
495
496 // Algorithm control
497 vtkTypeBool BoundaryCapping; // produce boundary surfaces or not
498 vtkTypeBool MergePoints; // merge near coincident points or not
499 vtkTypeBool Smoothing; // enable built-in smoothing process
500
501 // Internal classes Related to point location and smoothing control
503 vtkSmartPointer<vtkStaticPointLocator> Locator; // locator for finding proximal points
504
505 // Control the type of output mesh. Triangles by default.
507
508 // Related to internal Voronoi methods
509 double Padding; // amount to pad out input points bounding box
510 vtkTypeBool Validate; // Choose to validate and repair output
511 vtkTypeBool GenerateSmoothingStencils; // Produce smoothing stencils
512 unsigned char SmoothingConstraints[4]; // Specify which smoothing constraints are active
513 vtkIdType PointOfInterest; // specify a single input point to process
515 vtkIdType MaximumNumberOfHullClips; // limit the number of hull clips
516 double PruneTolerance; // prune tiny faces
517 unsigned int BatchSize; // process data in batches of specified size
518 int NumberOfThreads; // report on the number of threads used during processing
519 int NumberOfPrunes; // If spoke pruning is enabled, report number of pruning operations
520
521 // Satisfy pipeline-related API
523 int FillInputPortInformation(int port, vtkInformation* info) override;
524
525private:
527 void operator=(const vtkGeneralizedSurfaceNets3D&) = delete;
528};
529
530VTK_ABI_NAMESPACE_END
531#endif
adjust point positions using constrained smoothing
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
void GenerateValues(int numContours, double range[2])
Generate numLabels equally spaced labels between the specified range.
int GetNumberOfPrunes()
Return the number of hull prunes performed during execution.
void SetValue(int i, double value)
Set a particular label value at label number i.
MeshType
This enum is used to control the type of the output polygonal mesh.
vtkSmartPointer< vtkContourValues > Labels
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, type information, and printing.
void GenerateLabels(int numLabels, double range[2])
Generate numLabels equally spaced labels between the specified range.
double GetRelaxationFactor()
Convenience methods that delegate to the internal smoothing filter follow below.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
double GetLabel(int i)
Get the ith label value.
int GetNumberOfIterations()
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.
void GetValues(double *contourValues)
Fill a supplied list with label values.
vtkIdType GetNumberOfLabels()
Get the number of labels in the list of label values.
vtkMTimeType GetMTime() override
The modified time is also a function of the built in locator, smoothing filter, and label values.
void SetNumberOfContours(int number)
Set the number of labels to place into the list.
void AllSmoothingConstraintsOff()
If GenerateSmoothingStencils is on, this provides some control over each point's stencil creation.
void EdgeSmoothingConstraintOff()
If GenerateSmoothingStencils is on, this provides some control over each point's stencil creation.
double GetValue(int i)
Get the ith label value.
vtkSmartPointer< vtkIdTypeArray > PointsOfInterest
void SetNumberOfIterations(int n)
Convenience methods that delegate to the internal smoothing filter follow below.
void SetOutputMeshTypeToDefault()
Control the type of output mesh.
int GetNumberOfThreads()
Return the number of threads actually used during execution.
void SetOutputMeshTypeToPolygons()
Control the type of output mesh.
vtkSmartPointer< vtkConstrainedSmoothingFilter > Smoother
double * GetLabels()
Get a pointer to an array of labels.
virtual void SetSmoothingConstraints(unsigned char, unsigned char, unsigned char, unsigned char)
If GenerateSmoothingStencils is on, this provides some control over each point's stencil creation.
double * GetValues()
Get a pointer to an array of labels.
~vtkGeneralizedSurfaceNets3D() override=default
void SetRelaxationFactor(double f)
Convenience methods that delegate to the internal smoothing filter follow below.
static vtkGeneralizedSurfaceNets3D * New()
Standard methods for instantiation, type information, and printing.
void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
Generate numLabels equally spaced labels between the specified range.
vtkIdType GetNumberOfContours()
Get the number of labels in the list of label values.
vtkStaticPointLocator * GetLocator()
Retrieve the internal locator to manually configure it, for example specifying the number of points p...
void SetConstraintDistance(double d)
Convenience methods that delegate to the internal smoothing filter follow below.
void SetNumberOfLabels(int number)
Set the number of labels to place into the list.
void GetLabels(double *contourValues)
Fill a supplied list with label values.
void SetOutputMeshTypeToTriangles()
Control the type of output mesh.
vtkSmartPointer< vtkStaticPointLocator > Locator
virtual void SetOutputMeshType(int)
Control the type of output mesh.
void AllSmoothingConstraintsOn()
If GenerateSmoothingStencils is on, this provides some control over each point's stencil creation.
void GenerateValues(int numContours, double rangeStart, double rangeEnd)
Generate numLabels equally spaced labels between the specified range.
double GetConstraintDistance()
Convenience methods that delegate to the internal smoothing filter follow below.
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Hold a reference to a vtkObjectBase instance.
quickly locate points in 3-space
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:363
#define VTK_ID_MAX
Definition vtkType.h:367
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:318
#define VTK_INT_MAX
Definition vtkType.h:192
#define VTK_MARSHALAUTO