Defines a three-dimensional 4-node tetrahedral finite element using P1
interpolation.
More...


Public Member Functions | |
Tetra4 () | |
Default Constructor. | |
Tetra4 (const Element *el) | |
Constructor when data of Element el are given. | |
~Tetra4 () | |
Destructor. | |
void | set (const Element *el) |
Choose element by giving its pointer. | |
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 x , y and z partial derivatives of shape function associated to node i . More... | |
real_t | getVolume () const |
Return volume of element. | |
Point< real_t > | getRefCoord (const Point< real_t > &x) const |
Return reference coordinates of a point x in element. | |
bool | isIn (const Point< real_t > &x) |
Check whether point x is in current tetrahedron or not. | |
real_t | getInterpolate (const Point< real_t > &x, const LocalVect< real_t, 4 > &v) |
Return interpolated value at point of coordinate x | |
Point< real_t > | EdgeSh (size_t k, Point< real_t > s) |
Return edge shape function. More... | |
Point< real_t > | CurlEdgeSh (size_t k) |
Return curl of edge shape function. More... | |
real_t | getMaxEdgeLength () const |
Return maximal edge length of tetrahedron. | |
real_t | getMinEdgeLength () const |
Return minimal edge length of tetrahedron. | |
real_t | Sh (size_t i) const |
Return shape function of node i at given point. | |
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 4-node tetrahedral finite element using P1
interpolation.
The reference element is the right tetrahedron with four unit edges interpolation.
Member Function Documentation
Calculate shape function of node i
at a given point s
.
s
is a point in the reference tetrahedron.
Return x
, y
and z
partial derivatives of shape function associated to node i
.
Note that these are constant in element.
|
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.