VTK  9.6.20260615
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
226
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 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 VTK_MARSHALAUTO 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
322
324
334 vtkSetSmartPointerMacro(Transform, vtkAbstractTransform);
335 vtkGetSmartPointerMacro(Transform, vtkAbstractTransform);
337
339
348 vtkGetMacro(ProjectionPlaneMode, int);
350
358
360
369
371
382 vtkBooleanMacro(UseHilbertSorter, vtkTypeBool);
384
385protected:
387
389
390 double Alpha;
391 double Tolerance;
393 double Offset;
396
397 // Transform input points (if necessary)
399
400 int ProjectionPlaneMode; // selects the plane in 3D where the Delaunay triangulation will be
401 // computed.
402
403private:
404 vtkSmartPointer<vtkPolyData> Mesh; // the created mesh
405
406 // the raw points in double precision, and methods to access them
407 double* Points;
408 void SetPoint(vtkIdType id, double* x)
409 {
410 vtkIdType idx = 3 * id;
411 this->Points[idx] = x[0];
412 this->Points[idx + 1] = x[1];
413 this->Points[idx + 2] = x[2];
414 }
415 void GetPoint(vtkIdType id, double x[3])
416 {
417 double* ptr = this->Points + 3 * id;
418 x[0] = *ptr++;
419 x[1] = *ptr++;
420 x[2] = *ptr;
421 }
422
423 // Keep track of the bounding radius of all points (including the eight bounding points).
424 // This is used occasionally for numerical sanity checks to determine whether a point is
425 // within a circumcircle.
426 double BoundingRadius2;
427
428 int NumberOfDuplicatePoints;
429 int NumberOfDegeneracies;
430
431 // Various methods to support the Delaunay algorithm
432 int* RecoverBoundary(vtkPolyData* source);
433 int RecoverEdge(vtkPolyData* source, vtkIdType p1, vtkIdType p2);
434 void FillPolygons(vtkCellArray* polys, int* triUse);
435
436 int InCircle(double x[3], double x1[3], double x2[3], double x3[3]);
437
438 // CheckEdge() determines if the edge (p1,p2) of triangle tri satisfies the
439 // Delaunay criterion; if not it swaps the edge diagonal.
440 bool CheckEdge(vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2, vtkIdType tri);
441
442 int FillInputPortInformation(int, vtkInformation*) override;
443
444 vtkDelaunay2D(const vtkDelaunay2D&) = delete;
445 void operator=(const vtkDelaunay2D&) = delete;
446};
447
448VTK_ABI_NAMESPACE_END
449#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...
vtkTypeBool UseHilbertSorter
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.
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
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:363
#define VTK_DOUBLE_MAX
Definition vtkType.h:202
#define VTK_MARSHALAUTO