VTK  9.5.2
vtkPolygon.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
115#ifndef vtkPolygon_h
116#define vtkPolygon_h
117
118#include "vtkCell.h"
119#include "vtkCommonDataModelModule.h" // For export macro
120
121VTK_ABI_NAMESPACE_BEGIN
122class vtkDoubleArray;
123class vtkIdTypeArray;
124class vtkLine;
125class vtkPoints;
126class vtkQuad;
127class vtkTriangle;
129
130class VTKCOMMONDATAMODEL_EXPORT vtkPolygon : public vtkCell
131{
132public:
133 static vtkPolygon* New();
134 vtkTypeMacro(vtkPolygon, vtkCell);
135 void PrintSelf(ostream& os, vtkIndent indent) override;
136
138
141 int GetCellType() override { return VTK_POLYGON; }
142 int GetCellDimension() override { return 2; }
143 int GetNumberOfEdges() override { return this->GetNumberOfPoints(); }
144 int GetNumberOfFaces() override { return 0; }
145 vtkCell* GetEdge(int edgeId) override;
146 vtkCell* GetFace(int) override { return nullptr; }
147 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
148 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
149 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
150 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
151 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
152 vtkCellArray* tris, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
153 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
154 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
155 double& dist2, double weights[]) override;
156 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
157 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
158 double pcoords[3], int& subId) override;
159 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
161 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
162 int IsPrimaryCell() VTK_FUTURE_CONST override { return 0; }
164
171 double ComputeArea();
172
182 void InterpolateFunctions(const double x[3], double* sf) override;
183
185
189 static void ComputeNormal(vtkPoints* p, int numPts, const vtkIdType* pts, double n[3]);
190 static void ComputeNormal(vtkPoints* p, double n[3]);
191 static void ComputeNormal(vtkIdTypeArray* ids, vtkPoints* pts, double n[3]);
193
198 static void ComputeNormal(int numPts, double* pts, double n[3]);
199
206 bool IsConvex();
207
209
213 static bool IsConvex(vtkPoints* p, int numPts, const vtkIdType* pts);
214 static bool IsConvex(vtkIdTypeArray* ids, vtkPoints* p);
215 static bool IsConvex(vtkPoints* p);
217
219
243 static bool ComputeCentroid(
244 vtkPoints* p, int numPts, const vtkIdType* pts, double centroid[3], double tolerance);
245 static bool ComputeCentroid(vtkPoints* p, int numPts, const vtkIdType* pts, double centroid[3]);
246 static bool ComputeCentroid(vtkIdTypeArray* ids, vtkPoints* pts, double centroid[3]);
248
257 static double ComputeArea(vtkPoints* p, vtkIdType numPts, const vtkIdType* pts, double normal[3]);
258
267 double p0[3], double p10[3], double& l10, double p20[3], double& l20, double n[3]);
268
280 static int PointInPolygon(double x[3], int numPts, double* pts, double bounds[6], double n[3]);
281
282 // Needed to remove warning "member function does not override any
283 // base class virtual member function"
284 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override
285 {
286 return vtkCell::Triangulate(index, ptIds, pts);
287 }
288
294
302 int BoundedTriangulate(vtkIdList* outTris, double tol);
303
309 static double DistanceToPolygon(
310 double x[3], int numPts, double* pts, double bounds[6], double closest[3]);
311
320 static int IntersectPolygonWithPolygon(int npts, double* pts, double bounds[6], int npts2,
321 double* pts2, double bounds2[6], double tol, double x[3]);
322
335 vtkCell* cell1, vtkCell* cell2, double tol, double p0[3], double p1[3]);
336
338
344 vtkGetMacro(UseMVCInterpolation, bool);
345 vtkSetMacro(UseMVCInterpolation, bool);
347
349
357 vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
358 vtkGetMacro(Tolerance, double);
360
361protected:
363 ~vtkPolygon() override;
364
365 // Compute the interpolation functions using Mean Value Coordinate.
366 void InterpolateFunctionsUsingMVC(const double x[3], double* weights);
367
368 // variables used by instances of this class
369 double Tolerance; // Intersection tolerance set by public API
370 double Tol; // Internal tolerance set by ComputeBounds()
371 void ComputeTolerance(); // Compute the internal tolerance Tol
372
373 int SuccessfulTriangulation; // Stops recursive triangulation if necessary
374 vtkIdList* Tris; // Output triangulation placed here
375
376 // These are used for internal computation.
381
382 // Parameter indicating whether to use Mean Value Coordinate algorithm
383 // for interpolation. The parameter is false by default.
385
386 // Helper methods for triangulation------------------------------
387 // Made public for external access
388public:
389 // Ear cut triangulation options. The order in which vertices are
390 // removed are controlled by different measures. Changing this can
391 // make subtle differences in some cases. Historically the
392 // PERIMETER2_TO_AREA_RATIO has been used.
394 {
395 PERIMETER2_TO_AREA_RATIO = 0,
396 DOT_PRODUCT = 1,
397 BEST_QUALITY = 2
398 };
399
401
409 int EarCutTriangulation(int measure = PERIMETER2_TO_AREA_RATIO);
410 int EarCutTriangulation(vtkIdList* outTris, int measure = PERIMETER2_TO_AREA_RATIO);
412
414
421 int UnbiasedEarCutTriangulation(int seed, int measure = PERIMETER2_TO_AREA_RATIO);
423 int seed, vtkIdList* outTris, int measure = PERIMETER2_TO_AREA_RATIO);
425
426private:
427 vtkPolygon(const vtkPolygon&) = delete;
428 void operator=(const vtkPolygon&) = delete;
429};
430
431VTK_ABI_NAMESPACE_END
432#endif
object to represent cell connectivity
represent and manipulate cell attribute data
abstract class to specify cell behavior
Definition vtkCell.h:130
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)
Generate simplices of proper dimension.
vtkIdType GetNumberOfPoints() const
Return the number of points in the cell.
Definition vtkCell.h:213
abstract superclass for arrays of numeric data
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:133
dynamic, self-adjusting array of vtkIdType
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:108
cell represents a 1D line
Definition vtkLine.h:132
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:139
a cell that represents an n-sided polygon
Definition vtkPolygon.h:131
static int PointInPolygon(double x[3], int numPts, double *pts, double bounds[6], double n[3])
Determine whether a point is inside the specified polygon.
double ComputeArea()
Compute the area of a polygon.
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition vtkPolygon.h:141
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tris, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
See the vtkCell API for descriptions of these methods.
vtkTriangle * Triangle
Definition vtkPolygon.h:377
static bool IsConvex(vtkPoints *p, int numPts, const vtkIdType *pts)
Determine whether or not a polygon is convex.
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
See the vtkCell API for descriptions of these methods.
static void ComputeNormal(int numPts, double *pts, double n[3])
Compute the polygon normal from an array of points.
int UnbiasedEarCutTriangulation(int seed, int measure=PERIMETER2_TO_AREA_RATIO)
A fast triangulation method.
static int IntersectPolygonWithPolygon(int npts, double *pts, double bounds[6], int npts2, double *pts2, double bounds2[6], double tol, double x[3])
Method intersects two polygons.
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition vtkPolygon.h:143
static bool IsConvex(vtkPoints *p)
Determine whether or not a polygon is convex.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
Definition vtkPolygon.h:284
void ComputeTolerance()
bool UseMVCInterpolation
Definition vtkPolygon.h:384
static void ComputeNormal(vtkPoints *p, double n[3])
Computes the unit normal to the polygon.
double Tolerance
Definition vtkPolygon.h:369
void InterpolateFunctionsUsingMVC(const double x[3], double *weights)
int IsPrimaryCell() VTK_FUTURE_CONST override
See the vtkCell API for descriptions of these methods.
Definition vtkPolygon.h:162
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
Definition vtkPolygon.h:144
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
See the vtkCell API for descriptions of these methods.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int NonDegenerateTriangulate(vtkIdList *outTris)
Same as Triangulate(vtkIdList *outTris) but with a first pass to split the polygon into non-degenerat...
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition vtkPolygon.h:142
static double DistanceToPolygon(double x[3], int numPts, double *pts, double bounds[6], double closest[3])
Compute the distance of a point to a polygon.
int SuccessfulTriangulation
Definition vtkPolygon.h:373
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
See the vtkCell API for descriptions of these methods.
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
See the vtkCell API for descriptions of these methods.
vtkQuad * Quad
Definition vtkPolygon.h:378
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
See the vtkCell API for descriptions of these methods.
vtkIdList * Tris
Definition vtkPolygon.h:374
void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override
See the vtkCell API for descriptions of these methods.
~vtkPolygon() override
static bool IsConvex(vtkIdTypeArray *ids, vtkPoints *p)
Determine whether or not a polygon is convex.
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
See the vtkCell API for descriptions of these methods.
int ParameterizePolygon(double p0[3], double p10[3], double &l10, double p20[3], double &l20, double n[3])
Create a local s-t coordinate system for a polygon.
int EarCutTriangulation(vtkIdList *outTris, int measure=PERIMETER2_TO_AREA_RATIO)
A fast triangulation method.
vtkDoubleArray * TriScalars
Definition vtkPolygon.h:379
static void ComputeNormal(vtkIdTypeArray *ids, vtkPoints *pts, double n[3])
Computes the unit normal to the polygon.
int BoundedTriangulate(vtkIdList *outTris, double tol)
Triangulate polygon and enforce that the ratio of the smallest triangle area to the polygon area is g...
static bool ComputeCentroid(vtkPoints *p, int numPts, const vtkIdType *pts, double centroid[3], double tolerance)
Compute the centroid of a set of points.
int EarCutTriangulation(int measure=PERIMETER2_TO_AREA_RATIO)
A fast triangulation method.
vtkCell * GetFace(int) override
See the vtkCell API for descriptions of these methods.
Definition vtkPolygon.h:146
static bool ComputeCentroid(vtkPoints *p, int numPts, const vtkIdType *pts, double centroid[3])
Compute the centroid of a set of points.
static double ComputeArea(vtkPoints *p, vtkIdType numPts, const vtkIdType *pts, double normal[3])
Compute the area of a polygon in 3D.
void InterpolateFunctions(const double x[3], double *sf) override
Compute the interpolation functions/derivatives.
bool IsConvex()
Determine whether or not a polygon is convex.
vtkLine * Line
Definition vtkPolygon.h:380
static bool ComputeCentroid(vtkIdTypeArray *ids, vtkPoints *pts, double centroid[3])
Compute the centroid of a set of points.
int UnbiasedEarCutTriangulation(int seed, vtkIdList *outTris, int measure=PERIMETER2_TO_AREA_RATIO)
A fast triangulation method.
static vtkPolygon * New()
double Tol
Definition vtkPolygon.h:370
static void ComputeNormal(vtkPoints *p, int numPts, const vtkIdType *pts, double n[3])
Computes the unit normal to the polygon.
static int IntersectConvex2DCells(vtkCell *cell1, vtkCell *cell2, double tol, double p0[3], double p1[3])
Intersect two convex 2D polygons to produce a line segment as output.
a cell that represents a 2D quadrilateral
Definition vtkQuad.h:87
a cell that represents a triangle
@ VTK_POLYGON
Definition vtkCellType.h:44
int vtkIdType
Definition vtkType.h:332