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=========================================================================*/
124
125#ifndef vtkDelaunay2D_h
126#define vtkDelaunay2D_h
127
128#include "vtkFiltersCoreModule.h" // For export macro
129#include "vtkPolyDataAlgorithm.h"
130
132class vtkCellArray;
133class vtkIdList;
134class vtkPointSet;
135
136#define VTK_DELAUNAY_XY_PLANE 0
137#define VTK_SET_TRANSFORM_PLANE 1
138#define VTK_BEST_FITTING_PLANE 2
139
140class VTKFILTERSCORE_EXPORT vtkDelaunay2D : public vtkPolyDataAlgorithm
141{
142public:
144 void PrintSelf(ostream& os, vtkIndent indent) override;
145
151
162
172
177
179
185 vtkSetClampMacro(Alpha, double, 0.0, VTK_DOUBLE_MAX);
186 vtkGetMacro(Alpha, double);
188
190
195 vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
196 vtkGetMacro(Tolerance, double);
198
200
204 vtkSetClampMacro(Offset, double, 0.75, VTK_DOUBLE_MAX);
205 vtkGetMacro(Offset, double);
207
209
219
221
232 vtkGetObjectMacro(Transform, vtkAbstractTransform);
234
236
245 vtkGetMacro(ProjectionPlaneMode, int);
247
255
256protected:
258 ~vtkDelaunay2D() override;
259
261
262 double Alpha;
263 double Tolerance;
265 double Offset;
266
268
269 int ProjectionPlaneMode; // selects the plane in 3D where the Delaunay triangulation will be
270 // computed.
271
272private:
273 vtkPolyData* Mesh; // the created mesh
274 double* Points; // the raw points in double precision
275 void SetPoint(vtkIdType id, double* x)
276 {
277 vtkIdType idx = 3 * id;
278 this->Points[idx] = x[0];
279 this->Points[idx + 1] = x[1];
280 this->Points[idx + 2] = x[2];
281 }
282
283 void GetPoint(vtkIdType id, double x[3])
284 {
285 double* ptr = this->Points + 3 * id;
286 x[0] = *ptr++;
287 x[1] = *ptr++;
288 x[2] = *ptr;
289 }
290
291 int NumberOfDuplicatePoints;
292 int NumberOfDegeneracies;
293
294 int* RecoverBoundary(vtkPolyData* source);
295 int RecoverEdge(vtkPolyData* source, vtkIdType p1, vtkIdType p2);
296 void FillPolygons(vtkCellArray* polys, int* triUse);
297
298 int InCircle(double x[3], double x1[3], double x2[3], double x3[3]);
299 vtkIdType FindTriangle(double x[3], vtkIdType ptIds[3], vtkIdType tri, double tol,
300 vtkIdType nei[3], vtkIdList* neighbors);
301 bool CheckEdge(
302 vtkIdType ptId, double x[3], vtkIdType p1, vtkIdType p2, vtkIdType tri, bool recursive);
303
304 int FillInputPortInformation(int, vtkInformation*) override;
305
306private:
307 vtkDelaunay2D(const vtkDelaunay2D&) = delete;
308 void operator=(const vtkDelaunay2D&) = delete;
309};
310
311#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
~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:31
a simple class to control print indentation
Definition vtkIndent.h:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
concrete class for storing a set of points
Definition vtkPointSet.h:67
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:85
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