LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Cubic Lagrangian finite element on a quadrilateral reference element. More...
#include <lf/uscalfe/uscalfe.h>
Public Member Functions | |
FeLagrangeO3Quad (const FeLagrangeO3Quad &)=default | |
FeLagrangeO3Quad (FeLagrangeO3Quad &&) noexcept=default | |
FeLagrangeO3Quad & | operator= (const FeLagrangeO3Quad &)=default |
FeLagrangeO3Quad & | operator= (FeLagrangeO3Quad &&) noexcept=default |
FeLagrangeO3Quad ()=default | |
~FeLagrangeO3Quad () override=default | |
base::RefEl | RefEl () const override |
Tells the type of reference cell underlying the parametric finite element. | |
unsigned | Degree () const override |
Request the maximal polynomial degree of the basis functions in this finite element. | |
size_type | NumRefShapeFunctions () const override |
Total number of reference shape functions associated with the reference cell. | |
size_type | NumRefShapeFunctions (dim_t codim) const override |
One shape function attached to each node and two to each edge of the quadrilateral. Four interior shape functions. | |
size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const override |
One shape function attached to each node and two to each edge of the quadrilateral. Four interior shape functions. | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | EvalReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override |
Evaluation of all reference shape functions in a number of points. | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | GradientsReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override |
Computation of the gradients of all reference shape functions in a number of points. | |
Eigen::MatrixXd | EvaluationNodes () const override |
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom, nodal interpolation in the simplest case. | |
size_type | NumEvaluationNodes () const override |
Sixteen evaluation nodes. | |
![]() | |
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. | |
Static Private Attributes | |
static const FeLagrangeO3Segment< SCALAR > | krsf_segment_ |
static const Eigen::Matrix< int, 16, 2 > | ksegment_to_quad_mapping_ |
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. | |
Cubic Lagrangian finite element on a quadrilateral reference element.
This is a specialization of ScalarReferenceFiniteElement. Refer to its documentation
The shape functions are computed by a tensor product construction: The shape function associated with the interpolation node \( \mathbf{p} =(x_j,y_l)\) is computed as
\[ \hat{b}^{\mathbf{p}}(x_1,x_2)=\hat{b}^j(x_1)*\hat{b}^l(x_2) \]
Where \( \hat{b}^j\) and \( \hat{b}^l \) are the cubic Lagrangian shape functions on the reference line segment associated to the interpolation nodes \( x_j \) and \( y_l \), see Figure 119.
The first four shape functions are associated with the vertices, the next eight = 4*2 with edges of the quadrilateral. The last four shape functions are interior shape functions.
|
default |
|
defaultnoexcept |
|
default |
|
overridedefault |
|
inlineoverridevirtual |
Request the maximal polynomial degree of the basis functions in this finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Evaluation of all reference shape functions in a number of points.
refcoords | coordinates of N points in the reference cell passed as columns of a matrix of size dim x N, where dim is the dimension of the reference element, that is =0 for points, =1 for edges, =2 for triangles and quadrilaterals |
NumRefShapeFunctions() x refcoords.cols()
which contains the shape functions evaluated at every quadrature point.Concerning the numbering of local shape functions, please consult the documentation of lf::assemble::DofHandler or the documentation of the class.
There are three reference shape functions \(\hat{b}^1,\hat{b}^2,\hat{b}^3\) associated with the vertices of the reference triangle. Let us assume that the refcoords
argument is a 2x2 matrix \([\mathbf{x}_1\;\mathbf{x}_2]\), which corresponds to passing the coordinates of two points in the reference triangle. Then this method will return a 3x2
matrix:
\[ \begin{pmatrix}\hat{b}^1(\mathbf{x}_1) & \hat{b}^1(\mathbf{x}_2) \\ \hat{b}^2(\mathbf{x}_1) & \hat{b}^2(\mathbf{x}_2) \\ \hat{b}^3(\mathbf{x}_1)\ & \hat{b}^3(\mathbf{x}_2) \end{pmatrix} \]
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 1332 of file lagr_fe.h.
References lf::uscalfe::FeLagrangeO3Quad< SCALAR >::krsf_segment_, and lf::uscalfe::FeLagrangeO3Quad< SCALAR >::ksegment_to_quad_mapping_.
|
inlineoverridevirtual |
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom, nodal interpolation in the simplest case.
Every parametric scalar finite element implicitly defines a local interpolation operator by duality with the reference shape functions. This interpolation operator can be realized through evaluations at certain evaluation nodes, which are provided by this method.
The evaluation points must satisfy the following requirement: If the values of a function belonging to the span of the reference shape functions are known in the evaluation nodes, then this function is uniquely determined. This entails that the number of evaluation nodes must be at least as big as the number of reference shape functions.
For triangular Lagrangian finite elements of degree p the evaluation nodes, usually called "interpolation nodes" in this context, can be chosen as \(\left(\frac{j}{p},\frac{k}{p}\right),\; 0\leq j,k \leq p, j+k\leq p\).
For some finite element spaces the interpolation functional may be defined based on integrals over edges. In this case the evaluation nodes will be quadrature nodes for the approximate evaluation of these integrals.
The quadrature rule must be exact for the polynomials contained in the local finite element spaces.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Computation of the gradients of all reference shape functions in a number of points.
refcoords | coordinates of N points in the reference cell passed as columns of a matrix of size dim x N. |
NumRefShapeFunctions() x (dim * refcoords.cols())
where dim
is the dimension of the reference finite element. The gradients are returned in chunks of rows of this matrix.Concerning the numbering of local shape functions, please consult the documentation of lf::assemble::DofHandler.
There are three reference shape functions \(\hat{b}^1,\hat{b}^2,\hat{b}^3\) associated with the vertices of the reference triangle. Let us assume that the refcoords
argument is a 2x2 matrix \([\mathbf{x}_1\;\mathbf{x}_2]\), which corresponds to passing the coordinates of two points in the reference triangle. Then this method will return a 3x4
matrix:
\[ \begin{pmatrix} \grad^T\hat{b}^1(\mathbf{x}_1) & \grad^T\hat{b}^1(\mathbf{x}_2) \\ \grad^T\hat{b}^2(\mathbf{x}_1) & \grad^T\hat{b}^2(\mathbf{x}_2) \\ \grad^T\hat{b}^3(\mathbf{x}_1) & \grad^T\hat{b}^3(\mathbf{x}_2) \end{pmatrix} \]
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 1356 of file lagr_fe.h.
References lf::uscalfe::FeLagrangeO3Quad< SCALAR >::krsf_segment_, and lf::uscalfe::FeLagrangeO3Quad< SCALAR >::ksegment_to_quad_mapping_.
|
inlineoverridevirtual |
Sixteen evaluation nodes.
Tell the number of evaluation (interpolation) nodes.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 1426 of file lagr_fe.h.
References lf::uscalfe::FeLagrangeO3Quad< SCALAR >::NumRefShapeFunctions().
|
inlineoverridevirtual |
Total number of reference shape functions associated with the reference cell.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 1295 of file lagr_fe.h.
Referenced by lf::uscalfe::FeLagrangeO3Quad< SCALAR >::NumEvaluationNodes(), and lf::uscalfe::FeLagrangeO3Quad< SCALAR >::NumRefShapeFunctions().
|
inlineoverridevirtual |
One shape function attached to each node and two to each edge of the quadrilateral. Four interior shape functions.
The number of interior reference shape functions for sub-entities of a particular co-dimension.
codim | co-dimension of the subentity |
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
One shape function attached to each node and two to each edge of the quadrilateral. Four interior shape functions.
The number of interior reference shape functions for every sub-entity.
codim | do-dimension of the subentity |
subidx | local index of the sub-entity |
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 1323 of file lagr_fe.h.
References lf::uscalfe::FeLagrangeO3Quad< SCALAR >::NumRefShapeFunctions().
|
default |
|
defaultnoexcept |
|
inlineoverridevirtual |
Tells the type of reference cell underlying the parametric finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 1289 of file lagr_fe.h.
References lf::base::RefEl::kQuad().
|
staticprivate |
description of the shape functions on the reference line segment \( [0,1] \) that are used in the tensor product construction of the shape functions
Definition at line 1434 of file lagr_fe.h.
Referenced by lf::uscalfe::FeLagrangeO3Quad< SCALAR >::EvalReferenceShapeFunctions(), and lf::uscalfe::FeLagrangeO3Quad< SCALAR >::GradientsReferenceShapeFunctions().
|
staticprivate |
The shape functions are computed by a tensor product construction: The shape function associated with the interpolation node \( \mathbf{p} =(x_j,y_l)\) is computed as
\[ \hat{b}^{\mathbf{p}}(x_1,x_2)=\hat{b}^j(x_1)*\hat{b}^l(x_2) \]
Where \( \hat{b}^j\) and \( \hat{b}^l \) are the cubic Lagrangian shape functions on the reference line segment associated to the interpolation nodes \( x_j \) and \( y_l \).
This matrix of size NumRefShapeFunctions x 2 contains in each row the indices of the "segment" shape functions, that define the reference shape function via the above formula. if \(p\) was the \(k\)-th interpolation node, the \(k\)-th row would contain the indices \(j\) and \(l\).
Definition at line 1453 of file lagr_fe.h.
Referenced by lf::uscalfe::FeLagrangeO3Quad< SCALAR >::EvalReferenceShapeFunctions(), and lf::uscalfe::FeLagrangeO3Quad< SCALAR >::GradientsReferenceShapeFunctions().