Quad4 Class Reference
[Shape Function]
#include <Quad4.h>

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.
Public Member Functions | |
| double | check () const |
| Check element jacobian and number of nodes. | |
| Point< double > | DSh (size_t i) const |
| Return derivatives of shape function of node i at a given point. | |
| Point< double > | getCenter () const |
| Return coordinates of center of element. | |
| double | getDet () const |
| Return determinant of jacobian. | |
| Point< double > | getLocalPoint (const Point< double > &s) const |
| Localize a point in the element. | |
| Point< double > | getLocalPoint () const |
| Localize a point in the element. | |
| double | getMaxEdgeLength () const |
| Return maximal edge length of quadrilateral. | |
| double | getMinEdgeLength () const |
| Return minimal edge length of quadrilateral. | |
| Quad4 (const Side *side) | |
| Constructor when data of Side sd are given. | |
| Quad4 (const Element *element) | |
| Constructor when data of Element el are given. | |
| Quad4 () | |
| Default Constructor. | |
| void | setLocal (const Point< double > &s) |
| Initialize local point coordinates in element. | |
| double | Sh (size_t i, Point< double > s) const |
| Calculate shape function of node i at a given point s. | |
| double | Sh (size_t i) const |
| Return shape function of node i at given point. | |
| ~Quad4 () | |
| Destructor. | |
Member Function Documentation
| double check | ( | ) | const |
Check element jacobian and number of nodes.
- Return values:
-
m - > 0 : m is the length
- = 0 : zero length (=> Error)
| Point<double> DSh | ( | size_t | i | ) | const |
Return derivatives of shape function of node i at a given point.
Member function setLocal() must have been called before in order to calcuate relevant quantities.
Reimplemented from FEShape.
| double getDet | ( | ) | const |
Return determinant of jacobian.
Member function setLocal() must have been called before in order to calcuate relevant quantities.
Reimplemented from FEShape.
Localize a point in the element.
Return actual coordinates where s are coordinates in the reference element.
References FEShape::Sh().
| Point<double> getLocalPoint | ( | ) | const [inherited] |
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.
Reimplemented in Line3.
| void setLocal | ( | const Point< double > & | s | ) |
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.
| double Sh | ( | size_t | i, | |
| Point< double > | s | |||
| ) | const [inherited] |