VTK
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=========================================================================*/
131#ifndef vtkDelaunay2D_h
132#define vtkDelaunay2D_h
133
134#include "vtkFiltersCoreModule.h" // For export macro
135#include "vtkPolyDataAlgorithm.h"
136
138class vtkCellArray;
139class vtkIdList;
140class vtkPointSet;
141
142#define VTK_DELAUNAY_XY_PLANE 0
143#define VTK_SET_TRANSFORM_PLANE 1
144#define VTK_BEST_FITTING_PLANE 2
145
146class VTKFILTERSCORE_EXPORT vtkDelaunay2D : public vtkPolyDataAlgorithm
147{
148public:
150 void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
151
157
168
178
183
185
191 vtkSetClampMacro(Alpha,double,0.0,VTK_DOUBLE_MAX);
192 vtkGetMacro(Alpha,double);
194
196
201 vtkSetClampMacro(Tolerance,double,0.0,1.0);
202 vtkGetMacro(Tolerance,double);
204
206
210 vtkSetClampMacro(Offset,double,0.75,VTK_DOUBLE_MAX);
211 vtkGetMacro(Offset,double);
213
215
221 vtkSetMacro(BoundingTriangulation,int);
222 vtkGetMacro(BoundingTriangulation,int);
223 vtkBooleanMacro(BoundingTriangulation,int);
225
227
238 vtkGetObjectMacro(Transform, vtkAbstractTransform);
240
242
245 vtkSetClampMacro(ProjectionPlaneMode,int,
247 vtkGetMacro(ProjectionPlaneMode,int);
249
250protected:
252 ~vtkDelaunay2D() VTK_OVERRIDE;
253
254 int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) VTK_OVERRIDE;
255
256 vtkAbstractTransform * ComputeBestFittingPlane(vtkPointSet *input);
257
258 double Alpha;
259 double Tolerance;
260 int BoundingTriangulation;
261 double Offset;
262
264
265 int ProjectionPlaneMode; //selects the plane in 3D where the Delaunay triangulation will be computed.
266
267private:
268 vtkPolyData *Mesh; //the created mesh
269 double *Points; //the raw points in double precision
270 void SetPoint(vtkIdType id, double *x)
271 {vtkIdType idx=3*id;
272 this->Points[idx] = x[0];
273 this->Points[idx+1] = x[1];
274 this->Points[idx+2] = x[2];
275 }
276
277 void GetPoint(vtkIdType id, double x[3])
278 {double *ptr = this->Points + 3*id;
279 x[0] = *ptr++;
280 x[1] = *ptr++;
281 x[2] = *ptr;
282 }
283
284 int NumberOfDuplicatePoints;
285 int NumberOfDegeneracies;
286
287 int *RecoverBoundary(vtkPolyData *source);
288 int RecoverEdge(vtkPolyData* source, vtkIdType p1, vtkIdType p2);
289 void FillPolygons(vtkCellArray *polys, int *triUse);
290
291 int InCircle (double x[3], double x1[3], double x2[3], double x3[3]);
292 vtkIdType FindTriangle(double x[3], vtkIdType ptIds[3], vtkIdType tri,
293 double tol, vtkIdType nei[3], vtkIdList *neighbors);
294 void CheckEdge(vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2,
295 vtkIdType tri, bool recursive);
296
297 int FillInputPortInformation(int, vtkInformation*) VTK_OVERRIDE;
298
299private:
300 vtkDelaunay2D(const vtkDelaunay2D&) VTK_DELETE_FUNCTION;
301 void operator=(const vtkDelaunay2D&) VTK_DELETE_FUNCTION;
302};
303
304#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:51
create 2D Delaunay triangulation of input points
~vtkDelaunay2D() override
vtkPolyData * GetSource()
Get a pointer to the source object.
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.
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:37
a simple class to control print indentation
Definition: vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
abstract class for specifying dataset behavior
Definition: vtkPointSet.h:43
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:86
Transform
Definition: ADIOSDefs.h:40
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
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
int vtkIdType
Definition: vtkType.h:287
#define VTK_DOUBLE_MAX
Definition: vtkType.h:163