![]() |
LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Hierarchic Finite Elements of arbitrary degree on segments. More...
#include <lf/fe/fe.h>
Public Member Functions | |
| FeHierarchicSegment (const FeHierarchicSegment &)=default | |
| FeHierarchicSegment (FeHierarchicSegment &&) noexcept=default | |
| FeHierarchicSegment & | operator= (const FeHierarchicSegment &)=default |
| FeHierarchicSegment & | operator= (FeHierarchicSegment &&) noexcept=default |
| ~FeHierarchicSegment () override=default | |
| FeHierarchicSegment (unsigned degree, const quad::QuadRuleCache &qr_cache) | |
| lf::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. | |
| lf::base::size_type | NumRefShapeFunctions () const override |
| The local shape functions. | |
| lf::base::size_type | NumRefShapeFunctions (dim_t codim) const override |
| One shape function for each vertex, p-1 shape functions for the segment. | |
| lf::base::size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const override |
| One shape function for each vertex, p-1 shape functions for the segment. | |
| 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 |
| Evaluation nodes are the endpoints of the segment and the Chebyshev nodes of degree p-1 on the segment. | |
| lf::base::size_type | NumEvaluationNodes () const override |
| Tell the number of evaluation (interpolation) nodes. | |
| Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > | NodalValuesToDofs (const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > &nodevals) const override |
| Maps function evaluations to basis function coefficients. | |
Public Member Functions inherited from lf::fe::ScalarReferenceFiniteElement< SCALAR > | |
| virtual | ~ScalarReferenceFiniteElement ()=default |
| dim_t | Dimension () const |
| Returns the spatial dimension of the reference cell. | |
Private Attributes | |
| unsigned | degree_ |
| const lf::quad::QuadRule * | qr_dual_ |
Additional Inherited Members | |
Public Types inherited from lf::fe::ScalarReferenceFiniteElement< SCALAR > | |
| using | Scalar = SCALAR |
| The underlying scalar type. | |
Protected Member Functions inherited from lf::fe::ScalarReferenceFiniteElement< SCALAR > | |
| ScalarReferenceFiniteElement ()=default | |
| ScalarReferenceFiniteElement (const ScalarReferenceFiniteElement &)=default | |
| ScalarReferenceFiniteElement (ScalarReferenceFiniteElement &&) noexcept=default | |
| ScalarReferenceFiniteElement & | operator= (const ScalarReferenceFiniteElement &)=default |
| ScalarReferenceFiniteElement & | operator= (ScalarReferenceFiniteElement &&) noexcept=default |
Related Symbols inherited from lf::fe::ScalarReferenceFiniteElement< SCALAR > | |
| 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. | |
Hierarchic Finite Elements of arbitrary degree on segments.
The shape functions associated with the vertices are given by
\[ \begin{align*} \widehat{b^\cdot}^1(x) &:= 1 - x \\ \widehat{b^\cdot}^2(x) &:= x \end{align*} \]
and the interior basis functions associated with the segment itself are given by the integrated shifted Legendre polynomials
\[ \widehat{b^-}^i(x) = L_i(x) = \int_0^{x}\! P_{i-1}(\xi) \,\mathrm{d}\xi \quad\mbox{ for }\quad i= 2, \cdots, p \]
where \(P_i : [0, 1] \to \mathbb{R}\) is the shifted Legendre polynomial of degree \(i\).
To compute the basis function coefficients from point evaluations of a function, we make use of the dual basis given by
\[ \lambda_i^-[f] = \begin{cases} f(0) &\mbox{ for } i = 0 \\ f(1) &\mbox{ for } i = 1 \\ \frac{1}{2i - 1} \left [P_{i-1}(1)f(1) - P_{i-1}(0)f(0) - \int_0^1\! P_{i-1}'(x)f(x) \,\mathrm{d}x \right] &\mbox{ for } i \geq 2 \end{cases} \]
lf::mesh::Orientation::negative in order to keep a global ordering of DOFs on the cell interfaces.A complete description of the basis functions and dual basis can be found here.
Definition at line 242 of file hierarchic_fe.h.
|
default |
References FeHierarchicSegment().
Referenced by FeHierarchicSegment(), FeHierarchicSegment(), FeHierarchicSegment(), operator=(), operator=(), and ~FeHierarchicSegment().
|
defaultnoexcept |
References FeHierarchicSegment().
|
overridedefault |
References FeHierarchicSegment().
|
inline |
Definition at line 250 of file hierarchic_fe.h.
References degree_, FeHierarchicSegment(), qr_dual_, RefEl(), and lf::fe::ScalarReferenceFiniteElement< SCALAR >::ScalarReferenceFiniteElement().
|
inlinenodiscardoverridevirtual |
Request the maximal polynomial degree of the basis functions in this finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 259 of file hierarchic_fe.h.
References degree_.
|
inlinenodiscardoverridevirtual |
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 295 of file hierarchic_fe.h.
References degree_, lf::fe::ILegendre(), and NumRefShapeFunctions().
|
inlinenodiscardoverridevirtual |
Evaluation nodes are the endpoints of the segment and the Chebyshev nodes of degree p-1 on the segment.
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 >.
Definition at line 335 of file hierarchic_fe.h.
References NumEvaluationNodes(), and qr_dual_.
|
inlinenodiscardoverridevirtual |
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 312 of file hierarchic_fe.h.
References degree_, lf::fe::ILegendreDx(), and NumRefShapeFunctions().
|
inlinenodiscardoverridevirtual |
Maps function evaluations to basis function coefficients.
| nodevals | The value of the function at the evaluation nodes |
This function computes the basis function coefficients using the dual basis of the basis functions on the segment. It is given by
\[ \lambda_i^-[f] = \begin{cases} f(0) &\mbox{ for } i = 0 \\ f(1) &\mbox{ for } i = 1 \\ \frac{1}{2i - 1} \left [P_{i-1}(1)f(1) - P_{i-1}(0)f(0) - \int_0^1\! P_{i-1}'(x)f(x) \,\mathrm{d}x \right] &\mbox{ for } i \geq 2 \end{cases} \]
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 367 of file hierarchic_fe.h.
References lf::fe::ILegendreDx(), lf::fe::LegendreDx(), NumRefShapeFunctions(), qr_dual_, and lf::fe::transpose().
|
inlinenodiscardoverridevirtual |
Tell the number of evaluation (interpolation) nodes.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 348 of file hierarchic_fe.h.
References qr_dual_.
Referenced by EvaluationNodes().
|
inlinenodiscardoverridevirtual |
The local shape functions.
Total number of reference shape functions associated with the reference cell.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 265 of file hierarchic_fe.h.
References degree_.
Referenced by EvalReferenceShapeFunctions(), GradientsReferenceShapeFunctions(), NodalValuesToDofs(), and NumRefShapeFunctions().
|
inlinenodiscardoverridevirtual |
One shape function for each vertex, p-1 shape functions for the segment.
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 >.
Definition at line 274 of file hierarchic_fe.h.
References degree_.
|
inlinenodiscardoverridevirtual |
One shape function for each vertex, p-1 shape functions for the segment.
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 287 of file hierarchic_fe.h.
References NumRefShapeFunctions().
|
default |
References FeHierarchicSegment().
|
defaultnoexcept |
References FeHierarchicSegment().
|
inlinenodiscardoverridevirtual |
Tells the type of reference cell underlying the parametric finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 255 of file hierarchic_fe.h.
References lf::base::RefEl::kSegment().
Referenced by FeHierarchicSegment().
|
private |
Definition at line 392 of file hierarchic_fe.h.
Referenced by Degree(), EvalReferenceShapeFunctions(), FeHierarchicSegment(), GradientsReferenceShapeFunctions(), NumRefShapeFunctions(), and NumRefShapeFunctions().
|
private |
Definition at line 393 of file hierarchic_fe.h.
Referenced by EvaluationNodes(), FeHierarchicSegment(), NodalValuesToDofs(), and NumEvaluationNodes().