VTK  9.2.6
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=========================================================================*/
39#ifndef vtkImageMarchingCubes_h
40#define vtkImageMarchingCubes_h
41
42#include "vtkFiltersGeneralModule.h" // For export macro
44
45#include "vtkContourValues.h" // Needed for direct access to ContourValues
46
47class vtkCellArray;
48class vtkFloatArray;
49class vtkImageData;
50class vtkPoints;
51
52class VTKFILTERSGENERAL_EXPORT vtkImageMarchingCubes : public vtkPolyDataAlgorithm
53{
54public:
57 void PrintSelf(ostream& os, vtkIndent indent) override;
58
60
63 void SetValue(int i, double value);
64 double GetValue(int i);
65 double* GetValues();
66 void GetValues(double* contourValues);
67 void SetNumberOfContours(int number);
68 vtkIdType GetNumberOfContours();
69 void GenerateValues(int numContours, double range[2]);
70 void GenerateValues(int numContours, double rangeStart, double rangeEnd);
72
77
79
82 vtkSetMacro(ComputeScalars, vtkTypeBool);
83 vtkGetMacro(ComputeScalars, vtkTypeBool);
84 vtkBooleanMacro(ComputeScalars, vtkTypeBool);
86
88
93 vtkSetMacro(ComputeNormals, vtkTypeBool);
94 vtkGetMacro(ComputeNormals, vtkTypeBool);
95 vtkBooleanMacro(ComputeNormals, vtkTypeBool);
97
99
106 vtkSetMacro(ComputeGradients, vtkTypeBool);
107 vtkGetMacro(ComputeGradients, vtkTypeBool);
108 vtkBooleanMacro(ComputeGradients, vtkTypeBool);
110
111 // Should be protected, but the templated functions need these
116
122
123 vtkIdType GetLocatorPoint(int cellX, int cellY, int edge);
124 void AddLocatorPoint(int cellX, int cellY, int edge, vtkIdType ptId);
126
128
133 vtkSetMacro(InputMemoryLimit, vtkIdType);
134 vtkGetMacro(InputMemoryLimit, vtkIdType);
136
137protected:
140
143
145
151
154 int FillInputPortInformation(int port, vtkInformation* info) override;
155
156 void March(vtkImageData* inData, int chunkMin, int chunkMax, int numContours, double* values);
157 void InitializeLocator(int min0, int max0, int min1, int max1);
159 vtkIdType* GetLocatorPointer(int cellX, int cellY, int edge);
160
161private:
163 void operator=(const vtkImageMarchingCubes&) = delete;
164};
165
170inline void vtkImageMarchingCubes::SetValue(int i, double value)
171{
172 this->ContourValues->SetValue(i, value);
173}
174
179{
180 return this->ContourValues->GetValue(i);
181}
182
188{
189 return this->ContourValues->GetValues();
190}
191
197inline void vtkImageMarchingCubes::GetValues(double* contourValues)
198{
199 this->ContourValues->GetValues(contourValues);
200}
201
208{
209 this->ContourValues->SetNumberOfContours(number);
210}
211
216{
217 return this->ContourValues->GetNumberOfContours();
218}
219
224inline void vtkImageMarchingCubes::GenerateValues(int numContours, double range[2])
225{
226 this->ContourValues->GenerateValues(numContours, range);
227}
228
234 int numContours, double rangeStart, double rangeEnd)
235{
236 this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
237}
238
239#endif
object to represent cell connectivity
Definition: vtkCellArray.h:187
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 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.
dynamic, self-adjusting array of float
Definition: vtkFloatArray.h:42
topologically and geometrically regular array of data
Definition: vtkImageData.h:54
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:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate 3D points
Definition: vtkPoints.h:40
Superclass for algorithms that produce only polydata as output.
int vtkTypeBool
Definition: vtkABI.h:69
int vtkIdType
Definition: vtkType.h:332
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:287