VTK  9.4.20241218
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
150#ifndef vtkSurfaceNets3D_h
151#define vtkSurfaceNets3D_h
152
153#include "vtkConstrainedSmoothingFilter.h" // Perform mesh smoothing
154#include "vtkContourValues.h" // Needed for direct access to ContourValues
155#include "vtkFiltersCoreModule.h" // For export macro
156#include "vtkPolyData.h" // To support data caching
157#include "vtkPolyDataAlgorithm.h"
158
159#include <vector> // For selected seeds
160
161VTK_ABI_NAMESPACE_BEGIN
162
163class vtkImageData;
164
165class VTKFILTERSCORE_EXPORT vtkSurfaceNets3D : public vtkPolyDataAlgorithm
166{
167public:
169
175 void PrintSelf(ostream& os, vtkIndent indent) override;
177
183
185
198 void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
199 void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
201
203
206 double GetValue(int i) { return this->Labels->GetValue(i); }
207 double GetLabel(int i) { return this->Labels->GetValue(i); }
209
211
215 double* GetValues() { return this->Labels->GetValues(); }
216 double* GetLabels() { return this->Labels->GetValues(); }
218
220
225 void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
226 void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
228
230
237 void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
238 void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
240
242
245 vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
246 vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
248
250
254 void GenerateLabels(int numLabels, double range[2])
255 {
256 this->Labels->GenerateValues(numLabels, range);
257 }
258 void GenerateValues(int numContours, double range[2])
259 {
260 this->Labels->GenerateValues(numContours, range);
261 }
262 void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
263 {
264 this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
265 }
266 void GenerateValues(int numContours, double rangeStart, double rangeEnd)
267 {
268 this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
269 }
271
273
283 vtkSetMacro(BackgroundLabel, double);
284 vtkGetMacro(BackgroundLabel, double);
286
288
292 vtkSetMacro(ArrayComponent, int);
293 vtkGetMacro(ArrayComponent, int);
295
300 {
301 MESH_TYPE_DEFAULT = 0,
303 MESH_TYPE_QUADS
304 };
305
307
317 vtkSetClampMacro(OutputMeshType, int, MESH_TYPE_DEFAULT, MESH_TYPE_QUADS);
318 vtkGetMacro(OutputMeshType, int);
319 void SetOutputMeshTypeToDefault() { this->SetOutputMeshType(MESH_TYPE_DEFAULT); }
320 void SetOutputMeshTypeToTriangles() { this->SetOutputMeshType(MESH_TYPE_TRIANGLES); }
321 void SetOutputMeshTypeToQuads() { this->SetOutputMeshType(MESH_TYPE_QUADS); }
323
324 // The following code is used to control the smoothing process. Internally
325 // there is a vtkConstrainedSmoothingFilter that can be directly
326 // manipulated. In addition, methods that delegate to this filter are also
327 // provided.
328
330
335 vtkSetMacro(Smoothing, bool);
336 vtkGetMacro(Smoothing, bool);
337 vtkBooleanMacro(Smoothing, bool);
339
341
346 void SetNumberOfIterations(int n) { this->Smoother->SetNumberOfIterations(n); }
347 int GetNumberOfIterations() { return this->Smoother->GetNumberOfIterations(); }
348 void SetRelaxationFactor(double f) { this->Smoother->SetRelaxationFactor(f); }
349 double GetRelaxationFactor() { return this->Smoother->GetRelaxationFactor(); }
350 void SetConstraintDistance(double d) { this->Smoother->SetConstraintDistance(d); }
351 double GetConstraintDistance() { return this->Smoother->GetConstraintDistance(); }
352 void SetConstraintBox(double sx, double sy, double sz)
353 {
354 this->Smoother->SetConstraintBox(sx, sy, sz);
355 }
356 void SetConstraintBox(double s[3]) { this->Smoother->SetConstraintBox(s); }
357 double* GetConstraintBox() VTK_SIZEHINT(3) { return this->Smoother->GetConstraintBox(); }
358 void GetConstraintBox(double s[3]) { this->Smoother->GetConstraintBox(s); }
360 {
361 this->Smoother->SetConstraintStrategyToConstraintDistance();
362 }
364 {
365 this->Smoother->SetConstraintStrategyToConstraintBox();
366 }
367 int GetConstraintStrategy() { return this->Smoother->GetConstraintStrategy(); }
369
371
383 vtkSetMacro(AutomaticSmoothingConstraints, bool);
384 vtkGetMacro(AutomaticSmoothingConstraints, bool);
385 vtkBooleanMacro(AutomaticSmoothingConstraints, bool);
386 vtkSetClampMacro(ConstraintScale, double, 0, 100);
387 vtkGetMacro(ConstraintScale, double);
389
391
397 vtkSetMacro(OptimizedSmoothingStencils, bool);
398 vtkGetMacro(OptimizedSmoothingStencils, bool);
399 vtkBooleanMacro(OptimizedSmoothingStencils, bool);
401
403
412
413 // The following code is used to control what is produced for output.
414
430 {
431 OUTPUT_STYLE_DEFAULT = 0,
433 OUTPUT_STYLE_SELECTED
434 };
435
437
449 vtkSetClampMacro(OutputStyle, int, OUTPUT_STYLE_DEFAULT, OUTPUT_STYLE_SELECTED);
450 vtkGetMacro(OutputStyle, int);
451 void SetOutputStyleToDefault() { this->SetOutputStyle(OUTPUT_STYLE_DEFAULT); }
452 void SetOutputStyleToBoundary() { this->SetOutputStyle(OUTPUT_STYLE_BOUNDARY); }
453 void SetOutputStyleToSelected() { this->SetOutputStyle(OUTPUT_STYLE_SELECTED); }
455
457
462 void AddSelectedLabel(double label);
463 void DeleteSelectedLabel(double label);
465 double GetSelectedLabel(vtkIdType ithLabel);
467
472 {
473 TRIANGULATION_GREEDY = 0,
475 TRIANGULATION_MIN_AREA
476 };
477
479
489 vtkSetClampMacro(TriangulationStrategy, int, TRIANGULATION_GREEDY, TRIANGULATION_MIN_AREA);
490 vtkGetMacro(TriangulationStrategy, int);
491 void SetTriangulationStrategyToGreedy() { this->SetTriangulationStrategy(TRIANGULATION_GREEDY); }
493 {
494 this->SetTriangulationStrategy(TRIANGULATION_MIN_EDGE);
495 }
497 {
498 this->SetTriangulationStrategy(TRIANGULATION_MIN_AREA);
499 }
501
503
513 vtkSetMacro(DataCaching, bool);
514 vtkGetMacro(DataCaching, bool);
515 vtkBooleanMacro(DataCaching, bool);
517
518protected:
520 ~vtkSurfaceNets3D() override = default;
521
522 // Support visualization pipeline operations.
524 int FillInputPortInformation(int port, vtkInformation* info) override;
525
526 // Support the contouring operation.
532
533 // Support smoothing.
539
540 // Support data caching of the extracted surface nets. This is used to
541 // avoid repeated surface extraction when only smoothing filter
542 // parameters are modified.
549
550 // Support output style
552 std::vector<double> SelectedLabels;
554
555 // Support triangulation strategy
557
558private:
559 vtkSurfaceNets3D(const vtkSurfaceNets3D&) = delete;
560 void operator=(const vtkSurfaceNets3D&) = delete;
561};
562
563VTK_ABI_NAMESPACE_END
564#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.
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Hold a reference to a vtkObjectBase instance.
generate smoothed isocontours from segmented 3D image data (i.e., "label maps")
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_...
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.
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_...
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.
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.
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
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
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270
#define VTK_SIZEHINT(...)