VTK
9.5.20250625
|
A 3D cell defined by a set of polygonal faces. More...
#include <vtkPolyhedron.h>
Public Types | |
using | vtkPointIdMap = std::map< vtkIdType, vtkIdType > |
using | Status = vtkCellStatus |
Adopt vtkCellStatus to describe degenerate polyhedral cells. | |
![]() | |
typedef vtkCell | Superclass |
![]() | |
typedef vtkObject | Superclass |
Public Member Functions | |
double * | GetParametricCoords () override |
See vtkCell3D API for description of this method. | |
int | GetCellType () override |
See the vtkCell API for descriptions of these methods. | |
int | RequiresInitialization () override |
This cell requires that it be initialized prior to access. | |
void | Initialize () override |
The Initialize method builds up internal structures of vtkPolyhedron. | |
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. | |
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. | |
int | EvaluatePosition (const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override |
Satisfy the vtkCell API. | |
void | EvaluateLocation (int &subId, const double pcoords[3], double x[3], double *weights) override |
The inverse of EvaluatePosition. | |
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. | |
int | TriangulateLocalIds (int index, vtkIdList *ptIds) override |
Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh. | |
int | TriangulateFaces (vtkIdList *newFaces) |
Triangulate each face of the polyhedron. | |
int | TriangulateFaces (vtkCellArray *newFaces) |
Triangulate each face of the polyhedron. | |
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. | |
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). | |
int | GetParametricCenter (double pcoords[3]) override |
Return the center of the cell in parametric coordinates. | |
int | IsPrimaryCell () VTK_FUTURE_CONST override |
A polyhedron is a full-fledged primary cell. | |
int | RequiresExplicitFaceRepresentation () VTK_FUTURE_CONST override |
Satisfy the vtkCell API. | |
void | SetFaces (vtkIdType *faces) |
Set the faces of the polyhedron. | |
vtkIdType * | GetFaces () |
Get the faces of the polyhedron. | |
int | SetCellFaces (vtkCellArray *faces) |
Set the faces of the polyhedron. | |
int | IsInside (const double x[3], double tolerance) |
A method particular to vtkPolyhedron. | |
bool | IsConvex () |
Determine whether or not a polyhedron is convex. | |
Status | IsConvex (double planarThreshold) |
vtkPolyData * | GetPolyData () |
Construct polydata if no one exist, then return this->PolyData. | |
void | ShallowCopy (vtkCell *c) override |
Shallow copy of a polyhedron. | |
void | DeepCopy (vtkCell *c) override |
Deep copy of a polyhedron. | |
void | GetEdgePoints (vtkIdType edgeId, const vtkIdType *&pts) override |
See vtkCell3D API for description of these methods. | |
vtkIdType | GetFacePoints (vtkIdType faceId, const vtkIdType *&pts) override |
See vtkCell3D API for description of these methods. | |
void | GetEdgeToAdjacentFaces (vtkIdType edgeId, const vtkIdType *&pts) override |
See vtkCell3D API for description of these methods. | |
vtkIdType | GetFaceToAdjacentFaces (vtkIdType faceId, const vtkIdType *&faceIds) override |
See vtkCell3D API for description of these methods. | |
vtkIdType | GetPointToIncidentEdges (vtkIdType pointId, const vtkIdType *&edgeIds) override |
See vtkCell3D API for description of these methods. | |
vtkIdType | GetPointToIncidentFaces (vtkIdType pointId, const vtkIdType *&faceIds) override |
See vtkCell3D API for description of these methods. | |
vtkIdType | GetPointToOneRingPoints (vtkIdType pointId, const vtkIdType *&pts) override |
See vtkCell3D API for description of these methods. | |
bool | GetCentroid (double centroid[3]) const override |
See vtkCell3D API for description of these methods. | |
int | GetNumberOfEdges () override |
A polyhedron is represented internally by a set of polygonal faces. | |
vtkCell * | GetEdge (int) override |
A polyhedron is represented internally by a set of polygonal faces. | |
int | GetNumberOfFaces () override |
A polyhedron is represented internally by a set of polygonal faces. | |
vtkCell * | GetFace (int faceId) override |
A polyhedron is represented internally by a set of polygonal faces. | |
void | InterpolateFunctions (const double x[3], double *sf) override |
Compute the interpolation functions/derivatives (aka shape functions/derivatives). | |
void | InterpolateDerivs (const double x[3], double *derivs) override |
Compute the interpolation functions/derivatives (aka shape functions/derivatives). | |
vtkCellArray * | GetCellFaces () |
Get the faces of the polyhedron. | |
void | GetCellFaces (vtkCellArray *faces) |
Get the faces of the polyhedron. | |
![]() | |
virtual vtkTypeBool | IsA (const char *type) |
Return 1 if this class is the same type of (or a subclass of) the named class. | |
vtkCell3D * | NewInstance () const |
void | PrintSelf (ostream &os, vtkIndent indent) override |
Methods invoked by print to print information about the object including superclasses. | |
virtual void | GetEdgePoints (vtkIdType edgeId, const vtkIdType *&pts)=0 |
Get the pair of vertices that define an edge. | |
virtual vtkIdType | GetFacePoints (vtkIdType faceId, const vtkIdType *&pts)=0 |
Get the list of vertices that define a face. | |
virtual void | GetEdgeToAdjacentFaces (vtkIdType edgeId, const vtkIdType *&faceIds)=0 |
Get the ids of the two adjacent faces to edge of id edgeId. | |
virtual vtkIdType | GetFaceToAdjacentFaces (vtkIdType faceId, const vtkIdType *&faceIds)=0 |
Get the ids of the adjacent faces to face of id faceId. | |
virtual vtkIdType | GetPointToIncidentEdges (vtkIdType pointId, const vtkIdType *&edgeIds)=0 |
Get the ids of the incident edges to point of id pointId. | |
virtual vtkIdType | GetPointToIncidentFaces (vtkIdType pointId, const vtkIdType *&faceIds)=0 |
Get the ids of the incident faces point of id pointId. | |
virtual vtkIdType | GetPointToOneRingPoints (vtkIdType pointId, const vtkIdType *&pts)=0 |
Get the ids of a one-ring surrounding point of id pointId. | |
virtual bool | IsInsideOut () |
Returns true if the normals of the vtkCell3D point inside the cell. | |
virtual bool | GetCentroid (double centroid[3]) const =0 |
Computes the centroid of the cell. | |
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. | |
void | Clip (double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, 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. | |
int | GetCellDimension () override |
The topological dimension of the cell. | |
int | Inflate (double dist) override |
Inflates the cell. | |
virtual void | SetMergeTolerance (double) |
Set the tolerance for merging clip intersection points that are near the vertices of cells. | |
virtual double | GetMergeTolerance () |
Set the tolerance for merging clip intersection points that are near the vertices of cells. | |
![]() | |
virtual vtkTypeBool | IsA (const char *type) |
Return 1 if this class is the same type of (or a subclass of) the named class. | |
vtkCell * | NewInstance () 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 () VTK_FUTURE_CONST |
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 () VTK_FUTURE_CONST |
Explicit cells require additional representational information beyond the usual cell type and connectivity list information. | |
virtual int | RequiresExplicitFaceRepresentation () VTK_FUTURE_CONST |
Determine whether the cell requires explicit face representation, and methods for setting and getting the faces (see vtkPolyhedron for example usage of these methods). | |
vtkPoints * | GetPoints () |
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. | |
vtkIdList * | GetPointIds () |
Return the list of point ids defining the cell. | |
vtkIdType | GetPointId (int ptId) |
For cell point i, return the actual point id. | |
virtual vtkCell * | GetEdge (int edgeId)=0 |
Return the edge cell from the edgeId of the cell. | |
virtual vtkCell * | GetFace (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) |
Generate simplices of proper dimension. | |
virtual int | TriangulateIds (int index, vtkIdList *ptIds) |
Generate simplices of proper dimension. | |
virtual int | TriangulateLocalIds (int index, vtkIdList *ptIds)=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 () VTK_FUTURE_CONST |
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 pcoords[3], double *weight) |
Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this level. | |
virtual void | InterpolateDerivs (const double pcoords[3], double *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. | |
![]() | |
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) |
vtkTypeBool | InvokeEvent (unsigned long event) |
vtkTypeBool | InvokeEvent (const char *event) |
std::string | GetObjectDescription () const override |
The object description printed in messages and PrintSelf output. | |
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. | |
vtkCommand * | GetCommand (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. | |
vtkTypeBool | InvokeEvent (unsigned long event, void *callData) |
This method invokes an event and return whether the event was aborted or not. | |
vtkTypeBool | InvokeEvent (const char *event, void *callData) |
This method invokes an event and return whether the event was aborted or not. | |
virtual void | SetObjectName (const std::string &objectName) |
Set/get the name of this object for reporting purposes. | |
virtual std::string | GetObjectName () const |
Set/get the name of this object for reporting purposes. | |
![]() | |
const char * | GetClassName () const |
Return the class name as a string. | |
virtual std::string | GetObjectDescription () const |
The object description printed in messages and PrintSelf output. | |
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. | |
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. | |
virtual bool | UsesGarbageCollector () const |
Indicate whether the class uses vtkGarbageCollector or not. | |
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 | ConstructPolyData () |
void | ConstructLocator () |
![]() | |
virtual vtkObjectBase * | NewInstanceInternal () const |
vtkCell3D () | |
~vtkCell3D () override | |
![]() | |
virtual vtkObjectBase * | NewInstanceInternal () const |
vtkCell () | |
~vtkCell () override | |
![]() | |
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. | |
![]() | |
vtkObjectBase () | |
virtual | ~vtkObjectBase () |
virtual void | RegisterInternal (vtkObjectBase *, vtkTypeBool check) |
virtual void | UnRegisterInternal (vtkObjectBase *, vtkTypeBool check) |
virtual void | ReportReferences (vtkGarbageCollector *) |
virtual void | ObjectFinalize () |
vtkObjectBase (const vtkObjectBase &) | |
void | operator= (const vtkObjectBase &) |
Friends | |
class | vtkPolyhedronUtilities |
typedef vtkCell3D | Superclass |
Standard new methods. | |
static vtkPolyhedron * | New () |
Standard new methods. | |
static vtkTypeBool | IsTypeOf (const char *type) |
Standard new methods. | |
static vtkPolyhedron * | SafeDownCast (vtkObjectBase *o) |
Standard new methods. | |
virtual vtkTypeBool | IsA (const char *type) |
Standard new methods. | |
vtkPolyhedron * | NewInstance () const |
Standard new methods. | |
void | PrintSelf (ostream &os, vtkIndent indent) override |
Standard new methods. | |
virtual vtkObjectBase * | NewInstanceInternal () const |
Standard new methods. | |
Additional Inherited Members | |
![]() | |
static vtkTypeBool | IsTypeOf (const char *type) |
static vtkCell3D * | SafeDownCast (vtkObjectBase *o) |
![]() | |
static vtkTypeBool | IsTypeOf (const char *type) |
static vtkCell * | SafeDownCast (vtkObjectBase *o) |
![]() | |
static vtkObject * | New () |
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 (vtkTypeBool 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 vtkTypeBool | GetGlobalWarningDisplay () |
This is a global flag that controls whether any debug, warning or error messages are displayed. | |
![]() | |
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 vtkObjectBase * | New () |
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. | |
![]() | |
vtkPoints * | Points |
vtkIdList * | PointIds |
![]() | |
static vtkMallocingFunction | GetCurrentMallocFunction () |
static vtkReallocingFunction | GetCurrentReallocFunction () |
static vtkFreeingFunction | GetCurrentFreeFunction () |
static vtkFreeingFunction | GetAlternateFreeFunction () |
A 3D cell defined by a set of polygonal faces.
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:
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.
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:
Interpolation functions and weights are defined / computed using the method of Mean Value Coordinates (MVC). See the VTK class vtkMeanValueCoordinatesInterpolator for more information.
Definition at line 166 of file vtkPolyhedron.h.
using vtkPolyhedron::vtkPointIdMap = std::map<vtkIdType, vtkIdType> |
Definition at line 169 of file vtkPolyhedron.h.
using vtkPolyhedron::Status = vtkCellStatus |
Adopt vtkCellStatus to describe degenerate polyhedral cells.
Definition at line 172 of file vtkPolyhedron.h.
typedef vtkCell3D vtkPolyhedron::Superclass |
Standard new methods.
Definition at line 179 of file vtkPolyhedron.h.
|
protected |
|
overrideprotected |
|
static |
Standard new methods.
|
static |
Standard new methods.
|
virtual |
Standard new methods.
Reimplemented from vtkCell3D.
|
static |
Standard new methods.
|
protectedvirtual |
Standard new methods.
Reimplemented from vtkCell3D.
vtkPolyhedron * vtkPolyhedron::NewInstance | ( | ) | const |
Standard new methods.
|
overridevirtual |
Standard new methods.
Reimplemented from vtkCell3D.
|
inlineoverridevirtual |
See vtkCell3D API for description of these methods.
Implements vtkCell3D.
Definition at line 188 of file vtkPolyhedron.h.
|
inlineoverridevirtual |
See vtkCell3D API for description of these methods.
Implements vtkCell3D.
Definition at line 192 of file vtkPolyhedron.h.
|
inlineoverridevirtual |
See vtkCell3D API for description of these methods.
Implements vtkCell3D.
Definition at line 197 of file vtkPolyhedron.h.
|
inlineoverridevirtual |
See vtkCell3D API for description of these methods.
Implements vtkCell3D.
Definition at line 202 of file vtkPolyhedron.h.
|
inlineoverridevirtual |
See vtkCell3D API for description of these methods.
Implements vtkCell3D.
Definition at line 208 of file vtkPolyhedron.h.
|
overridevirtual |
See vtkCell3D API for description of these methods.
Implements vtkCell3D.
|
inlineoverridevirtual |
See vtkCell3D API for description of these methods.
Implements vtkCell3D.
Definition at line 215 of file vtkPolyhedron.h.
|
inlineoverridevirtual |
See vtkCell3D API for description of these methods.
Implements vtkCell3D.
Definition at line 221 of file vtkPolyhedron.h.
|
overridevirtual |
|
inlineoverridevirtual |
See the vtkCell API for descriptions of these methods.
Implements vtkCell.
Definition at line 236 of file vtkPolyhedron.h.
|
inlineoverridevirtual |
This cell requires that it be initialized prior to access.
Reimplemented from vtkCell.
Definition at line 241 of file vtkPolyhedron.h.
|
overridevirtual |
The Initialize method builds up internal structures of vtkPolyhedron.
Reimplemented from vtkCell.
|
overridevirtual |
A polyhedron is represented internally by a set of polygonal faces.
These faces can be processed to explicitly determine edges.
Implements vtkCell.
|
overridevirtual |
A polyhedron is represented internally by a set of polygonal faces.
These faces can be processed to explicitly determine edges.
Implements vtkCell.
|
overridevirtual |
A polyhedron is represented internally by a set of polygonal faces.
These faces can be processed to explicitly determine edges.
Implements vtkCell.
|
overridevirtual |
A polyhedron is represented internally by a set of polygonal faces.
These faces can be processed to explicitly determine edges.
Implements vtkCell.
|
overridevirtual |
|
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.
Reimplemented from vtkCell3D.
|
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.
|
overridevirtual |
|
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.
|
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.
Implements vtkCell.
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).
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).
|
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.
|
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.
|
overridevirtual |
Return the center of the cell in parametric coordinates.
In this cell, the parametric location (within its bounds) of the centroid of its points is returned.
Reimplemented from vtkCell.
|
inlineoverridevirtual |
A polyhedron is a full-fledged primary cell.
Reimplemented from vtkCell.
Definition at line 366 of file vtkPolyhedron.h.
|
overridevirtual |
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
Here we use the MVC calculation process to compute the interpolation functions.
Reimplemented from vtkCell.
|
overridevirtual |
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
Here we use the MVC calculation process to compute the interpolation functions.
Reimplemented from vtkCell.
|
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 383 of file vtkPolyhedron.h.
void vtkPolyhedron::SetFaces | ( | vtkIdType * | faces | ) |
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.
faces | 1-dimensional array with the following structure : [ NbOfFaces,
NbOfPtsFace1, face1Pt1, face1Pt2, …, face1PtNbOfPtsFace1,
NbOfPtsFace2, face2Pt1, face2Pt2, …, face2PtNbOfPtsFace2,
…,
NbOfPtsFaceN, faceNPt1, faceNPt2, …, faceNPtNbOfPtsFaceN ]
|
vtkIdType * vtkPolyhedron::GetFaces | ( | ) |
Get the faces of the polyhedron.
Face are expressed as sequences of global point IDs .
int vtkPolyhedron::SetCellFaces | ( | vtkCellArray * | faces | ) |
Set the faces of the polyhedron.
Symmetric method to GetCellFaces
faces | vtkCellArray that stores a contiguous list of polygonal faces with their corresponding global point IDs defining a polyhedron. |
vtkCellArray * vtkPolyhedron::GetCellFaces | ( | ) |
Get the faces of the polyhedron.
Face are expressed as sequences of global point IDs .
faces | vtkCellArray that stores the list of polygonal faces with their corresponding global point IDs |
void vtkPolyhedron::GetCellFaces | ( | vtkCellArray * | faces | ) |
Get the faces of the polyhedron.
Face are expressed as sequences of global point IDs .
faces | vtkCellArray that stores the list of polygonal faces with their corresponding global point IDs |
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.
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.
Some variants have a return value which indicates a particular reason for the cell to be marked non-convex (currently non-planar face(s) or a concave arrangement of neighboring faces). These variants accept a planarThreshold for the degree of non-planarity allowed. The planarThreshold is the maximum ratio of the distance out of the plane of the face compared to the size in the plane of the polygon. If a face exceeds planarThreshold, then IsConvex() will return false. Pass planarThreshold < 0 to ignore non-planar faces. The default planarThreshold is 0.1.
Status vtkPolyhedron::IsConvex | ( | double | planarThreshold | ) |
vtkPolyData * vtkPolyhedron::GetPolyData | ( | ) |
Construct polydata if no one exist, then return this->PolyData.
|
overridevirtual |
Shallow copy of a polyhedron.
Reimplemented from vtkCell.
|
overridevirtual |
Deep copy of a polyhedron.
Reimplemented from vtkCell.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
friend |
Definition at line 542 of file vtkPolyhedron.h.
Definition at line 492 of file vtkPolyhedron.h.
|
protected |
Definition at line 493 of file vtkPolyhedron.h.
Definition at line 494 of file vtkPolyhedron.h.
|
protected |
Definition at line 495 of file vtkPolyhedron.h.
Definition at line 496 of file vtkPolyhedron.h.
|
protected |
Definition at line 500 of file vtkPolyhedron.h.
|
protected |
Definition at line 504 of file vtkPolyhedron.h.
|
protected |
Definition at line 507 of file vtkPolyhedron.h.
|
protected |
Definition at line 508 of file vtkPolyhedron.h.
|
protected |
Definition at line 509 of file vtkPolyhedron.h.
|
protected |
Definition at line 510 of file vtkPolyhedron.h.
|
protected |
Definition at line 519 of file vtkPolyhedron.h.
|
protected |
Definition at line 520 of file vtkPolyhedron.h.
|
protected |
Definition at line 523 of file vtkPolyhedron.h.
|
protected |
Definition at line 529 of file vtkPolyhedron.h.
|
protected |
Definition at line 530 of file vtkPolyhedron.h.
|
protected |
Definition at line 532 of file vtkPolyhedron.h.
|
protected |
Definition at line 533 of file vtkPolyhedron.h.
Definition at line 535 of file vtkPolyhedron.h.
|
protected |
Definition at line 536 of file vtkPolyhedron.h.