VTK
vtkHyperOctreeCutter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkHyperOctreeCutter.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 =========================================================================*/
49 #ifndef vtkHyperOctreeCutter_h
50 #define vtkHyperOctreeCutter_h
51 
52 #include "vtkFiltersHyperTreeModule.h" // For export macro
53 #include "vtkPolyDataAlgorithm.h"
54 
55 #include "vtkContourValues.h" // Needed for inline methods
56 
57 #include "vtkCutter.h" // for VTK_SORT_BY_VALUE and VTK_SORT_BY_CELL
58 
59 //#define VTK_SORT_BY_VALUE 0
60 //#define VTK_SORT_BY_CELL 1
61 // This does not really belong here, ut it is for a temporary
62 // fix until this filter can be converted to geernate unstructured grids.
63 //#define VTK_NUMBER_OF_CELL_TYPES 68
64 
67 class vtkHyperOctree;
70 class vtkTetra;
73 
74 class VTKFILTERSHYPERTREE_EXPORT vtkHyperOctreeCutter : public vtkPolyDataAlgorithm
75 {
76 public:
78  void PrintSelf(ostream& os, vtkIndent indent);
79 
84  static vtkHyperOctreeCutter *New();
85 
90  void SetValue(int i, double value)
91  {this->ContourValues->SetValue(i,value);}
92 
96  double GetValue(int i)
97  {return this->ContourValues->GetValue(i);}
98 
103  double *GetValues()
104  {return this->ContourValues->GetValues();}
105 
111  void GetValues(double *contourValues)
112  {this->ContourValues->GetValues(contourValues);}
113 
119  void SetNumberOfContours(int number)
120  {this->ContourValues->SetNumberOfContours(number);}
121 
126  {return this->ContourValues->GetNumberOfContours();}
127 
132  void GenerateValues(int numContours, double range[2])
133  {this->ContourValues->GenerateValues(numContours, range);}
134 
139  void GenerateValues(int numContours, double rangeStart, double rangeEnd)
140  {this->ContourValues->GenerateValues(numContours, rangeStart, rangeEnd);}
141 
147 
149 
152  virtual void SetCutFunction(vtkImplicitFunction*);
153  vtkGetObjectMacro(CutFunction,vtkImplicitFunction);
155 
157 
162  vtkSetMacro(GenerateCutScalars,int);
163  vtkGetMacro(GenerateCutScalars,int);
164  vtkBooleanMacro(GenerateCutScalars,int);
166 
168 
172  void SetLocator(vtkIncrementalPointLocator *locator);
173  vtkGetObjectMacro(Locator,vtkIncrementalPointLocator);
175 
177 
192  vtkSetClampMacro(SortBy,int,VTK_SORT_BY_VALUE,VTK_SORT_BY_CELL);
193  vtkGetMacro(SortBy,int);
195  {this->SetSortBy(VTK_SORT_BY_VALUE);}
197  {this->SetSortBy(VTK_SORT_BY_CELL);}
199 
201 
204  const char *GetSortByAsString()
205  {
206  if ( this->SortBy == VTK_SORT_BY_VALUE )
207  {
208  return "SortByValue";
209  }
210  else
211  {
212  return "SortByCell";
213  }
214  }
216 
221  void CreateDefaultLocator();
222 
223 protected:
226 
230 
236  void CutNode(vtkHyperOctreeCursor *cursor,
237  int level,
238  double bounds[6]);
239 
241 
242 
244  int SortBy;
247 
250 
251 
255 
260  vtkHyperOctreeCursor *Sibling; // to avoid allocation in the loop
261 
262  int Iter; // iterate over contour values in case of VTK_SORT_BY_CELL
263 
264 
268 
271 
272  vtkIdType CellTypeCounter[65536]; // up-to-65536 points per octant
274  vtkIdType TemplateCounter; // record the number of octants that succceed
275  // to use the template triangulator
276 
277  // in VTK_SORT_BY_VALUE case, rejection test need to combine all values.
278  int *AllLess;
281 
282 private:
283  vtkHyperOctreeCutter(const vtkHyperOctreeCutter&) VTK_DELETE_FUNCTION;
284  void operator=(const vtkHyperOctreeCutter&) VTK_DELETE_FUNCTION;
285 };
286 
287 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:39
vtkHyperOctreeCutter::OutPD
vtkPointData * OutPD
Definition: vtkHyperOctreeCutter.h:258
vtkHyperOctreeCutter::NewPolys
vtkCellArray * NewPolys
Definition: vtkHyperOctreeCutter.h:254
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:37
vtkX3D::value
@ value
Definition: vtkX3D.h:220
vtkIdType
int vtkIdType
Definition: vtkType.h:287
vtkPolygon
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:45
vtkHyperOctreeCutter::NewVerts
vtkCellArray * NewVerts
Definition: vtkHyperOctreeCutter.h:252
vtkDataSetAttributes
represent and manipulate attribute data in a dataset
Definition: vtkDataSetAttributes.h:58
vtkInformationVector
Store zero or more vtkInformation instances.
Definition: vtkInformationVector.h:41
vtkHyperOctreeCutter::Tetra
vtkTetra * Tetra
Definition: vtkHyperOctreeCutter.h:266
VTK_SORT_BY_VALUE
#define VTK_SORT_BY_VALUE
Definition: vtkCutter.h:59
vtkHyperOctreeCutter::TetScalars
vtkDoubleArray * TetScalars
Definition: vtkHyperOctreeCutter.h:267
vtkHyperOctreeCutter::Output
vtkPolyData * Output
Definition: vtkHyperOctreeCutter.h:249
vtkHyperOctreeCutter::Grabber
vtkHyperOctreeClipCutPointsGrabber * Grabber
Definition: vtkHyperOctreeCutter.h:280
vtkPolyDataAlgorithm::RequestUpdateExtent
virtual int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
vtkHyperOctreeCursor
Objects that can traverse hyperoctree nodes.
Definition: vtkHyperOctreeCursor.h:51
vtkX3D::range
@ range
Definition: vtkX3D.h:238
vtkHyperOctreeCutter::CutFunction
vtkImplicitFunction * CutFunction
Definition: vtkHyperOctreeCutter.h:240
vtkPolyDataAlgorithm::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkHyperOctreeClipCutPointsGrabber
A concrete implementation of vtkHyperOctreePointsGrabber used by vtkClipHyperOctree and vtkHyperOctre...
Definition: vtkHyperOctreeClipCutPointsGrabber.h:36
vtkOrderedTriangulator
helper class to generate triangulations
Definition: vtkOrderedTriangulator.h:117
vtkHyperOctreeCutter::SetSortByToSortByValue
void SetSortByToSortByValue()
Definition: vtkHyperOctreeCutter.h:194
vtkMTimeType
vtkTypeUInt64 vtkMTimeType
Definition: vtkType.h:248
vtkPolyDataAlgorithm.h
vtkHyperOctreeCutter::Locator
vtkIncrementalPointLocator * Locator
Definition: vtkHyperOctreeCutter.h:243
vtkImplicitFunction
abstract interface for implicit functions
Definition: vtkImplicitFunction.h:58
vtkX3D::level
@ level
Definition: vtkX3D.h:395
vtkX3D::port
@ port
Definition: vtkX3D.h:447
vtkHyperOctreeCutter::SetNumberOfContours
void SetNumberOfContours(int number)
Set the number of contours to place into the list.
Definition: vtkHyperOctreeCutter.h:119
vtkHyperOctreeCutter::GetValue
double GetValue(int i)
Get the ith contour value.
Definition: vtkHyperOctreeCutter.h:96
vtkHyperOctreeCutter::GetSortByAsString
const char * GetSortByAsString()
Return the sorting procedure as a descriptive character string.
Definition: vtkHyperOctreeCutter.h:204
vtkPolyDataAlgorithm::RequestData
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
vtkHyperOctreeCutter::InCD
vtkDataSetAttributes * InCD
Definition: vtkHyperOctreeCutter.h:256
vtkObject::GetMTime
virtual vtkMTimeType GetMTime()
Return this object's modified time.
vtkHyperOctreeCutter::Sibling
vtkHyperOctreeCursor * Sibling
Definition: vtkHyperOctreeCutter.h:260
vtkHyperOctreeCutter::Iter
int Iter
Definition: vtkHyperOctreeCutter.h:262
vtkContourValues
helper object to manage setting and generating contour values
Definition: vtkContourValues.h:35
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:38
vtkHyperOctreeCutter::ContourValues
vtkContourValues * ContourValues
Definition: vtkHyperOctreeCutter.h:245
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:39
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:50
vtkHyperOctreeCutter::NewLines
vtkCellArray * NewLines
Definition: vtkHyperOctreeCutter.h:253
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:51
vtkPolyDataAlgorithm::FillInputPortInformation
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkHyperOctreeCutter::CellScalars
vtkDoubleArray * CellScalars
Definition: vtkHyperOctreeCutter.h:265
vtkHyperOctreeCutter::GenerateValues
void GenerateValues(int numContours, double rangeStart, double rangeEnd)
Generate numContours equally spaced contour values between specified range.
Definition: vtkHyperOctreeCutter.h:139
vtkHyperOctreeCutter::Polygon
vtkPolygon * Polygon
Definition: vtkHyperOctreeCutter.h:270
vtkHyperOctreeCutter::Triangulator
vtkOrderedTriangulator * Triangulator
Definition: vtkHyperOctreeCutter.h:259
vtkHyperOctreeCutter::Input
vtkHyperOctree * Input
Definition: vtkHyperOctreeCutter.h:248
vtkContourValues.h
vtkHyperOctreeCutter
Cut vtkHyperOctree with user-specified implicit function.
Definition: vtkHyperOctreeCutter.h:74
vtkBooleanMacro
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkInformation
Store vtkAlgorithm input/output information.
Definition: vtkInformation.h:86
vtkSetMacro
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkHyperOctreeCutter::AllGreater
int * AllGreater
Definition: vtkHyperOctreeCutter.h:279
vtkX3D::info
@ info
Definition: vtkX3D.h:376
vtkHyperOctreeCutter::GenerateCutScalars
int GenerateCutScalars
Definition: vtkHyperOctreeCutter.h:246
vtkPolyData
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:85
vtkHyperOctreeCutter::TotalCounter
vtkIdType TotalCounter
Definition: vtkHyperOctreeCutter.h:273
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:41
vtkHyperOctreeCutter::GetNumberOfContours
int GetNumberOfContours()
Get the number of contours in the list of contour values.
Definition: vtkHyperOctreeCutter.h:125
vtkCutter.h
vtkHyperOctreeCutter::GetValues
void GetValues(double *contourValues)
Fill a supplied list with contour values.
Definition: vtkHyperOctreeCutter.h:111
vtkHyperOctree
A dataset structured as a tree where each node has exactly 2^n children.
Definition: vtkHyperOctree.h:143
vtkHyperOctreeCutter::GetValues
double * GetValues()
Get a pointer to an array of contour values.
Definition: vtkHyperOctreeCutter.h:103
vtkTetra
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:47
vtkHyperOctreeCutter::Pts
vtkPoints * Pts
Definition: vtkHyperOctreeCutter.h:269
VTK_SORT_BY_CELL
#define VTK_SORT_BY_CELL
Definition: vtkCutter.h:60
vtkHyperOctreeCutter::SetSortByToSortByCell
void SetSortByToSortByCell()
Definition: vtkHyperOctreeCutter.h:196
vtkHyperOctreeCutter::SortBy
int SortBy
Definition: vtkHyperOctreeCutter.h:244
vtkHyperOctreeCutter::AllLess
int * AllLess
Definition: vtkHyperOctreeCutter.h:278
vtkHyperOctreeCutter::OutCD
vtkCellData * OutCD
Definition: vtkHyperOctreeCutter.h:257
vtkPolyDataAlgorithm::New
static vtkPolyDataAlgorithm * New()
vtkHyperOctreeCutter::SetValue
void SetValue(int i, double value)
Set a particular contour value at contour number i.
Definition: vtkHyperOctreeCutter.h:90
vtkHyperOctreeCutter::TemplateCounter
vtkIdType TemplateCounter
Definition: vtkHyperOctreeCutter.h:274
vtkHyperOctreeCutter::GenerateValues
void GenerateValues(int numContours, double range[2])
Generate numContours equally spaced contour values between specified range.
Definition: vtkHyperOctreeCutter.h:132
vtkPolyDataAlgorithm
Superclass for algorithms that produce only polydata as output.
Definition: vtkPolyDataAlgorithm.h:44