LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Quadratic Lagrangian finite element on a triangular reference element. More...
#include <lf/uscalfe/uscalfe.h>
Public Member Functions | |
FeLagrangeO2Tria (const FeLagrangeO2Tria &)=default | |
FeLagrangeO2Tria (FeLagrangeO2Tria &&) noexcept=default | |
FeLagrangeO2Tria & | operator= (const FeLagrangeO2Tria &)=default |
FeLagrangeO2Tria & | operator= (FeLagrangeO2Tria &&) noexcept=default |
FeLagrangeO2Tria ()=default | |
~FeLagrangeO2Tria () override=default | |
base::RefEl | RefEl () const override |
Tells the type of reference cell underlying the parametric finite element. | |
unsigned | Degree () const override |
Quadratic Lagrangian finite elements rely on polynomials of degree 2. | |
size_type | NumRefShapeFunctions () const override |
Six local shape functions are associated with a triangular cell in the case of quadratic Lagrangian finite elements. | |
size_type | NumRefShapeFunctions (dim_t codim) const override |
One shape function attached to each node and one to each edge of the triangle. There are no interior shape functions. | |
size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t) const override |
One shape function attached to each node and one to each edge of the triangle. There are no interior shape functions. | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | EvalReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override |
Point evaluation of reference shape functions. | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | GradientsReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override |
Point evaluations of gradients of reference shape functions. | |
Eigen::MatrixXd | EvaluationNodes () const override |
Evaluation nodes are the three vertices of the triangle and the three midpoints of its edges. | |
size_type | NumEvaluationNodes () const override |
Six evaluation nodes, the same number as local shape functions. | |
![]() | |
virtual | ~ScalarReferenceFiniteElement ()=default |
dim_t | Dimension () const |
Returns the spatial dimension of the reference cell. | |
virtual Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > | NodalValuesToDofs (const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > &nodvals) const |
Computes the linear combination of reference shape functions matching function values at evaluation nodes. | |
Additional Inherited Members | |
![]() | |
using | Scalar = SCALAR |
The underlying scalar type. | |
![]() | |
ScalarReferenceFiniteElement ()=default | |
ScalarReferenceFiniteElement (const ScalarReferenceFiniteElement &)=default | |
ScalarReferenceFiniteElement (ScalarReferenceFiniteElement &&) noexcept=default | |
ScalarReferenceFiniteElement & | operator= (const ScalarReferenceFiniteElement &)=default |
ScalarReferenceFiniteElement & | operator= (ScalarReferenceFiniteElement &&) noexcept=default |
![]() | |
template<class SCALAR > | |
void | PrintInfo (std::ostream &o, const ScalarReferenceFiniteElement< SCALAR > &srfe, unsigned int ctrl=0) |
Print information about a ScalarReferenceFiniteElement to the given stream. | |
template<typename SCALAR > | |
std::ostream & | operator<< (std::ostream &o, const ScalarReferenceFiniteElement< SCALAR > &fe_desc) |
Stream output operator: just calls the ScalarReferenceFiniteElement::print() method. | |
Quadratic Lagrangian finite element on a triangular reference element.
This is a specialization of fe::ScalarReferenceFiniteElement.
The reference shape functions are combinations of the barycentric coordinate functions on the reference triangle, see Example 2.6.1.2.
The first three shape functions are associated to the vertices, the other three shape functions to the edges of the triangle. There are no interior shape functions, see Equation 2.6.1.6:
\begin{eqnarray} \hat{b}^1({\bf x}) &=& 2(1-x_1-x_2)(\frac12 -x_1 -x_2)\;,\\ \hat{b}^2({\bf x}) &=& 2x_1(x_1-\frac12)\;,\\ \hat{b}^3({\bf x}) &=& 2x_2(x_2-\frac12)\;,\\ \hat{b}^4({\bf x}) &=& 4(1-x_1-x_2)x_1\;,\\ \hat{b}^5({\bf x}) &=& 4x_1x_2\;,\\ \hat{b}^6({\bf x}) &=& 4(1-x_1-x_2)x_2\;. \end{eqnarray}
The numbering convention for the local shape functions is explained in Example 2.7.4.12.
|
default |
|
defaultnoexcept |
|
default |
|
overridedefault |
|
inlineoverridevirtual |
Quadratic Lagrangian finite elements rely on polynomials of degree 2.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Point evaluation of reference shape functions.
The formulas of the reference shape functions are given in the class documentation of FeLagrangeO2Tria.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Evaluation nodes are the three vertices of the triangle and the three midpoints of its edges.
The reference shape functions as implemented by this class are a cardinal basis of the space of quadratic two-variate polynomials with respect to these evaluation nodes.
The numbering of evluation nodes is fixed by the numbering of the local shape functions.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Point evaluations of gradients of reference shape functions.
Refer to Example 2.7.5.7 for related explanations
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Six evaluation nodes, the same number as local shape functions.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 543 of file lagr_fe.h.
References lf::uscalfe::FeLagrangeO2Tria< SCALAR >::NumRefShapeFunctions().
|
inlineoverridevirtual |
Six local shape functions are associated with a triangular cell in the case of quadratic Lagrangian finite elements.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 426 of file lagr_fe.h.
Referenced by lf::uscalfe::FeLagrangeO2Tria< SCALAR >::NumEvaluationNodes().
|
inlineoverridevirtual |
One shape function attached to each node and one to each edge of the triangle. There are no interior shape functions.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
One shape function attached to each node and one to each edge of the triangle. There are no interior shape functions.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
default |
|
defaultnoexcept |
|
inlineoverridevirtual |
Tells the type of reference cell underlying the parametric finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 412 of file lagr_fe.h.
References lf::base::RefEl::kTria().