VTK  9.0.3
vtkImageMarchingCubes.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkImageMarchingCubes.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
36 #ifndef vtkImageMarchingCubes_h
37 #define vtkImageMarchingCubes_h
38 
39 #include "vtkFiltersGeneralModule.h" // For export macro
40 #include "vtkPolyDataAlgorithm.h"
41 
42 #include "vtkContourValues.h" // Needed for direct access to ContourValues
43 
44 class vtkCellArray;
45 class vtkFloatArray;
46 class vtkImageData;
47 class vtkPoints;
48 
49 class VTKFILTERSGENERAL_EXPORT vtkImageMarchingCubes : public vtkPolyDataAlgorithm
50 {
51 public:
54  void PrintSelf(ostream& os, vtkIndent indent) override;
55 
57 
60  void SetValue(int i, double value);
61  double GetValue(int i);
62  double* GetValues();
63  void GetValues(double* contourValues);
64  void SetNumberOfContours(int number);
65  vtkIdType GetNumberOfContours();
66  void GenerateValues(int numContours, double range[2]);
67  void GenerateValues(int numContours, double rangeStart, double rangeEnd);
69 
73  vtkMTimeType GetMTime() override;
74 
76 
79  vtkSetMacro(ComputeScalars, vtkTypeBool);
80  vtkGetMacro(ComputeScalars, vtkTypeBool);
81  vtkBooleanMacro(ComputeScalars, vtkTypeBool);
83 
85 
90  vtkSetMacro(ComputeNormals, vtkTypeBool);
91  vtkGetMacro(ComputeNormals, vtkTypeBool);
92  vtkBooleanMacro(ComputeNormals, vtkTypeBool);
94 
96 
103  vtkSetMacro(ComputeGradients, vtkTypeBool);
104  vtkGetMacro(ComputeGradients, vtkTypeBool);
105  vtkBooleanMacro(ComputeGradients, vtkTypeBool);
107 
108  // Should be protected, but the templated functions need these
113 
119 
120  vtkIdType GetLocatorPoint(int cellX, int cellY, int edge);
121  void AddLocatorPoint(int cellX, int cellY, int edge, vtkIdType ptId);
123 
125 
130  vtkSetMacro(InputMemoryLimit, vtkIdType);
131  vtkGetMacro(InputMemoryLimit, vtkIdType);
133 
134 protected:
137 
140 
142 
148 
152 
153  void March(vtkImageData* inData, int chunkMin, int chunkMax, int numContours, double* values);
154  void InitializeLocator(int min0, int max0, int min1, int max1);
156  vtkIdType* GetLocatorPointer(int cellX, int cellY, int edge);
157 
158 private:
160  void operator=(const vtkImageMarchingCubes&) = delete;
161 };
162 
167 inline void vtkImageMarchingCubes::SetValue(int i, double value)
168 {
169  this->ContourValues->SetValue(i, value);
170 }
171 
176 {
177  return this->ContourValues->GetValue(i);
178 }
179 
185 {
186  return this->ContourValues->GetValues();
187 }
188 
194 inline void vtkImageMarchingCubes::GetValues(double* contourValues)
195 {
196  this->ContourValues->GetValues(contourValues);
197 }
198 
205 {
206  this->ContourValues->SetNumberOfContours(number);
207 }
208 
213 {
214  return this->ContourValues->GetNumberOfContours();
215 }
216 
221 inline void vtkImageMarchingCubes::GenerateValues(int numContours, double range[2])
222 {
223  this->ContourValues->GenerateValues(numContours, range);
224 }
225 
231  int numContours, double rangeStart, double rangeEnd)
232 {
233  this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
234 }
235 
236 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:180
helper object to manage setting and generating contour values
int GetNumberOfContours()
Return the number of contours in the.
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
void SetNumberOfContours(const int number)
Set the number of contours to place into the list.
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
Definition: vtkFloatArray.h:36
topologically and geometrically regular array of data
Definition: vtkImageData.h:42
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:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition: vtkPoints.h:34
Superclass for algorithms that produce only polydata as output.
@ info
Definition: vtkX3D.h:382
@ value
Definition: vtkX3D.h:226
@ port
Definition: vtkX3D.h:453
@ range
Definition: vtkX3D.h:244
int vtkTypeBool
Definition: vtkABI.h:69
int vtkIdType
Definition: vtkType.h:338
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:293