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

cell represents a biquadratic, 24-node isoparametric hexahedron More...

#include <vtkBiQuadraticQuadraticHexahedron.h>

Inheritance diagram for vtkBiQuadraticQuadraticHexahedron:
[legend]
Collaboration diagram for vtkBiQuadraticQuadraticHexahedron:
[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...
 
vtkBiQuadraticQuadraticHexahedronNewInstance () 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 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. More...
 
int EvaluatePosition (const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
 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 EvaluateLocation (int &subId, const double pcoords[3], double x[3], double *weights) override
 Determine global coordinate (x[3]) from subId and parametric coordinates. More...
 
int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts) override
 Generate simplices of proper dimension. More...
 
void Derivatives (int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
 Compute derivatives given cell subId and parametric coordinates. More...
 
double * GetParametricCoords () override
 Return a contiguous array of parametric coordinates of the points defining this 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 biquadratic hexahedron 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...
 
void JacobianInverse (const double pcoords[3], double **inverse, double derivs[72])
 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) override
 Implement the vtkCell API. More...
 
vtkCellGetFace (int) override
 Implement the vtkCell API. More...
 
void InterpolateFunctions (const double pcoords[3], double weights[24]) override
 Compute the interpolation functions/derivatives (aka shape functions/derivatives) More...
 
void InterpolateDerivs (const double pcoords[3], double derivs[72]) override
 Compute the interpolation functions/derivatives (aka shape functions/derivatives) More...
 
- Public Member Functions inherited from vtkNonLinearCell
virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class. More...
 
vtkNonLinearCellNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. More...
 
int IsLinear () override
 Non-linear cells require special treatment (tessellation) when converting to graphics primitives (during mapping). More...
 
- Public Member Functions inherited from vtkCell
virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class. More...
 
vtkCellNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. More...
 
void Initialize (int npts, const 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 GetCellType ()=0
 Return the type of cell. More...
 
virtual int GetCellDimension ()=0
 Return the topological dimensional of the cell (0,1,2, or 3). More...
 
virtual int IsLinear ()
 Non-linear cells require special treatment beyond the usual cell type and connectivity list information. 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 () const
 Return the number of points in the cell. More...
 
virtual int GetNumberOfEdges ()=0
 Return the number of edges in the cell. More...
 
virtual int GetNumberOfFaces ()=0
 Return the number of faces 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 vtkCellGetEdge (int edgeId)=0
 Return the edge cell from the edgeId of the cell. More...
 
virtual vtkCellGetFace (int faceId)=0
 Return the face cell from the faceId of the cell. More...
 
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. 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...
 
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. More...
 
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. More...
 
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. More...
 
virtual int Inflate (double dist)
 Inflates the cell. More...
 
virtual double ComputeBoundingSphere (double center[3]) const
 Computes the bounding sphere of the cell. More...
 
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. More...
 
virtual int Triangulate (int index, vtkIdList *ptIds, vtkPoints *pts)=0
 Generate simplices of proper dimension. More...
 
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. More...
 
void GetBounds (double bounds[6])
 Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax). More...
 
double * GetBounds ()
 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 int GetParametricCenter (double pcoords[3])
 Return center of the cell in parametric coordinates. 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 double * GetParametricCoords ()
 Return a contiguous array of parametric coordinates of the points defining this cell. 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))
 
virtual int IntersectWithCell (vtkCell *other, double tol=0.0)
 Intersects with an other cell. More...
 
virtual int IntersectWithCell (vtkCell *other, const vtkBoundingBox &boudingBox, const vtkBoundingBox &otherBoundingBox, double tol=0.0)
 Intersects with an other cell. More...
 
- 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 PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. 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 vtkTypeBool IsA (const char *name)
 Return 1 if this class is the same type of (or a subclass of) the named class. More...
 
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). 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...
 
bool GetIsInMemkind () const
 A local state flag that remembers whether this object lives in the normal or extended memory space. 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 vtkBiQuadraticQuadraticHexahedronNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkBiQuadraticQuadraticHexahedronSafeDownCast (vtkObjectBase *o)
 
static void InterpolationFunctions (const double pcoords[3], double weights[24])
 
static void InterpolationDerivs (const double pcoords[3], double derivs[72])
 
static const vtkIdTypeGetEdgeArray (vtkIdType edgeId)
 Return the ids of the vertices defining edge/face (edgeId/‘faceId’). More...
 
static const vtkIdTypeGetFaceArray (vtkIdType 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 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). More...
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 
static void SetMemkindDirectory (const char *directoryname)
 The name of a directory, ideally mounted -o dax, to memory map an extended memory space within. More...
 
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. More...
 

Protected Member Functions

virtual vtkObjectBaseNewInstanceInternal () const
 
 vtkBiQuadraticQuadraticHexahedron ()
 
 ~vtkBiQuadraticQuadraticHexahedron () override
 
void Subdivide (vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
 
- Protected Member Functions inherited from vtkNonLinearCell
virtual vtkObjectBaseNewInstanceInternal () const
 
 vtkNonLinearCell ()
 
 ~vtkNonLinearCell () override=default
 
- Protected Member Functions inherited from vtkCell
virtual vtkObjectBaseNewInstanceInternal () const
 
 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 RegisterInternal (vtkObjectBase *, vtkTypeBool check)
 
virtual void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check)
 
virtual void ReportReferences (vtkGarbageCollector *)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

vtkQuadraticEdgeEdge
 
vtkQuadraticQuadFace
 
vtkBiQuadraticQuadBiQuadFace
 
vtkHexahedronHex
 
vtkPointDataPointData
 
vtkCellDataCellData
 
vtkDoubleArrayCellScalars
 
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
std::atomic< int32_t > ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Additional Inherited Members

- Public Attributes inherited from vtkCell
vtkPointsPoints
 
vtkIdListPointIds
 
- Static Protected Member Functions inherited from vtkObjectBase
static vtkMallocingFunction GetCurrentMallocFunction ()
 
static vtkReallocingFunction GetCurrentReallocFunction ()
 
static vtkFreeingFunction GetCurrentFreeFunction ()
 
static vtkFreeingFunction GetAlternateFreeFunction ()
 

Detailed Description

cell represents a biquadratic, 24-node isoparametric hexahedron

vtkBiQuadraticQuadraticHexahedron is a concrete implementation of vtkNonLinearCell to represent a three-dimensional, 24-node isoparametric biquadratic hexahedron. The interpolation is the standard finite element, biquadratic-quadratic isoparametric shape function. The cell includes mid-edge and center-face nodes. The ordering of the 24 points defining the cell is point ids (0-7,8-19, 20-23) where point ids 0-7 are the eight corner vertices of the cube; followed by twelve midedge nodes (8-19), nodes 20-23 are the center-face nodes. Note that these midedge nodes correspond lie on the edges defined by (0,1), (1,2), (2,3), (3,0), (4,5), (5,6), (6,7), (7,4), (0,4), (1,5), (2,6), (3,7). The center face nodes laying in quad 22-(0,1,5,4), 21-(1,2,6,5), 23-(2,3,7,6) and 22-(3,0,4,7)

*
* top
*  7--14--6
*  |      |
* 15      13
*  |      |
*  4--12--5
*
*  middle
* 19--23--18
*  |      |
* 20      21
*  |      |
* 16--22--17
*
* bottom
*  3--10--2
*  |      |
* 11      9
*  |      |
*  0-- 8--1
*
* 
See also
vtkQuadraticEdge vtkQuadraticTriangle vtkQuadraticTetra vtkQuadraticQuad vtkQuadraticPyramid vtkQuadraticWedge
Thanks:
Thanks to Soeren Gebbert who developed this class and integrated it into VTK 5.0.
Online Examples:

Definition at line 99 of file vtkBiQuadraticQuadraticHexahedron.h.

Member Typedef Documentation

◆ Superclass

Definition at line 103 of file vtkBiQuadraticQuadraticHexahedron.h.

Constructor & Destructor Documentation

◆ vtkBiQuadraticQuadraticHexahedron()

vtkBiQuadraticQuadraticHexahedron::vtkBiQuadraticQuadraticHexahedron ( )
protected

◆ ~vtkBiQuadraticQuadraticHexahedron()

vtkBiQuadraticQuadraticHexahedron::~vtkBiQuadraticQuadraticHexahedron ( )
overrideprotected

Member Function Documentation

◆ New()

static vtkBiQuadraticQuadraticHexahedron * vtkBiQuadraticQuadraticHexahedron::New ( )
static

◆ IsTypeOf()

static vtkTypeBool vtkBiQuadraticQuadraticHexahedron::IsTypeOf ( const char *  type)
static

◆ IsA()

virtual vtkTypeBool vtkBiQuadraticQuadraticHexahedron::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.

◆ SafeDownCast()

static vtkBiQuadraticQuadraticHexahedron * vtkBiQuadraticQuadraticHexahedron::SafeDownCast ( vtkObjectBase o)
static

◆ NewInstanceInternal()

virtual vtkObjectBase * vtkBiQuadraticQuadraticHexahedron::NewInstanceInternal ( ) const
protectedvirtual

Reimplemented from vtkNonLinearCell.

◆ NewInstance()

vtkBiQuadraticQuadraticHexahedron * vtkBiQuadraticQuadraticHexahedron::NewInstance ( ) const

◆ PrintSelf()

void vtkBiQuadraticQuadraticHexahedron::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 vtkObject.

◆ GetCellType()

int vtkBiQuadraticQuadraticHexahedron::GetCellType ( )
inlineoverridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 111 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ GetCellDimension()

int vtkBiQuadraticQuadraticHexahedron::GetCellDimension ( )
inlineoverridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 112 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ GetNumberOfEdges()

int vtkBiQuadraticQuadraticHexahedron::GetNumberOfEdges ( )
inlineoverridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 113 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ GetNumberOfFaces()

int vtkBiQuadraticQuadraticHexahedron::GetNumberOfFaces ( )
inlineoverridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

Definition at line 114 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ GetEdge()

vtkCell * vtkBiQuadraticQuadraticHexahedron::GetEdge ( int  )
overridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

◆ GetFace()

vtkCell * vtkBiQuadraticQuadraticHexahedron::GetFace ( int  )
overridevirtual

Implement the vtkCell API.

See the vtkCell API for descriptions of these methods.

Implements vtkCell.

◆ CellBoundary()

int vtkBiQuadraticQuadraticHexahedron::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.

◆ Contour()

void vtkBiQuadraticQuadraticHexahedron::Contour ( double  value,
vtkDataArray cellScalars,
vtkIncrementalPointLocator locator,
vtkCellArray verts,
vtkCellArray lines,
vtkCellArray polys,
vtkPointData inPd,
vtkPointData outPd,
vtkCellData inCd,
vtkIdType  cellId,
vtkCellData outCd 
)
overridevirtual

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 nullptr, then no interpolation is performed. Also, if the output cell data is non-nullptr, 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.)

Implements vtkCell.

◆ EvaluatePosition()

int vtkBiQuadraticQuadraticHexahedron::EvaluatePosition ( const double  x[3],
double  closestPoint[3],
int &  subId,
double  pcoords[3],
double &  dist2,
double  weights[] 
)
overridevirtual

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.

Implements vtkCell.

◆ EvaluateLocation()

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

Determine global coordinate (x[3]) from subId and parametric coordinates.

Also returns interpolation weights. (The number of weights is equal to the number of points in the cell.)

Implements vtkCell.

◆ Triangulate()

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

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.

Implements vtkCell.

◆ Derivatives()

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

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

Implements vtkCell.

◆ GetParametricCoords()

double * vtkBiQuadraticQuadraticHexahedron::GetParametricCoords ( )
overridevirtual

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-nullptr 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 from vtkCell.

◆ Clip()

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

Clip this biquadratic hexahedron using scalar value provided.

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

Implements vtkCell.

◆ IntersectWithLine()

int vtkBiQuadraticQuadraticHexahedron::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.

◆ InterpolationFunctions()

static void vtkBiQuadraticQuadraticHexahedron::InterpolationFunctions ( const double  pcoords[3],
double  weights[24] 
)
static

◆ InterpolationDerivs()

static void vtkBiQuadraticQuadraticHexahedron::InterpolationDerivs ( const double  pcoords[3],
double  derivs[72] 
)
static

◆ InterpolateFunctions()

void vtkBiQuadraticQuadraticHexahedron::InterpolateFunctions ( const double  pcoords[3],
double  weights[24] 
)
inlineoverride

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

Definition at line 154 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ InterpolateDerivs()

void vtkBiQuadraticQuadraticHexahedron::InterpolateDerivs ( const double  pcoords[3],
double  derivs[72] 
)
inlineoverride

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

Definition at line 158 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ GetEdgeArray()

static const vtkIdType * vtkBiQuadraticQuadraticHexahedron::GetEdgeArray ( vtkIdType  edgeId)
static

Return the ids of the vertices defining edge/face (edgeId/‘faceId’).

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

Note
The return type changed. It used to be int*, it is now const vtkIdType*. This is so ids are unified between vtkCell and vtkPoints.

◆ GetFaceArray()

static const vtkIdType * vtkBiQuadraticQuadraticHexahedron::GetFaceArray ( vtkIdType  faceId)
static

Return the ids of the vertices defining edge/face (edgeId/‘faceId’).

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

Note
The return type changed. It used to be int*, it is now const vtkIdType*. This is so ids are unified between vtkCell and vtkPoints.

◆ JacobianInverse()

void vtkBiQuadraticQuadraticHexahedron::JacobianInverse ( const double  pcoords[3],
double **  inverse,
double  derivs[72] 
)

Given parametric coordinates compute inverse Jacobian transformation matrix.

Returns 9 elements of 3x3 inverse Jacobian plus interpolation function derivatives.

◆ Subdivide()

void vtkBiQuadraticQuadraticHexahedron::Subdivide ( vtkPointData inPd,
vtkCellData inCd,
vtkIdType  cellId,
vtkDataArray cellScalars 
)
protected

Member Data Documentation

◆ Edge

vtkQuadraticEdge* vtkBiQuadraticQuadraticHexahedron::Edge
protected

Definition at line 186 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ Face

vtkQuadraticQuad* vtkBiQuadraticQuadraticHexahedron::Face
protected

Definition at line 187 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ BiQuadFace

vtkBiQuadraticQuad* vtkBiQuadraticQuadraticHexahedron::BiQuadFace
protected

Definition at line 188 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ Hex

vtkHexahedron* vtkBiQuadraticQuadraticHexahedron::Hex
protected

Definition at line 189 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ PointData

vtkPointData* vtkBiQuadraticQuadraticHexahedron::PointData
protected

Definition at line 190 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ CellData

vtkCellData* vtkBiQuadraticQuadraticHexahedron::CellData
protected

Definition at line 191 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ CellScalars

vtkDoubleArray* vtkBiQuadraticQuadraticHexahedron::CellScalars
protected

Definition at line 192 of file vtkBiQuadraticQuadraticHexahedron.h.

◆ Scalars

vtkDoubleArray* vtkBiQuadraticQuadraticHexahedron::Scalars
protected

Definition at line 193 of file vtkBiQuadraticQuadraticHexahedron.h.


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