VTK
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkPolygon Class Reference

a cell that represents an n-sided polygon More...

#include <vtkPolygon.h>

Inheritance diagram for vtkPolygon:
[legend]
Collaboration diagram for vtkPolygon:
[legend]

Public Types

typedef vtkCell Superclass
 
- Public Types inherited from vtkCell
typedef vtkObject Superclass
 

Public Member Functions

virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class. More...
 
vtkPolygonNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. More...
 
int GetCellType () override
 See the vtkCell API for descriptions of these methods. More...
 
int GetCellDimension () override
 Return the topological dimensional of the cell (0,1,2, or 3). More...
 
int GetNumberOfEdges () override
 Return the number of edges in the cell. More...
 
int GetNumberOfFaces () override
 Return the number of faces in the cell. More...
 
vtkCellGetEdge (int edgeId) override
 Return the edge cell from the edgeId of the cell. More...
 
vtkCellGetFace (int) override
 Return the face cell from the faceId of the cell. More...
 
int CellBoundary (int subId, const double pcoords[3], vtkIdList *pts) override
 Given parametric coordinates of a point, return the closest cell boundary, and whether the point is inside or outside of the cell. More...
 
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
 Generate contouring primitives. More...
 
void Clip (double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tris, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
 Cut (or clip) the cell based on the input cellScalars and the specified value. More...
 
int EvaluatePosition (const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
 Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; evaluate parametric coordinates, sub-cell id (!=0 only if cell is composite), distance squared of point x[3] to cell (in particular, the sub-cell indicated), closest point on cell to x[3] (unless closestPoint is null, in which case, the closest point and dist2 are not found), and interpolation weights in cell. More...
 
void EvaluateLocation (int &subId, const double pcoords[3], double x[3], double *weights) override
 Determine global coordinate (x[3]) from subId and parametric coordinates. More...
 
int IntersectWithLine (const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
 Intersect with a ray. More...
 
int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts) override
 Generate simplices of proper dimension. More...
 
void Derivatives (int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
 Compute derivatives given cell subId and parametric coordinates. More...
 
int IsPrimaryCell () override
 Return whether this cell type has a fixed topology or whether the topology varies depending on the data (e.g., vtkConvexPointSet). More...
 
double ComputeArea ()
 Compute the area of a polygon. More...
 
void InterpolateFunctions (const double x[3], double *sf) override
 Compute the interpolation functions/derivatives. More...
 
bool IsConvex ()
 Determine whether or not a polygon is convex. More...
 
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. More...
 
int Triangulate (vtkIdList *outTris)
 Triangulate this polygon. More...
 
int NonDegenerateTriangulate (vtkIdList *outTris)
 Same as Triangulate(vtkIdList *outTris) but with a first pass to split the polygon into non-degenerate polygons. More...
 
int BoundedTriangulate (vtkIdList *outTris, double tol)
 Triangulate polygon and enforce that the ratio of the smallest triangle area to the polygon area is greater than a user-defined tolerance. More...
 
virtual bool GetUseMVCInterpolation ()
 Set/Get the flag indicating whether to use Mean Value Coordinate for the interpolation. More...
 
virtual void SetUseMVCInterpolation (bool)
 
- Public Member Functions inherited from vtkCell
vtkCellNewInstance () const
 
void Initialize (int npts, vtkIdType *pts, vtkPoints *p)
 Initialize cell from outside with point ids and point coordinates specified. More...
 
void Initialize (int npts, vtkPoints *p)
 Initialize the cell with point coordinates specified. More...
 
virtual void ShallowCopy (vtkCell *c)
 Copy this cell by reference counting the internal data structures. More...
 
virtual void DeepCopy (vtkCell *c)
 Copy this cell by completely copying internal data structures. More...
 
virtual int IsLinear ()
 Non-linear cells require special treatment beyond the usual cell type and connectivity list information. More...
 
virtual int RequiresInitialization ()
 Some cells require initialization prior to access. More...
 
virtual void Initialize ()
 
virtual int IsExplicitCell ()
 Explicit cells require additional representational information beyond the usual cell type and connectivity list information. More...
 
virtual int RequiresExplicitFaceRepresentation ()
 Determine whether the cell requires explicit face representation, and methods for setting and getting the faces (see vtkPolyhedron for example usage of these methods). More...
 
virtual void SetFaces (vtkIdType *vtkNotUsed(faces))
 
virtual vtkIdTypeGetFaces ()
 
vtkPointsGetPoints ()
 Get the point coordinates for the cell. More...
 
vtkIdType GetNumberOfPoints ()
 Return the number of points in the cell. More...
 
vtkIdListGetPointIds ()
 Return the list of point ids defining the cell. More...
 
vtkIdType GetPointId (int ptId)
 For cell point i, return the actual point id. More...
 
void GetBounds (double bounds[6])
 Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). More...
 
doubleGetBounds ()
 Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). More...
 
double GetLength2 ()
 Compute Length squared of cell (i.e., bounding box diagonal squared). More...
 
virtual int GetParametricCenter (double pcoords[3])
 Return center of the cell in parametric coordinates. More...
 
virtual double GetParametricDistance (const double pcoords[3])
 Return the distance of the parametric coordinate provided to the cell. More...
 
virtual doubleGetParametricCoords ())
 Return a contiguous array of parametric coordinates of the points defining this cell. More...
 
virtual void InterpolateFunctions (const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(weight))
 Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this level. More...
 
virtual void InterpolateDerivs (const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 Turn debugging output on. More...
 
virtual void DebugOff ()
 Turn debugging output off. More...
 
bool GetDebug ()
 Get the value of the debug flag. More...
 
void SetDebug (bool debugFlag)
 Set the value of the debug flag. More...
 
virtual void Modified ()
 Update the modification time for this object. More...
 
virtual vtkMTimeType GetMTime ()
 Return this object's modified time. More...
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 
vtkCommandGetCommand (unsigned long tag)
 
void RemoveObserver (vtkCommand *)
 
void RemoveObservers (unsigned long event, vtkCommand *)
 
void RemoveObservers (const char *event, vtkCommand *)
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 
void RemoveObserver (unsigned long tag)
 
void RemoveObservers (unsigned long event)
 
void RemoveObservers (const char *event)
 
void RemoveAllObservers ()
 
vtkTypeBool HasObserver (unsigned long event)
 
vtkTypeBool HasObserver (const char *event)
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Allow user to set the AbortFlagOn() with the return value of the callback method. More...
 
int InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
int InvokeEvent (const char *event, void *callData)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string. More...
 
virtual void Delete ()
 Delete a VTK object. More...
 
virtual void FastDelete ()
 Delete a reference to this object. More...
 
void InitializeObjectBase ()
 
void Print (ostream &os)
 Print an object to an ostream. More...
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 
virtual void Register (vtkObjectBase *o)
 Increase the reference count (mark as used by another object). More...
 
virtual void UnRegister (vtkObjectBase *o)
 Decrease the reference count (release by another object). More...
 
int GetReferenceCount ()
 Return the current reference count of this object. More...
 
void SetReferenceCount (int)
 Sets the reference count. More...
 
void PrintRevisions (ostream &)
 Legacy. More...
 

Static Public Member Functions

static vtkPolygonNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkPolygonSafeDownCast (vtkObjectBase *o)
 
static void ComputeNormal (vtkPoints *p, int numPts, vtkIdType *pts, double n[3])
 Computes the unit normal to the polygon. More...
 
static void ComputeNormal (vtkPoints *p, double n[3])
 
static void ComputeNormal (vtkIdTypeArray *ids, vtkPoints *pts, double n[3])
 
static void ComputeNormal (int numPts, double *pts, double n[3])
 Compute the polygon normal from an array of points. More...
 
static bool IsConvex (vtkPoints *p, int numPts, vtkIdType *pts)
 Determine whether or not a polygon is convex. More...
 
static bool IsConvex (vtkIdTypeArray *ids, vtkPoints *p)
 
static bool IsConvex (vtkPoints *p)
 
static bool ComputeCentroid (vtkPoints *p, int numPts, vtkIdType *pts, double centroid[3])
 Compute the centroid of a set of points. More...
 
static bool ComputeCentroid (vtkIdTypeArray *ids, vtkPoints *pts, double centroid[3])
 
static double ComputeArea (vtkPoints *p, vtkIdType numPts, vtkIdType *pts, double normal[3])
 Compute the area of a polygon in 3D. More...
 
static int PointInPolygon (double x[3], int numPts, double *pts, double bounds[6], double n[3])
 Determine whether point is inside polygon. More...
 
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. More...
 
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. More...
 
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. More...
 
- Static Public Member Functions inherited from vtkCell
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkCellSafeDownCast (vtkObjectBase *o)
 
- Static Public Member Functions inherited from vtkObject
static vtkObjectNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 
static void BreakOnError ()
 This method is called when vtkErrorMacro executes. More...
 
static void SetGlobalWarningDisplay (int val)
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOn ()
 
static void GlobalWarningDisplayOff ()
 
static int GetGlobalWarningDisplay ()
 
- Static Public Member Functions inherited from vtkObjectBase
static vtkTypeBool IsTypeOf (const char *name)
 Return 1 if this class type is the same type of (or a subclass of) the named class. More...
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 

Protected Member Functions

virtual vtkObjectBaseNewInstanceInternal () const
 
 vtkPolygon ()
 
 ~vtkPolygon () override
 
void InterpolateFunctionsUsingMVC (const double x[3], double *weights)
 
int EarCutTriangulation ()
 A fast triangulation method. More...
 
int UnbiasedEarCutTriangulation (int seed)
 A fast triangulation method. More...
 
- Protected Member Functions inherited from vtkCell
 vtkCell ()
 
 ~vtkCell () override
 
- Protected Member Functions inherited from vtkObject
 vtkObject ()
 
 ~vtkObject () override
 
void RegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=nullptr)
 These methods allow a command to exclusively grab all events. More...
 
void InternalReleaseFocus ()
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void CollectRevisions (ostream &)
 
virtual void ReportReferences (vtkGarbageCollector *)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

double Tolerance
 
int SuccessfulTriangulation
 
double Normal [3]
 
vtkIdListTris
 
vtkTriangleTriangle
 
vtkQuadQuad
 
vtkDoubleArrayTriScalars
 
vtkLineLine
 
bool UseMVCInterpolation
 
- Protected Attributes inherited from vtkCell
double Bounds [6]
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
- Protected Attributes inherited from vtkObjectBase
vtkAtomicInt32 ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Additional Inherited Members

- Public Attributes inherited from vtkCell
vtkPointsPoints
 
vtkIdListPointIds
 

Detailed Description

a cell that represents an n-sided polygon

vtkPolygon is a concrete implementation of vtkCell to represent a 2D n-sided polygon. The polygons cannot have any internal holes, and cannot self-intersect. Define the polygon with n-points ordered in the counter- clockwise direction; do not repeat the last point.

Examples:
vtkPolygon (Examples)
Tests:
vtkPolygon (Tests)

Definition at line 45 of file vtkPolygon.h.

Member Typedef Documentation

◆ Superclass

Definition at line 49 of file vtkPolygon.h.

Constructor & Destructor Documentation

◆ vtkPolygon()

vtkPolygon::vtkPolygon ( )
protected

◆ ~vtkPolygon()

vtkPolygon::~vtkPolygon ( )
overrideprotected

Member Function Documentation

◆ New()

static vtkPolygon* vtkPolygon::New ( )
static

◆ IsTypeOf()

static vtkTypeBool vtkPolygon::IsTypeOf ( const char *  type)
static

◆ IsA()

virtual vtkTypeBool vtkPolygon::IsA ( const char *  name)
virtual

Return 1 if this class is the same type of (or a subclass of) the named class.

Returns 0 otherwise. This method works in combination with vtkTypeMacro found in vtkSetGet.h.

Reimplemented from vtkCell.

◆ SafeDownCast()

static vtkPolygon* vtkPolygon::SafeDownCast ( vtkObjectBase o)
static

◆ NewInstanceInternal()

virtual vtkObjectBase* vtkPolygon::NewInstanceInternal ( ) const
protectedvirtual

Reimplemented from vtkCell.

◆ NewInstance()

vtkPolygon* vtkPolygon::NewInstance ( ) const

◆ PrintSelf()

void vtkPolygon::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
overridevirtual

Methods invoked by print to print information about the object including superclasses.

Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from vtkCell.

◆ GetCellType()

int vtkPolygon::GetCellType ( )
inlineoverridevirtual

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 56 of file vtkPolygon.h.

◆ GetCellDimension()

int vtkPolygon::GetCellDimension ( )
inlineoverridevirtual

Return the topological dimensional of the cell (0,1,2, or 3).

Implements vtkCell.

Definition at line 57 of file vtkPolygon.h.

◆ GetNumberOfEdges()

int vtkPolygon::GetNumberOfEdges ( )
inlineoverridevirtual

Return the number of edges in the cell.

Implements vtkCell.

Definition at line 58 of file vtkPolygon.h.

◆ GetNumberOfFaces()

int vtkPolygon::GetNumberOfFaces ( )
inlineoverridevirtual

Return the number of faces in the cell.

Implements vtkCell.

Definition at line 59 of file vtkPolygon.h.

◆ GetEdge()

vtkCell* vtkPolygon::GetEdge ( int  edgeId)
overridevirtual

Return the edge cell from the edgeId of the cell.

Implements vtkCell.

◆ GetFace()

vtkCell* vtkPolygon::GetFace ( int  faceId)
inlineoverridevirtual

Return the face cell from the faceId of the cell.

Implements vtkCell.

Definition at line 61 of file vtkPolygon.h.

◆ CellBoundary()

int vtkPolygon::CellBoundary ( int  subId,
const double  pcoords[3],
vtkIdList pts 
)
overridevirtual

Given parametric coordinates of a point, return the closest cell boundary, and whether the point is inside or outside of the cell.

The cell boundary is defined by a list of points (pts) that specify a face (3D cell), edge (2D cell), or vertex (1D cell). If the return value of the method is != 0, then the point is inside the cell.

Implements vtkCell.

◆ Contour()

void vtkPolygon::Contour ( double  value,
vtkDataArray cellScalars,
vtkIncrementalPointLocator locator,
vtkCellArray verts,
vtkCellArray lines,
vtkCellArray polys,
vtkPointData inPd,
vtkPointData outPd,
vtkCellData inCd,
vtkIdType  cellId,
vtkCellData outCd 
)
overridevirtual

Generate contouring primitives.

The scalar list cellScalars are scalar values at each cell point. The point locator is essentially a points list that merges points as they are inserted (i.e., prevents duplicates). Contouring primitives can be vertices, lines, or polygons. It is possible to interpolate point data along the edge by providing input and output point data - if outPd is nullptr, then no interpolation is performed. Also, if the output cell data is non-nullptr, the cell data from the contoured cell is passed to the generated contouring primitives. (Note: the CopyAllocate() method must be invoked on both the output cell and point data. The cellId refers to the cell from which the cell data is copied.)

Implements vtkCell.

◆ Clip()

void vtkPolygon::Clip ( double  value,
vtkDataArray cellScalars,
vtkIncrementalPointLocator locator,
vtkCellArray connectivity,
vtkPointData inPd,
vtkPointData outPd,
vtkCellData inCd,
vtkIdType  cellId,
vtkCellData outCd,
int  insideOut 
)
overridevirtual

Cut (or clip) the cell based on the input cellScalars and the specified value.

The output of the clip operation will be one or more cells of the same topological dimension as the original cell. The flag insideOut controls what part of the cell is considered inside - normally cell points whose scalar value is greater than "value" are considered inside. If insideOut is on, this is reversed. Also, if the output cell data is non-nullptr, the cell data from the clipped cell is passed to the generated contouring primitives. (Note: the CopyAllocate() method must be invoked on both the output cell and point data. The cellId refers to the cell from which the cell data is copied.)

Implements vtkCell.

◆ EvaluatePosition()

int vtkPolygon::EvaluatePosition ( const double  x[3],
double  closestPoint[3],
int subId,
double  pcoords[3],
double dist2,
double  weights[] 
)
overridevirtual

Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; evaluate parametric coordinates, sub-cell id (!=0 only if cell is composite), distance squared of point x[3] to cell (in particular, the sub-cell indicated), closest point on cell to x[3] (unless closestPoint is null, in which case, the closest point and dist2 are not found), and interpolation weights in cell.

(The number of weights is equal to the number of points defining the cell). Note: on rare occasions a -1 is returned from the method. This means that numerical error has occurred and all data returned from this method should be ignored. Also, inside/outside is determine parametrically. That is, a point is inside if it satisfies parametric limits. This can cause problems for cells of topological dimension 2 or less, since a point in 3D can project onto the cell within parametric limits but be "far" from the cell. Thus the value dist2 may be checked to determine true in/out.

Implements vtkCell.

◆ EvaluateLocation()

void vtkPolygon::EvaluateLocation ( int subId,
const double  pcoords[3],
double  x[3],
double weights 
)
overridevirtual

Determine global coordinate (x[3]) from subId and parametric coordinates.

Also returns interpolation weights. (The number of weights is equal to the number of points in the cell.)

Implements vtkCell.

◆ IntersectWithLine()

int vtkPolygon::IntersectWithLine ( const double  p1[3],
const double  p2[3],
double  tol,
double t,
double  x[3],
double  pcoords[3],
int subId 
)
overridevirtual

Intersect with a ray.

Return parametric coordinates (both line and cell) and global intersection coordinates, given ray definition p1[3], p2[3] and tolerance tol. The method returns non-zero value if intersection occurs. A parametric distance t between 0 and 1 along the ray representing the intersection point, the point coordinates x[3] in data coordinates and also pcoords[3] in parametric coordinates. subId is the index within the cell if a composed cell like a triangle strip.

Implements vtkCell.

◆ Triangulate() [1/2]

int vtkPolygon::Triangulate ( int  index,
vtkIdList ptIds,
vtkPoints pts 
)
overridevirtual

Generate simplices of proper dimension.

If cell is 3D, tetrahedron are generated; if 2D triangles; if 1D lines; if 0D points. The form of the output is a sequence of points, each n+1 points (where n is topological cell dimension) defining a simplex. The index is a parameter that controls which triangulation to use (if more than one is possible). If numerical degeneracy encountered, 0 is returned, otherwise 1 is returned. This method does not insert new points: all the points that define the simplices are the points that define the cell.

Implements vtkCell.

◆ Derivatives()

void vtkPolygon::Derivatives ( int  subId,
const double  pcoords[3],
const double values,
int  dim,
double derivs 
)
overridevirtual

Compute derivatives given cell subId and parametric coordinates.

The values array is a series of data value(s) at the cell points. There is a one-to-one correspondence between cell point and data value(s). Dim is the number of data values per cell point. Derivs are derivatives in the x-y-z coordinate directions for each data value. Thus, if computing derivatives for a scalar function in a hexahedron, dim=1, 8 values are supplied, and 3 deriv values are returned (i.e., derivatives in x-y-z directions). On the other hand, if computing derivatives of velocity (vx,vy,vz) dim=3, 24 values are supplied ((vx,vy,vz)1, (vx,vy,vz)2, ....()8), and 9 deriv values are returned ((d(vx)/dx),(d(vx)/dy),(d(vx)/dz), (d(vy)/dx),(d(vy)/dy), (d(vy)/dz), (d(vz)/dx),(d(vz)/dy),(d(vz)/dz)).

Implements vtkCell.

◆ IsPrimaryCell()

int vtkPolygon::IsPrimaryCell ( )
inlineoverridevirtual

Return whether this cell type has a fixed topology or whether the topology varies depending on the data (e.g., vtkConvexPointSet).

This compares to composite cells that are typically composed of primary cells (e.g., a triangle strip composite cell is made up of triangle primary cells).

Reimplemented from vtkCell.

Definition at line 83 of file vtkPolygon.h.

◆ ComputeArea() [1/2]

double vtkPolygon::ComputeArea ( )

Compute the area of a polygon.

This is a convenience function which simply calls static double ComputeArea(vtkPoints *p, vtkIdType numPts, vtkIdType *pts, double normal[3]); with the appropriate parameters from the instantiated vtkPolygon.

◆ InterpolateFunctions()

void vtkPolygon::InterpolateFunctions ( const double  x[3],
double sf 
)
override

Compute the interpolation functions/derivatives.

(aka shape functions/derivatives) Two interpolation algorithms are available: 1/r^2 and Mean Value Coordinate. The former is used by default. To use the second algorithm, set UseMVCInterpolation to be true. The function assumes the input point lies on the polygon plane without checking that.

◆ ComputeNormal() [1/4]

static void vtkPolygon::ComputeNormal ( vtkPoints p,
int  numPts,
vtkIdType pts,
double  n[3] 
)
static

Computes the unit normal to the polygon.

If pts=nullptr, point indexing is assumed to be {0, 1, ..., numPts-1}.

◆ ComputeNormal() [2/4]

static void vtkPolygon::ComputeNormal ( vtkPoints p,
double  n[3] 
)
static

◆ ComputeNormal() [3/4]

static void vtkPolygon::ComputeNormal ( vtkIdTypeArray ids,
vtkPoints pts,
double  n[3] 
)
static

◆ ComputeNormal() [4/4]

static void vtkPolygon::ComputeNormal ( int  numPts,
double pts,
double  n[3] 
)
static

Compute the polygon normal from an array of points.

This version assumes that the polygon is convex, and looks for the first valid normal.

◆ IsConvex() [1/4]

bool vtkPolygon::IsConvex ( )

Determine whether or not a polygon is convex.

This is a convenience function that simply calls static bool IsConvex(int numPts, vtkIdType *pts, vtkPoints *p) with the appropriate parameters from the instantiated vtkPolygon.

◆ IsConvex() [2/4]

static bool vtkPolygon::IsConvex ( vtkPoints p,
int  numPts,
vtkIdType pts 
)
static

Determine whether or not a polygon is convex.

If pts=nullptr, point indexing is assumed to be {0, 1, ..., numPts-1}.

◆ IsConvex() [3/4]

static bool vtkPolygon::IsConvex ( vtkIdTypeArray ids,
vtkPoints p 
)
static

◆ IsConvex() [4/4]

static bool vtkPolygon::IsConvex ( vtkPoints p)
static

◆ ComputeCentroid() [1/2]

static bool vtkPolygon::ComputeCentroid ( vtkPoints p,
int  numPts,
vtkIdType pts,
double  centroid[3] 
)
static

Compute the centroid of a set of points.

Returns false if the computation is invalid (this occurs when numPts=0 or when ids is empty).

◆ ComputeCentroid() [2/2]

static bool vtkPolygon::ComputeCentroid ( vtkIdTypeArray ids,
vtkPoints pts,
double  centroid[3] 
)
static

◆ ComputeArea() [2/2]

static double vtkPolygon::ComputeArea ( vtkPoints p,
vtkIdType  numPts,
vtkIdType pts,
double  normal[3] 
)
static

Compute the area of a polygon in 3D.

The area is returned, as well as the normal (a side effect of using this method). If you desire to compute the area of a triangle, use vtkTriangleArea which is faster. If pts==nullptr, point indexing is supposed to be {0, 1, ..., numPts-1}. If you already have a vtkPolygon instantiated, a convenience function, ComputeArea() is provided.

◆ ParameterizePolygon()

int vtkPolygon::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.

The point p0 is the origin of the local system, p10 is s-axis vector, and p20 is the t-axis vector. (These are expressed in the modeling coordinate system and are vectors of dimension [3].) The values l20 and l20 are the lengths of the vectors p10 and p20, and n is the polygon normal.

◆ PointInPolygon()

static int vtkPolygon::PointInPolygon ( double  x[3],
int  numPts,
double pts,
double  bounds[6],
double  n[3] 
)
static

Determine whether point is inside polygon.

Function uses ray-casting to determine if point is inside polygon. Works for arbitrary polygon shape (e.g., non-convex). Returns 0 if point is not in polygon; 1 if it is. Can also return -1 to indicate degenerate polygon.

◆ Triangulate() [2/2]

int vtkPolygon::Triangulate ( vtkIdList outTris)

Triangulate this polygon.

The user must provide the vtkIdList outTris. On output, the outTris list contains the ids of the points defining the triangulation. The ids are ordered into groups of three: each three-group defines one triangle.

◆ NonDegenerateTriangulate()

int vtkPolygon::NonDegenerateTriangulate ( vtkIdList outTris)

Same as Triangulate(vtkIdList *outTris) but with a first pass to split the polygon into non-degenerate polygons.

◆ BoundedTriangulate()

int vtkPolygon::BoundedTriangulate ( vtkIdList outTris,
double  tol 
)

Triangulate polygon and enforce that the ratio of the smallest triangle area to the polygon area is greater than a user-defined tolerance.

The user must provide the vtkIdList outTris. On output, the outTris list contains the ids of the points defining the triangulation. The ids are ordered into groups of three: each three-group defines one triangle.

◆ DistanceToPolygon()

static double vtkPolygon::DistanceToPolygon ( double  x[3],
int  numPts,
double pts,
double  bounds[6],
double  closest[3] 
)
static

Compute the distance of a point to a polygon.

The closest point on the polygon is also returned. The bounds should be provided to accelerate the computation.

◆ IntersectPolygonWithPolygon()

static int vtkPolygon::IntersectPolygonWithPolygon ( int  npts,
double pts,
double  bounds[6],
int  npts2,
double pts2,
double  bounds2[3],
double  tol,
double  x[3] 
)
static

Method intersects two polygons.

You must supply the number of points and point coordinates (npts, *pts) and the bounding box (bounds) of the two polygons. Also supply a tolerance squared for controlling error. The method returns 1 if there is an intersection, and 0 if not. A single point of intersection x[3] is also returned if there is an intersection.

◆ IntersectConvex2DCells()

static int vtkPolygon::IntersectConvex2DCells ( vtkCell cell1,
vtkCell cell2,
double  tol,
double  p0[3],
double  p1[3] 
)
static

Intersect two convex 2D polygons to produce a line segment as output.

The return status of the methods indicated no intersection (returns 0); a single point of intersection (returns 1); or a line segment (i.e., two points of intersection, returns 2). The points of intersection are returned in the arrays p0 and p1. If less than two points of intersection are generated then p1 and/or p0 may be indeterminiate. Finally, if the two convex polygons are parallel, then "0" is returned (i.e., no intersection) even if the triangles lie on one another.

◆ GetUseMVCInterpolation()

virtual bool vtkPolygon::GetUseMVCInterpolation ( )
virtual

Set/Get the flag indicating whether to use Mean Value Coordinate for the interpolation.

If true, InterpolateFunctions() uses the Mean Value Coordinate to compute weights. Otherwise, the conventional 1/r^2 method is used. The UseMVCInterpolation parameter is set to false by default.

◆ SetUseMVCInterpolation()

virtual void vtkPolygon::SetUseMVCInterpolation ( bool  )
virtual

◆ InterpolateFunctionsUsingMVC()

void vtkPolygon::InterpolateFunctionsUsingMVC ( const double  x[3],
double weights 
)
protected

◆ EarCutTriangulation()

int vtkPolygon::EarCutTriangulation ( )
protected

A fast triangulation method.

Uses recursive divide and conquer based on plane splitting to reduce loop into triangles. The cell (e.g., triangle) is presumed properly initialized (i.e., Points and PointIds).

◆ UnbiasedEarCutTriangulation()

int vtkPolygon::UnbiasedEarCutTriangulation ( int  seed)
protected

A fast triangulation method.

Uses recursive divide and conquer based on plane splitting to reduce loop into triangles. The cell (e.g., triangle) is presumed properly initialized (i.e., Points and PointIds). Unlike EarCutTriangulation(), vertices are visited sequentially without preference to angle.

Member Data Documentation

◆ Tolerance

double vtkPolygon::Tolerance
protected

Definition at line 258 of file vtkPolygon.h.

◆ SuccessfulTriangulation

int vtkPolygon::SuccessfulTriangulation
protected

Definition at line 259 of file vtkPolygon.h.

◆ Normal

double vtkPolygon::Normal[3]
protected

Definition at line 260 of file vtkPolygon.h.

◆ Tris

vtkIdList* vtkPolygon::Tris
protected

Definition at line 261 of file vtkPolygon.h.

◆ Triangle

vtkTriangle* vtkPolygon::Triangle
protected

Definition at line 262 of file vtkPolygon.h.

◆ Quad

vtkQuad* vtkPolygon::Quad
protected

Definition at line 263 of file vtkPolygon.h.

◆ TriScalars

vtkDoubleArray* vtkPolygon::TriScalars
protected

Definition at line 264 of file vtkPolygon.h.

◆ Line

vtkLine* vtkPolygon::Line
protected

Definition at line 265 of file vtkPolygon.h.

◆ UseMVCInterpolation

bool vtkPolygon::UseMVCInterpolation
protected

Definition at line 269 of file vtkPolygon.h.


The documentation for this class was generated from the following file: