VTK
|
a 3D cell defined by a set of polygonal faces More...
#include <vtkPolyhedron.h>
a 3D cell defined by a set of polygonal faces
vtkPolyhedron is a concrete implementation that represents a 3D cell defined by a set of polygonal faces. The polyhedron should be watertight, non-self-intersecting and manifold (each edge is used twice).
Interpolation functions and weights are defined / computed using the method of Mean Value Coordinates (MVC). See the VTK class vtkMeanValueCoordinatesInterpolator for more information.
The class assumes that the polyhedron is non-convex. However, the polygonal faces should be planar. Non-planar polygonal faces will definitely cause problems, especially in severely warped situations.
Definition at line 59 of file vtkPolyhedron.h.
typedef vtkCell3D vtkPolyhedron::Superclass |
vtkPolyhedron::vtkPolyhedron | ( | ) | [protected] |
vtkPolyhedron::~vtkPolyhedron | ( | ) | [protected] |
static vtkPolyhedron* vtkPolyhedron::New | ( | ) | [static] |
Standard new methods.
Reimplemented from vtkObject.
static int vtkPolyhedron::IsTypeOf | ( | const char * | type | ) | [static] |
Standard new methods.
Reimplemented from vtkCell3D.
virtual int vtkPolyhedron::IsA | ( | const char * | type | ) | [virtual] |
Standard new methods.
Reimplemented from vtkCell3D.
static vtkPolyhedron* vtkPolyhedron::SafeDownCast | ( | vtkObjectBase * | o | ) | [static] |
Standard new methods.
Reimplemented from vtkCell3D.
virtual vtkObjectBase* vtkPolyhedron::NewInstanceInternal | ( | ) | const [protected, virtual] |
Standard new methods.
Reimplemented from vtkCell3D.
vtkPolyhedron* vtkPolyhedron::NewInstance | ( | ) | const |
Standard new methods.
Reimplemented from vtkCell3D.
void vtkPolyhedron::PrintSelf | ( | ostream & | os, |
vtkIndent | indent | ||
) | [virtual] |
Standard new methods.
Reimplemented from vtkCell3D.
virtual void vtkPolyhedron::GetEdgePoints | ( | int | vtkNotUsededgeId, |
int *& | vtkNotUsedpts | ||
) | [inline, virtual] |
See vtkCell3D API for description of these methods.
Definition at line 71 of file vtkPolyhedron.h.
virtual void vtkPolyhedron::GetFacePoints | ( | int | vtkNotUsedfaceId, |
int *& | vtkNotUsedpts | ||
) | [inline, virtual] |
See vtkCell3D API for description of these methods.
Definition at line 72 of file vtkPolyhedron.h.
virtual double* vtkPolyhedron::GetParametricCoords | ( | ) | [virtual] |
virtual int vtkPolyhedron::GetCellType | ( | ) | [inline, virtual] |
See the vtkCell API for descriptions of these methods.
Implements vtkCell.
Definition at line 77 of file vtkPolyhedron.h.
virtual int vtkPolyhedron::RequiresInitialization | ( | ) | [inline, virtual] |
This cell requires that it be initialized prior to access.
Reimplemented from vtkCell.
Definition at line 81 of file vtkPolyhedron.h.
virtual void vtkPolyhedron::Initialize | ( | ) | [virtual] |
This cell requires that it be initialized prior to access.
Reimplemented from vtkCell.
virtual int vtkPolyhedron::GetNumberOfEdges | ( | ) | [virtual] |
A polyhedron is represented internally by a set of polygonal faces. These faces can be processed to explicitly determine edges.
Implements vtkCell.
virtual vtkCell* vtkPolyhedron::GetEdge | ( | int | ) | [virtual] |
A polyhedron is represented internally by a set of polygonal faces. These faces can be processed to explicitly determine edges.
Implements vtkCell.
virtual int vtkPolyhedron::GetNumberOfFaces | ( | ) | [virtual] |
A polyhedron is represented internally by a set of polygonal faces. These faces can be processed to explicitly determine edges.
Implements vtkCell.
virtual vtkCell* vtkPolyhedron::GetFace | ( | int | faceId | ) | [virtual] |
A polyhedron is represented internally by a set of polygonal faces. These faces can be processed to explicitly determine edges.
Implements vtkCell.
virtual 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 | ||
) | [virtual] |
virtual void vtkPolyhedron::Clip | ( | double | value, |
vtkDataArray * | scalars, | ||
vtkIncrementalPointLocator * | locator, | ||
vtkCellArray * | connectivity, | ||
vtkPointData * | inPd, | ||
vtkPointData * | outPd, | ||
vtkCellData * | inCd, | ||
vtkIdType | cellId, | ||
vtkCellData * | outCd, | ||
int | insideOut | ||
) | [virtual] |
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. Note: the algorithm assumes water-tight polyhedron cells.
Reimplemented from vtkCell3D.
virtual int vtkPolyhedron::EvaluatePosition | ( | double | x[3], |
double * | closestPoint, | ||
int & | subId, | ||
double | pcoords[3], | ||
double & | dist2, | ||
double * | weights | ||
) | [virtual] |
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.
virtual void vtkPolyhedron::EvaluateLocation | ( | int & | subId, |
double | pcoords[3], | ||
double | x[3], | ||
double * | weights | ||
) | [virtual] |
The inverse of EvaluatePosition. Note the weights should be the MVC weights.
Implements vtkCell.
virtual int vtkPolyhedron::IntersectWithLine | ( | double | p1[3], |
double | p2[3], | ||
double | tol, | ||
double & | t, | ||
double | x[3], | ||
double | pcoords[3], | ||
int & | subId | ||
) | [virtual] |
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 the number of intersection points.
Implements vtkCell.
virtual int vtkPolyhedron::Triangulate | ( | int | index, |
vtkIdList * | ptIds, | ||
vtkPoints * | pts | ||
) | [virtual] |
Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh. This method works well for a convex polyhedron but may return wrong result in a concave case. 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.
virtual void vtkPolyhedron::Derivatives | ( | int | subId, |
double | pcoords[3], | ||
double * | values, | ||
int | dim, | ||
double * | derivs | ||
) | [virtual] |
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.
virtual int vtkPolyhedron::CellBoundary | ( | int | subId, |
double | pcoords[3], | ||
vtkIdList * | pts | ||
) | [virtual] |
Find the boundary face closest to the point defined by the pcoords[3] and subId of the cell (subId can be ignored).
Implements vtkCell.
int vtkPolyhedron::GetParametricCenter | ( | double | pcoords[3] | ) | [inline, virtual] |
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 301 of file vtkPolyhedron.h.
int vtkPolyhedron::IsPrimaryCell | ( | ) | [inline, virtual] |
A polyhedron is a full-fledged primary cell.
Reimplemented from vtkCell.
Definition at line 184 of file vtkPolyhedron.h.
virtual void vtkPolyhedron::InterpolateFunctions | ( | double | x[3], |
double * | sf | ||
) | [virtual] |
Compute the interpolation functions/derivatives (aka shape functions/derivatives). Here we use the MVC calculation process to compute the interpolation functions.
virtual void vtkPolyhedron::InterpolateDerivs | ( | double | x[3], |
double * | derivs | ||
) | [virtual] |
Compute the interpolation functions/derivatives (aka shape functions/derivatives). Here we use the MVC calculation process to compute the interpolation functions.
virtual int vtkPolyhedron::RequiresExplicitFaceRepresentation | ( | ) | [inline, virtual] |
Methods supporting the definition of faces. Note that the GetFaces() returns a list of faces in vtkCellArray form; use the method GetNumberOfFaces() to determine the number of faces in the list. The SetFaces() method is also in vtkCellArray form, except that it begins with a leading count indicating the total number of faces in the list.
Reimplemented from vtkCell.
Definition at line 201 of file vtkPolyhedron.h.
virtual void vtkPolyhedron::SetFaces | ( | vtkIdType * | faces | ) | [virtual] |
Methods supporting the definition of faces. Note that the GetFaces() returns a list of faces in vtkCellArray form; use the method GetNumberOfFaces() to determine the number of faces in the list. The SetFaces() method is also in vtkCellArray form, except that it begins with a leading count indicating the total number of faces in the list.
virtual vtkIdType* vtkPolyhedron::GetFaces | ( | ) | [virtual] |
Methods supporting the definition of faces. Note that the GetFaces() returns a list of faces in vtkCellArray form; use the method GetNumberOfFaces() to determine the number of faces in the list. The SetFaces() method is also in vtkCellArray form, except that it begins with a leading count indicating the total number of faces in the list.
Reimplemented from vtkCell.
int vtkPolyhedron::IsInside | ( | double | x[3], |
double | tolerance | ||
) |
Construct polydata if no one exist, then return this->PolyData
int vtkPolyhedron::GenerateEdges | ( | ) | [protected] |
void vtkPolyhedron::GenerateFaces | ( | ) | [protected] |
void vtkPolyhedron::ComputeBounds | ( | ) | [protected] |
void vtkPolyhedron::ComputeParametricCoordinate | ( | double | x[3], |
double | pc[3] | ||
) | [protected] |
void vtkPolyhedron::ComputePositionFromParametricCoordinate | ( | double | pc[3], |
double | x[3] | ||
) | [protected] |
void vtkPolyhedron::ConstructPolyData | ( | ) | [protected] |
void vtkPolyhedron::ConstructLocator | ( | ) | [protected] |
int vtkPolyhedron::InternalContour | ( | double | value, |
int | insideOut, | ||
vtkIncrementalPointLocator * | locator, | ||
vtkDataArray * | inScalars, | ||
vtkDataArray * | outScalars, | ||
vtkPointData * | inPd, | ||
vtkPointData * | outPd, | ||
vtkCellArray * | contourPolys, | ||
vtkIdToIdVectorMapType & | faceToPointsMap, | ||
vtkIdToIdVectorMapType & | pointToFacesMap, | ||
vtkIdToIdMapType & | pointIdMap | ||
) | [protected] |
int vtkPolyhedron::IntersectWithContour | ( | double | value, |
int | insideOut, | ||
vtkDataArray * | inScalars | ||
) | [protected] |
vtkLine* vtkPolyhedron::Line [protected] |
Definition at line 221 of file vtkPolyhedron.h.
vtkTriangle* vtkPolyhedron::Triangle [protected] |
Definition at line 222 of file vtkPolyhedron.h.
vtkQuad* vtkPolyhedron::Quad [protected] |
Definition at line 223 of file vtkPolyhedron.h.
vtkPolygon* vtkPolyhedron::Polygon [protected] |
Definition at line 224 of file vtkPolyhedron.h.
vtkTetra* vtkPolyhedron::Tetra [protected] |
Definition at line 225 of file vtkPolyhedron.h.
vtkIdTypeArray* vtkPolyhedron::GlobalFaces [protected] |
Definition at line 226 of file vtkPolyhedron.h.
vtkIdTypeArray* vtkPolyhedron::FaceLocations [protected] |
Definition at line 227 of file vtkPolyhedron.h.
vtkPointIdMap* vtkPolyhedron::PointIdMap [protected] |
Definition at line 234 of file vtkPolyhedron.h.
int vtkPolyhedron::EdgesGenerated [protected] |
Definition at line 238 of file vtkPolyhedron.h.
vtkEdgeTable* vtkPolyhedron::EdgeTable [protected] |
Definition at line 239 of file vtkPolyhedron.h.
vtkIdTypeArray* vtkPolyhedron::Edges [protected] |
Definition at line 240 of file vtkPolyhedron.h.
vtkIdTypeArray* vtkPolyhedron::Faces [protected] |
Definition at line 246 of file vtkPolyhedron.h.
int vtkPolyhedron::FacesGenerated [protected] |
Definition at line 247 of file vtkPolyhedron.h.
int vtkPolyhedron::BoundsComputed [protected] |
Definition at line 251 of file vtkPolyhedron.h.
int vtkPolyhedron::PolyDataConstructed [protected] |
Definition at line 257 of file vtkPolyhedron.h.
vtkPolyData* vtkPolyhedron::PolyData [protected] |
Definition at line 258 of file vtkPolyhedron.h.
vtkCellArray* vtkPolyhedron::Polys [protected] |
Definition at line 259 of file vtkPolyhedron.h.
vtkIdTypeArray* vtkPolyhedron::PolyConnectivity [protected] |
Definition at line 260 of file vtkPolyhedron.h.
int vtkPolyhedron::LocatorConstructed [protected] |
Definition at line 262 of file vtkPolyhedron.h.
vtkCellLocator* vtkPolyhedron::CellLocator [protected] |
Definition at line 263 of file vtkPolyhedron.h.
vtkIdList* vtkPolyhedron::CellIds [protected] |
Definition at line 265 of file vtkPolyhedron.h.
vtkGenericCell* vtkPolyhedron::Cell [protected] |
Definition at line 266 of file vtkPolyhedron.h.