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 =========================================================================*/
130 #ifndef vtkDelaunay2D_h
131 #define vtkDelaunay2D_h
132 
133 #include "vtkFiltersCoreModule.h" // For export macro
134 #include "vtkPolyDataAlgorithm.h"
135 
137 class vtkCellArray;
138 class vtkIdList;
139 class vtkPointSet;
140 
141 #define VTK_DELAUNAY_XY_PLANE 0
142 #define VTK_SET_TRANSFORM_PLANE 1
143 #define VTK_BEST_FITTING_PLANE 2
144 
146 {
147 public:
149  void PrintSelf(ostream& os, vtkIndent indent);
150 
153  static vtkDelaunay2D *New();
154 
162  void SetSourceData(vtkPolyData *);
163 
170  void SetSourceConnection(vtkAlgorithmOutput *algOutput);
171 
173  vtkPolyData *GetSource();
174 
176 
180  vtkSetClampMacro(Alpha,double,0.0,VTK_DOUBLE_MAX);
181  vtkGetMacro(Alpha,double);
183 
185 
188  vtkSetClampMacro(Tolerance,double,0.0,1.0);
189  vtkGetMacro(Tolerance,double);
191 
193 
195  vtkSetClampMacro(Offset,double,0.75,VTK_DOUBLE_MAX);
196  vtkGetMacro(Offset,double);
198 
200 
204  vtkSetMacro(BoundingTriangulation,int);
205  vtkGetMacro(BoundingTriangulation,int);
206  vtkBooleanMacro(BoundingTriangulation,int);
208 
210 
218  virtual void SetTransform(vtkAbstractTransform*);
219  vtkGetObjectMacro(Transform, vtkAbstractTransform);
221 
223 
224  vtkSetClampMacro(ProjectionPlaneMode,int,
226  vtkGetMacro(ProjectionPlaneMode,int);
228 
229 protected:
230  vtkDelaunay2D();
231  ~vtkDelaunay2D();
232 
234 
235  vtkAbstractTransform * ComputeBestFittingPlane(vtkPointSet *input);
236 
237  double Alpha;
238  double Tolerance;
240  double Offset;
241 
243 
244  int ProjectionPlaneMode; //selects the plane in 3D where the Delaunay triangulation will be computed.
245 
246 private:
247  vtkPolyData *Mesh; //the created mesh
248  double *Points; //the raw points in double precision
249  void SetPoint(vtkIdType id, double *x)
250  {vtkIdType idx=3*id;
251  this->Points[idx] = x[0];
252  this->Points[idx+1] = x[1];
253  this->Points[idx+2] = x[2];
254  }
255 
256  void GetPoint(vtkIdType id, double x[3])
257  {double *ptr = this->Points + 3*id;
258  x[0] = *ptr++;
259  x[1] = *ptr++;
260  x[2] = *ptr;
261  }
262 
263  int NumberOfDuplicatePoints;
264  int NumberOfDegeneracies;
265 
266  int *RecoverBoundary(vtkPolyData *source);
267  int RecoverEdge(vtkIdType p1, vtkIdType p2);
268  void FillPolygons(vtkCellArray *polys, int *triUse);
269 
270  int InCircle (double x[3], double x1[3], double x2[3], double x3[3]);
271  vtkIdType FindTriangle(double x[3], vtkIdType ptIds[3], vtkIdType tri,
272  double tol, vtkIdType nei[3], vtkIdList *neighbors);
273  void CheckEdge(vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2,
274  vtkIdType tri, bool recursive);
275 
276  virtual int FillInputPortInformation(int, vtkInformation*);
277 
278 private:
279  vtkDelaunay2D(const vtkDelaunay2D&); // Not implemented.
280  void operator=(const vtkDelaunay2D&); // Not implemented.
281 };
282 
283 #endif
virtual int FillInputPortInformation(int port, vtkInformation *info)
#define VTK_DOUBLE_MAX
Definition: vtkType.h:142
Store vtkAlgorithm input/output information.
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
#define VTKFILTERSCORE_EXPORT
abstract class for specifying dataset behavior
Definition: vtkPointSet.h:44
int vtkIdType
Definition: vtkType.h:275
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:83
Proxy object to connect input/output ports.
static vtkPolyDataAlgorithm * New()
void PrintSelf(ostream &os, vtkIndent indent)
Superclass for algorithms that produce only polydata as output.
a simple class to control print indentation
Definition: vtkIndent.h:38
list of point or cell ids
Definition: vtkIdList.h:35
superclass for all geometric transformations
vtkAbstractTransform * Transform
#define VTK_BEST_FITTING_PLANE
boost::graph_traits< vtkGraph * >::vertex_descriptor source(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
object to represent cell connectivity
Definition: vtkCellArray.h:49
create 2D Delaunay triangulation of input points
Store zero or more vtkInformation instances.
#define VTK_DELAUNAY_XY_PLANE
int BoundingTriangulation