VTK  9.1.0
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

enum  EarCutMeasureTypes { PERIMETER2_TO_AREA_RATIO = 0 , DOT_PRODUCT = 1 , BEST_QUALITY = 2 }
 
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.
 
vtkPolygonNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses.
 
double ComputeArea ()
 Compute the area of a polygon.
 
void InterpolateFunctions (const double x[3], double *sf) override
 Compute the interpolation functions/derivatives.
 
bool IsConvex ()
 Determine whether or not a polygon is convex.
 
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 Triangulate (vtkIdList *outTris)
 Triangulate this polygon.
 
int NonDegenerateTriangulate (vtkIdList *outTris)
 Same as Triangulate(vtkIdList *outTris) but with a first pass to split the polygon into non-degenerate polygons.
 
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.
 
virtual void SetTolerance (double)
 Specify an internal tolerance for operations requiring polygon triangulation.
 
virtual double GetTolerance ()
 
int GetCellType () override
 See the vtkCell API for descriptions of these methods.
 
int GetCellDimension () override
 See the vtkCell API for descriptions of these methods.
 
int GetNumberOfEdges () override
 See the vtkCell API for descriptions of these methods.
 
int GetNumberOfFaces () override
 See the vtkCell API for descriptions of these methods.
 
vtkCellGetEdge (int edgeId) override
 See the vtkCell API for descriptions of these methods.
 
vtkCellGetFace (int) override
 See the vtkCell API for descriptions of these methods.
 
int CellBoundary (int subId, const double pcoords[3], vtkIdList *pts) override
 See the vtkCell API for descriptions of these methods.
 
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.
 
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.
 
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 EvaluateLocation (int &subId, const double pcoords[3], double x[3], double *weights) override
 See the vtkCell API for descriptions of these methods.
 
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 Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts) override
 See the vtkCell API for descriptions of these methods.
 
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.
 
int IsPrimaryCell () override
 See the vtkCell API for descriptions of these methods.
 
virtual bool GetUseMVCInterpolation ()
 Set/Get the flag indicating whether to use Mean Value Coordinate for the interpolation.
 
virtual void SetUseMVCInterpolation (bool)
 Set/Get the flag indicating whether to use Mean Value Coordinate for the interpolation.
 
int EarCutTriangulation (int measure=PERIMETER2_TO_AREA_RATIO)
 A fast triangulation method.
 
int EarCutTriangulation (vtkIdList *outTris, int measure=PERIMETER2_TO_AREA_RATIO)
 A fast triangulation method.
 
int UnbiasedEarCutTriangulation (int seed, int measure=PERIMETER2_TO_AREA_RATIO)
 A fast triangulation method.
 
int UnbiasedEarCutTriangulation (int seed, vtkIdList *outTris, int measure=PERIMETER2_TO_AREA_RATIO)
 A fast triangulation method.
 
- Public Member Functions inherited from vtkCell
virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class.
 
vtkCellNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses.
 
void Initialize (int npts, const vtkIdType *pts, vtkPoints *p)
 Initialize cell from outside with point ids and point coordinates specified.
 
void Initialize (int npts, vtkPoints *p)
 Initialize the cell with point coordinates specified.
 
virtual void ShallowCopy (vtkCell *c)
 Copy this cell by reference counting the internal data structures.
 
virtual void DeepCopy (vtkCell *c)
 Copy this cell by completely copying internal data structures.
 
virtual int GetCellType ()=0
 Return the type of cell.
 
virtual int GetCellDimension ()=0
 Return the topological dimensional of the cell (0,1,2, or 3).
 
virtual int IsLinear ()
 Non-linear cells require special treatment beyond the usual cell type and connectivity list information.
 
virtual int RequiresInitialization ()
 Some cells require initialization prior to access.
 
virtual void Initialize ()
 
virtual int IsExplicitCell ()
 Explicit cells require additional representational information beyond the usual cell type and connectivity list information.
 
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).
 
virtual void SetFaces (vtkIdType *vtkNotUsed(faces))
 
virtual vtkIdTypeGetFaces ()
 
vtkPointsGetPoints ()
 Get the point coordinates for the cell.
 
vtkIdType GetNumberOfPoints () const
 Return the number of points in the cell.
 
virtual int GetNumberOfEdges ()=0
 Return the number of edges in the cell.
 
virtual int GetNumberOfFaces ()=0
 Return the number of faces in the cell.
 
vtkIdListGetPointIds ()
 Return the list of point ids defining the cell.
 
vtkIdType GetPointId (int ptId)
 For cell point i, return the actual point id.
 
virtual vtkCellGetEdge (int edgeId)=0
 Return the edge cell from the edgeId of the cell.
 
virtual vtkCellGetFace (int faceId)=0
 Return the face cell from the faceId of the cell.
 
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 inside or outside of the cell.
 
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; 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.
 
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.
 
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 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.
 
virtual int Inflate (double dist)
 Inflates the cell.
 
virtual double ComputeBoundingSphere (double center[3]) const
 Computes the bounding sphere of the cell.
 
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.
 
virtual int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts)=0
 Generate simplices of proper dimension.
 
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.
 
void GetBounds (double bounds[6])
 Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax).
 
double * GetBounds ()
 Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax).
 
double GetLength2 ()
 Compute Length squared of cell (i.e., bounding box diagonal squared).
 
virtual int GetParametricCenter (double pcoords[3])
 Return center of the cell in parametric coordinates.
 
virtual double GetParametricDistance (const double pcoords[3])
 Return the distance of the parametric coordinate provided to the cell.
 
virtual int IsPrimaryCell ()
 Return whether this cell type has a fixed topology or whether the topology varies depending on the data (e.g., vtkConvexPointSet).
 
virtual double * GetParametricCoords ()
 Return a contiguous array of parametric coordinates of the points defining this cell.
 
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.
 
virtual void InterpolateDerivs (const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
 
virtual int IntersectWithCell (vtkCell *other, double tol=0.0)
 Intersects with an other cell.
 
virtual int IntersectWithCell (vtkCell *other, const vtkBoundingBox &boudingBox, const vtkBoundingBox &otherBoundingBox, double tol=0.0)
 Intersects with an other cell.
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 Turn debugging output on.
 
virtual void DebugOff ()
 Turn debugging output off.
 
bool GetDebug ()
 Get the value of the debug flag.
 
void SetDebug (bool debugFlag)
 Set the value of the debug flag.
 
virtual void Modified ()
 Update the modification time for this object.
 
virtual vtkMTimeType GetMTime ()
 Return this object's modified time.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses.
 
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)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkCommandGetCommand (unsigned long tag)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObserver (vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObservers (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObservers (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
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.
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events.
 
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.
 
int InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
int InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string.
 
virtual vtkTypeBool IsA (const char *name)
 Return 1 if this class is the same type of (or a subclass of) the named class.
 
virtual vtkIdType GetNumberOfGenerationsFromBase (const char *name)
 Given the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class).
 
virtual void Delete ()
 Delete a VTK object.
 
virtual void FastDelete ()
 Delete a reference to this object.
 
void InitializeObjectBase ()
 
void Print (ostream &os)
 Print an object to an ostream.
 
virtual void Register (vtkObjectBase *o)
 Increase the reference count (mark as used by another object).
 
virtual void UnRegister (vtkObjectBase *o)
 Decrease the reference count (release by another object).
 
int GetReferenceCount ()
 Return the current reference count of this object.
 
void SetReferenceCount (int)
 Sets the reference count.
 
bool GetIsInMemkind () const
 A local state flag that remembers whether this object lives in the normal or extended memory space.
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses.
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses.
 

Static Public Member Functions

static vtkPolygonNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkPolygonSafeDownCast (vtkObjectBase *o)
 
static void ComputeNormal (int numPts, double *pts, double n[3])
 Compute the polygon normal from an array of points.
 
static double ComputeArea (vtkPoints *p, vtkIdType numPts, const vtkIdType *pts, double normal[3])
 Compute the area of a polygon in 3D.
 
static int PointInPolygon (double x[3], int numPts, double *pts, double bounds[6], double n[3])
 Determine whether point is inside polygon.
 
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.
 
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.
 
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.
 
static void ComputeNormal (vtkPoints *p, int numPts, const vtkIdType *pts, double n[3])
 Computes the unit normal to the polygon.
 
static void ComputeNormal (vtkPoints *p, double n[3])
 Computes the unit normal to the polygon.
 
static void ComputeNormal (vtkIdTypeArray *ids, vtkPoints *pts, double n[3])
 Computes the unit normal to the polygon.
 
static bool IsConvex (vtkPoints *p, int numPts, vtkIdType *pts)
 Determine whether or not a polygon is convex.
 
static bool IsConvex (vtkIdTypeArray *ids, vtkPoints *p)
 Determine whether or not a polygon is convex.
 
static bool IsConvex (vtkPoints *p)
 Determine whether or not a polygon is convex.
 
static bool ComputeCentroid (vtkPoints *p, int numPts, const vtkIdType *pts, double centroid[3])
 Compute the centroid of a set of points.
 
static bool ComputeCentroid (vtkIdTypeArray *ids, vtkPoints *pts, double centroid[3])
 Compute the centroid of a set of points.
 
- 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.
 
static void BreakOnError ()
 This method is called when vtkErrorMacro executes.
 
static void SetGlobalWarningDisplay (int val)
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static void GlobalWarningDisplayOn ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static void GlobalWarningDisplayOff ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static int GetGlobalWarningDisplay ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
- 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.
 
static vtkIdType GetNumberOfGenerationsFromBaseType (const char *name)
 Given a the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class).
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
 
static void SetMemkindDirectory (const char *directoryname)
 The name of a directory, ideally mounted -o dax, to memory map an extended memory space within.
 
static bool GetUsingMemkind ()
 A global state flag that controls whether vtkObjects are constructed in the usual way (the default) or within the extended memory space.
 

Protected Member Functions

virtual vtkObjectBaseNewInstanceInternal () const
 
 vtkPolygon ()
 
 ~vtkPolygon () override
 
void InterpolateFunctionsUsingMVC (const double x[3], double *weights)
 
void ComputeTolerance ()
 
- Protected Member Functions inherited from vtkCell
virtual vtkObjectBaseNewInstanceInternal () const
 
 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.
 
void InternalReleaseFocus ()
 These methods allow a command to exclusively grab all events.
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void RegisterInternal (vtkObjectBase *, vtkTypeBool check)
 
virtual void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check)
 
virtual void ReportReferences (vtkGarbageCollector *)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

double Tolerance
 
double Tol
 
int SuccessfulTriangulation
 
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
std::atomic< int32_t > ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Additional Inherited Members

- Public Attributes inherited from vtkCell
vtkPointsPoints
 
vtkIdListPointIds
 
- Static Protected Member Functions inherited from vtkObjectBase
static vtkMallocingFunction GetCurrentMallocFunction ()
 
static vtkReallocingFunction GetCurrentReallocFunction ()
 
static vtkFreeingFunction GetCurrentFreeFunction ()
 
static vtkFreeingFunction GetAlternateFreeFunction ()
 

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.

Online Examples:

Definition at line 148 of file vtkPolygon.h.

Member Typedef Documentation

◆ Superclass

Definition at line 152 of file vtkPolygon.h.

Member Enumeration Documentation

◆ EarCutMeasureTypes

Enumerator
PERIMETER2_TO_AREA_RATIO 
DOT_PRODUCT 
BEST_QUALITY 

Definition at line 387 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 159 of file vtkPolygon.h.

◆ GetCellDimension()

int vtkPolygon::GetCellDimension ( )
inlineoverridevirtual

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 160 of file vtkPolygon.h.

◆ GetNumberOfEdges()

int vtkPolygon::GetNumberOfEdges ( )
inlineoverridevirtual

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 161 of file vtkPolygon.h.

◆ GetNumberOfFaces()

int vtkPolygon::GetNumberOfFaces ( )
inlineoverridevirtual

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 162 of file vtkPolygon.h.

◆ GetEdge()

vtkCell * vtkPolygon::GetEdge ( int  edgeId)
overridevirtual

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

◆ GetFace()

vtkCell * vtkPolygon::GetFace ( int  )
inlineoverridevirtual

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 164 of file vtkPolygon.h.

◆ CellBoundary()

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

See the vtkCell API for descriptions of these methods.

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

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

◆ Clip()

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

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

◆ EvaluatePosition()

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

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

◆ EvaluateLocation()

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

See the vtkCell API for descriptions of these methods.

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

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

◆ Triangulate() [1/2]

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

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

◆ Derivatives()

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

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

◆ IsPrimaryCell()

int vtkPolygon::IsPrimaryCell ( )
inlineoverridevirtual

See the vtkCell API for descriptions of these methods.

Reimplemented from vtkCell.

Definition at line 180 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,
const 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

Computes the unit normal to the polygon.

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

◆ ComputeNormal() [3/4]

static void vtkPolygon::ComputeNormal ( vtkIdTypeArray ids,
vtkPoints 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() [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

Determine whether or not a polygon is convex.

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

◆ IsConvex() [4/4]

static bool vtkPolygon::IsConvex ( vtkPoints p)
static

Determine whether or not a polygon is convex.

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

◆ ComputeCentroid() [1/2]

static bool vtkPolygon::ComputeCentroid ( vtkPoints p,
int  numPts,
const 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

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).

◆ ComputeArea() [2/2]

static double vtkPolygon::ComputeArea ( vtkPoints p,
vtkIdType  numPts,
const 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 (i.e., not the associated polygon->PointIds, rather the index into the polygon->Points array). The ids are ordered into groups of three: each three-group defines one triangle. The method returns non-zero if the triangulation is successful.

◆ 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

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.

◆ SetTolerance()

virtual void vtkPolygon::SetTolerance ( double  )
virtual

Specify an internal tolerance for operations requiring polygon triangulation.

(For example, clipping and contouring operations proceed by first triangulating the polygon, and then clipping/contouring the resulting triangles.) This is a normalized tolerance value multiplied by the diagonal length of the polygon bounding box. Is it used to determine whether potential triangulation edges intersect one another.

◆ GetTolerance()

virtual double vtkPolygon::GetTolerance ( )
virtual

◆ InterpolateFunctionsUsingMVC()

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

◆ ComputeTolerance()

void vtkPolygon::ComputeTolerance ( )
protected

◆ EarCutTriangulation() [1/2]

int vtkPolygon::EarCutTriangulation ( int  measure = PERIMETER2_TO_AREA_RATIO)

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). Ears can be removed using different measures (the measures indicate convexity plus characterize the local geometry around each vertex).

◆ EarCutTriangulation() [2/2]

int vtkPolygon::EarCutTriangulation ( vtkIdList outTris,
int  measure = PERIMETER2_TO_AREA_RATIO 
)

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). Ears can be removed using different measures (the measures indicate convexity plus characterize the local geometry around each vertex).

◆ UnbiasedEarCutTriangulation() [1/2]

int vtkPolygon::UnbiasedEarCutTriangulation ( int  seed,
int  measure = PERIMETER2_TO_AREA_RATIO 
)

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.

◆ UnbiasedEarCutTriangulation() [2/2]

int vtkPolygon::UnbiasedEarCutTriangulation ( int  seed,
vtkIdList outTris,
int  measure = PERIMETER2_TO_AREA_RATIO 
)

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 363 of file vtkPolygon.h.

◆ Tol

double vtkPolygon::Tol
protected

Definition at line 364 of file vtkPolygon.h.

◆ SuccessfulTriangulation

int vtkPolygon::SuccessfulTriangulation
protected

Definition at line 367 of file vtkPolygon.h.

◆ Tris

vtkIdList* vtkPolygon::Tris
protected

Definition at line 368 of file vtkPolygon.h.

◆ Triangle

vtkTriangle* vtkPolygon::Triangle
protected

Definition at line 371 of file vtkPolygon.h.

◆ Quad

vtkQuad* vtkPolygon::Quad
protected

Definition at line 372 of file vtkPolygon.h.

◆ TriScalars

vtkDoubleArray* vtkPolygon::TriScalars
protected

Definition at line 373 of file vtkPolygon.h.

◆ Line

vtkLine* vtkPolygon::Line
protected

Definition at line 374 of file vtkPolygon.h.

◆ UseMVCInterpolation

bool vtkPolygon::UseMVCInterpolation
protected

Definition at line 378 of file vtkPolygon.h.


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