VTK  9.0.20210518
vtkPolygon.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkPolygon.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 =========================================================================*/
25 #ifndef vtkPolygon_h
26 #define vtkPolygon_h
27 
28 #include "vtkCell.h"
29 #include "vtkCommonDataModelModule.h" // For export macro
30 
31 class vtkDoubleArray;
32 class vtkIdTypeArray;
33 class vtkLine;
34 class vtkPoints;
35 class vtkQuad;
36 class vtkTriangle;
38 
39 class VTKCOMMONDATAMODEL_EXPORT vtkPolygon : public vtkCell
40 {
41 public:
42  static vtkPolygon* New();
43  vtkTypeMacro(vtkPolygon, vtkCell);
44  void PrintSelf(ostream& os, vtkIndent indent) override;
45 
47 
50  int GetCellType() override { return VTK_POLYGON; }
51  int GetCellDimension() override { return 2; }
52  int GetNumberOfEdges() override { return this->GetNumberOfPoints(); }
53  int GetNumberOfFaces() override { return 0; }
54  vtkCell* GetEdge(int edgeId) override;
55  vtkCell* GetFace(int) override { return nullptr; }
56  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
57  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
58  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
59  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
60  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
61  vtkCellArray* tris, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
62  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
63  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
64  double& dist2, double weights[]) override;
65  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
66  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
67  double pcoords[3], int& subId) override;
68  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
70  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
71  int IsPrimaryCell() override { return 0; }
73 
80  double ComputeArea();
81 
91  void InterpolateFunctions(const double x[3], double* sf) override;
92 
94 
98  static void ComputeNormal(vtkPoints* p, int numPts, const vtkIdType* pts, double n[3]);
99  static void ComputeNormal(vtkPoints* p, double n[3]);
100  static void ComputeNormal(vtkIdTypeArray* ids, vtkPoints* pts, double n[3]);
102 
107  static void ComputeNormal(int numPts, double* pts, double n[3]);
108 
115  bool IsConvex();
116 
118 
122  static bool IsConvex(vtkPoints* p, int numPts, vtkIdType* pts);
123  static bool IsConvex(vtkIdTypeArray* ids, vtkPoints* p);
124  static bool IsConvex(vtkPoints* p);
126 
128 
132  static bool ComputeCentroid(vtkPoints* p, int numPts, const vtkIdType* pts, double centroid[3]);
133  static bool ComputeCentroid(vtkIdTypeArray* ids, vtkPoints* pts, double centroid[3]);
135 
144  static double ComputeArea(vtkPoints* p, vtkIdType numPts, const vtkIdType* pts, double normal[3]);
145 
154  double p0[3], double p10[3], double& l10, double p20[3], double& l20, double n[3]);
155 
162  static int PointInPolygon(double x[3], int numPts, double* pts, double bounds[6], double n[3]);
163 
172  int Triangulate(vtkIdList* outTris);
173 
179 
187  int BoundedTriangulate(vtkIdList* outTris, double tol);
188 
194  static double DistanceToPolygon(
195  double x[3], int numPts, double* pts, double bounds[6], double closest[3]);
196 
205  static int IntersectPolygonWithPolygon(int npts, double* pts, double bounds[6], int npts2,
206  double* pts2, double bounds2[3], double tol, double x[3]);
207 
220  vtkCell* cell1, vtkCell* cell2, double tol, double p0[3], double p1[3]);
221 
223 
229  vtkGetMacro(UseMVCInterpolation, bool);
230  vtkSetMacro(UseMVCInterpolation, bool);
232 
234 
242  vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
243  vtkGetMacro(Tolerance, double);
245 
246 protected:
248  ~vtkPolygon() override;
249 
250  // Compute the interpolation functions using Mean Value Coordinate.
251  void InterpolateFunctionsUsingMVC(const double x[3], double* weights);
252 
253  // variables used by instances of this class
254  double Tolerance; // Intersection tolerance
255  double Tol;
256  int SuccessfulTriangulation; // Stops recursive tri. if necessary
257  double Normal[3]; // polygon normal
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------------------------------
276 
285 
286 private:
287  vtkPolygon(const vtkPolygon&) = delete;
288  void operator=(const vtkPolygon&) = delete;
289 };
290 
291 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:34
vtkPolygon::IsConvex
static bool IsConvex(vtkPoints *p)
Determine whether or not a polygon is convex.
vtkPolygon::NonDegenerateTriangulate
int NonDegenerateTriangulate(vtkIdList *outTris)
Same as Triangulate(vtkIdList *outTris) but with a first pass to split the polygon into non-degenerat...
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:33
vtkPolygon::GetNumberOfFaces
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:53
vtkX3D::value
@ value
Definition: vtkX3D.h:226
vtkIdType
int vtkIdType
Definition: vtkType.h:332
vtkPolygon::ComputeArea
double ComputeArea()
Compute the area of a polygon.
vtkPolygon::GetFace
vtkCell * GetFace(int) override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:55
vtkPolygon::ComputeCentroid
static bool ComputeCentroid(vtkIdTypeArray *ids, vtkPoints *pts, double centroid[3])
Compute the centroid of a set of points.
vtkPolygon
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:40
vtkPolygon::Triangulate
int Triangulate(vtkIdList *outTris)
Triangulate this polygon.
vtkPolygon::ComputeNormal
static void ComputeNormal(vtkIdTypeArray *ids, vtkPoints *pts, double n[3])
Computes the unit normal to the polygon.
vtkX3D::Normal
@ Normal
Definition: vtkX3D.h:51
vtkPolygon::Contour
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::Tolerance
double Tolerance
Definition: vtkPolygon.h:254
vtkPolygon::IntersectConvex2DCells
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.
vtkPolygon::BoundedTriangulate
int BoundedTriangulate(vtkIdList *outTris, double tol)
Triangulate polygon and enforce that the ratio of the smallest triangle area to the polygon area is g...
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
vtkPolygon::IsConvex
static bool IsConvex(vtkPoints *p, int numPts, vtkIdType *pts)
Determine whether or not a polygon is convex.
vtkPolygon::Triangle
vtkTriangle * Triangle
Definition: vtkPolygon.h:259
vtkPolygon::GetCellType
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:50
vtkPolygon::InterpolateFunctionsUsingMVC
void InterpolateFunctionsUsingMVC(const double x[3], double *weights)
vtkPolygon::Clip
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.
vtkLine
cell represents a 1D line
Definition: vtkLine.h:31
vtkPolygon::IsConvex
static bool IsConvex(vtkIdTypeArray *ids, vtkPoints *p)
Determine whether or not a polygon is convex.
vtkPolygon::CellBoundary
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
See the vtkCell API for descriptions of these methods.
vtkCell.h
vtkPolygon::Tol
double Tol
Definition: vtkPolygon.h:255
vtkPolygon::SuccessfulTriangulation
int SuccessfulTriangulation
Definition: vtkPolygon.h:256
vtkPolygon::UnbiasedEarCutTriangulation
int UnbiasedEarCutTriangulation(int seed)
A fast triangulation method.
vtkPolygon::IntersectWithLine
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.
vtkPolygon::~vtkPolygon
~vtkPolygon() override
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:58
vtkPolygon::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkPolygon::PointInPolygon
static int PointInPolygon(double x[3], int numPts, double *pts, double bounds[6], double n[3])
Determine whether point is inside polygon.
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:34
vtkPolygon::DistanceToPolygon
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.
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:181
vtkPolygon::EvaluatePosition
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.
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:52
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:31
vtkPolygon::TriScalars
vtkDoubleArray * TriScalars
Definition: vtkPolygon.h:261
vtkPolygon::GetEdge
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
vtkPolygon::Tris
vtkIdList * Tris
Definition: vtkPolygon.h:258
vtkPolygon::GetNumberOfEdges
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:52
vtkPolygon::EvaluateLocation
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
See the vtkCell API for descriptions of these methods.
vtkPolygon::UseMVCInterpolation
bool UseMVCInterpolation
Definition: vtkPolygon.h:266
vtkPolygon::EarCutTriangulation
int EarCutTriangulation()
A fast triangulation method.
vtkTriangle
a cell that represents a triangle
Definition: vtkTriangle.h:36
vtkPolygon::InterpolateFunctions
void InterpolateFunctions(const double x[3], double *sf) override
Compute the interpolation functions/derivatives.
vtkPolygon::vtkPolygon
vtkPolygon()
vtkPolygon::ParameterizePolygon
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.
vtkPolygon::IntersectPolygonWithPolygon
static int IntersectPolygonWithPolygon(int npts, double *pts, double bounds[6], int npts2, double *pts2, double bounds2[3], double tol, double x[3])
Method intersects two polygons.
vtkPolygon::Derivatives
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.
vtkIdTypeArray
dynamic, self-adjusting array of vtkIdType
Definition: vtkIdTypeArray.h:36
vtkPolygon::ComputeNormal
static void ComputeNormal(int numPts, double *pts, double n[3])
Compute the polygon normal from an array of points.
VTK_POLYGON
@ VTK_POLYGON
Definition: vtkCellType.h:53
vtkPolygon::IsPrimaryCell
int IsPrimaryCell() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:71
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:36
vtkPolygon::ComputeCentroid
static bool ComputeCentroid(vtkPoints *p, int numPts, const vtkIdType *pts, double centroid[3])
Compute the centroid of a set of points.
vtkPolygon::GetCellDimension
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPolygon.h:51
vtkPolygon::Line
vtkLine * Line
Definition: vtkPolygon.h:262
vtkPolygon::ComputeNormal
static void ComputeNormal(vtkPoints *p, int numPts, const vtkIdType *pts, double n[3])
Computes the unit normal to the polygon.
vtkCell::GetNumberOfPoints
vtkIdType GetNumberOfPoints() const
Return the number of points in the cell.
Definition: vtkCell.h:138
vtkPolygon::ComputeArea
static double ComputeArea(vtkPoints *p, vtkIdType numPts, const vtkIdType *pts, double normal[3])
Compute the area of a polygon in 3D.
vtkPolygon::Triangulate
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
See the vtkCell API for descriptions of these methods.
vtkPolygon::New
static vtkPolygon * New()
vtkX3D::index
@ index
Definition: vtkX3D.h:252
vtkPolygon::ComputeNormal
static void ComputeNormal(vtkPoints *p, double n[3])
Computes the unit normal to the polygon.
vtkQuad
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:36
vtkPolygon::IsConvex
bool IsConvex()
Determine whether or not a polygon is convex.
vtkPolygon::Quad
vtkQuad * Quad
Definition: vtkPolygon.h:260