VTK  9.0.20200814
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 =========================================================================*/
94 #ifndef vtkDelaunay3D_h
95 #define vtkDelaunay3D_h
96 
97 #include "vtkFiltersCoreModule.h" // For export macro
99 
100 class vtkIdList;
101 class vtkPointLocator;
102 class vtkPointSet;
103 class vtkPoints;
104 class vtkTetraArray;
106 
107 class VTKFILTERSCORE_EXPORT vtkDelaunay3D : public vtkUnstructuredGridAlgorithm
108 {
109 public:
111  void PrintSelf(ostream& os, vtkIndent indent) override;
112 
117  static vtkDelaunay3D* New();
118 
120 
129  vtkSetClampMacro(Alpha, double, 0.0, VTK_DOUBLE_MAX);
130  vtkGetMacro(Alpha, double);
132 
134 
137  vtkSetMacro(AlphaTets, vtkTypeBool);
138  vtkGetMacro(AlphaTets, vtkTypeBool);
139  vtkBooleanMacro(AlphaTets, vtkTypeBool);
141 
143 
146  vtkSetMacro(AlphaTris, vtkTypeBool);
147  vtkGetMacro(AlphaTris, vtkTypeBool);
148  vtkBooleanMacro(AlphaTris, vtkTypeBool);
150 
152 
155  vtkSetMacro(AlphaLines, vtkTypeBool);
156  vtkGetMacro(AlphaLines, vtkTypeBool);
157  vtkBooleanMacro(AlphaLines, vtkTypeBool);
159 
161 
164  vtkSetMacro(AlphaVerts, vtkTypeBool);
165  vtkGetMacro(AlphaVerts, vtkTypeBool);
166  vtkBooleanMacro(AlphaVerts, vtkTypeBool);
168 
170 
175  vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
176  vtkGetMacro(Tolerance, double);
178 
180 
184  vtkSetClampMacro(Offset, double, 2.5, VTK_DOUBLE_MAX);
185  vtkGetMacro(Offset, double);
187 
189 
195  vtkSetMacro(BoundingTriangulation, vtkTypeBool);
196  vtkGetMacro(BoundingTriangulation, vtkTypeBool);
197  vtkBooleanMacro(BoundingTriangulation, vtkTypeBool);
199 
201 
205  void SetLocator(vtkIncrementalPointLocator* locator);
206  vtkGetObjectMacro(Locator, vtkIncrementalPointLocator);
208 
213  void CreateDefaultLocator();
214 
227  vtkUnstructuredGrid* InitPointInsertion(
228  double center[3], double length, vtkIdType numPts, vtkPoints*& points);
229 
240  void InsertPoint(
241  vtkUnstructuredGrid* Mesh, vtkPoints* points, vtkIdType id, double x[3], vtkIdList* holeTetras);
242 
249  void EndPointInsertion();
250 
254  vtkMTimeType GetMTime() override;
255 
257 
262  vtkSetMacro(OutputPointsPrecision, int);
263  vtkGetMacro(OutputPointsPrecision, int);
265 
266 protected:
267  vtkDelaunay3D();
268  ~vtkDelaunay3D() override;
269 
271 
272  double Alpha;
277  double Tolerance;
279  double Offset;
281 
282  vtkIncrementalPointLocator* Locator; // help locate points faster
283 
284  vtkTetraArray* TetraArray; // used to keep track of circumspheres/neighbors
285  int FindTetra(vtkUnstructuredGrid* Mesh, double x[3], vtkIdType tetraId, int depth);
286  int InSphere(double x[3], vtkIdType tetraId);
287  void InsertTetra(vtkUnstructuredGrid* Mesh, vtkPoints* points, vtkIdType tetraId);
288 
289  int NumberOfDuplicatePoints; // keep track of bad data
291 
292  // Keep track of number of references to points to avoid new/delete calls
294 
295  vtkIdType FindEnclosingFaces(double x[3], vtkUnstructuredGrid* Mesh, vtkIdList* tetras,
296  vtkIdList* faces, vtkIncrementalPointLocator* Locator);
297 
298  int FillInputPortInformation(int, vtkInformation*) override;
299 
300 private: // members added for performance
301  vtkIdList* Tetras; // used in InsertPoint
302  vtkIdList* Faces; // used in InsertPoint
303  vtkIdList* CheckedTetras; // used by InsertPoint
304 
305 private:
306  vtkDelaunay3D(const vtkDelaunay3D&) = delete;
307  void operator=(const vtkDelaunay3D&) = delete;
308 };
309 
310 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:33
vtkDelaunay3D::AlphaLines
vtkTypeBool AlphaLines
Definition: vtkDelaunay3D.h:275
vtkDelaunay3D::AlphaTets
vtkTypeBool AlphaTets
Definition: vtkDelaunay3D.h:273
vtkDelaunay3D::Tolerance
double Tolerance
Definition: vtkDelaunay3D.h:277
vtkDelaunay3D::Alpha
double Alpha
Definition: vtkDelaunay3D.h:272
vtkIdType
int vtkIdType
Definition: vtkType.h:330
vtkInformationVector
Store zero or more vtkInformation instances.
Definition: vtkInformationVector.h:35
vtkPointLocator
quickly locate points in 3-space
Definition: vtkPointLocator.h:50
vtkDelaunay3D::TetraArray
vtkTetraArray * TetraArray
Definition: vtkDelaunay3D.h:284
vtkX3D::length
Definition: vtkX3D.h:399
vtkX3D::center
Definition: vtkX3D.h:236
vtkDelaunay3D::AlphaTris
vtkTypeBool AlphaTris
Definition: vtkDelaunay3D.h:274
vtkDelaunay3D::OutputPointsPrecision
int OutputPointsPrecision
Definition: vtkDelaunay3D.h:280
vtkDelaunay3D::Locator
vtkIncrementalPointLocator * Locator
Definition: vtkDelaunay3D.h:282
vtkX3D::points
Definition: vtkX3D.h:452
vtkObject::GetMTime
virtual vtkMTimeType GetMTime()
Return this object's modified time.
vtkUnstructuredGridAlgorithm::RequestData
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
vtkDelaunay3D::NumberOfDuplicatePoints
int NumberOfDuplicatePoints
Definition: vtkDelaunay3D.h:289
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:51
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:30
vtkUnstructuredGridAlgorithm::FillInputPortInformation
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkInformation
Store vtkAlgorithm input/output information.
Definition: vtkInformation.h:73
vtkUnstructuredGridAlgorithm::New
static vtkUnstructuredGridAlgorithm * New()
vtkPointSet
concrete class for storing a set of points
Definition: vtkPointSet.h:66
vtkDelaunay3D::References
int * References
Definition: vtkDelaunay3D.h:293
vtkDelaunay3D
create 3D Delaunay triangulation of input points
Definition: vtkDelaunay3D.h:107
vtkUnstructuredGridAlgorithm::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkUnstructuredGridAlgorithm.h
vtkUnstructuredGridAlgorithm
Superclass for algorithms that produce only unstructured grid as output.
Definition: vtkUnstructuredGridAlgorithm.h:40
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:93
vtkDelaunay3D::Offset
double Offset
Definition: vtkDelaunay3D.h:279
vtkDelaunay3D::NumberOfDegeneracies
int NumberOfDegeneracies
Definition: vtkDelaunay3D.h:290
VTK_DOUBLE_MAX
#define VTK_DOUBLE_MAX
Definition: vtkType.h:157
vtkDelaunay3D::AlphaVerts
vtkTypeBool AlphaVerts
Definition: vtkDelaunay3D.h:276
vtkTypeBool
int vtkTypeBool
Definition: vtkABI.h:69
vtkMTimeType
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:285
vtkDelaunay3D::BoundingTriangulation
vtkTypeBool BoundingTriangulation
Definition: vtkDelaunay3D.h:278