VTK  9.4.20241218
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
227#ifndef vtkDelaunay2D_h
228#define vtkDelaunay2D_h
229
230#include "vtkAbstractTransform.h" // For point transformation
231#include "vtkFiltersCoreModule.h" // For export macro
232#include "vtkPolyDataAlgorithm.h"
233
234VTK_ABI_NAMESPACE_BEGIN
235class vtkCellArray;
236class vtkIdList;
237class vtkPointSet;
238
239#define VTK_DELAUNAY_XY_PLANE 0
240#define VTK_SET_TRANSFORM_PLANE 1
241#define VTK_BEST_FITTING_PLANE 2
242
243class VTKFILTERSCORE_EXPORT vtkDelaunay2D : public vtkPolyDataAlgorithm
244{
245public:
247 void PrintSelf(ostream& os, vtkIndent indent) override;
248
254
265
275
280
282
288 vtkSetClampMacro(Alpha, double, 0.0, VTK_DOUBLE_MAX);
289 vtkGetMacro(Alpha, double);
291
293
298 vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
299 vtkGetMacro(Tolerance, double);
301
303
307 vtkSetClampMacro(Offset, double, 0.75, VTK_DOUBLE_MAX);
308 vtkGetMacro(Offset, double);
310
312
318 vtkSetMacro(BoundingTriangulation, vtkTypeBool);
319 vtkGetMacro(BoundingTriangulation, vtkTypeBool);
320 vtkBooleanMacro(BoundingTriangulation, vtkTypeBool);
322
324
337
339
347 vtkSetClampMacro(ProjectionPlaneMode, int, VTK_DELAUNAY_XY_PLANE, VTK_BEST_FITTING_PLANE);
348 vtkGetMacro(ProjectionPlaneMode, int);
350
358
360
365 vtkSetMacro(RandomPointInsertion, vtkTypeBool);
366 vtkGetMacro(RandomPointInsertion, vtkTypeBool);
367 vtkBooleanMacro(RandomPointInsertion, vtkTypeBool);
369
370protected:
372
374
375 double Alpha;
376 double Tolerance;
378 double Offset;
380
381 // Transform input points (if necessary)
383
384 int ProjectionPlaneMode; // selects the plane in 3D where the Delaunay triangulation will be
385 // computed.
386
387private:
388 vtkSmartPointer<vtkPolyData> Mesh; // the created mesh
389
390 // the raw points in double precision, and methods to access them
391 double* Points;
392 void SetPoint(vtkIdType id, double* x)
393 {
394 vtkIdType idx = 3 * id;
395 this->Points[idx] = x[0];
396 this->Points[idx + 1] = x[1];
397 this->Points[idx + 2] = x[2];
398 }
399 void GetPoint(vtkIdType id, double x[3])
400 {
401 double* ptr = this->Points + 3 * id;
402 x[0] = *ptr++;
403 x[1] = *ptr++;
404 x[2] = *ptr;
405 }
406
407 // Keep track of the bounding radius of all points (including the eight bounding points).
408 // This is used occasionally for numerical sanity checks to determine whether a point is
409 // within a circumcircle.
410 double BoundingRadius2;
411
412 int NumberOfDuplicatePoints;
413 int NumberOfDegeneracies;
414
415 // Various methods to support the Delaunay algorithm
416 int* RecoverBoundary(vtkPolyData* source);
417 int RecoverEdge(vtkPolyData* source, vtkIdType p1, vtkIdType p2);
418 void FillPolygons(vtkCellArray* polys, int* triUse);
419
420 int InCircle(double x[3], double x1[3], double x2[3], double x3[3]);
421 vtkIdType FindTriangle(double x[3], vtkIdType ptIds[3], vtkIdType tri, double tol,
422 vtkIdType nei[3], vtkIdList* neighbors);
423
424 // CheckEdge() is a recursive function to determine if triangles satisfy the Delaunay
425 // criterion. To prevent segfaults due to excessive recursion, recursion depth is limited.
426 bool CheckEdge(vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2, vtkIdType tri,
427 bool recursive, unsigned int depth);
428
429 int FillInputPortInformation(int, vtkInformation*) override;
430
431 vtkDelaunay2D(const vtkDelaunay2D&) = delete;
432 void operator=(const vtkDelaunay2D&) = delete;
433};
434
435VTK_ABI_NAMESPACE_END
436#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
create 2D Delaunay triangulation of input points
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
vtkSetSmartPointerMacro(Transform, vtkAbstractTransform)
Set / get the transform which is applied to points to generate a 2D problem.
void SetSourceData(vtkPolyData *)
Specify the source object used to specify constrained edges and loops.
vtkGetSmartPointerMacro(Transform, vtkAbstractTransform)
Set / get the transform which is applied to points to generate a 2D problem.
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:133
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
concrete class for storing a set of points
Definition vtkPointSet.h:98
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Hold a reference to a vtkObjectBase instance.
int vtkTypeBool
Definition vtkABI.h:64
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:315
#define VTK_DOUBLE_MAX
Definition vtkType.h:154