VTK  9.5.20250905
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#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
159
160#include <vector> // For selected seeds
161
162VTK_ABI_NAMESPACE_BEGIN
163
164class vtkImageData;
165
166class VTKFILTERSCORE_EXPORT VTK_MARSHALAUTO vtkSurfaceNets3D : public vtkPolyDataAlgorithm
167{
168public:
170
176 void PrintSelf(ostream& os, vtkIndent indent) override;
178
184
186
199 void SetValue(int i, double value) { this->Labels->SetValue(i, value); }
200 void SetLabel(int i, double value) { this->Labels->SetValue(i, value); }
202
204
207 double GetValue(int i) { return this->Labels->GetValue(i); }
208 double GetLabel(int i) { return this->Labels->GetValue(i); }
210
212
216 double* GetValues() { return this->Labels->GetValues(); }
217 double* GetLabels() { return this->Labels->GetValues(); }
219
221
226 void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); }
227 void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); }
229
231
238 void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); }
239 void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); }
241
243
246 vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); }
247 vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); }
249
251
255 void GenerateLabels(int numLabels, double range[2])
256 {
257 this->Labels->GenerateValues(numLabels, range);
258 }
259 void GenerateValues(int numContours, double range[2])
260 {
261 this->Labels->GenerateValues(numContours, range);
262 }
263 void GenerateLabels(int numLabels, double rangeStart, double rangeEnd)
264 {
265 this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd);
266 }
267 void GenerateValues(int numContours, double rangeStart, double rangeEnd)
268 {
269 this->Labels->GenerateValues(numContours, rangeStart, rangeEnd);
270 }
272
274
284 vtkSetMacro(BackgroundLabel, double);
285 vtkGetMacro(BackgroundLabel, double);
287
289
293 vtkSetMacro(ArrayComponent, int);
294 vtkGetMacro(ArrayComponent, int);
296
301 {
302 MESH_TYPE_DEFAULT = 0,
304 MESH_TYPE_QUADS
305 };
306
308
318 vtkSetClampMacro(OutputMeshType, int, MESH_TYPE_DEFAULT, MESH_TYPE_QUADS);
319 vtkGetMacro(OutputMeshType, int);
320 void SetOutputMeshTypeToDefault() { this->SetOutputMeshType(MESH_TYPE_DEFAULT); }
321 void SetOutputMeshTypeToTriangles() { this->SetOutputMeshType(MESH_TYPE_TRIANGLES); }
322 void SetOutputMeshTypeToQuads() { this->SetOutputMeshType(MESH_TYPE_QUADS); }
324
325 // The following code is used to control the smoothing process. Internally
326 // there is a vtkConstrainedSmoothingFilter that can be directly
327 // manipulated. In addition, methods that delegate to this filter are also
328 // provided.
329
331
336 vtkSetMacro(Smoothing, bool);
337 vtkGetMacro(Smoothing, bool);
338 vtkBooleanMacro(Smoothing, bool);
340
342
347 void SetNumberOfIterations(int n) { this->Smoother->SetNumberOfIterations(n); }
348 int GetNumberOfIterations() { return this->Smoother->GetNumberOfIterations(); }
349 void SetRelaxationFactor(double f) { this->Smoother->SetRelaxationFactor(f); }
350 double GetRelaxationFactor() { return this->Smoother->GetRelaxationFactor(); }
351 void SetConstraintDistance(double d) { this->Smoother->SetConstraintDistance(d); }
352 double GetConstraintDistance() { return this->Smoother->GetConstraintDistance(); }
353 void SetConstraintBox(double sx, double sy, double sz)
354 {
355 this->Smoother->SetConstraintBox(sx, sy, sz);
356 }
357 void SetConstraintBox(double s[3]) { this->Smoother->SetConstraintBox(s); }
358 double* GetConstraintBox() VTK_SIZEHINT(3) { return this->Smoother->GetConstraintBox(); }
359 void GetConstraintBox(double s[3]) { this->Smoother->GetConstraintBox(s); }
361 {
362 this->Smoother->SetConstraintStrategyToConstraintDistance();
363 }
365 {
366 this->Smoother->SetConstraintStrategyToConstraintBox();
367 }
368 int GetConstraintStrategy() { return this->Smoother->GetConstraintStrategy(); }
370
372
384 vtkSetMacro(AutomaticSmoothingConstraints, bool);
385 vtkGetMacro(AutomaticSmoothingConstraints, bool);
386 vtkBooleanMacro(AutomaticSmoothingConstraints, bool);
387 vtkSetClampMacro(ConstraintScale, double, 0, 100);
388 vtkGetMacro(ConstraintScale, double);
390
392
398 vtkSetMacro(OptimizedSmoothingStencils, bool);
399 vtkGetMacro(OptimizedSmoothingStencils, bool);
400 vtkBooleanMacro(OptimizedSmoothingStencils, bool);
402
404
411 vtkGetSmartPointerMacro(Smoother, vtkConstrainedSmoothingFilter);
413
414 // The following code is used to control what is produced for output.
415
431 {
432 OUTPUT_STYLE_DEFAULT = 0,
434 OUTPUT_STYLE_SELECTED
435 };
436
438
450 vtkSetClampMacro(OutputStyle, int, OUTPUT_STYLE_DEFAULT, OUTPUT_STYLE_SELECTED);
451 vtkGetMacro(OutputStyle, int);
452 void SetOutputStyleToDefault() { this->SetOutputStyle(OUTPUT_STYLE_DEFAULT); }
453 void SetOutputStyleToBoundary() { this->SetOutputStyle(OUTPUT_STYLE_BOUNDARY); }
454 void SetOutputStyleToSelected() { this->SetOutputStyle(OUTPUT_STYLE_SELECTED); }
456
458
463 void AddSelectedLabel(double label);
464 void DeleteSelectedLabel(double label);
466 double GetSelectedLabel(vtkIdType ithLabel);
468
473 {
474 TRIANGULATION_GREEDY = 0,
476 TRIANGULATION_MIN_AREA
477 };
478
480
490 vtkSetClampMacro(TriangulationStrategy, int, TRIANGULATION_GREEDY, TRIANGULATION_MIN_AREA);
491 vtkGetMacro(TriangulationStrategy, int);
492 void SetTriangulationStrategyToGreedy() { this->SetTriangulationStrategy(TRIANGULATION_GREEDY); }
494 {
495 this->SetTriangulationStrategy(TRIANGULATION_MIN_EDGE);
496 }
498 {
499 this->SetTriangulationStrategy(TRIANGULATION_MIN_AREA);
500 }
502
504
514 vtkSetMacro(DataCaching, bool);
515 vtkGetMacro(DataCaching, bool);
516 vtkBooleanMacro(DataCaching, bool);
518
519protected:
521 ~vtkSurfaceNets3D() override = default;
522
523 // Support visualization pipeline operations.
525 int FillInputPortInformation(int port, vtkInformation* info) override;
526
527 // Support the contouring operation.
533
534 // Support smoothing.
540
541 // Support data caching of the extracted surface nets. This is used to
542 // avoid repeated surface extraction when only smoothing filter
543 // parameters are modified.
550
551 // Support output style
553 std::vector<double> SelectedLabels;
555
556 // Support triangulation strategy
558
559private:
560 vtkSurfaceNets3D(const vtkSurfaceNets3D&) = delete;
561 void operator=(const vtkSurfaceNets3D&) = delete;
562};
563
564VTK_ABI_NAMESPACE_END
565#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_...
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:332
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287
#define VTK_SIZEHINT(...)
#define VTK_MARSHALAUTO