VTK  9.1.0
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 =========================================================================*/
55 #ifndef vtkImageMarchingCubes_h
56 #define vtkImageMarchingCubes_h
57 
58 #include "vtkFiltersGeneralModule.h" // For export macro
59 #include "vtkPolyDataAlgorithm.h"
60 
61 #include "vtkContourValues.h" // Needed for direct access to ContourValues
62 
63 class vtkCellArray;
64 class vtkFloatArray;
65 class vtkImageData;
66 class vtkPoints;
67 
68 class VTKFILTERSGENERAL_EXPORT vtkImageMarchingCubes : public vtkPolyDataAlgorithm
69 {
70 public:
73  void PrintSelf(ostream& os, vtkIndent indent) override;
74 
76 
79  void SetValue(int i, double value);
80  double GetValue(int i);
81  double* GetValues();
82  void GetValues(double* contourValues);
83  void SetNumberOfContours(int number);
84  vtkIdType GetNumberOfContours();
85  void GenerateValues(int numContours, double range[2]);
86  void GenerateValues(int numContours, double rangeStart, double rangeEnd);
88 
92  vtkMTimeType GetMTime() override;
93 
95 
98  vtkSetMacro(ComputeScalars, vtkTypeBool);
99  vtkGetMacro(ComputeScalars, vtkTypeBool);
100  vtkBooleanMacro(ComputeScalars, vtkTypeBool);
102 
104 
109  vtkSetMacro(ComputeNormals, vtkTypeBool);
110  vtkGetMacro(ComputeNormals, vtkTypeBool);
111  vtkBooleanMacro(ComputeNormals, vtkTypeBool);
113 
115 
122  vtkSetMacro(ComputeGradients, vtkTypeBool);
123  vtkGetMacro(ComputeGradients, vtkTypeBool);
124  vtkBooleanMacro(ComputeGradients, vtkTypeBool);
126 
127  // Should be protected, but the templated functions need these
132 
138 
139  vtkIdType GetLocatorPoint(int cellX, int cellY, int edge);
140  void AddLocatorPoint(int cellX, int cellY, int edge, vtkIdType ptId);
142 
144 
149  vtkSetMacro(InputMemoryLimit, vtkIdType);
150  vtkGetMacro(InputMemoryLimit, vtkIdType);
152 
153 protected:
156 
159 
161 
167 
171 
172  void March(vtkImageData* inData, int chunkMin, int chunkMax, int numContours, double* values);
173  void InitializeLocator(int min0, int max0, int min1, int max1);
175  vtkIdType* GetLocatorPointer(int cellX, int cellY, int edge);
176 
177 private:
179  void operator=(const vtkImageMarchingCubes&) = delete;
180 };
181 
186 inline void vtkImageMarchingCubes::SetValue(int i, double value)
187 {
188  this->ContourValues->SetValue(i, value);
189 }
190 
195 {
196  return this->ContourValues->GetValue(i);
197 }
198 
204 {
205  return this->ContourValues->GetValues();
206 }
207 
213 inline void vtkImageMarchingCubes::GetValues(double* contourValues)
214 {
215  this->ContourValues->GetValues(contourValues);
216 }
217 
224 {
225  this->ContourValues->SetNumberOfContours(number);
226 }
227 
232 {
233  return this->ContourValues->GetNumberOfContours();
234 }
235 
240 inline void vtkImageMarchingCubes::GenerateValues(int numContours, double range[2])
241 {
242  this->ContourValues->GenerateValues(numContours, range);
243 }
244 
250  int numContours, double rangeStart, double rangeEnd)
251 {
252  this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
253 }
254 
255 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:290
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
topologically and geometrically regular array of data
Definition: vtkImageData.h:157
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:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition: vtkPoints.h:143
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:332
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:287