VTK  9.3.0
vtkImageMarchingCubes.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
27#ifndef vtkImageMarchingCubes_h
28#define vtkImageMarchingCubes_h
29
30#include "vtkFiltersGeneralModule.h" // For export macro
32
33#include "vtkContourValues.h" // Needed for direct access to ContourValues
34
35VTK_ABI_NAMESPACE_BEGIN
36class vtkCellArray;
37class vtkFloatArray;
38class vtkImageData;
39class vtkPoints;
40
41class VTKFILTERSGENERAL_EXPORT vtkImageMarchingCubes : public vtkPolyDataAlgorithm
42{
43public:
46 void PrintSelf(ostream& os, vtkIndent indent) override;
47
49
52 void SetValue(int i, double value);
53 double GetValue(int i);
54 double* GetValues();
55 void GetValues(double* contourValues);
56 void SetNumberOfContours(int number);
57 vtkIdType GetNumberOfContours();
58 void GenerateValues(int numContours, double range[2]);
59 void GenerateValues(int numContours, double rangeStart, double rangeEnd);
61
66
68
71 vtkSetMacro(ComputeScalars, vtkTypeBool);
72 vtkGetMacro(ComputeScalars, vtkTypeBool);
73 vtkBooleanMacro(ComputeScalars, vtkTypeBool);
75
77
82 vtkSetMacro(ComputeNormals, vtkTypeBool);
83 vtkGetMacro(ComputeNormals, vtkTypeBool);
84 vtkBooleanMacro(ComputeNormals, vtkTypeBool);
86
88
95 vtkSetMacro(ComputeGradients, vtkTypeBool);
96 vtkGetMacro(ComputeGradients, vtkTypeBool);
97 vtkBooleanMacro(ComputeGradients, vtkTypeBool);
99
100 // Should be protected, but the templated functions need these
105
111
112 vtkIdType GetLocatorPoint(int cellX, int cellY, int edge);
113 void AddLocatorPoint(int cellX, int cellY, int edge, vtkIdType ptId);
115
117
122 vtkSetMacro(InputMemoryLimit, vtkIdType);
123 vtkGetMacro(InputMemoryLimit, vtkIdType);
125
126protected:
129
132
134
140
143 int FillInputPortInformation(int port, vtkInformation* info) override;
144
145 void March(vtkImageData* inData, int chunkMin, int chunkMax, int numContours, double* values);
146 void InitializeLocator(int min0, int max0, int min1, int max1);
148 vtkIdType* GetLocatorPointer(int cellX, int cellY, int edge);
149
150private:
152 void operator=(const vtkImageMarchingCubes&) = delete;
153};
154
159inline void vtkImageMarchingCubes::SetValue(int i, double value)
160{
161 this->ContourValues->SetValue(i, value);
162}
163
168{
169 return this->ContourValues->GetValue(i);
170}
171
177{
178 return this->ContourValues->GetValues();
179}
180
186inline void vtkImageMarchingCubes::GetValues(double* contourValues)
187{
188 this->ContourValues->GetValues(contourValues);
189}
190
197{
198 this->ContourValues->SetNumberOfContours(number);
199}
200
205{
206 return this->ContourValues->GetNumberOfContours();
207}
208
213inline void vtkImageMarchingCubes::GenerateValues(int numContours, double range[2])
214{
215 this->ContourValues->GenerateValues(numContours, range);
216}
217
223 int numContours, double rangeStart, double rangeEnd)
224{
225 this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
226}
227
228VTK_ABI_NAMESPACE_END
229#endif
object to represent cell connectivity
helper object to manage setting and generating contour values
double * GetValues()
Return a pointer to a list of contour values.
int GetNumberOfContours()
Return the number of contours in the.
void SetNumberOfContours(int number)
Set the number of contours to place into the list.
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
void SetValue(int i, double value)
Set the ith contour value.
double GetValue(int i)
Get the ith contour value.
dynamic, self-adjusting array of float
topologically and geometrically regular array of data
generate isosurface(s) from volume/images
vtkIdType GetNumberOfContours()
Get the number of contours in the list of contour values.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void March(vtkImageData *inData, int chunkMin, int chunkMax, int numContours, double *values)
vtkContourValues * ContourValues
static vtkImageMarchingCubes * New()
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void AddLocatorPoint(int cellX, int cellY, int edge, vtkIdType ptId)
vtkIdType * GetLocatorPointer(int cellX, int cellY, int edge)
double GetValue(int i)
Get the ith contour value.
vtkMTimeType GetMTime() override
Because we delegate to vtkContourValues & refer to vtkImplicitFunction.
double * GetValues()
Get a pointer to an array of contour values.
~vtkImageMarchingCubes() override
void InitializeLocator(int min0, int max0, int min1, int max1)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void SetValue(int i, double value)
Methods to set contour values.
void SetNumberOfContours(int number)
Set the number of contours to place into the list.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
vtkIdType GetLocatorPoint(int cellX, int cellY, int edge)
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition vtkPoints.h:29
Superclass for algorithms that produce only polydata as output.
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:270