VTK
vtkDelaunay3D.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkDelaunay3D.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 =========================================================================*/
100 #ifndef vtkDelaunay3D_h
101 #define vtkDelaunay3D_h
102 
103 #include "vtkFiltersCoreModule.h" // For export macro
105 
106 class vtkIdList;
107 class vtkPointLocator;
108 class vtkPointSet;
109 class vtkPoints;
110 class vtkTetraArray;
112 
113 class VTKFILTERSCORE_EXPORT vtkDelaunay3D : public vtkUnstructuredGridAlgorithm
114 {
115 public:
117  void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
118 
123  static vtkDelaunay3D *New();
124 
126 
135  vtkSetClampMacro(Alpha,double,0.0,VTK_DOUBLE_MAX);
136  vtkGetMacro(Alpha,double);
138 
140 
143  vtkSetMacro(AlphaTets,int);
144  vtkGetMacro(AlphaTets,int);
145  vtkBooleanMacro(AlphaTets,int);
147 
149 
152  vtkSetMacro(AlphaTris,int);
153  vtkGetMacro(AlphaTris,int);
154  vtkBooleanMacro(AlphaTris,int);
156 
158 
161  vtkSetMacro(AlphaLines,int);
162  vtkGetMacro(AlphaLines,int);
163  vtkBooleanMacro(AlphaLines,int);
165 
167 
170  vtkSetMacro(AlphaVerts,int);
171  vtkGetMacro(AlphaVerts,int);
172  vtkBooleanMacro(AlphaVerts,int);
174 
176 
181  vtkSetClampMacro(Tolerance,double,0.0,1.0);
182  vtkGetMacro(Tolerance,double);
184 
186 
190  vtkSetClampMacro(Offset,double,2.5,VTK_DOUBLE_MAX);
191  vtkGetMacro(Offset,double);
193 
195 
201  vtkSetMacro(BoundingTriangulation,int);
202  vtkGetMacro(BoundingTriangulation,int);
203  vtkBooleanMacro(BoundingTriangulation,int);
205 
207 
212  vtkGetObjectMacro(Locator,vtkIncrementalPointLocator);
214 
220 
234  vtkIdType numPts, vtkPoints* &pts);
235 
247  vtkIdType id, double x[3], vtkIdList *holeTetras);
248 
256 
260  vtkMTimeType GetMTime() VTK_OVERRIDE;
261 
263 
268  vtkSetMacro(OutputPointsPrecision,int);
269  vtkGetMacro(OutputPointsPrecision,int);
271 
272 protected:
274  ~vtkDelaunay3D() VTK_OVERRIDE;
275 
276  int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) VTK_OVERRIDE;
277 
278  double Alpha;
279  int AlphaTets;
280  int AlphaTris;
281  int AlphaLines;
282  int AlphaVerts;
283  double Tolerance;
284  int BoundingTriangulation;
285  double Offset;
286  int OutputPointsPrecision;
287 
288  vtkIncrementalPointLocator *Locator; //help locate points faster
289 
290  vtkTetraArray *TetraArray; //used to keep track of circumspheres/neighbors
291  int FindTetra(vtkUnstructuredGrid *Mesh, double x[3], vtkIdType tetId,
292  int depth);
293  int InSphere(double x[3], vtkIdType tetraId);
294  void InsertTetra(vtkUnstructuredGrid *Mesh, vtkPoints *pts,
295  vtkIdType tetraId);
296 
297  int NumberOfDuplicatePoints; //keep track of bad data
298  int NumberOfDegeneracies;
299 
300  // Keep track of number of references to points to avoid new/delete calls
301  int *References;
302 
303  vtkIdType FindEnclosingFaces(double x[3], vtkUnstructuredGrid *Mesh,
304  vtkIdList *tetras, vtkIdList *faces,
305  vtkIncrementalPointLocator *Locator);
306 
307  int FillInputPortInformation(int, vtkInformation*) VTK_OVERRIDE;
308 private: //members added for performance
309  vtkIdList *Tetras; //used in InsertPoint
310  vtkIdList *Faces; //used in InsertPoint
311  vtkIdList *CheckedTetras; //used by InsertPoint
312 
313 private:
314  vtkDelaunay3D(const vtkDelaunay3D&) VTK_DELETE_FUNCTION;
315  void operator=(const vtkDelaunay3D&) VTK_DELETE_FUNCTION;
316 };
317 
318 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:40
vtkDelaunay3D::EndPointInsertion
void EndPointInsertion()
Invoke this method after all points have been inserted.
vtkIdType
int vtkIdType
Definition: vtkType.h:287
vtkInformationVector
Store zero or more vtkInformation instances.
Definition: vtkInformationVector.h:42
vtkPointLocator
quickly locate points in 3-space
Definition: vtkPointLocator.h:54
vtkX3D::length
@ length
Definition: vtkX3D.h:393
vtkX3D::center
@ center
Definition: vtkX3D.h:230
vtkMTimeType
vtkTypeUInt64 vtkMTimeType
Definition: vtkType.h:248
vtkX3D::points
@ points
Definition: vtkX3D.h:446
vtkDelaunay3D::New
static vtkDelaunay3D * New()
Construct object with Alpha = 0.0; Tolerance = 0.001; Offset = 2.5; BoundingTriangulation turned off.
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:40
vtkDelaunay3D::CreateDefaultLocator
void CreateDefaultLocator()
Create default locator.
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:52
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:37
vtkDelaunay3D::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkBooleanMacro
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkInformation
Store vtkAlgorithm input/output information.
Definition: vtkInformation.h:87
vtkSetMacro
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkDelaunay3D::InsertPoint
void InsertPoint(vtkUnstructuredGrid *Mesh, vtkPoints *points, vtkIdType id, double x[3], vtkIdList *holeTetras)
This is a helper method used with InitPointInsertion() to create tetrahedronalizations of points.
vtkPointSet
abstract class for specifying dataset behavior
Definition: vtkPointSet.h:43
vtkDelaunay3D::SetLocator
void SetLocator(vtkIncrementalPointLocator *locator)
Set / get a spatial locator for merging points.
vtkDelaunay3D
create 3D Delaunay triangulation of input points
Definition: vtkDelaunay3D.h:114
vtkUnstructuredGridAlgorithm.h
vtkUnstructuredGridAlgorithm
Superclass for algorithms that produce only unstructured grid as output.
Definition: vtkUnstructuredGridAlgorithm.h:41
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:83
vtkDelaunay3D::InitPointInsertion
vtkUnstructuredGrid * InitPointInsertion(double center[3], double length, vtkIdType numPts, vtkPoints *&pts)
This is a helper method used with InsertPoint() to create tetrahedronalizations of points.
vtkDelaunay3D::GetMTime
vtkMTimeType GetMTime() override
Return the MTime also considering the locator.
VTK_DOUBLE_MAX
#define VTK_DOUBLE_MAX
Definition: vtkType.h:163