VTK  9.3.20240318
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
vtkPolyhedron Class Reference

A 3D cell defined by a set of polygonal faces. More...

#include <vtkPolyhedron.h>

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

Public Types

typedef std::map< vtkIdType, vtkIdTypevtkPointIdMap
 
- Public Types inherited from vtkCell3D
typedef vtkCell Superclass
 
- Public Types inherited from vtkCell
typedef vtkObject Superclass
 

Public Member Functions

double * GetParametricCoords () override
 See vtkCell3D API for description of this method. More...
 
int GetCellType () override
 See the vtkCell API for descriptions of these methods. More...
 
int RequiresInitialization () override
 This cell requires that it be initialized prior to access. More...
 
void Initialize () override
 The Initialize method builds up internal structures of vtkPolyhedron. More...
 
void Contour (double value, vtkDataArray *scalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override
 Satisfy the vtkCell API. More...
 
void Clip (double value, vtkDataArray *scalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
 Satisfy the vtkCell API. More...
 
int EvaluatePosition (const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
 Satisfy the vtkCell API. More...
 
void EvaluateLocation (int &subId, const double pcoords[3], double x[3], double *weights) override
 The inverse of EvaluatePosition. 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 the line (p1,p2) with a given tolerance tol to determine a point of intersection x[3] with parametric coordinate t along the line. More...
 
int TriangulateLocalIds (int index, vtkIdList *ptIds) override
 Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh. More...
 
int TriangulateFaces (vtkIdList *newFaces)
 Triangulate each face of the polyhedron. More...
 
int TriangulateFaces (vtkCellArray *newFaces)
 Triangulate each face of the polyhedron. More...
 
void Derivatives (int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
 Computes derivatives at the point specified by the parameter coordinate. More...
 
int CellBoundary (int subId, const double pcoords[3], vtkIdList *pts) override
 Find the boundary face closest to the point defined by the pcoords[3] and subId of the cell (subId can be ignored). More...
 
int GetParametricCenter (double pcoords[3]) override
 Return the center of the cell in parametric coordinates. More...
 
int IsPrimaryCell () override
 A polyhedron is a full-fledged primary cell. More...
 
int RequiresExplicitFaceRepresentation () override
 Satisfy the vtkCell API. More...
 
void SetFaces (vtkIdType *faces) override
 Set the faces of the polyhedron. More...
 
vtkIdTypeGetFaces () override
 Get the faces of the polyhedron. More...
 
int SetCellFaces (vtkCellArray *faces)
 Set the faces of the polyhedron. More...
 
int IsInside (const double x[3], double tolerance)
 A method particular to vtkPolyhedron. More...
 
bool IsConvex ()
 Determine whether or not a polyhedron is convex. More...
 
vtkPolyDataGetPolyData ()
 Construct polydata if no one exist, then return this->PolyData. More...
 
void GetEdgePoints (vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
 See vtkCell3D API for description of these methods. More...
 
vtkIdType GetFacePoints (vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(pts)) override
 See vtkCell3D API for description of these methods. More...
 
void GetEdgeToAdjacentFaces (vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
 See vtkCell3D API for description of these methods. More...
 
vtkIdType GetFaceToAdjacentFaces (vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(faceIds)) override
 See vtkCell3D API for description of these methods. More...
 
vtkIdType GetPointToIncidentEdges (vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(edgeIds)) override
 See vtkCell3D API for description of these methods. More...
 
vtkIdType GetPointToIncidentFaces (vtkIdType pointId, const vtkIdType *&faceIds) override
 See vtkCell3D API for description of these methods. More...
 
vtkIdType GetPointToOneRingPoints (vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(pts)) override
 See vtkCell3D API for description of these methods. More...
 
bool GetCentroid (double vtkNotUsed(centroid)[3]) const override
 See vtkCell3D API for description of these methods. More...
 
int GetNumberOfEdges () override
 A polyhedron is represented internally by a set of polygonal faces. More...
 
vtkCellGetEdge (int) override
 A polyhedron is represented internally by a set of polygonal faces. More...
 
int GetNumberOfFaces () override
 A polyhedron is represented internally by a set of polygonal faces. More...
 
vtkCellGetFace (int faceId) override
 A polyhedron is represented internally by a set of polygonal faces. More...
 
void InterpolateFunctions (const double x[3], double *sf) override
 Compute the interpolation functions/derivatives (aka shape functions/derivatives). More...
 
void InterpolateDerivs (const double x[3], double *derivs) override
 Compute the interpolation functions/derivatives (aka shape functions/derivatives). More...
 
vtkCellArrayGetCellFaces ()
 Get the faces of the polyhedron. More...
 
void GetCellFaces (vtkCellArray *faces)
 Get the faces of the polyhedron. More...
 
- Public Member Functions inherited from vtkCell3D
vtkCell3DNewInstance () const
 
virtual void GetEdgePoints (vtkIdType edgeId, const vtkIdType *&pts)=0
 Get the pair of vertices that define an edge. More...
 
virtual vtkIdType GetFacePoints (vtkIdType faceId, const vtkIdType *&pts)=0
 Get the list of vertices that define a face. More...
 
virtual void GetEdgeToAdjacentFaces (vtkIdType edgeId, const vtkIdType *&faceIds)=0
 Get the ids of the two adjacent faces to edge of id edgeId. More...
 
virtual vtkIdType GetFaceToAdjacentFaces (vtkIdType faceId, const vtkIdType *&faceIds)=0
 Get the ids of the adjacent faces to face of id faceId. More...
 
virtual vtkIdType GetPointToIncidentEdges (vtkIdType pointId, const vtkIdType *&edgeIds)=0
 Get the ids of the incident edges to point of id pointId. More...
 
virtual vtkIdType GetPointToOneRingPoints (vtkIdType pointId, const vtkIdType *&pts)=0
 Get the ids of a one-ring surrounding point of id pointId. More...
 
virtual bool IsInsideOut ()
 Returns true if the normals of the vtkCell3D point inside the cell. More...
 
virtual bool GetCentroid (double centroid[3]) const =0
 Computes the centroid of the cell. More...
 
int GetCellDimension () override
 The topological dimension of the cell. More...
 
int Inflate (double dist) override
 Inflates the cell. More...
 
virtual void SetMergeTolerance (double)
 Set the tolerance for merging clip intersection points that are near the vertices of cells. More...
 
virtual double GetMergeTolerance ()
 Set the tolerance for merging clip intersection points that are near the vertices of cells. More...
 
- Public Member Functions inherited from vtkCell
vtkCellNewInstance () const
 
void Initialize (int npts, const 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 IsExplicitCell ()
 Explicit cells require additional representational information beyond the usual cell type and connectivity list information. More...
 
virtual void SetFaces (vtkIdType *vtkNotUsed(faces))
 
vtkPointsGetPoints ()
 Get the point coordinates for the cell. More...
 
vtkIdType GetNumberOfPoints () const
 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...
 
virtual double ComputeBoundingSphere (double center[3]) const
 Computes the bounding sphere of the cell. More...
 
virtual int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts)
 Generate simplices of proper dimension. More...
 
virtual int TriangulateIds (int index, vtkIdList *ptIds)
 Generate simplices of proper dimension. More...
 
void GetBounds (double bounds[6])
 Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). More...
 
double * GetBounds ()
 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 double GetParametricDistance (const double pcoords[3])
 Return the distance of the parametric coordinate provided to the 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))
 
virtual int IntersectWithCell (vtkCell *other, double tol=0.0)
 Intersects with an other cell. More...
 
virtual int IntersectWithCell (vtkCell *other, const vtkBoundingBox &boudingBox, const vtkBoundingBox &otherBoundingBox, double tol=0.0)
 Intersects with an other cell. More...
 
- 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...
 
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)
 
vtkTypeBool InvokeEvent (unsigned long event)
 
vtkTypeBool InvokeEvent (const char *event)
 
std::string GetObjectDescription () const override
 The object description printed in messages and PrintSelf output. 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)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkCommandGetCommand (unsigned long tag)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObserver (vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
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)
 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, 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...
 
vtkTypeBool InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
vtkTypeBool InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
virtual void SetObjectName (const std::string &objectName)
 Set/get the name of this object for reporting purposes. More...
 
virtual std::string GetObjectName () const
 Set/get the name of this object for reporting purposes. More...
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string. More...
 
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). 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...
 
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...
 
bool GetIsInMemkind () const
 A local state flag that remembers whether this object lives in the normal or extended memory space. More...
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 
virtual bool UsesGarbageCollector () const
 Indicate whether the class uses vtkGarbageCollector or not. More...
 

Protected Member Functions

 vtkPolyhedron ()
 
 ~vtkPolyhedron () override
 
int GenerateEdges ()
 
void GenerateFaces ()
 
void ComputeBounds ()
 
void ComputeParametricCoordinate (const double x[3], double pc[3])
 
void ComputePositionFromParametricCoordinate (const double pc[3], double x[3])
 
void GeneratePointToIncidentFacesAndValenceAtPoint ()
 
void ConstructPolyData ()
 
void ConstructLocator ()
 
- Protected Member Functions inherited from vtkCell3D
 vtkCell3D ()
 
 ~vtkCell3D () override
 
- 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 ()
 These methods allow a command to exclusively grab all events. More...
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void ReportReferences (vtkGarbageCollector *)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

vtkLineLine
 
vtkTriangleTriangle
 
vtkQuadQuad
 
vtkPolygonPolygon
 
vtkTetraTetra
 
vtkCellArrayGlobalFaces
 
vtkIdTypeArrayLegacyGlobalFaces
 
vtkPointIdMapPointIdMap
 
int EdgesGenerated
 
vtkEdgeTableEdgeTable
 
vtkIdTypeArrayEdges
 
vtkIdTypeArrayEdgeFaces
 
vtkCellArrayFaces
 
int FacesGenerated
 
int BoundsComputed
 
int PolyDataConstructed
 
vtkPolyDataPolyData
 
int LocatorConstructed
 
vtkCellLocatorCellLocator
 
vtkIdListCellIds
 
vtkGenericCellCell
 
vtkIdType ** PointToIncidentFaces
 
vtkIdTypeValenceAtPoint
 
- Protected Attributes inherited from vtkCell3D
vtkOrderedTriangulatorTriangulator
 
double MergeTolerance
 
vtkTetraClipTetra
 
vtkDoubleArrayClipScalars
 
- Protected Attributes inherited from vtkCell
double Bounds [6]
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
std::string ObjectName
 
- Protected Attributes inherited from vtkObjectBase
std::atomic< int32_t > ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Friends

class vtkPolyhedronUtilities
 
typedef vtkCell3D Superclass
 Standard new methods. More...
 
static vtkPolyhedronNew ()
 Standard new methods. More...
 
static vtkTypeBool IsTypeOf (const char *type)
 Standard new methods. More...
 
static vtkPolyhedronSafeDownCast (vtkObjectBase *o)
 Standard new methods. More...
 
virtual vtkTypeBool IsA (const char *type)
 Standard new methods. More...
 
vtkPolyhedronNewInstance () const
 Standard new methods. More...
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Standard new methods. More...
 
virtual vtkObjectBaseNewInstanceInternal () const
 Standard new methods. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from vtkCell3D
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkCell3DSafeDownCast (vtkObjectBase *o)
 
- 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 (vtkTypeBool val)
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOn ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOff ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static vtkTypeBool GetGlobalWarningDisplay ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
- 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 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). More...
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 
static void SetMemkindDirectory (const char *directoryname)
 The name of a directory, ideally mounted -o dax, to memory map an extended memory space within. More...
 
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. More...
 
- 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 3D cell defined by a set of polygonal faces.

Instantiation

vtkPolyhedron is a concrete implementation of vtkCell that represents a 3D cell defined by a set of polygonal faces.

To instantiate a vtkPolyhedron, like any vtkCell, one needs to define the following structures:

Note that the ordering of points coordinates or IDs is not important. However, it MUST be consistent between the two lists.

Unlike other kinds of cells (e.g. vtkVoxel), the topology is not directly deduced from points coordinates or point IDs ordering; it must be explicitly defined by providing a list of faces (see the SetFaces() method). Each face is represented as a sequence of global point Ids.

Once point coordinates, point IDs and faces are defined, the Initialize() method should be called in order to setup the internal structures and finalize the construction of the polyhedron.

Here is an example of vtkPolyhedron instantiation:

// 9 +------+.11
// |`. | `.
// |13`+--+---+ 15
// | | | |
// 8 +---+--+.10|
// `. | `.|
// 12`+------+ 14
//
// (Global IDs are arbitrarily chosen between 8 and 15)
// Insert point coordinates
polyhedron->GetPoints()->SetNumberOfPoints(8);
polyhedron->GetPoints()->SetPoint(0, 0., 0., 0.); // 8
polyhedron->GetPoints()->SetPoint(1, 0., 1., 0.); // 9
polyhedron->GetPoints()->SetPoint(2, 1., 0., 0.); // 10
polyhedron->GetPoints()->SetPoint(3, 1., 1., 0.); // 11
polyhedron->GetPoints()->SetPoint(4, 0., 0., 1.); // 12
polyhedron->GetPoints()->SetPoint(5, 0., 1., 1.); // 13
polyhedron->GetPoints()->SetPoint(6, 1., 0., 1.); // 14
polyhedron->GetPoints()->SetPoint(7, 1., 1., 1.); // 15
// Insert point IDs (global IDs)
// Note that the canonical ordering (0, 1, ..., 8) is used
// to correlate point Ids and coordinates
polyhedron->GetPointIds()->Allocate(8);
for (int i = 8; i < 16; ++i)
{
polyhedron->GetPointIds()->InsertNextId(i);
}
// Describe faces, indexed on global IDs
vtkIdType faces[31] = { 6, // Number of faces
4, 9 , 11, 10, 8 , // Number of points in the face + point IDs
4, 11, 15, 14, 10, // Faces are described counter-clockwise
4, 15, 13, 12, 14, // looking from the "outside" of the cell
4, 13, 9 , 8 , 12,
4, 14, 12, 8 , 10,
4, 13, 15, 11, 9 };
polyhedron->SetFaces(faces);
// Initialize the polyhedron
// This will build internal structures and should be done before the proper
// use of the cell.
polyhedron->Initialize();
int vtkIdType
Definition: vtkType.h:315

Specifications

Polyhedrons described by this class must conform to some criteria in order to avoid errors and guarantee good results in terms of visualization and processing.

These specifications are described as follows. Polyhedrons must:

In a more global perspective, polyhedrons must be watertight and manifold. In particular, each edge of the polyhedron must be adjacent to exactly two faces. Several algorithms like contour, clip or slice will assume that each edge of the polyhedron is adjacent to exactly two faces and will definitely lead to bad results (and generate numerous warnings) if this criterion is not fulfilled.

Limitations

The class does not require the polyhedron to be convex. However, the support of concave polyhedrons is currently limited. Concavity can lead to bad results with some filters, including:

Other details

Interpolation functions and weights are defined / computed using the method of Mean Value Coordinates (MVC). See the VTK class vtkMeanValueCoordinatesInterpolator for more information.

See also
vtkCell3D vtkConvexPointSet vtkMeanValueCoordinatesInterpolator vtkPolyhedronUtilities
Online Examples:

Tests:
vtkPolyhedron (Tests)

Definition at line 164 of file vtkPolyhedron.h.

Member Typedef Documentation

◆ vtkPointIdMap

Definition at line 167 of file vtkPolyhedron.h.

◆ Superclass

Standard new methods.

Definition at line 174 of file vtkPolyhedron.h.

Constructor & Destructor Documentation

◆ vtkPolyhedron()

vtkPolyhedron::vtkPolyhedron ( )
protected

◆ ~vtkPolyhedron()

vtkPolyhedron::~vtkPolyhedron ( )
overrideprotected

Member Function Documentation

◆ New()

static vtkPolyhedron* vtkPolyhedron::New ( )
static

Standard new methods.

◆ IsTypeOf()

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

Standard new methods.

◆ IsA()

virtual vtkTypeBool vtkPolyhedron::IsA ( const char *  type)
virtual

Standard new methods.

Reimplemented from vtkCell3D.

◆ SafeDownCast()

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

Standard new methods.

◆ NewInstanceInternal()

virtual vtkObjectBase* vtkPolyhedron::NewInstanceInternal ( ) const
protectedvirtual

Standard new methods.

Reimplemented from vtkCell3D.

◆ NewInstance()

vtkPolyhedron* vtkPolyhedron::NewInstance ( ) const

Standard new methods.

◆ PrintSelf()

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

Standard new methods.

Reimplemented from vtkCell3D.

◆ GetEdgePoints()

void vtkPolyhedron::GetEdgePoints ( vtkIdType   vtkNotUsededgeId,
const vtkIdType *&  vtkNotUsedpts 
)
inlineoverride

See vtkCell3D API for description of these methods.

Warning
These method are unimplemented in vtkPolyhedron

Definition at line 183 of file vtkPolyhedron.h.

◆ GetFacePoints()

vtkIdType vtkPolyhedron::GetFacePoints ( vtkIdType   vtkNotUsedfaceId,
const vtkIdType *&  vtkNotUsedpts 
)
inlineoverride

See vtkCell3D API for description of these methods.

Warning
These method are unimplemented in vtkPolyhedron

Definition at line 187 of file vtkPolyhedron.h.

◆ GetEdgeToAdjacentFaces()

void vtkPolyhedron::GetEdgeToAdjacentFaces ( vtkIdType   vtkNotUsededgeId,
const vtkIdType *&  vtkNotUsedpts 
)
inlineoverride

See vtkCell3D API for description of these methods.

Warning
These method are unimplemented in vtkPolyhedron

Definition at line 192 of file vtkPolyhedron.h.

◆ GetFaceToAdjacentFaces()

vtkIdType vtkPolyhedron::GetFaceToAdjacentFaces ( vtkIdType   vtkNotUsedfaceId,
const vtkIdType *&  vtkNotUsedfaceIds 
)
inlineoverride

See vtkCell3D API for description of these methods.

Warning
These method are unimplemented in vtkPolyhedron

Definition at line 197 of file vtkPolyhedron.h.

◆ GetPointToIncidentEdges()

vtkIdType vtkPolyhedron::GetPointToIncidentEdges ( vtkIdType   vtkNotUsedpointId,
const vtkIdType *&  vtkNotUsededgeIds 
)
inlineoverride

See vtkCell3D API for description of these methods.

Warning
These method are unimplemented in vtkPolyhedron

Definition at line 203 of file vtkPolyhedron.h.

◆ GetPointToIncidentFaces()

vtkIdType vtkPolyhedron::GetPointToIncidentFaces ( vtkIdType  pointId,
const vtkIdType *&  faceIds 
)
overridevirtual

See vtkCell3D API for description of these methods.

Warning
These method are unimplemented in vtkPolyhedron

Implements vtkCell3D.

◆ GetPointToOneRingPoints()

vtkIdType vtkPolyhedron::GetPointToOneRingPoints ( vtkIdType   vtkNotUsedpointId,
const vtkIdType *&  vtkNotUsedpts 
)
inlineoverride

See vtkCell3D API for description of these methods.

Warning
These method are unimplemented in vtkPolyhedron

Definition at line 210 of file vtkPolyhedron.h.

◆ GetCentroid()

bool vtkPolyhedron::GetCentroid ( double   vtkNotUsed(centroid)[3]) const
inlineoverride

See vtkCell3D API for description of these methods.

Warning
These method are unimplemented in vtkPolyhedron

Definition at line 216 of file vtkPolyhedron.h.

◆ GetParametricCoords()

double* vtkPolyhedron::GetParametricCoords ( )
overridevirtual

See vtkCell3D API for description of this method.

Reimplemented from vtkCell.

◆ GetCellType()

int vtkPolyhedron::GetCellType ( )
inlineoverridevirtual

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 231 of file vtkPolyhedron.h.

◆ RequiresInitialization()

int vtkPolyhedron::RequiresInitialization ( )
inlineoverridevirtual

This cell requires that it be initialized prior to access.

Reimplemented from vtkCell.

Definition at line 236 of file vtkPolyhedron.h.

◆ Initialize()

void vtkPolyhedron::Initialize ( )
overridevirtual

The Initialize method builds up internal structures of vtkPolyhedron.

Warning
This method should be called after setting the point coordinates, point IDs and faces of the polyhedron in order to finalize its construction.

Reimplemented from vtkCell.

◆ GetNumberOfEdges()

int vtkPolyhedron::GetNumberOfEdges ( )
overridevirtual

A polyhedron is represented internally by a set of polygonal faces.

These faces can be processed to explicitly determine edges.

Implements vtkCell.

◆ GetEdge()

vtkCell* vtkPolyhedron::GetEdge ( int  )
overridevirtual

A polyhedron is represented internally by a set of polygonal faces.

These faces can be processed to explicitly determine edges.

Implements vtkCell.

◆ GetNumberOfFaces()

int vtkPolyhedron::GetNumberOfFaces ( )
overridevirtual

A polyhedron is represented internally by a set of polygonal faces.

These faces can be processed to explicitly determine edges.

Implements vtkCell.

◆ GetFace()

vtkCell* vtkPolyhedron::GetFace ( int  faceId)
overridevirtual

A polyhedron is represented internally by a set of polygonal faces.

These faces can be processed to explicitly determine edges.

Implements vtkCell.

◆ Contour()

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

Satisfy the vtkCell API.

This method contours the input polyhedron and outputs a polygon. When the result polygon is not planar, it will be triangulated.

Warning
The current implementation assumes water-tight and manifold polyhedron cells.

Reimplemented from vtkCell3D.

◆ Clip()

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

Satisfy the vtkCell API.

This method clips the input polyhedron and outputs a new polyhedron. The face information of the output polyhedron is encoded in the output vtkCellArray using a special format: CellLength [nCellFaces, nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...]. Use the static method vtkUnstructuredGrid::DecomposePolyhedronCellArray to convert it into a standard format.

Warning
The current implementation assumes water-tight and manifold polyhedron cells.

Reimplemented from vtkCell3D.

◆ EvaluatePosition()

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

Satisfy the vtkCell API.

The subId is ignored and zero is always returned. The parametric coordinates pcoords are normalized values in the bounding box of the polyhedron. The weights are determined by evaluating the MVC coordinates. The dist is always zero if the point x[3] is inside the polyhedron; otherwise it's the distance to the surface.

Implements vtkCell.

◆ EvaluateLocation()

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

The inverse of EvaluatePosition.

Note the weights should be the MVC weights.

Implements vtkCell.

◆ IntersectWithLine()

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

Intersect the line (p1,p2) with a given tolerance tol to determine a point of intersection x[3] with parametric coordinate t along the line.

The parametric coordinates are returned as well (subId can be ignored). Returns true if the line intersects a face.

Implements vtkCell.

◆ TriangulateLocalIds()

int vtkPolyhedron::TriangulateLocalIds ( int  index,
vtkIdList ptIds 
)
overridevirtual

Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh.

Once triangulation has been performed, the results are saved in ptIds and pts. The ptIds is a vtkIdList with 4xn number of ids (n is the number of result tetrahedrons). The first 4 represent the point ids of the first tetrahedron, the second 4 represents the point ids of the second tetrahedron and so on. The point ids represent global dataset ids. The points of result tetrahedons are stored in pts. Note that there are 4xm output points (m is the number of points in the original polyhedron). A point may be stored multiple times when it is shared by more than one tetrahedrons. The points stored in pts are ordered the same as they are listed in ptIds.

Warning
This method works well for a convex polyhedron but may return wrong result in a concave case.

Implements vtkCell.

◆ TriangulateFaces() [1/2]

int vtkPolyhedron::TriangulateFaces ( vtkIdList newFaces)

Triangulate each face of the polyhedron.

This method internally use the vtkCell::Triangulate method on each face (so the triangulation method vary depending on the 2D cell type corresponding to the face).

Warning
Can lead to bad results with non-planar faces.

◆ TriangulateFaces() [2/2]

int vtkPolyhedron::TriangulateFaces ( vtkCellArray newFaces)

Triangulate each face of the polyhedron.

This method internally use the vtkCell::Triangulate method on each face (so the triangulation method vary depending on the 2D cell type corresponding to the face).

Warning
Can lead to bad results with non-planar faces.

◆ Derivatives()

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

Computes derivatives at the point specified by the parameter coordinate.

Current implementation uses all vertices and subId is not used. To accelerate the speed, the future implementation can triangulate and extract the local tetrahedron from subId and pcoords, then evaluate derivatives on the local tetrahedron.

Implements vtkCell.

◆ CellBoundary()

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

Find the boundary face closest to the point defined by the pcoords[3] and subId of the cell (subId can be ignored).

Implements vtkCell.

◆ GetParametricCenter()

int vtkPolyhedron::GetParametricCenter ( double  pcoords[3])
inlineoverridevirtual

Return the center of the cell in parametric coordinates.

In this cell, the center of the bounding box is returned.

Reimplemented from vtkCell.

Definition at line 530 of file vtkPolyhedron.h.

◆ IsPrimaryCell()

int vtkPolyhedron::IsPrimaryCell ( )
inlineoverridevirtual

A polyhedron is a full-fledged primary cell.

Reimplemented from vtkCell.

Definition at line 361 of file vtkPolyhedron.h.

◆ InterpolateFunctions()

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

Compute the interpolation functions/derivatives (aka shape functions/derivatives).

Here we use the MVC calculation process to compute the interpolation functions.

◆ InterpolateDerivs()

void vtkPolyhedron::InterpolateDerivs ( const double  x[3],
double *  derivs 
)
override

Compute the interpolation functions/derivatives (aka shape functions/derivatives).

Here we use the MVC calculation process to compute the interpolation functions.

◆ RequiresExplicitFaceRepresentation()

int vtkPolyhedron::RequiresExplicitFaceRepresentation ( )
inlineoverridevirtual

Satisfy the vtkCell API.

Always return true, because vtkPolyhedron needs explicit faces definition in order to describe the topology of the cell.

Reimplemented from vtkCell.

Definition at line 378 of file vtkPolyhedron.h.

◆ SetFaces()

void vtkPolyhedron::SetFaces ( vtkIdType faces)
override

Set the faces of the polyhedron.

Face are expressed as sequences of global point IDs . The SetFaces method will require a copy from internal unstructured grid layout.

Parameters
faces1-dimensional array with the following structure :
[ NbOfFaces,
NbOfPtsFace1, face1Pt1, face1Pt2, …, face1PtNbOfPtsFace1,
NbOfPtsFace2, face2Pt1, face2Pt2, …, face2PtNbOfPtsFace2,
…,
NbOfPtsFaceN, faceNPt1, faceNPt2, …, faceNPtNbOfPtsFaceN ]
This ordering corresponds to the legacy vtkCellArray form, with in addition a leading count indicating the total number of faces in the list.

◆ GetFaces()

vtkIdType* vtkPolyhedron::GetFaces ( )
overridevirtual

Get the faces of the polyhedron.

Face are expressed as sequences of global point IDs .

Returns
A 1-dimentional array with the following structure :
[ NbOfPtsFace1, face1Pt1, face1Pt2, …, face1PtNbOfPtsFace1,
NbOfPtsFace2, face2Pt1, face2Pt2, …, face2PtNbOfPtsFace2,
…,
NbOfPtsFaceN, faceNPt1, faceNPt2, …, faceNPtNbOfPtsFaceN ]
This ordering corresponds to the legacy vtkCellArray form. Note that unlike the SetFaces method, the total faces number leading count is missing. In order to get the number of faces, please use the vtkPolyhedron::GetNumberOfFaces() method.

Reimplemented from vtkCell.

◆ SetCellFaces()

int vtkPolyhedron::SetCellFaces ( vtkCellArray faces)

Set the faces of the polyhedron.

Symmetric method to GetCellFaces

Parameters
facesvtkCellArray that stores a contiguous list of polygonal faces with their corresponding global point IDs defining a polyhedron.

◆ GetCellFaces() [1/2]

vtkCellArray* vtkPolyhedron::GetCellFaces ( )

Get the faces of the polyhedron.

Face are expressed as sequences of global point IDs .

Parameters
facesvtkCellArray that stores the list of polygonal faces with their corresponding global point IDs

◆ GetCellFaces() [2/2]

void vtkPolyhedron::GetCellFaces ( vtkCellArray faces)

Get the faces of the polyhedron.

Face are expressed as sequences of global point IDs .

Parameters
facesvtkCellArray that stores the list of polygonal faces with their corresponding global point IDs

◆ IsInside()

int vtkPolyhedron::IsInside ( const double  x[3],
double  tolerance 
)

A method particular to vtkPolyhedron.

It determines whether a point x[3] is inside the polyhedron or not (returns 1 is the point is inside, 0 otherwise). The tolerance is expressed in normalized space; i.e., a fraction of the size of the bounding box.

◆ IsConvex()

bool vtkPolyhedron::IsConvex ( )

Determine whether or not a polyhedron is convex.

This method is adapted from Devillers et al., "Checking the Convexity of Polytopes and the Planarity of Subdivisions", Computational Geometry, Volume 11, Issues 3 - 4, December 1998, Pages 187 - 208.

◆ GetPolyData()

vtkPolyData* vtkPolyhedron::GetPolyData ( )

Construct polydata if no one exist, then return this->PolyData.

◆ GenerateEdges()

int vtkPolyhedron::GenerateEdges ( )
protected

◆ GenerateFaces()

void vtkPolyhedron::GenerateFaces ( )
protected

◆ ComputeBounds()

void vtkPolyhedron::ComputeBounds ( )
protected

◆ ComputeParametricCoordinate()

void vtkPolyhedron::ComputeParametricCoordinate ( const double  x[3],
double  pc[3] 
)
protected

◆ ComputePositionFromParametricCoordinate()

void vtkPolyhedron::ComputePositionFromParametricCoordinate ( const double  pc[3],
double  x[3] 
)
protected

◆ GeneratePointToIncidentFacesAndValenceAtPoint()

void vtkPolyhedron::GeneratePointToIncidentFacesAndValenceAtPoint ( )
protected

◆ ConstructPolyData()

void vtkPolyhedron::ConstructPolyData ( )
protected

◆ ConstructLocator()

void vtkPolyhedron::ConstructLocator ( )
protected

Friends And Related Function Documentation

◆ vtkPolyhedronUtilities

friend class vtkPolyhedronUtilities
friend

Definition at line 526 of file vtkPolyhedron.h.

Member Data Documentation

◆ Line

vtkLine* vtkPolyhedron::Line
protected

Definition at line 464 of file vtkPolyhedron.h.

◆ Triangle

vtkTriangle* vtkPolyhedron::Triangle
protected

Definition at line 465 of file vtkPolyhedron.h.

◆ Quad

vtkQuad* vtkPolyhedron::Quad
protected

Definition at line 466 of file vtkPolyhedron.h.

◆ Polygon

vtkPolygon* vtkPolyhedron::Polygon
protected

Definition at line 467 of file vtkPolyhedron.h.

◆ Tetra

vtkTetra* vtkPolyhedron::Tetra
protected

Definition at line 468 of file vtkPolyhedron.h.

◆ GlobalFaces

vtkCellArray* vtkPolyhedron::GlobalFaces
protected

Definition at line 472 of file vtkPolyhedron.h.

◆ LegacyGlobalFaces

vtkIdTypeArray* vtkPolyhedron::LegacyGlobalFaces
protected

Definition at line 475 of file vtkPolyhedron.h.

◆ PointIdMap

vtkPointIdMap* vtkPolyhedron::PointIdMap
protected

Definition at line 482 of file vtkPolyhedron.h.

◆ EdgesGenerated

int vtkPolyhedron::EdgesGenerated
protected

Definition at line 485 of file vtkPolyhedron.h.

◆ EdgeTable

vtkEdgeTable* vtkPolyhedron::EdgeTable
protected

Definition at line 486 of file vtkPolyhedron.h.

◆ Edges

vtkIdTypeArray* vtkPolyhedron::Edges
protected

Definition at line 487 of file vtkPolyhedron.h.

◆ EdgeFaces

vtkIdTypeArray* vtkPolyhedron::EdgeFaces
protected

Definition at line 488 of file vtkPolyhedron.h.

◆ Faces

vtkCellArray* vtkPolyhedron::Faces
protected

Definition at line 497 of file vtkPolyhedron.h.

◆ FacesGenerated

int vtkPolyhedron::FacesGenerated
protected

Definition at line 498 of file vtkPolyhedron.h.

◆ BoundsComputed

int vtkPolyhedron::BoundsComputed
protected

Definition at line 501 of file vtkPolyhedron.h.

◆ PolyDataConstructed

int vtkPolyhedron::PolyDataConstructed
protected

Definition at line 509 of file vtkPolyhedron.h.

◆ PolyData

vtkPolyData* vtkPolyhedron::PolyData
protected

Definition at line 510 of file vtkPolyhedron.h.

◆ LocatorConstructed

int vtkPolyhedron::LocatorConstructed
protected

Definition at line 512 of file vtkPolyhedron.h.

◆ CellLocator

vtkCellLocator* vtkPolyhedron::CellLocator
protected

Definition at line 513 of file vtkPolyhedron.h.

◆ CellIds

vtkIdList* vtkPolyhedron::CellIds
protected

Definition at line 515 of file vtkPolyhedron.h.

◆ Cell

vtkGenericCell* vtkPolyhedron::Cell
protected

Definition at line 516 of file vtkPolyhedron.h.

◆ PointToIncidentFaces

vtkIdType** vtkPolyhedron::PointToIncidentFaces
protected

Definition at line 519 of file vtkPolyhedron.h.

◆ ValenceAtPoint

vtkIdType* vtkPolyhedron::ValenceAtPoint
protected

Definition at line 520 of file vtkPolyhedron.h.


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