VTK  9.1.0
vtkDelaunay2D.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkDelaunay2D.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=========================================================================*/
234#ifndef vtkDelaunay2D_h
235#define vtkDelaunay2D_h
236
237#include "vtkFiltersCoreModule.h" // For export macro
238#include "vtkPolyDataAlgorithm.h"
239
241class vtkCellArray;
242class vtkIdList;
243class vtkPointSet;
244
245#define VTK_DELAUNAY_XY_PLANE 0
246#define VTK_SET_TRANSFORM_PLANE 1
247#define VTK_BEST_FITTING_PLANE 2
248
249class VTKFILTERSCORE_EXPORT vtkDelaunay2D : public vtkPolyDataAlgorithm
250{
251public:
253 void PrintSelf(ostream& os, vtkIndent indent) override;
254
260
271
281
286
288
294 vtkSetClampMacro(Alpha, double, 0.0, VTK_DOUBLE_MAX);
295 vtkGetMacro(Alpha, double);
297
299
304 vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
305 vtkGetMacro(Tolerance, double);
307
309
313 vtkSetClampMacro(Offset, double, 0.75, VTK_DOUBLE_MAX);
314 vtkGetMacro(Offset, double);
316
318
324 vtkSetMacro(BoundingTriangulation, vtkTypeBool);
325 vtkGetMacro(BoundingTriangulation, vtkTypeBool);
326 vtkBooleanMacro(BoundingTriangulation, vtkTypeBool);
328
330
341 vtkGetObjectMacro(Transform, vtkAbstractTransform);
343
345
353 vtkSetClampMacro(ProjectionPlaneMode, int, VTK_DELAUNAY_XY_PLANE, VTK_BEST_FITTING_PLANE);
354 vtkGetMacro(ProjectionPlaneMode, int);
356
364
365protected:
367 ~vtkDelaunay2D() override;
368
370
371 double Alpha;
372 double Tolerance;
374 double Offset;
375
377
378 int ProjectionPlaneMode; // selects the plane in 3D where the Delaunay triangulation will be
379 // computed.
380
381private:
382 vtkPolyData* Mesh; // the created mesh
383 double* Points; // the raw points in double precision
384 void SetPoint(vtkIdType id, double* x)
385 {
386 vtkIdType idx = 3 * id;
387 this->Points[idx] = x[0];
388 this->Points[idx + 1] = x[1];
389 this->Points[idx + 2] = x[2];
390 }
391
392 void GetPoint(vtkIdType id, double x[3])
393 {
394 double* ptr = this->Points + 3 * id;
395 x[0] = *ptr++;
396 x[1] = *ptr++;
397 x[2] = *ptr;
398 }
399
400 int NumberOfDuplicatePoints;
401 int NumberOfDegeneracies;
402
403 int* RecoverBoundary(vtkPolyData* source);
404 int RecoverEdge(vtkPolyData* source, vtkIdType p1, vtkIdType p2);
405 void FillPolygons(vtkCellArray* polys, int* triUse);
406
407 int InCircle(double x[3], double x1[3], double x2[3], double x3[3]);
408 vtkIdType FindTriangle(double x[3], vtkIdType ptIds[3], vtkIdType tri, double tol,
409 vtkIdType nei[3], vtkIdList* neighbors);
410 void CheckEdge(
411 vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2, vtkIdType tri, bool recursive);
412
413 int FillInputPortInformation(int, vtkInformation*) override;
414
415private:
416 vtkDelaunay2D(const vtkDelaunay2D&) = delete;
417 void operator=(const vtkDelaunay2D&) = delete;
418};
419
420#endif
void GetPoint(const int i, const int j, const int k, double pnt[3])
superclass for all geometric transformations
Proxy object to connect input/output ports.
object to represent cell connectivity
Definition: vtkCellArray.h:290
create 2D Delaunay triangulation of input points
~vtkDelaunay2D() override
vtkPolyData * GetSource()
Get a pointer to the source object.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
virtual void SetTransform(vtkAbstractTransform *)
Set / get the transform which is applied to points to generate a 2D problem.
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.
static vtkAbstractTransform * ComputeBestFittingPlane(vtkPointSet *input)
This method computes the best fit plane to a set of points represented by a vtkPointSet.
vtkTypeBool BoundingTriangulation
vtkAbstractTransform * Transform
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:140
a simple class to control print indentation
Definition: vtkIndent.h:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
concrete class for storing a set of points
Definition: vtkPointSet.h:106
Superclass for algorithms that produce only polydata as output.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:195
int vtkTypeBool
Definition: vtkABI.h:69
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:165