(9)

Defines a 6-node pentahedral finite element using P1 interpolation in local coordinates (s.x,s.y) and Q1 isoparametric interpolation in local coordinates (s.x,s.z) and (s.y,s.z). More...

Inheritance diagram for Penta6:
FEShape

Public Member Functions

 Penta6 ()
 Default Constructor.
 
 Penta6 (const Element *element)
 Constructor when data of Element el are given. More...
 
 ~Penta6 ()
 Destructor.
 
void set (const Element *el)
 Choose element by giving its pointer.
 
void setLocal (const Point< real_t > &s)
 Initialize local point coordinates in element. More...
 
vector< Point< real_t > > DSh () const
 Return partial derivatives of shape functions of element nodes. More...
 
real_t getMaxEdgeLength () const
 Return Maximum length of pentahedron edges.
 
real_t getMinEdgeLength () const
 Return Mimimum length of pentahedron edges.
 
- Public Member Functions inherited from FEShape
 FEShape ()
 Default Constructor.
 
 FEShape (const Element *el)
 Constructor for an element. More...
 
 FEShape (const Side *sd)
 Constructor for a side. More...
 
virtual ~FEShape ()
 Destructor.
 
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...
 
real_t getDet () const
 Return determinant of jacobian. More...
 
Point< real_tgetCenter () const
 Return coordinates of center of element.
 
Point< real_tgetLocalPoint () const
 Localize a point in the element. More...
 
Point< real_tgetLocalPoint (const Point< real_t > &s) const
 Localize a point in the element. More...
 

Detailed Description

Defines a 6-node pentahedral finite element using P1 interpolation in local coordinates (s.x,s.y) and Q1 isoparametric interpolation in local coordinates (s.x,s.z) and (s.y,s.z).

The reference element is the cartesian product of the standard reference triangle with the line [-1,1]. The nodes are ordered as follows: Node 1 in reference element is at s=(1,0,0) Node 2 in reference element is at s=(0,1,0) Node 3 in reference element is at s=(0,0,0) Node 4 in reference element is at s=(1,0,1) Node 5 in reference element is at s=(0,1,1) Node 6 in reference element is at s=(0,0,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 setLocal() must be invoked.

Author
Rachid Touzani

Constructor & Destructor Documentation

◆ Penta6()

Penta6 ( const Element element)

Constructor when data of Element el are given.

Parameters
[in]elementPointer to Element

Member Function Documentation

◆ DSh()

vector<Point<real_t> > DSh ( ) const

Return partial derivatives of shape functions of element nodes.

Returns
LocalVect instance of partial derivatives of shape functions e.g. dsh(i).x, dsh(i).y, are partial derivatives of the i-th shape function.
Note
The local point at which the derivatives are computed must be chosen before by using the member function setLocal

◆ setLocal()

void setLocal ( const Point< real_t > &  s)

Initialize local point coordinates in element.

Parameters
[in]sPoint in the reference element This function computes jacobian, shape functions and their partial derivatives at s. Other member functions only return these values.