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 Types inherited from vtkCell
typedef vtkObject Superclass
 

Public Member Functions

virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class. More...
 
vtkQuadraticLinearWedgeNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. More...
 
int CellBoundary (int subId, const double pcoords[3], vtkIdList *pts) override
 Given parametric coordinates of a point, return the closest cell boundary, and whether the point is inside or outside of the cell. More...
 
void Clip (double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tetras, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
 Clip this quadratic linear wedge using scalar value provided. More...
 
int IntersectWithLine (const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
 Line-edge intersection. More...
 
int GetParametricCenter (double pcoords[3]) override
 Return the center of the quadratic linear wedge in parametric coordinates. More...
 
void JacobianInverse (const double pcoords[3], double **inverse, double derivs[45])
 Given parametric coordinates compute inverse Jacobian transformation matrix. More...
 
int GetCellType () override
 Implement the vtkCell API. More...
 
int GetCellDimension () override
 Implement the vtkCell API. More...
 
int GetNumberOfEdges () override
 Implement the vtkCell API. More...
 
int GetNumberOfFaces () override
 Implement the vtkCell API. More...
 
vtkCellGetEdge (int edgeId) override
 Implement the vtkCell API. More...
 
vtkCellGetFace (int faceId) override
 Implement the vtkCell API. More...
 
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
 The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value. More...
 
int EvaluatePosition (const double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights) override
 The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value. More...
 
void EvaluateLocation (int &subId, const double pcoords[3], double x[3], double *weights) override
 The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value. More...
 
int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts) override
 The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value. More...
 
void Derivatives (int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
 The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value. More...
 
doubleGetParametricCoords () override
 The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value. More...
 
void InterpolateFunctions (const double pcoords[3], double weights[15]) override
 Compute the interpolation functions/derivatives (aka shape functions/derivatives) More...
 
void InterpolateDerivs (const double pcoords[3], double derivs[45]) override
 Compute the interpolation functions/derivatives (aka shape functions/derivatives) More...
 
- Public Member Functions inherited from vtkNonLinearCell
vtkNonLinearCellNewInstance () const
 
int IsLinear () override
 Non-linear cells require special treatment (tessellation) when converting to graphics primitives (during mapping). More...
 
- Public Member Functions inherited from vtkCell
vtkCellNewInstance () const
 
void Initialize (int npts, 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 RequiresInitialization ()
 Some cells require initialization prior to access. More...
 
virtual void Initialize ()
 
virtual int IsExplicitCell ()
 Explicit cells require additional representational information beyond the usual cell type and connectivity list information. More...
 
virtual int RequiresExplicitFaceRepresentation ()
 Determine whether the cell requires explicit face representation, and methods for setting and getting the faces (see vtkPolyhedron for example usage of these methods). More...
 
virtual void SetFaces (vtkIdType *vtkNotUsed(faces))
 
virtual vtkIdTypeGetFaces ()
 
vtkPointsGetPoints ()
 Get the point coordinates for the cell. More...
 
vtkIdType GetNumberOfPoints ()
 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 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. More...
 
void GetBounds (double bounds[6])
 Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). More...
 
doubleGetBounds ()
 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 int IsPrimaryCell ()
 Return whether this cell type has a fixed topology or whether the topology varies depending on the data (e.g., vtkConvexPointSet). 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))
 
- 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)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
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...
 
int InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
int InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string. 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...
 
virtual 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...
 
void PrintRevisions (ostream &)
 Legacy. 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...
 

Static Public Member Functions

static vtkQuadraticLinearWedgeNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkQuadraticLinearWedgeSafeDownCast (vtkObjectBase *o)
 
static void InterpolationFunctions (const double pcoords[3], double weights[15])
 
static void InterpolationDerivs (const double pcoords[3], double derivs[45])
 
static intGetEdgeArray (int edgeId)
 Return the ids of the vertices defining edge/face (edgeId/`faceId'). More...
 
static intGetFaceArray (int faceId)
 Return the ids of the vertices defining edge/face (edgeId/`faceId'). More...
 
- Static Public Member Functions inherited from vtkNonLinearCell
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkNonLinearCellSafeDownCast (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 (int 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 int 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 vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 

Protected Member Functions

virtual vtkObjectBaseNewInstanceInternal () const
 
 vtkQuadraticLinearWedge ()
 
 ~vtkQuadraticLinearWedge () override
 
- Protected Member Functions inherited from vtkNonLinearCell
 vtkNonLinearCell ()
 
 ~vtkNonLinearCell () 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 CollectRevisions (ostream &)
 
virtual void ReportReferences (vtkGarbageCollector *)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

vtkQuadraticEdgeQuadEdge
 
vtkLineEdge
 
vtkQuadraticTriangleTriangleFace
 
vtkQuadraticLinearQuadFace
 
vtkWedgeWedge
 
vtkDoubleArrayScalars
 
- Protected Attributes inherited from vtkCell
double Bounds [6]
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
- Protected Attributes inherited from vtkObjectBase
vtkAtomicInt32 ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Additional Inherited Members

- Public Attributes inherited from vtkCell
vtkPointsPoints
 
vtkIdListPointIds
 

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) don't 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 55 of file vtkQuadraticLinearWedge.h.

Member Typedef Documentation

Definition at line 59 of file vtkQuadraticLinearWedge.h.

Constructor & Destructor Documentation

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

Member Function Documentation

static vtkQuadraticLinearWedge* vtkQuadraticLinearWedge::New ( )
static
static vtkTypeBool vtkQuadraticLinearWedge::IsTypeOf ( const char *  type)
static
virtual vtkTypeBool vtkQuadraticLinearWedge::IsA ( const char *  name)
virtual

Return 1 if this class is the same type of (or a subclass of) the named class.

Returns 0 otherwise. This method works in combination with vtkTypeMacro found in vtkSetGet.h.

Reimplemented from vtkNonLinearCell.

static vtkQuadraticLinearWedge* vtkQuadraticLinearWedge::SafeDownCast ( vtkObjectBase o)
static
virtual vtkObjectBase* vtkQuadraticLinearWedge::NewInstanceInternal ( ) const
protectedvirtual

Reimplemented from vtkNonLinearCell.

vtkQuadraticLinearWedge* vtkQuadraticLinearWedge::NewInstance ( ) const
void vtkQuadraticLinearWedge::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
overridevirtual

Methods invoked by print to print information about the object including superclasses.

Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from vtkNonLinearCell.

int vtkQuadraticLinearWedge::GetCellType ( )
inlineoverridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 67 of file vtkQuadraticLinearWedge.h.

int vtkQuadraticLinearWedge::GetCellDimension ( )
inlineoverridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 68 of file vtkQuadraticLinearWedge.h.

int vtkQuadraticLinearWedge::GetNumberOfEdges ( )
inlineoverridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 69 of file vtkQuadraticLinearWedge.h.

int vtkQuadraticLinearWedge::GetNumberOfFaces ( )
inlineoverridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 70 of file vtkQuadraticLinearWedge.h.

vtkCell* vtkQuadraticLinearWedge::GetEdge ( int  edgeId)
overridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

vtkCell* vtkQuadraticLinearWedge::GetFace ( int  faceId)
overridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

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

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.

Implements vtkCell.

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 
)
overridevirtual

The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value.

Implements vtkCell.

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

The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value.

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

The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value.

Implements vtkCell.

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

The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value.

Implements vtkCell.

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

The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value.

Implements vtkCell.

double* vtkQuadraticLinearWedge::GetParametricCoords ( )
overridevirtual

The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided scalar value.

Reimplemented from vtkCell.

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

Clip this quadratic linear wedge using scalar value provided.

Like contouring, except that it cuts the hex to produce linear tetrahedron.

Implements vtkCell.

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

Line-edge intersection.

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

Implements vtkCell.

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

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

Reimplemented from vtkCell.

Definition at line 175 of file vtkQuadraticLinearWedge.h.

static void vtkQuadraticLinearWedge::InterpolationFunctions ( const double  pcoords[3],
double  weights[15] 
)
static
static void vtkQuadraticLinearWedge::InterpolationDerivs ( const double  pcoords[3],
double  derivs[45] 
)
static
void vtkQuadraticLinearWedge::InterpolateFunctions ( const double  pcoords[3],
double  weights[15] 
)
inlineoverride

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

Definition at line 133 of file vtkQuadraticLinearWedge.h.

void vtkQuadraticLinearWedge::InterpolateDerivs ( const double  pcoords[3],
double  derivs[45] 
)
inlineoverride

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

Definition at line 137 of file vtkQuadraticLinearWedge.h.

static int* vtkQuadraticLinearWedge::GetEdgeArray ( int  edgeId)
static

Return the ids of the vertices defining edge/face (edgeId/`faceId').

Ids are related to the cell, not to the dataset.

static int* vtkQuadraticLinearWedge::GetFaceArray ( int  faceId)
static

Return the ids of the vertices defining edge/face (edgeId/`faceId').

Ids are related to the cell, not to the dataset.

void vtkQuadraticLinearWedge::JacobianInverse ( const 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 162 of file vtkQuadraticLinearWedge.h.

vtkLine* vtkQuadraticLinearWedge::Edge
protected

Definition at line 163 of file vtkQuadraticLinearWedge.h.

vtkQuadraticTriangle* vtkQuadraticLinearWedge::TriangleFace
protected

Definition at line 164 of file vtkQuadraticLinearWedge.h.

vtkQuadraticLinearQuad* vtkQuadraticLinearWedge::Face
protected

Definition at line 165 of file vtkQuadraticLinearWedge.h.

vtkWedge* vtkQuadraticLinearWedge::Wedge
protected

Definition at line 166 of file vtkQuadraticLinearWedge.h.

vtkDoubleArray* vtkQuadraticLinearWedge::Scalars
protected

Definition at line 167 of file vtkQuadraticLinearWedge.h.


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