Defines a 4-node quadrilateral finite element using Q1
isoparametric interpolation.
More...


Public Member Functions | |
Quad4 () | |
Default Constructor. | |
Quad4 (const Element *element) | |
Constructor when data of Element el are given. More... | |
Quad4 (const Side *side) | |
Constructor when data of Side sd are given. More... | |
~Quad4 () | |
Destructor. | |
void | set (const Element *el) |
Choose element by giving its pointer. | |
void | set (const Side *sd) |
Choose side by giving its pointer. | |
void | setLocal (const Point< real_t > &s) |
Initialize local point coordinates in element. More... | |
Point< real_t > | DSh (size_t i) const |
Return derivatives of shape function of node i at a given point. More... | |
Point< real_t > | Grad (const LocalVect< real_t, 4 > &u, const Point< real_t > &s) |
Return gradient of a function defined at element nodes. More... | |
real_t | getMaxEdgeLength () const |
Return maximal edge length of quadrilateral. | |
real_t | getMinEdgeLength () const |
Return minimal edge length of quadrilateral. | |
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_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 4-node quadrilateral finite element using Q1
isoparametric interpolation.
The reference element is the square [-1,1]x[-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 setLocal() must be invoked.
Constructor & Destructor Documentation
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 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.
Return gradient of a function defined at element nodes.
- Parameters
-
[in] u Vector of values at nodes [in] s Local coordinates (in [-1,1]
) of point where the gradient is evaluated
- Returns
- Value of gradient
- Note
- If the derivatives of shape functions were not computed before calling this function (by calling setLocal), this function will compute them
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
|
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.