VTK
vtkHyperOctreeContourFilter.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkHyperOctreeContourFilter.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=========================================================================*/
51#ifndef vtkHyperOctreeContourFilter_h
52#define vtkHyperOctreeContourFilter_h
53
54#include "vtkFiltersHyperTreeModule.h" // For export macro
56
57#include "vtkContourValues.h" // Needed for inline methods
58#include "vtkCutter.h" // for VTK_SORT_BY_VALUE
59
61class vtkHyperOctree;
63class vtkTetra;
65
68class vtkIdTypeArray;
69class vtkHyperOctreeContourPointsGrabber;
70class vtkBitArray;
71
72class VTKFILTERSHYPERTREE_EXPORT vtkHyperOctreeContourFilter : public vtkPolyDataAlgorithm
73{
74public:
76 void PrintSelf(ostream& os, vtkIndent indent);
77
83
92 void SetValue(int i, double value)
93 {
94 this->ContourValues->SetValue(i,value);
95 }
96
100 double GetValue(int i)
101 {
102 return this->ContourValues->GetValue(i);
103 }
104
109 double *GetValues()
110 {
111 return this->ContourValues->GetValues();
112 }
113
119 void GetValues(double *contourValues)
120 {
121 this->ContourValues->GetValues(contourValues);
122 }
123
129 void SetNumberOfContours(int number)
130 {
131 this->ContourValues->SetNumberOfContours(number);
132 }
133
138 {
139 return this->ContourValues->GetNumberOfContours();
140 }
141
146 void GenerateValues(int numContours, double range[2])
147 {
148 this->ContourValues->GenerateValues(numContours, range);
149 }
150
155 void GenerateValues(int numContours, double
156 rangeStart, double rangeEnd)
157 {
158 this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);
159 }
160
165
167
172 vtkGetObjectMacro(Locator,vtkIncrementalPointLocator);
174
180
181protected:
184
185 virtual int RequestData(vtkInformation* request,
186 vtkInformationVector** inputVector,
187 vtkInformationVector* outputVector);
192
197
201 double ComputePointValue(int ptIndices[3]);
202
204
207
208 vtkIdList *CellPts; // for 2D case
209
212
216
222
223 vtkHyperOctreeCursor *Sibling; // to avoid allocation in the loop
224
225
229
231
234
235 vtkIdType CellTypeCounter[65536]; // up-to-65536 points per octant
237 vtkIdType TemplateCounter; // record the number of octants that succceed
238 // to use the template triangulator
239
241 vtkHyperOctreeContourPointsGrabber *Grabber;
242
245 int Iter; // iterate over contour values in case of VTK_SORT_BY_CELL
246
248 double LeftValue;
249 double LeftCoord;
250
251 friend class vtkHyperOctreeContourPointsGrabber;
252
253private:
255 void operator=(const vtkHyperOctreeContourFilter&) VTK_DELETE_FUNCTION;
256};
257#endif
dynamic, self-adjusting array of bits
Definition: vtkBitArray.h:37
object to represent cell connectivity
Definition: vtkCellArray.h:51
represent and manipulate cell attribute data
Definition: vtkCellData.h:39
helper object to manage setting and generating contour values
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
dynamic, self-adjusting array of double
generate isosurfaces/isolines from scalar values
vtkIncrementalPointLocator * Locator
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
void GetValues(double *contourValues)
Fill a supplied list with contour values.
static vtkHyperOctreeContourFilter * New()
Construct object with initial range (0,1) and single contour value of 0.0.
double ComputePointValue(int ptIndices[3])
(i,j,k) are point coordinates at last level
void SetLocator(vtkIncrementalPointLocator *locator)
Set / get a spatial locator for merging points.
double GetValue(int i)
Get the ith contour value.
void SetNumberOfContours(int number)
Set the number of contours to place into the list.
virtual int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
void CreateDefaultLocator()
Create default locator.
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
void ContourNode()
Do the recursive contour of the node pointed by Cursor.
vtkHyperOctreeContourPointsGrabber * Grabber
double * GetValues()
Get a pointer to an array of contour values.
int GetNumberOfContours()
Get the number of contours in the list of contour values.
virtual int FillInputPortInformation(int port, vtkInformation *info)
Fill the input port information objects for this algorithm.
void SetValue(int i, double value)
Methods to set / get contour values.
void GenerateValues(int numContours, double rangeStart, double rangeEnd)
Generate numContours equally spaced contour values between specified range.
vtkMTimeType GetMTime()
Modified GetMTime Because we delegate to vtkContourValues.
void PrintSelf(ostream &os, vtkIndent indent)
Methods invoked by print to print information about the object including superclasses.
Objects that can traverse hyperoctree nodes.
A dataset structured as a tree where each node has exactly 2^n children.
list of point or cell ids
Definition: vtkIdList.h:37
dynamic, self-adjusting array of vtkIdType
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
cell represents a 1D line
Definition: vtkLine.h:36
helper class to generate triangulations
represent and manipulate point attribute data
Definition: vtkPointData.h:38
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:86
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:46
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:48
dynamic, self-adjusting array of unsigned char
dataset represents arbitrary combinations of all possible cell types
@ info
Definition: vtkX3D.h:376
@ value
Definition: vtkX3D.h:220
@ port
Definition: vtkX3D.h:447
@ range
Definition: vtkX3D.h:238
int vtkIdType
Definition: vtkType.h:287
vtkTypeUInt64 vtkMTimeType
Definition: vtkType.h:248