Defines a three-dimensional 8-node hexahedral finite element using Q1-isoparametric interpolation. More...


Public Member Functions | |
Hexa8 () | |
Default Constructor. | |
Hexa8 (const Element *el) | |
Constructor when data of Element el are given. | |
~Hexa8 () | |
Destructor. | |
void | setLocal (const Point< real_t > &s) |
Initialize local point coordinates in element. More... | |
Point< real_t > | DSh (size_t i) |
Return x , y and z partial derivatives of shape function of node i at a given point. More... | |
void | atGauss1 (LocalVect< Point< real_t >, 8 > &dsh, real_t &w) |
Calculate shape function derivatives and integration weights for 1-point Gauss rule. | |
void | atGauss2 (LocalMatrix< Point< real_t >, 8, 8 > &dsh, LocalVect< real_t, 8 > &w) |
Calculate shape function derivatives and integration weights for 2x2x2-point Gauss rule. | |
real_t | getMaxEdgeLength () const |
Return maximal edge length. | |
real_t | getMinEdgeLength () const |
Return minimal edge length. | |
real_t | Sh (size_t i) const |
Return shape function of node i at given point. | |
real_t | Sh (size_t i, Point< real_t > s) const |
Calculate shape function of node i at a given point s . More... | |
Point< real_t > | DSh (size_t i) const |
Return derivatives of shape function of node i at a given point. More... | |
real_t | getDet () const |
Return determinant of jacobian. More... | |
Point< real_t > | getCenter () const |
Return coordinates of center of element. | |
Point< real_t > | getLocalPoint () const |
Localize a point in the element. More... | |
Point< real_t > | getLocalPoint (const Point< real_t > &s) const |
Localize a point in the element. More... | |
Detailed Description
Defines a three-dimensional 8-node hexahedral finite element using Q1-isoparametric interpolation.
The reference element is the cube [-1,1]*[-1,1]*[-1,1]
. The user must take care to the fact that determinant of jacobian and other quantities depend on the point in the reference element where they are calculated. For this, before any utilization of shape functions or jacobian, function getLocal(s) must be invoked.
Member Function Documentation
Initialize local point coordinates in element.
- Parameters
-
[in] s Point in the reference element This function computes jacobian, shape functions and their partial derivatives at s
. Other member functions only return these values.
Return x
, y
and z
partial derivatives of shape function of node i
at a given point.
Member function setLocal must have been called before in order to calculate relevant quantities.
Calculate shape function of node i
at a given point s
.
- Parameters
-
[in] i Local node label [in] s Point in the reference triangle where the shape function is evaluated
Return derivatives of shape function of node i
at a given point.
If the transformation (Reference element -> Actual element) is not affine, member function setLocal()
must have been called before in order to calcuate relevant quantities.
- Parameters
-
[in] i Partial derivative index (1, 2 or 3)
|
inherited |
Return determinant of jacobian.
If the transformation (Reference element -> Actual element) is not affine, member function setLocal() must have been called before in order to calcuate relevant quantities.
Localize a point in the element.
Return actual coordinates in the reference element. If the transformation (Reference element -> Actual element) is not affine, member function setLocal() must have been called before in order to calcuate relevant quantities.