VTK
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes
vtkCell Class Reference

abstract class to specify cell behavior More...

#include <vtkCell.h>

Inheritance diagram for vtkCell:
Inheritance graph
[legend]
Collaboration diagram for vtkCell:
Collaboration graph
[legend]

List of all members.

Public Types

typedef vtkObject Superclass

Public Member Functions

virtual int IsA (const char *type)
vtkCellNewInstance () const
void PrintSelf (ostream &os, vtkIndent indent)
void Initialize (int npts, vtkIdType *pts, vtkPoints *p)
virtual void ShallowCopy (vtkCell *c)
virtual void DeepCopy (vtkCell *c)
virtual int GetCellType ()=0
virtual int GetCellDimension ()=0
virtual int IsLinear ()
virtual int IsExplicitCell ()
vtkPointsGetPoints ()
vtkIdType GetNumberOfPoints ()
virtual int GetNumberOfEdges ()=0
virtual int GetNumberOfFaces ()=0
vtkIdListGetPointIds ()
vtkIdType GetPointId (int ptId)
virtual vtkCellGetEdge (int edgeId)=0
virtual vtkCellGetFace (int faceId)=0
virtual int CellBoundary (int subId, double pcoords[3], vtkIdList *pts)=0
virtual int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts)=0
void GetBounds (double bounds[6])
doubleGetBounds ()
double GetLength2 ()
virtual int GetParametricCenter (double pcoords[3])
virtual double GetParametricDistance (double pcoords[3])
virtual int IsPrimaryCell ()
virtual doubleGetParametricCoords ()
virtual int RequiresInitialization ()
virtual void Initialize ()
virtual int RequiresExplicitFaceRepresentation ()
virtual void SetFaces (vtkIdType *vtkNotUsed(faces))
virtual vtkIdTypeGetFaces ()
virtual int EvaluatePosition (double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights)=0
virtual void EvaluateLocation (int &subId, double pcoords[3], double x[3], double *weights)=0
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
virtual void Clip (double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)=0
virtual int IntersectWithLine (double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
virtual void Derivatives (int subId, double pcoords[3], double *values, int dim, double *derivs)=0
virtual void InterpolateFunctions (double pcoords[3], double weights[3])
virtual void InterpolateDerivs (double pcoords[3], double derivs[3])

Static Public Member Functions

static int IsTypeOf (const char *type)
static vtkCellSafeDownCast (vtkObjectBase *o)

Public Attributes

vtkPointsPoints
vtkIdListPointIds

Protected Member Functions

virtual vtkObjectBaseNewInstanceInternal () const
 vtkCell ()
 ~vtkCell ()

Protected Attributes

double Bounds [6]

Detailed Description

abstract class to specify cell behavior

vtkCell is an abstract class that specifies the interfaces for data cells. Data cells are simple topological elements like points, lines, polygons, and tetrahedra of which visualization datasets are composed. In some cases visualization datasets may explicitly represent cells (e.g., vtkPolyData, vtkUnstructuredGrid), and in some cases, the datasets are implicitly composed of cells (e.g., vtkStructuredPoints).

Warning:
The #define VTK_CELL_SIZE is a parameter used to construct cells and provide a general guideline for controlling object execution. This parameter is not a hard boundary: you can create cells with more points.
See also:
vtkHexahedron vtkLine vtkPixel vtkPolyLine vtkPolyVertex vtkPolygon vtkQuad vtkTetra vtkTriangle vtkTriangleStrip vtkVertex vtkVoxel vtkWedge vtkPyramid
Tests:
vtkCell (Tests)

Definition at line 58 of file vtkCell.h.


Member Typedef Documentation


Constructor & Destructor Documentation

vtkCell::vtkCell ( ) [protected]
vtkCell::~vtkCell ( ) [protected]

Member Function Documentation

static int vtkCell::IsTypeOf ( const char *  name) [static]
virtual int vtkCell::IsA ( const char *  name) [virtual]
static vtkCell* vtkCell::SafeDownCast ( vtkObjectBase o) [static]
virtual vtkObjectBase* vtkCell::NewInstanceInternal ( ) const [protected, virtual]
void vtkCell::PrintSelf ( ostream &  os,
vtkIndent  indent 
) [virtual]
void vtkCell::Initialize ( int  npts,
vtkIdType pts,
vtkPoints p 
)

Initialize cell from outside with point ids and point coordinates specified.

virtual void vtkCell::ShallowCopy ( vtkCell c) [virtual]

Copy this cell by reference counting the internal data structures. This is safe if you want a "read-only" copy. If you modify the cell you might wish to use DeepCopy().

Reimplemented in vtkGenericCell.

virtual void vtkCell::DeepCopy ( vtkCell c) [virtual]

Copy this cell by completely copying internal data structures. This is slower but safer than ShallowCopy().

Reimplemented in vtkGenericCell.

virtual int vtkCell::GetCellType ( ) [pure virtual]
virtual int vtkCell::GetCellDimension ( ) [pure virtual]
virtual int vtkCell::IsLinear ( ) [inline, virtual]

Non-linear cells require special treatment beyond the usual cell type and connectivity list information. Most cells in VTK are implicit cells.

Reimplemented in vtkGenericCell, and vtkNonLinearCell.

Definition at line 86 of file vtkCell.h.

virtual int vtkCell::RequiresInitialization ( ) [inline, virtual]

Some cells require initialization prior to access. For example, they may have to triangulate themselves or set up internal data structures.

Reimplemented in vtkPolyhedron, vtkConvexPointSet, and vtkGenericCell.

Definition at line 92 of file vtkCell.h.

virtual void vtkCell::Initialize ( ) [inline, virtual]

Some cells require initialization prior to access. For example, they may have to triangulate themselves or set up internal data structures.

Reimplemented in vtkPolyhedron, vtkConvexPointSet, and vtkGenericCell.

Definition at line 93 of file vtkCell.h.

virtual int vtkCell::IsExplicitCell ( ) [inline, virtual]

Explicit cells require additional representational information beyond the usual cell type and connectivity list information. Most cells in VTK are implicit cells.

Reimplemented in vtkExplicitCell.

Definition at line 99 of file vtkCell.h.

virtual int vtkCell::RequiresExplicitFaceRepresentation ( ) [inline, virtual]

Determine whether the cell requires explicit face representation, and methods for setting and getting the faces (see vtkPolyhedron for example usage of these methods).

Reimplemented in vtkPolyhedron, and vtkGenericCell.

Definition at line 105 of file vtkCell.h.

virtual void vtkCell::SetFaces ( vtkIdType vtkNotUsedfaces) [inline, virtual]

Determine whether the cell requires explicit face representation, and methods for setting and getting the faces (see vtkPolyhedron for example usage of these methods).

Definition at line 106 of file vtkCell.h.

virtual vtkIdType* vtkCell::GetFaces ( ) [inline, virtual]

Determine whether the cell requires explicit face representation, and methods for setting and getting the faces (see vtkPolyhedron for example usage of these methods).

Reimplemented in vtkPolyhedron, and vtkGenericCell.

Definition at line 107 of file vtkCell.h.

Get the point coordinates for the cell.

Definition at line 111 of file vtkCell.h.

Return the number of points in the cell.

Definition at line 114 of file vtkCell.h.

virtual int vtkCell::GetNumberOfEdges ( ) [pure virtual]
virtual int vtkCell::GetNumberOfFaces ( ) [pure virtual]

Return the list of point ids defining the cell.

Definition at line 123 of file vtkCell.h.

vtkIdType vtkCell::GetPointId ( int  ptId) [inline]

For cell point i, return the actual point id.

Definition at line 126 of file vtkCell.h.

virtual vtkCell* vtkCell::GetEdge ( int  edgeId) [pure virtual]
virtual vtkCell* vtkCell::GetFace ( int  faceId) [pure virtual]
virtual int vtkCell::CellBoundary ( int  subId,
double  pcoords[3],
vtkIdList pts 
) [pure virtual]

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

Implemented in vtkPolyhedron, vtkConvexPointSet, vtkTriQuadraticHexahedron, vtkBiQuadraticQuadraticHexahedron, vtkTetra, vtkBiQuadraticQuadraticWedge, vtkQuadraticLinearWedge, vtkQuadraticPyramid, vtkBiQuadraticQuad, vtkVertex, vtkPentagonalPrism, vtkPyramid, vtkWedge, vtkBiQuadraticTriangle, vtkHexagonalPrism, vtkPolyLine, vtkQuadraticLinearQuad, vtkHexahedron, vtkQuadraticWedge, vtkQuadraticHexahedron, vtkVoxel, vtkQuadraticTetra, vtkGenericCell, vtkQuadraticQuad, vtkQuadraticTriangle, vtkQuadraticEdge, vtkCubicLine, vtkPolygon, vtkTriangle, vtkTriangleStrip, vtkQuad, vtkPixel, vtkPolyVertex, vtkLine, and vtkEmptyCell.

virtual int vtkCell::EvaluatePosition ( double  x[3],
double closestPoint,
int subId,
double  pcoords[3],
double dist2,
double weights 
) [pure virtual]

Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; evaluate parametric coordinates, sub-cell id (!=0 only if cell is composite), distance squared of point x[3] to cell (in particular, the sub-cell indicated), closest point on cell to x[3] (unless closestPoint is null, in which case, the closest point and dist2 are not found), and interpolation weights in cell. (The number of weights is equal to the number of points defining the cell). Note: on rare occasions a -1 is returned from the method. This means that numerical error has occurred and all data returned from this method should be ignored. Also, inside/outside is determine parametrically. That is, a point is inside if it satisfies parametric limits. This can cause problems for cells of topological dimension 2 or less, since a point in 3D can project onto the cell within parametric limits but be "far" from the cell. Thus the value dist2 may be checked to determine true in/out.

Implemented in vtkPolyhedron, vtkConvexPointSet, vtkTriQuadraticHexahedron, vtkBiQuadraticQuadraticHexahedron, vtkQuadraticLinearWedge, vtkBiQuadraticQuadraticWedge, vtkPolyLine, vtkQuadraticPyramid, vtkTetra, vtkPyramid, vtkWedge, vtkBiQuadraticTriangle, vtkHexahedron, vtkBiQuadraticQuad, vtkPentagonalPrism, vtkQuadraticLinearQuad, vtkQuadraticWedge, vtkHexagonalPrism, vtkQuadraticHexahedron, vtkVoxel, vtkPolygon, vtkQuadraticTetra, vtkTriangleStrip, vtkQuadraticQuad, vtkQuadraticTriangle, vtkQuadraticEdge, vtkCubicLine, vtkPixel, vtkGenericCell, vtkTriangle, vtkPolyVertex, vtkQuad, vtkEmptyCell, vtkVertex, and vtkLine.

virtual void vtkCell::EvaluateLocation ( int subId,
double  pcoords[3],
double  x[3],
double weights 
) [pure virtual]
virtual void vtkCell::Contour ( double  value,
vtkDataArray cellScalars,
vtkIncrementalPointLocator locator,
vtkCellArray verts,
vtkCellArray lines,
vtkCellArray polys,
vtkPointData inPd,
vtkPointData outPd,
vtkCellData inCd,
vtkIdType  cellId,
vtkCellData outCd 
) [pure virtual]

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

Implemented in vtkTriQuadraticHexahedron, vtkBiQuadraticQuadraticHexahedron, vtkPolyhedron, vtkConvexPointSet, vtkBiQuadraticQuad, vtkQuadraticLinearWedge, vtkVertex, vtkBiQuadraticQuadraticWedge, vtkQuadraticPyramid, vtkPyramid, vtkWedge, vtkBiQuadraticTriangle, vtkGenericCell, vtkPolyLine, vtkQuadraticLinearQuad, vtkHexahedron, vtkQuadraticWedge, vtkQuadraticHexahedron, vtkTetra, vtkVoxel, vtkQuadraticTetra, vtkQuadraticQuad, vtkQuadraticTriangle, vtkQuadraticEdge, vtkCubicLine, vtkPolygon, vtkTriangle, vtkTriangleStrip, vtkCell3D, vtkQuad, vtkPixel, vtkPolyVertex, vtkLine, and vtkEmptyCell.

virtual void vtkCell::Clip ( double  value,
vtkDataArray cellScalars,
vtkIncrementalPointLocator locator,
vtkCellArray connectivity,
vtkPointData inPd,
vtkPointData outPd,
vtkCellData inCd,
vtkIdType  cellId,
vtkCellData outCd,
int  insideOut 
) [pure virtual]

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

Implemented in vtkBiQuadraticQuadraticHexahedron, vtkTriQuadraticHexahedron, vtkPolyhedron, vtkConvexPointSet, vtkQuadraticLinearWedge, vtkBiQuadraticQuad, vtkBiQuadraticQuadraticWedge, vtkQuadraticPyramid, vtkBiQuadraticTriangle, vtkQuadraticWedge, vtkQuadraticHexahedron, vtkQuadraticLinearQuad, vtkCubicLine, vtkQuadraticTetra, vtkQuadraticTriangle, vtkQuadraticQuad, vtkTriangle, vtkQuad, vtkQuadraticEdge, vtkCell3D, vtkGenericCell, vtkPolyLine, vtkTetra, vtkLine, vtkPolygon, vtkTriangleStrip, vtkPixel, vtkPolyVertex, vtkEmptyCell, and vtkVertex.

virtual int vtkCell::IntersectWithLine ( double  p1[3],
double  p2[3],
double  tol,
double t,
double  x[3],
double  pcoords[3],
int subId 
) [pure virtual]
virtual int vtkCell::Triangulate ( int  index,
vtkIdList ptIds,
vtkPoints pts 
) [pure virtual]

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

Implemented in vtkPolyhedron, vtkConvexPointSet, vtkBiQuadraticQuadraticHexahedron, vtkTriQuadraticHexahedron, vtkVertex, vtkPolyLine, vtkQuadraticLinearWedge, vtkBiQuadraticQuadraticWedge, vtkQuadraticPyramid, vtkTetra, vtkPyramid, vtkWedge, vtkGenericCell, vtkHexahedron, vtkBiQuadraticTriangle, vtkPentagonalPrism, vtkHexagonalPrism, vtkVoxel, vtkBiQuadraticQuad, vtkPixel, vtkPolygon, vtkQuadraticLinearQuad, vtkQuadraticWedge, vtkTriangleStrip, vtkQuadraticHexahedron, vtkQuadraticTetra, vtkQuadraticQuad, vtkQuadraticTriangle, vtkQuadraticEdge, vtkCubicLine, vtkPolyVertex, vtkTriangle, vtkQuad, vtkEmptyCell, and vtkLine.

virtual void vtkCell::Derivatives ( int  subId,
double  pcoords[3],
double values,
int  dim,
double derivs 
) [pure virtual]

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

Implemented in vtkPolyhedron, vtkConvexPointSet, vtkBiQuadraticQuadraticHexahedron, vtkTriQuadraticHexahedron, vtkVertex, vtkPolyLine, vtkQuadraticLinearWedge, vtkBiQuadraticQuadraticWedge, vtkQuadraticPyramid, vtkTetra, vtkPyramid, vtkWedge, vtkGenericCell, vtkHexahedron, vtkBiQuadraticTriangle, vtkPentagonalPrism, vtkHexagonalPrism, vtkVoxel, vtkBiQuadraticQuad, vtkPixel, vtkPolygon, vtkQuadraticLinearQuad, vtkQuadraticWedge, vtkTriangleStrip, vtkQuadraticHexahedron, vtkQuadraticTetra, vtkQuadraticQuad, vtkQuadraticTriangle, vtkQuadraticEdge, vtkCubicLine, vtkPolyVertex, vtkTriangle, vtkQuad, vtkEmptyCell, and vtkLine.

void vtkCell::GetBounds ( double  bounds[6])

Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). Copy result into user provided array.

Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). Return pointer to array of six double values.

Compute Length squared of cell (i.e., bounding box diagonal squared).

virtual int vtkCell::GetParametricCenter ( double  pcoords[3]) [virtual]

Return center of the cell in parametric coordinates. Note that the parametric center is not always located at (0.5,0.5,0.5). The return value is the subId that the center is in (if a composite cell). If you want the center in x-y-z space, invoke the EvaluateLocation() method.

Reimplemented in vtkPolyhedron, vtkConvexPointSet, vtkTriangle, vtkQuadraticLinearWedge, vtkBiQuadraticQuadraticWedge, vtkQuadraticPyramid, vtkBiQuadraticQuad, vtkBiQuadraticTriangle, vtkQuadraticWedge, vtkQuadraticLinearQuad, vtkQuadraticTriangle, vtkQuadraticTetra, vtkQuadraticQuad, vtkTetra, vtkQuadraticEdge, vtkPolyLine, vtkCubicLine, vtkPyramid, vtkWedge, vtkPentagonalPrism, vtkVertex, vtkGenericCell, vtkHexagonalPrism, vtkTriangleStrip, vtkLine, vtkPolyVertex, vtkQuad, and vtkPixel.

virtual double vtkCell::GetParametricDistance ( double  pcoords[3]) [virtual]

Return the distance of the parametric coordinate provided to the cell. If inside the cell, a distance of zero is returned. This is used during picking to get the correct cell picked. (The tolerance will occasionally allow cells to be picked who are not really intersected "inside" the cell.)

Reimplemented in vtkTriangle, vtkBiQuadraticTriangle, vtkQuadraticTriangle, vtkQuadraticTetra, vtkTetra, and vtkCubicLine.

virtual int vtkCell::IsPrimaryCell ( ) [inline, virtual]

Return whether this cell type has a fixed topology or whether the topology varies depending on the data (e.g., vtkConvexPointSet). This compares to composite cells that are typically composed of primary cells (e.g., a triangle strip composite cell is made up of triangle primary cells).

Reimplemented in vtkPolyhedron, vtkConvexPointSet, vtkPolyLine, vtkGenericCell, vtkPolygon, vtkTriangleStrip, and vtkPolyVertex.

Definition at line 282 of file vtkCell.h.

virtual double* vtkCell::GetParametricCoords ( ) [virtual]

Return a contiguous array of parametric coordinates of the points defining this cell. In other words, (px,py,pz, px,py,pz, etc..) The coordinates are ordered consistent with the definition of the point ordering for the cell. This method returns a non-NULL pointer when the cell is a primary type (i.e., IsPrimaryCell() is true). Note that 3D parametric coordinates are returned no matter what the topological dimension of the cell.

Reimplemented in vtkBiQuadraticQuadraticHexahedron, vtkTriQuadraticHexahedron, vtkQuadraticLinearWedge, vtkBiQuadraticQuadraticWedge, vtkQuadraticPyramid, vtkTetra, vtkGenericCell, vtkPyramid, vtkWedge, vtkHexahedron, vtkBiQuadraticTriangle, vtkPentagonalPrism, vtkHexagonalPrism, vtkBiQuadraticQuad, vtkPixel, vtkQuadraticLinearQuad, vtkQuadraticWedge, vtkQuadraticHexahedron, vtkQuadraticTetra, vtkQuadraticQuad, vtkQuadraticTriangle, vtkQuadraticEdge, vtkCubicLine, vtkPolyhedron, vtkTriangle, vtkQuad, vtkLine, vtkVertex, vtkConvexPointSet, and vtkVoxel.

virtual void vtkCell::InterpolateFunctions ( double  pcoords[3],
double  weights[3] 
) [inline, virtual]

Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this level. Typically overridden in subclasses.

Reimplemented in vtkQuadraticEdge, and vtkTriangle.

Definition at line 298 of file vtkCell.h.

virtual void vtkCell::InterpolateDerivs ( double  pcoords[3],
double  derivs[3] 
) [inline, virtual]

Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this level. Typically overridden in subclasses.

Reimplemented in vtkVertex, and vtkQuadraticEdge.

Definition at line 303 of file vtkCell.h.


Member Data Documentation

Definition at line 311 of file vtkCell.h.

Definition at line 312 of file vtkCell.h.

double vtkCell::Bounds[6] [protected]

Definition at line 318 of file vtkCell.h.


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