VTK  9.6.1
vtkDelaunay2D.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
117
118#ifndef vtkDelaunay2D_h
119#define vtkDelaunay2D_h
120
121#include "vtkAbstractTransform.h" // For point transformation
122#include "vtkFiltersCoreModule.h" // For export macro
123#include "vtkPolyDataAlgorithm.h"
124#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
125
126VTK_ABI_NAMESPACE_BEGIN
127class vtkCellArray;
128class vtkIdList;
129class vtkPointSet;
130
131#define VTK_DELAUNAY_XY_PLANE 0
132#define VTK_SET_TRANSFORM_PLANE 1
133#define VTK_BEST_FITTING_PLANE 2
134
135class VTKFILTERSCORE_EXPORT VTK_MARSHALAUTO vtkDelaunay2D : public vtkPolyDataAlgorithm
136{
137public:
139 void PrintSelf(ostream& os, vtkIndent indent) override;
140
146
157
167
172
174
180 vtkSetClampMacro(Alpha, double, 0.0, VTK_DOUBLE_MAX);
181 vtkGetMacro(Alpha, double);
183
185
190 vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
191 vtkGetMacro(Tolerance, double);
193
195
199 vtkSetClampMacro(Offset, double, 0.75, VTK_DOUBLE_MAX);
200 vtkGetMacro(Offset, double);
202
204
214
216
226 vtkSetSmartPointerMacro(Transform, vtkAbstractTransform);
227 vtkGetSmartPointerMacro(Transform, vtkAbstractTransform);
229
231
240 vtkGetMacro(ProjectionPlaneMode, int);
242
250
252
261
262protected:
264
266
267 double Alpha;
268 double Tolerance;
270 double Offset;
272
273 // Transform input points (if necessary)
275
276 int ProjectionPlaneMode; // selects the plane in 3D where the Delaunay triangulation will be
277 // computed.
278
279private:
280 vtkSmartPointer<vtkPolyData> Mesh; // the created mesh
281
282 // the raw points in double precision, and methods to access them
283 double* Points;
284 void SetPoint(vtkIdType id, double* x)
285 {
286 vtkIdType idx = 3 * id;
287 this->Points[idx] = x[0];
288 this->Points[idx + 1] = x[1];
289 this->Points[idx + 2] = x[2];
290 }
291 void GetPoint(vtkIdType id, double x[3])
292 {
293 double* ptr = this->Points + 3 * id;
294 x[0] = *ptr++;
295 x[1] = *ptr++;
296 x[2] = *ptr;
297 }
298
299 // Keep track of the bounding radius of all points (including the eight bounding points).
300 // This is used occasionally for numerical sanity checks to determine whether a point is
301 // within a circumcircle.
302 double BoundingRadius2;
303
304 int NumberOfDuplicatePoints;
305 int NumberOfDegeneracies;
306
307 // Various methods to support the Delaunay algorithm
308 int* RecoverBoundary(vtkPolyData* source);
309 int RecoverEdge(vtkPolyData* source, vtkIdType p1, vtkIdType p2);
310 void FillPolygons(vtkCellArray* polys, int* triUse);
311
312 int InCircle(double x[3], double x1[3], double x2[3], double x3[3]);
313 vtkIdType FindTriangle(double x[3], vtkIdType ptIds[3], vtkIdType tri, double tol,
314 vtkIdType nei[3], vtkIdList* neighbors);
315
316 // CheckEdge() is a recursive function to determine if triangles satisfy the Delaunay
317 // criterion. To prevent segfaults due to excessive recursion, recursion depth is limited.
318 bool CheckEdge(vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2, vtkIdType tri,
319 bool recursive, unsigned int depth);
320
321 int FillInputPortInformation(int, vtkInformation*) override;
322
323 vtkDelaunay2D(const vtkDelaunay2D&) = delete;
324 void operator=(const vtkDelaunay2D&) = delete;
325};
326
327VTK_ABI_NAMESPACE_END
328#endif
void GetPoint(int i, int j, int k, double pnt[3])
superclass for all geometric transformations
Proxy object to connect input/output ports.
object to represent cell connectivity
vtkPolyData * GetSource()
Get a pointer to the source object.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
static vtkDelaunay2D * New()
Construct object with Alpha = 0.0; Tolerance = 0.001; Offset = 1.25; BoundingTriangulation turned off...
void SetSourceConnection(vtkAlgorithmOutput *algOutput)
Specify the source object used to specify constrained edges and loops.
vtkSmartPointer< vtkAbstractTransform > Transform
static vtkAbstractTransform * ComputeBestFittingPlane(vtkPointSet *input)
This method computes the best fit plane to a set of points represented by a vtkPointSet.
vtkTypeBool BoundingTriangulation
vtkTypeBool RandomPointInsertion
void SetSourceData(vtkPolyData *)
Specify the source object used to specify constrained edges and loops.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
list of point or cell ids
Definition vtkIdList.h:24
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
concrete class for storing a set of points
Definition vtkPointSet.h:59
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:72
Hold a reference to a vtkObjectBase instance.
int vtkTypeBool
Definition vtkABI.h:64
virtual int FillInputPortInformation(int port, vtkInformation *info)
Fill the input port information objects for this algorithm.
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
#define VTK_BEST_FITTING_PLANE
#define VTK_DELAUNAY_XY_PLANE
int vtkIdType
Definition vtkType.h:368
#define VTK_DOUBLE_MAX
Definition vtkType.h:207
#define VTK_MARSHALAUTO