VTK  9.3.1
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
16 #ifndef vtkPolygon_h
17 #define vtkPolygon_h
18 
19 #include "vtkCell.h"
20 #include "vtkCommonDataModelModule.h" // For export macro
21 
22 VTK_ABI_NAMESPACE_BEGIN
23 class vtkDoubleArray;
24 class vtkIdTypeArray;
25 class vtkLine;
26 class vtkPoints;
27 class vtkQuad;
28 class vtkTriangle;
30 
31 class VTKCOMMONDATAMODEL_EXPORT vtkPolygon : public vtkCell
32 {
33 public:
34  static vtkPolygon* New();
35  vtkTypeMacro(vtkPolygon, vtkCell);
36  void PrintSelf(ostream& os, vtkIndent indent) override;
37 
39 
42  int GetCellType() override { return VTK_POLYGON; }
43  int GetCellDimension() override { return 2; }
44  int GetNumberOfEdges() override { return this->GetNumberOfPoints(); }
45  int GetNumberOfFaces() override { return 0; }
46  vtkCell* GetEdge(int edgeId) override;
47  vtkCell* GetFace(int) override { return nullptr; }
48  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
49  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
50  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
51  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
52  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
53  vtkCellArray* tris, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
54  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
55  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
56  double& dist2, double weights[]) override;
57  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
58  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
59  double pcoords[3], int& subId) override;
60  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
61  void Derivatives(
62  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
63  int IsPrimaryCell() override { return 0; }
65 
72  double ComputeArea();
73 
83  void InterpolateFunctions(const double x[3], double* sf) override;
84 
86 
90  static void ComputeNormal(vtkPoints* p, int numPts, const vtkIdType* pts, double n[3]);
91  static void ComputeNormal(vtkPoints* p, double n[3]);
92  static void ComputeNormal(vtkIdTypeArray* ids, vtkPoints* pts, double n[3]);
94 
99  static void ComputeNormal(int numPts, double* pts, double n[3]);
100 
107  bool IsConvex();
108 
110 
114  static bool IsConvex(vtkPoints* p, int numPts, const vtkIdType* pts);
115  static bool IsConvex(vtkIdTypeArray* ids, vtkPoints* p);
116  static bool IsConvex(vtkPoints* p);
118 
120 
124  static bool ComputeCentroid(vtkPoints* p, int numPts, const vtkIdType* pts, double centroid[3]);
125  static bool ComputeCentroid(vtkIdTypeArray* ids, vtkPoints* pts, double centroid[3]);
127 
136  static double ComputeArea(vtkPoints* p, vtkIdType numPts, const vtkIdType* pts, double normal[3]);
137 
145  int ParameterizePolygon(
146  double p0[3], double p10[3], double& l10, double p20[3], double& l20, double n[3]);
147 
159  static int PointInPolygon(double x[3], int numPts, double* pts, double bounds[6], double n[3]);
160 
169  int Triangulate(vtkIdList* outTris);
170 
175  int NonDegenerateTriangulate(vtkIdList* outTris);
176 
184  int BoundedTriangulate(vtkIdList* outTris, double tol);
185 
191  static double DistanceToPolygon(
192  double x[3], int numPts, double* pts, double bounds[6], double closest[3]);
193 
202  static int IntersectPolygonWithPolygon(int npts, double* pts, double bounds[6], int npts2,
203  double* pts2, double bounds2[6], double tol, double x[3]);
204 
216  static int IntersectConvex2DCells(
217  vtkCell* cell1, vtkCell* cell2, double tol, double p0[3], double p1[3]);
218 
220 
226  vtkGetMacro(UseMVCInterpolation, bool);
227  vtkSetMacro(UseMVCInterpolation, bool);
229 
231 
239  vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
240  vtkGetMacro(Tolerance, double);
242 
243 protected:
244  vtkPolygon();
245  ~vtkPolygon() override;
246 
247  // Compute the interpolation functions using Mean Value Coordinate.
248  void InterpolateFunctionsUsingMVC(const double x[3], double* weights);
249 
250  // variables used by instances of this class
251  double Tolerance; // Intersection tolerance set by public API
252  double Tol; // Internal tolerance set by ComputeBounds()
253  void ComputeTolerance(); // Compute the internal tolerance Tol
254 
255  int SuccessfulTriangulation; // Stops recursive triangulation if necessary
256  vtkIdList* Tris; // Output triangulation placed here
257 
258  // These are used for internal computation.
263 
264  // Parameter indicating whether to use Mean Value Coordinate algorithm
265  // for interpolation. The parameter is false by default.
267 
268  // Helper methods for triangulation------------------------------
269  // Made public for external access
270 public:
271  // Ear cut triangulation options. The order in which vertices are
272  // removed are controlled by different measures. Changing this can
273  // make subtle differences in some cases. Historically the
274  // PERIMETER2_TO_AREA_RATIO has been used.
276  {
277  PERIMETER2_TO_AREA_RATIO = 0,
278  DOT_PRODUCT = 1,
279  BEST_QUALITY = 2
280  };
281 
283 
291  int EarCutTriangulation(int measure = PERIMETER2_TO_AREA_RATIO);
292  int EarCutTriangulation(vtkIdList* outTris, int measure = PERIMETER2_TO_AREA_RATIO);
294 
296 
303  int UnbiasedEarCutTriangulation(int seed, int measure = PERIMETER2_TO_AREA_RATIO);
304  int UnbiasedEarCutTriangulation(
305  int seed, vtkIdList* outTris, int measure = PERIMETER2_TO_AREA_RATIO);
307 
308 private:
309  vtkPolygon(const vtkPolygon&) = delete;
310  void operator=(const vtkPolygon&) = delete;
311 };
312 
313 VTK_ABI_NAMESPACE_END
314 #endif
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
represent and manipulate point attribute data
Definition: vtkPointData.h:29
vtkIdType GetNumberOfPoints() const
Return the number of points in the cell.
Definition: vtkCell.h:130
double Tolerance
Definition: vtkPolygon.h:251
represent and manipulate cell attribute data
Definition: vtkCellData.h:30
Abstract class in support of both point location and point insertion.
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
virtual void InterpolateFunctions(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(weight))
Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this leve...
Definition: vtkCell.h:380
double Tol
Definition: vtkPolygon.h:252
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:27
dynamic, self-adjusting array of vtkIdType
int vtkIdType
Definition: vtkType.h:315
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:44
virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
dynamic, self-adjusting array of double
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:45
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:42
cell represents a 1D line
Definition: vtkLine.h:22
abstract class to specify cell behavior
Definition: vtkCell.h:49
virtual void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkDoubleArray * TriScalars
Definition: vtkPolygon.h:261
a simple class to control print indentation
Definition: vtkIndent.h:28
list of point or cell ids
Definition: vtkIdList.h:22
int SuccessfulTriangulation
Definition: vtkPolygon.h:255
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:44
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:31
bool UseMVCInterpolation
Definition: vtkPolygon.h:266
virtual void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)=0
Cut (or clip) the cell based on the input cellScalars and the specified value.
vtkCell * GetFace(int) override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:47
virtual int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[])=0
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
object to represent cell connectivity
Definition: vtkCellArray.h:175
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
vtkLine * Line
Definition: vtkPolygon.h:262
int IsPrimaryCell() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:63
a cell that represents a triangle
Definition: vtkTriangle.h:27
vtkIdList * Tris
Definition: vtkPolygon.h:256
virtual void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd)=0
Generate contouring primitives.
virtual void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
vtkQuad * Quad
Definition: vtkPolygon.h:260
vtkTriangle * Triangle
Definition: vtkPolygon.h:259
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:43
represent and manipulate 3D points
Definition: vtkPoints.h:28