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 =========================================================================*/
128 #ifndef vtkDelaunay2D_h
129 #define vtkDelaunay2D_h
130 
131 #include "vtkFiltersCoreModule.h" // For export macro
132 #include "vtkPolyDataAlgorithm.h"
133 
135 class vtkCellArray;
136 class vtkIdList;
137 class vtkPointSet;
138 
139 #define VTK_DELAUNAY_XY_PLANE 0
140 #define VTK_SET_TRANSFORM_PLANE 1
141 #define VTK_BEST_FITTING_PLANE 2
142 
143 class VTKFILTERSCORE_EXPORT vtkDelaunay2D : public vtkPolyDataAlgorithm
144 {
145 public:
147  void PrintSelf(ostream& os, vtkIndent indent) override;
148 
153  static vtkDelaunay2D* New();
154 
164  void SetSourceData(vtkPolyData*);
165 
174  void SetSourceConnection(vtkAlgorithmOutput* algOutput);
175 
179  vtkPolyData* GetSource();
180 
182 
188  vtkSetClampMacro(Alpha, double, 0.0, VTK_DOUBLE_MAX);
189  vtkGetMacro(Alpha, double);
191 
193 
198  vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
199  vtkGetMacro(Tolerance, double);
201 
203 
207  vtkSetClampMacro(Offset, double, 0.75, VTK_DOUBLE_MAX);
208  vtkGetMacro(Offset, double);
210 
212 
218  vtkSetMacro(BoundingTriangulation, vtkTypeBool);
219  vtkGetMacro(BoundingTriangulation, vtkTypeBool);
220  vtkBooleanMacro(BoundingTriangulation, vtkTypeBool);
222 
224 
234  virtual void SetTransform(vtkAbstractTransform*);
235  vtkGetObjectMacro(Transform, vtkAbstractTransform);
237 
239 
247  vtkSetClampMacro(ProjectionPlaneMode, int, VTK_DELAUNAY_XY_PLANE, VTK_BEST_FITTING_PLANE);
248  vtkGetMacro(ProjectionPlaneMode, int);
250 
257  static vtkAbstractTransform* ComputeBestFittingPlane(vtkPointSet* input);
258 
259 protected:
260  vtkDelaunay2D();
261  ~vtkDelaunay2D() override;
262 
264 
265  double Alpha;
266  double Tolerance;
268  double Offset;
269 
271 
272  int ProjectionPlaneMode; // selects the plane in 3D where the Delaunay triangulation will be
273  // computed.
274 
275 private:
276  vtkPolyData* Mesh; // the created mesh
277  double* Points; // the raw points in double precision
278  void SetPoint(vtkIdType id, double* x)
279  {
280  vtkIdType idx = 3 * id;
281  this->Points[idx] = x[0];
282  this->Points[idx + 1] = x[1];
283  this->Points[idx + 2] = x[2];
284  }
285 
286  void GetPoint(vtkIdType id, double x[3])
287  {
288  double* ptr = this->Points + 3 * id;
289  x[0] = *ptr++;
290  x[1] = *ptr++;
291  x[2] = *ptr;
292  }
293 
294  int NumberOfDuplicatePoints;
295  int NumberOfDegeneracies;
296 
297  int* RecoverBoundary(vtkPolyData* source);
298  int RecoverEdge(vtkPolyData* source, vtkIdType p1, vtkIdType p2);
299  void FillPolygons(vtkCellArray* polys, int* triUse);
300 
301  int InCircle(double x[3], double x1[3], double x2[3], double x3[3]);
302  vtkIdType FindTriangle(double x[3], vtkIdType ptIds[3], vtkIdType tri, double tol,
303  vtkIdType nei[3], vtkIdList* neighbors);
304  bool CheckEdge(
305  vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2, vtkIdType tri, bool recursive);
306 
307  int FillInputPortInformation(int, vtkInformation*) override;
308 
309 private:
310  vtkDelaunay2D(const vtkDelaunay2D&) = delete;
311  void operator=(const vtkDelaunay2D&) = delete;
312 };
313 
314 #endif
#define VTK_DOUBLE_MAX
Definition: vtkType.h:165
Store vtkAlgorithm input/output information.
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
concrete class for storing a set of points
Definition: vtkPointSet.h:69
int vtkIdType
Definition: vtkType.h:332
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:90
void GetPoint(const int i, const int j, const int k, double pnt[3])
Proxy object to connect input/output ports.
static vtkPolyDataAlgorithm * New()
int vtkTypeBool
Definition: vtkABI.h:69
Superclass for algorithms that produce only polydata as output.
a simple class to control print indentation
Definition: vtkIndent.h:39
list of point or cell ids
Definition: vtkIdList.h:33
superclass for all geometric transformations
vtkAbstractTransform * Transform
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
#define VTK_BEST_FITTING_PLANE
object to represent cell connectivity
Definition: vtkCellArray.h:186
create 2D Delaunay triangulation of input points
vtkTypeBool BoundingTriangulation
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
Store zero or more vtkInformation instances.
#define VTK_DELAUNAY_XY_PLANE
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.