VTK  9.5.20250905
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#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
234
235VTK_ABI_NAMESPACE_BEGIN
236class vtkCellArray;
237class vtkIdList;
238class vtkPointSet;
239
240#define VTK_DELAUNAY_XY_PLANE 0
241#define VTK_SET_TRANSFORM_PLANE 1
242#define VTK_BEST_FITTING_PLANE 2
243
244class VTKFILTERSCORE_EXPORT VTK_MARSHALAUTO vtkDelaunay2D : public vtkPolyDataAlgorithm
245{
246public:
248 void PrintSelf(ostream& os, vtkIndent indent) override;
249
255
266
276
281
283
289 vtkSetClampMacro(Alpha, double, 0.0, VTK_DOUBLE_MAX);
290 vtkGetMacro(Alpha, double);
292
294
299 vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
300 vtkGetMacro(Tolerance, double);
302
304
308 vtkSetClampMacro(Offset, double, 0.75, VTK_DOUBLE_MAX);
309 vtkGetMacro(Offset, double);
311
313
319 vtkSetMacro(BoundingTriangulation, vtkTypeBool);
320 vtkGetMacro(BoundingTriangulation, vtkTypeBool);
321 vtkBooleanMacro(BoundingTriangulation, vtkTypeBool);
323
325
335 vtkSetSmartPointerMacro(Transform, vtkAbstractTransform);
336 vtkGetSmartPointerMacro(Transform, vtkAbstractTransform);
338
340
348 vtkSetClampMacro(ProjectionPlaneMode, int, VTK_DELAUNAY_XY_PLANE, VTK_BEST_FITTING_PLANE);
349 vtkGetMacro(ProjectionPlaneMode, int);
351
359
361
366 vtkSetMacro(RandomPointInsertion, vtkTypeBool);
367 vtkGetMacro(RandomPointInsertion, vtkTypeBool);
368 vtkBooleanMacro(RandomPointInsertion, vtkTypeBool);
370
371protected:
373
375
376 double Alpha;
377 double Tolerance;
379 double Offset;
381
382 // Transform input points (if necessary)
384
385 int ProjectionPlaneMode; // selects the plane in 3D where the Delaunay triangulation will be
386 // computed.
387
388private:
389 vtkSmartPointer<vtkPolyData> Mesh; // the created mesh
390
391 // the raw points in double precision, and methods to access them
392 double* Points;
393 void SetPoint(vtkIdType id, double* x)
394 {
395 vtkIdType idx = 3 * id;
396 this->Points[idx] = x[0];
397 this->Points[idx + 1] = x[1];
398 this->Points[idx + 2] = x[2];
399 }
400 void GetPoint(vtkIdType id, double x[3])
401 {
402 double* ptr = this->Points + 3 * id;
403 x[0] = *ptr++;
404 x[1] = *ptr++;
405 x[2] = *ptr;
406 }
407
408 // Keep track of the bounding radius of all points (including the eight bounding points).
409 // This is used occasionally for numerical sanity checks to determine whether a point is
410 // within a circumcircle.
411 double BoundingRadius2;
412
413 int NumberOfDuplicatePoints;
414 int NumberOfDegeneracies;
415
416 // Various methods to support the Delaunay algorithm
417 int* RecoverBoundary(vtkPolyData* source);
418 int RecoverEdge(vtkPolyData* source, vtkIdType p1, vtkIdType p2);
419 void FillPolygons(vtkCellArray* polys, int* triUse);
420
421 int InCircle(double x[3], double x1[3], double x2[3], double x3[3]);
422 vtkIdType FindTriangle(double x[3], vtkIdType ptIds[3], vtkIdType tri, double tol,
423 vtkIdType nei[3], vtkIdList* neighbors);
424
425 // CheckEdge() is a recursive function to determine if triangles satisfy the Delaunay
426 // criterion. To prevent segfaults due to excessive recursion, recursion depth is limited.
427 bool CheckEdge(vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2, vtkIdType tri,
428 bool recursive, unsigned int depth);
429
430 int FillInputPortInformation(int, vtkInformation*) override;
431
432 vtkDelaunay2D(const vtkDelaunay2D&) = delete;
433 void operator=(const vtkDelaunay2D&) = delete;
434};
435
436VTK_ABI_NAMESPACE_END
437#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
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: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:332
#define VTK_DOUBLE_MAX
Definition vtkType.h:171
#define VTK_MARSHALAUTO