Jive reference manual
List of all members | Static Public Member Functions
jive::geom::QuadraticTetrahedron Class Reference

Constructs shape objects that encapsulate the geometrical properties of quadratic tetrahedron elements. More...

Inheritance diagram for jive::geom::QuadraticTetrahedron:
Inheritance graph

Static Public Member Functions

static jem::Ref< InternalShapegetShape (int ipCount=11, int bndIpCount=7)
 Returns a shape with a specified number of integration points. More...
 
static jem::Ref< InternalShapegetShape (const Matrix &intPoints, const Vector &intWeights, const Matrix &bndIntPoints, const Vector &bndIntWeights)
 Returns a shape with a specific set of integration points. More...
 
static jem::Ref< InternalShapegetShape (const Matrix &intPoints, const Vector &intWeights, const ShapeBoundary &boundary)
 Returns a shape with a specific set of boundaries. More...
 

Additional Inherited Members

- Public Types inherited from jive::geom::TypeDefs
typedef jem::util::Properties Properties
 A type representing a set of properties. More...
 
typedef util::Vector Vector
 A vector type. More...
 
typedef util::Matrix Matrix
 A matrix type. More...
 

Detailed Description

#include <jive/geom/Tetrahedron.h>

The QuadraticTetrahedron class provides functions for constructing InternalShape objects that encapsulate the geometrical properties of tetrahedron elements with ten nodes and quadratic, isoparametric shape functions. The following figure shows how the nodes are numbered (left image) and where they are located in the local coordinate system (right image).

QuadraticTetrahedron.jpg

Note that the local node coordinates of a quadratic tetrahedron are identical to the node coordinates of a standard quadratic tetrahedron.

A quadratic tetrahedron has four boundaries. The table below indicates which boundaries are attached to which nodes.

Boundary index

Node indices

0 0, 5, 4, 3, 2, 1
1 0, 1, 2, 7, 9, 6
2 2, 3, 4, 8, 9, 7
3 0, 6, 9, 8, 4, 5

The normals on the boundaries point away from the interior of the tetrahedron. The boundaries are normally represented by BoundaryShape objects that are constructed by the QuadraticBoundaryTriangle class. However, you can also specify your own set of BoundaryShape objects to be used.

Note that the QuadraticTetrahedron class has only static members; it should be used as a mini name space.

See also
LinearTetrahedron

Member Function Documentation

static jem::Ref<InternalShape> jive::geom::QuadraticTetrahedron::getShape ( int  ipCount = 11,
int  bndIpCount = 7 
)
static

Returns an InternalShape object that encapsulates the geometrical properties of a quadratic tetrahedron element. This function is similar to LinearTetrahedron::getShape().

Parameters
ipCount- the number of integration points within the element.
bndIpCount- the number of integration points on each boundary.
Returns
An InternalShape object representing a quadratic tetrahedron element.
static jem::Ref<InternalShape> jive::geom::QuadraticTetrahedron::getShape ( const Matrix intPoints,
const Vector intWeights,
const Matrix bndIntPoints,
const Vector bndIntWeights 
)
static

Returns an InternalShape object that represents a quadratic tetrahedron element with a specific set of integration points. This function is similar to LinearTetrahedron::getShape().

Parameters
intPoints- a matrix containing the local coordinates of the integration points within the element.
intWeights- a vector containing the weights of the internal integration points in the local coordinate system.
bndIntPoints- a matrix containing the local coordinates of the integration points on each boundary.
bndIntWeights- a vector containing the weights of the boundary integration points in the local coordinate system of the boundaries.
Precondition
intPoints.size(0) == 3 &&
intPoints.size(1) == intWeights.size() &&
bndIntPoints.size(0) == 2 &&
bndIntPoints.size(1) == bndIntWeights.size()
Returns
An InternalShape object that encapsulates the geometrical properties of a quadratic tetrahedron element.
static jem::Ref<InternalShape> jive::geom::QuadraticTetrahedron::getShape ( const Matrix intPoints,
const Vector intWeights,
const ShapeBoundary boundary 
)
static

Returns an InternalShape object that represents a quadratic tetrahedron element with a specific set of integration points and boundaries. This function is similar to LinearTetrahedron::getShape().

Parameters
intPoints- a matrix containing the local coordinates of the integration points within the element.
intWeights- a vector containing the weights of the internal integration points in the local coordinate system.
boundary- the boundaries of the shape object that is returned by this function.
Precondition
intPoints.size(0) == 3 &&
intPoints.size(1) == intWeights.size()
Returns
An InternalShape object that encapsulates the geometrical properties of a quadratic tetrahedron element.