VTK
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkQuadraticLinearWedge Class Reference

cell represents a, 12-node isoparametric wedge More...

#include <vtkQuadraticLinearWedge.h>

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

Public Types

typedef vtkNonLinearCell Superclass
 
- Public Types inherited from vtkNonLinearCell
typedef vtkCell Superclass
 

Public Member Functions

virtual const char * GetClassName ()
 
virtual int IsA (const char *type)
 
void PrintSelf (ostream &os, vtkIndent indent)
 
int CellBoundary (int subId, double pcoords[3], vtkIdList *pts)
 
int GetParametricCenter (double pcoords[3])
 
void JacobianInverse (double pcoords[3], double **inverse, double derivs[45])
 
int GetCellType ()
 
int GetCellDimension ()
 
int GetNumberOfEdges ()
 
int GetNumberOfFaces ()
 
vtkCell * GetEdge (int edgeId)
 
vtkCell * GetFace (int faceId)
 
void Contour (double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd)
 
int EvaluatePosition (double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights)
 
void EvaluateLocation (int &subId, double pcoords[3], double x[3], double *weights)
 
int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts)
 
void Derivatives (int subId, double pcoords[3], double *values, int dim, double *derivs)
 
virtual double * GetParametricCoords ()
 
void Clip (double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tetras, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)
 
int IntersectWithLine (double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)
 
- Public Member Functions inherited from vtkNonLinearCell
void PrintSelf (ostream &os, vtkIndent indent)
 
virtual int IsLinear ()
 

Static Public Member Functions

static vtkQuadraticLinearWedgeNew ()
 
static int IsTypeOf (const char *type)
 
static vtkQuadraticLinearWedgeSafeDownCast (vtkObject *o)
 
- Static Public Member Functions inherited from vtkNonLinearCell
static int IsTypeOf (const char *type)
 
static vtkNonLinearCellSafeDownCast (vtkObject *o)
 

Protected Member Functions

 vtkQuadraticLinearWedge ()
 
 ~vtkQuadraticLinearWedge ()
 
- Protected Member Functions inherited from vtkNonLinearCell
 vtkNonLinearCell ()
 
 ~vtkNonLinearCell ()
 

Protected Attributes

vtkQuadraticEdge * QuadEdge
 
vtkLineEdge
 
vtkQuadraticTriangle * TriangleFace
 
vtkQuadraticLinearQuad * Face
 
vtkWedge * Wedge
 
vtkDoubleArrayScalars
 
static void InterpolationFunctions (double pcoords[3], double weights[15])
 
static void InterpolationDerivs (double pcoords[3], double derivs[45])
 
static int * GetEdgeArray (int edgeId)
 
static int * GetFaceArray (int faceId)
 
virtual void InterpolateFunctions (double pcoords[3], double weights[15])
 
virtual void InterpolateDerivs (double pcoords[3], double derivs[45])
 

Detailed Description

cell represents a, 12-node isoparametric wedge

vtkQuadraticLinearWedge is a concrete implementation of vtkNonLinearCell to represent a three-dimensional, 12-node isoparametric linear quadratic wedge. The interpolation is the standard finite element, quadratic isoparametric shape function in xy - layer and the linear functions in z - direction. The cell includes mid-edge node in the triangle edges. The ordering of the 12 points defining the cell is point ids (0-5,6-12) where point ids 0-5 are the six corner vertices of the wedge; followed by six midedge nodes (6-12). Note that these midedge nodes correspond lie on the edges defined by (0,1), (1,2), (2,0), (3,4), (4,5), (5,3). The Edges (0,3), (1,4), (2,5) dont have midedge nodes.

See Also
vtkQuadraticEdge vtkQuadraticTriangle vtkQuadraticTetra vtkQuadraticHexahedron vtkQuadraticQuad vtkQuadraticPyramid
Thanks:
Thanks to Soeren Gebbert who developed this class and integrated it into VTK 5.0.
Tests:
vtkQuadraticLinearWedge (Tests)

Definition at line 53 of file vtkQuadraticLinearWedge.h.

Member Typedef Documentation

Definition at line 57 of file vtkQuadraticLinearWedge.h.

Constructor & Destructor Documentation

vtkQuadraticLinearWedge::vtkQuadraticLinearWedge ( )
protected
vtkQuadraticLinearWedge::~vtkQuadraticLinearWedge ( )
protected

Member Function Documentation

static vtkQuadraticLinearWedge* vtkQuadraticLinearWedge::New ( )
static
virtual const char* vtkQuadraticLinearWedge::GetClassName ( )
virtual

Reimplemented from vtkNonLinearCell.

static int vtkQuadraticLinearWedge::IsTypeOf ( const char *  type)
static
virtual int vtkQuadraticLinearWedge::IsA ( const char *  type)
virtual

Reimplemented from vtkNonLinearCell.

static vtkQuadraticLinearWedge* vtkQuadraticLinearWedge::SafeDownCast ( vtkObject o)
static
void vtkQuadraticLinearWedge::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
int vtkQuadraticLinearWedge::GetCellType ( )
inline

Implement the vtkCell API. See the vtkCell API for descriptions of these methods.

Definition at line 63 of file vtkQuadraticLinearWedge.h.

int vtkQuadraticLinearWedge::GetCellDimension ( )
inline

Implement the vtkCell API. See the vtkCell API for descriptions of these methods.

Definition at line 64 of file vtkQuadraticLinearWedge.h.

int vtkQuadraticLinearWedge::GetNumberOfEdges ( )
inline

Implement the vtkCell API. See the vtkCell API for descriptions of these methods.

Definition at line 65 of file vtkQuadraticLinearWedge.h.

int vtkQuadraticLinearWedge::GetNumberOfFaces ( )
inline

Implement the vtkCell API. See the vtkCell API for descriptions of these methods.

Definition at line 66 of file vtkQuadraticLinearWedge.h.

vtkCell* vtkQuadraticLinearWedge::GetEdge ( int  edgeId)

Implement the vtkCell API. See the vtkCell API for descriptions of these methods.

vtkCell* vtkQuadraticLinearWedge::GetFace ( int  faceId)

Implement the vtkCell API. See the vtkCell API for descriptions of these methods.

int vtkQuadraticLinearWedge::CellBoundary ( int  subId,
double  pcoords[3],
vtkIdList pts 
)
void vtkQuadraticLinearWedge::Contour ( double  value,
vtkDataArray cellScalars,
vtkIncrementalPointLocator *  locator,
vtkCellArray verts,
vtkCellArray lines,
vtkCellArray polys,
vtkPointData inPd,
vtkPointData outPd,
vtkCellData *  inCd,
vtkIdType  cellId,
vtkCellData *  outCd 
)

The quadratic linear wege is splitted into 4 linear wedges, each of them is contoured by a provided scalar value

int vtkQuadraticLinearWedge::EvaluatePosition ( double  x[3],
double *  closestPoint,
int &  subId,
double  pcoords[3],
double &  dist2,
double *  weights 
)

The quadratic linear wege is splitted into 4 linear wedges, each of them is contoured by a provided scalar value

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

The quadratic linear wege is splitted into 4 linear wedges, each of them is contoured by a provided scalar value

int vtkQuadraticLinearWedge::Triangulate ( int  index,
vtkIdList ptIds,
vtkPoints pts 
)

The quadratic linear wege is splitted into 4 linear wedges, each of them is contoured by a provided scalar value

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

The quadratic linear wege is splitted into 4 linear wedges, each of them is contoured by a provided scalar value

virtual double* vtkQuadraticLinearWedge::GetParametricCoords ( )
virtual

The quadratic linear wege is splitted into 4 linear wedges, each of them is contoured by a provided scalar value

void vtkQuadraticLinearWedge::Clip ( double  value,
vtkDataArray cellScalars,
vtkIncrementalPointLocator *  locator,
vtkCellArray tetras,
vtkPointData inPd,
vtkPointData outPd,
vtkCellData *  inCd,
vtkIdType  cellId,
vtkCellData *  outCd,
int  insideOut 
)

Clip this quadratic linear wedge using scalar value provided. Like contouring, except that it cuts the hex to produce linear tetrahedron.

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

Line-edge intersection. Intersection has to occur within [0,1] parametric coordinates and with specified tolerance.

int vtkQuadraticLinearWedge::GetParametricCenter ( double  pcoords[3])
inline

Return the center of the quadratic linear wedge in parametric coordinates.

Definition at line 157 of file vtkQuadraticLinearWedge.h.

static void vtkQuadraticLinearWedge::InterpolationFunctions ( double  pcoords[3],
double  weights[15] 
)
static
static void vtkQuadraticLinearWedge::InterpolationDerivs ( double  pcoords[3],
double  derivs[45] 
)
static
virtual void vtkQuadraticLinearWedge::InterpolateFunctions ( double  pcoords[3],
double  weights[15] 
)
inlinevirtual
Deprecated:
Replaced by vtkQuadraticLinearWedge::InterpolateFunctions as of VTK 5.2

Definition at line 120 of file vtkQuadraticLinearWedge.h.

virtual void vtkQuadraticLinearWedge::InterpolateDerivs ( double  pcoords[3],
double  derivs[45] 
)
inlinevirtual
Deprecated:
Replaced by vtkQuadraticLinearWedge::InterpolateFunctions as of VTK 5.2

Definition at line 124 of file vtkQuadraticLinearWedge.h.

static int* vtkQuadraticLinearWedge::GetEdgeArray ( int  edgeId)
static
static int* vtkQuadraticLinearWedge::GetFaceArray ( int  faceId)
static
void vtkQuadraticLinearWedge::JacobianInverse ( double  pcoords[3],
double **  inverse,
double  derivs[45] 
)

Given parametric coordinates compute inverse Jacobian transformation matrix. Returns 9 elements of 3x3 inverse Jacobian plus interpolation function derivatives.

Member Data Documentation

vtkQuadraticEdge* vtkQuadraticLinearWedge::QuadEdge
protected

Definition at line 144 of file vtkQuadraticLinearWedge.h.

vtkLine* vtkQuadraticLinearWedge::Edge
protected

Definition at line 145 of file vtkQuadraticLinearWedge.h.

vtkQuadraticTriangle* vtkQuadraticLinearWedge::TriangleFace
protected

Definition at line 146 of file vtkQuadraticLinearWedge.h.

vtkQuadraticLinearQuad* vtkQuadraticLinearWedge::Face
protected

Definition at line 147 of file vtkQuadraticLinearWedge.h.

vtkWedge* vtkQuadraticLinearWedge::Wedge
protected

Definition at line 148 of file vtkQuadraticLinearWedge.h.

vtkDoubleArray* vtkQuadraticLinearWedge::Scalars
protected

Definition at line 149 of file vtkQuadraticLinearWedge.h.


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