VTK  9.2.6
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 =========================================================================*/
237 #ifndef vtkDelaunay2D_h
238 #define vtkDelaunay2D_h
239 
240 #include "vtkFiltersCoreModule.h" // For export macro
241 #include "vtkPolyDataAlgorithm.h"
242 
244 class vtkCellArray;
245 class vtkIdList;
246 class vtkPointSet;
247 
248 #define VTK_DELAUNAY_XY_PLANE 0
249 #define VTK_SET_TRANSFORM_PLANE 1
250 #define VTK_BEST_FITTING_PLANE 2
251 
252 class VTKFILTERSCORE_EXPORT vtkDelaunay2D : public vtkPolyDataAlgorithm
253 {
254 public:
256  void PrintSelf(ostream& os, vtkIndent indent) override;
257 
262  static vtkDelaunay2D* New();
263 
274 
284 
289 
291 
297  vtkSetClampMacro(Alpha, double, 0.0, VTK_DOUBLE_MAX);
298  vtkGetMacro(Alpha, double);
300 
302 
307  vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
308  vtkGetMacro(Tolerance, double);
310 
312 
316  vtkSetClampMacro(Offset, double, 0.75, VTK_DOUBLE_MAX);
317  vtkGetMacro(Offset, double);
319 
321 
327  vtkSetMacro(BoundingTriangulation, vtkTypeBool);
328  vtkGetMacro(BoundingTriangulation, vtkTypeBool);
329  vtkBooleanMacro(BoundingTriangulation, vtkTypeBool);
331 
333 
344  vtkGetObjectMacro(Transform, vtkAbstractTransform);
346 
348 
356  vtkSetClampMacro(ProjectionPlaneMode, int, VTK_DELAUNAY_XY_PLANE, VTK_BEST_FITTING_PLANE);
357  vtkGetMacro(ProjectionPlaneMode, int);
359 
367 
368 protected:
370  ~vtkDelaunay2D() override;
371 
373 
374  double Alpha;
375  double Tolerance;
377  double Offset;
378 
380 
381  int ProjectionPlaneMode; // selects the plane in 3D where the Delaunay triangulation will be
382  // computed.
383 
384 private:
385  vtkPolyData* Mesh; // the created mesh
386  double* Points; // the raw points in double precision
387  void SetPoint(vtkIdType id, double* x)
388  {
389  vtkIdType idx = 3 * id;
390  this->Points[idx] = x[0];
391  this->Points[idx + 1] = x[1];
392  this->Points[idx + 2] = x[2];
393  }
394 
395  void GetPoint(vtkIdType id, double x[3])
396  {
397  double* ptr = this->Points + 3 * id;
398  x[0] = *ptr++;
399  x[1] = *ptr++;
400  x[2] = *ptr;
401  }
402 
403  int NumberOfDuplicatePoints;
404  int NumberOfDegeneracies;
405 
406  int* RecoverBoundary(vtkPolyData* source);
407  int RecoverEdge(vtkPolyData* source, vtkIdType p1, vtkIdType p2);
408  void FillPolygons(vtkCellArray* polys, int* triUse);
409 
410  int InCircle(double x[3], double x1[3], double x2[3], double x3[3]);
411  vtkIdType FindTriangle(double x[3], vtkIdType ptIds[3], vtkIdType tri, double tol,
412  vtkIdType nei[3], vtkIdList* neighbors);
413  bool CheckEdge(
414  vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2, vtkIdType tri, bool recursive);
415 
416  int FillInputPortInformation(int, vtkInformation*) override;
417 
418 private:
419  vtkDelaunay2D(const vtkDelaunay2D&) = delete;
420  void operator=(const vtkDelaunay2D&) = delete;
421 };
422 
423 #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:296
create 2D Delaunay triangulation of input points
~vtkDelaunay2D() override
static vtkAbstractTransform * ComputeBestFittingPlane(vtkPointSet *input)
This method computes the best fit plane to a set of points represented by a vtkPointSet.
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...
virtual void SetTransform(vtkAbstractTransform *)
Set / get the transform which is applied to points to generate a 2D problem.
void SetSourceConnection(vtkAlgorithmOutput *algOutput)
Specify the source object used to specify constrained edges and loops.
vtkPolyData * GetSource()
Get a pointer to the source object.
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:143
a simple class to control print indentation
Definition: vtkIndent.h:119
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
concrete class for storing a set of points
Definition: vtkPointSet.h:109
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:200
@ Transform
Definition: vtkX3D.h:47
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