VTK  9.3.20240419
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
46 #ifndef vtkImageMarchingCubes_h
47 #define vtkImageMarchingCubes_h
48 
49 #include "vtkFiltersGeneralModule.h" // For export macro
50 #include "vtkPolyDataAlgorithm.h"
51 
52 #include "vtkContourValues.h" // Needed for direct access to ContourValues
53 
54 VTK_ABI_NAMESPACE_BEGIN
55 class vtkCellArray;
56 class vtkFloatArray;
57 class vtkImageData;
58 class vtkPoints;
59 
60 class VTKFILTERSGENERAL_EXPORT vtkImageMarchingCubes : public vtkPolyDataAlgorithm
61 {
62 public:
65  void PrintSelf(ostream& os, vtkIndent indent) override;
66 
68 
71  void SetValue(int i, double value);
72  double GetValue(int i);
73  double* GetValues();
74  void GetValues(double* contourValues);
75  void SetNumberOfContours(int number);
76  vtkIdType GetNumberOfContours();
77  void GenerateValues(int numContours, double range[2]);
78  void GenerateValues(int numContours, double rangeStart, double rangeEnd);
80 
84  vtkMTimeType GetMTime() override;
85 
87 
90  vtkSetMacro(ComputeScalars, vtkTypeBool);
91  vtkGetMacro(ComputeScalars, vtkTypeBool);
92  vtkBooleanMacro(ComputeScalars, vtkTypeBool);
94 
96 
101  vtkSetMacro(ComputeNormals, vtkTypeBool);
102  vtkGetMacro(ComputeNormals, vtkTypeBool);
103  vtkBooleanMacro(ComputeNormals, vtkTypeBool);
105 
107 
114  vtkSetMacro(ComputeGradients, vtkTypeBool);
115  vtkGetMacro(ComputeGradients, vtkTypeBool);
116  vtkBooleanMacro(ComputeGradients, vtkTypeBool);
118 
119  // Should be protected, but the templated functions need these
124 
130 
131  vtkIdType GetLocatorPoint(int cellX, int cellY, int edge);
132  void AddLocatorPoint(int cellX, int cellY, int edge, vtkIdType ptId);
134 
136 
141  vtkSetMacro(InputMemoryLimit, vtkIdType);
142  vtkGetMacro(InputMemoryLimit, vtkIdType);
144 
145 protected:
148 
151 
153 
159 
163 
164  void March(vtkImageData* inData, int chunkMin, int chunkMax, int numContours, double* values);
165  void InitializeLocator(int min0, int max0, int min1, int max1);
167  vtkIdType* GetLocatorPointer(int cellX, int cellY, int edge);
168 
169 private:
171  void operator=(const vtkImageMarchingCubes&) = delete;
172 };
173 
178 inline void vtkImageMarchingCubes::SetValue(int i, double value)
179 {
180  this->ContourValues->SetValue(i, value);
181 }
182 
187 {
188  return this->ContourValues->GetValue(i);
189 }
190 
196 {
197  return this->ContourValues->GetValues();
198 }
199 
205 inline void vtkImageMarchingCubes::GetValues(double* contourValues)
206 {
207  this->ContourValues->GetValues(contourValues);
208 }
209 
216 {
217  this->ContourValues->SetNumberOfContours(number);
218 }
219 
224 {
225  return this->ContourValues->GetNumberOfContours();
226 }
227 
232 inline void vtkImageMarchingCubes::GenerateValues(int numContours, double range[2])
233 {
234  this->ContourValues->GenerateValues(numContours, range);
235 }
236 
242  int numContours, double rangeStart, double rangeEnd)
243 {
244  this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
245 }
246 
247 VTK_ABI_NAMESPACE_END
248 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:286
helper object to manage setting and generating 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.
double * GetValues()
Return a pointer to a list of contour values.
dynamic, self-adjusting array of float
topologically and geometrically regular array of data
Definition: vtkImageData.h:156
generate isosurface(s) from volume/images
vtkIdType GetNumberOfContours()
Get the number of contours in the list of contour values.
static vtkImageMarchingCubes * New()
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void March(vtkImageData *inData, int chunkMin, int chunkMax, int numContours, double *values)
vtkIdType * GetLocatorPointer(int cellX, int cellY, int edge)
vtkContourValues * ContourValues
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void AddLocatorPoint(int cellX, int cellY, int edge, vtkIdType ptId)
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:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition: vtkPoints.h:139
Superclass for algorithms that produce only polydata as output.
@ info
Definition: vtkX3D.h:376
@ value
Definition: vtkX3D.h:220
@ port
Definition: vtkX3D.h:447
@ range
Definition: vtkX3D.h:238
int vtkTypeBool
Definition: vtkABI.h:64
int vtkIdType
Definition: vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270