LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Hierarchic Finite Elements of arbitrary degree on triangles. More...
#include <lf/fe/fe.h>
Public Member Functions | |
FeHierarchicTria (const FeHierarchicTria &)=default | |
FeHierarchicTria (FeHierarchicTria &&) noexcept=default | |
FeHierarchicTria & | operator= (const FeHierarchicTria &)=default |
FeHierarchicTria & | operator= (FeHierarchicTria &&) noexcept=default |
~FeHierarchicTria () override=default | |
FeHierarchicTria (unsigned interior_degree, std::array< unsigned, 3 > edge_degrees, const quad::QuadRuleCache &qr_cache, std::span< const lf::mesh::Orientation > rel_orient) | |
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 on the edges and max(0, (p-2)*(p-1)/2) shape functions on the triangle. | |
lf::base::size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const override |
One shape function for each vertex, p-1 shape functions on the edges and max(0, (p-2)*(p-1)/2) shape functions on the triangle. | |
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 vertices, the points of a Gauss quadrature rule for each edge and the points of a quadrature rule on the interior of the triangle. | |
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 |
Computes the linear combination of reference shape functions matching function values at evaluation nodes. | |
![]() | |
virtual | ~ScalarReferenceFiniteElement ()=default |
dim_t | Dimension () const |
Returns the spatial dimension of the reference cell. | |
Private Member Functions | |
Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > | NodalValuesToVertexDofs (const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > &nodevals) const |
Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > | NodalValuesToEdgeDofs (const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > &nodevals) const |
Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > | NodalValuesToFaceDofs (const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > &nodevals) const |
Private Attributes | |
unsigned | interior_degree_ |
std::array< unsigned, 3 > | edge_degrees_ |
std::array< const quad::QuadRule *, 3 > | qr_dual_edge_ |
const lf::quad::QuadRule * | qr_dual_tria_ |
std::span< const lf::mesh::Orientation > | rel_orient_ |
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. | |
Hierarchic Finite Elements of arbitrary degree on triangles.
The shape functions associated with the vertices are given by the barycentric coordinates
\[ \begin{align*} \widehat{b^{\cdot}}^1(\widehat{x}, \widehat{y}) &= \lambda_1 = 1 - \widehat{x} - \widehat{y} \\ \widehat{b^{\cdot}}^2(\widehat{x}, \widehat{y}) &= \lambda_2 = \widehat{x} \\ \widehat{b^{\cdot}}^3(\widehat{x}, \widehat{y}) &= \lambda_3 = \widehat{y} \end{align*} \]
The basis functions associated with the triangle edges are given by the homogenized integrated Legendre polynomials
\[ \begin{align*} \widehat{b^-}^i(\widehat{x}, \widehat{y}) &= (\lambda_1 + \lambda_2)^i L_i\left(\frac{\lambda_2}{\lambda_1+\lambda_2}\right) &\mbox{ for edge 1 and } i = 2, \cdots, p \\ \widehat{b^-}^i(\widehat{x}, \widehat{y}) &= (\lambda_2 + \lambda_3)^i L_i\left(\frac{\lambda_3}{\lambda_2+\lambda_3}\right) &\mbox{ for edge 2 and } i = 2, \cdots, p \\ \widehat{b^-}^i(\widehat{x}, \widehat{y}) &= (\lambda_3 + \lambda_1)^i L_i\left(\frac{\lambda_1}{\lambda_3+\lambda_1}\right) &\mbox{ for edge 3 and } i = 2, \cdots, p \end{align*} \]
Note that the basis function on a specific edge is always zero on the other two edges. This is needed to guarantee continuity of the function space.
The basis functions associated with the interior of the triangle are given by the edge basis functions multiplied with an integrated Jacobi polynomial to force the value of the basis function to be zero on all edges.
\[ \widehat{b^{\triangle}}^{ij}(\widehat{x}, \widehat{y}) = (\lambda_1 + \lambda_2)^i L_i\left(\frac{\lambda_2}{\lambda_1+\lambda_2}\right) L_j^{2i}(\lambda_3) \quad\mbox{ for } i \geq 2, j \geq 1, i+j = 3, \cdots, p \]
where \(L_i : [0, 1] \to \mathbb{R}\) and \(L_i^{\alpha} : [0, 1] \to \mathbb{R}\) are the integrated shifted Legendre and integrated shifted Jacobi polynomials respectively:
\[ \begin{align*} L_i(x) &= \int_0^x\! P_{i-1}(\xi) \,\mathrm{d}\xi \\ L_i^{\alpha}(x) &= \int_0^x\! P_{i-1}^{(\alpha, 0)}(\xi) \,\mathrm{d}\xi \end{align*} \]
To compute the basis function coefficients from point evaluations of a function, we make use of the dual basis. For the vertices it is simply given by
\[ \lambda_i^{\cdot}[f] = \begin{cases} f(0, 0) &\mbox{ for } i = 1 \\ f(1, 0) &\mbox{ for } i = 2 \\ f(0, 1) &\mbox{ for } i = 3 \end{cases} \]
For the dual basis on the edges, we simply apply the segment dual basis along the edges of the triangle
\[ \lambda_i^-[f] = \begin{cases} \frac{1}{2i-1} \left[ P_{i-1}(1)f(1, 0) - P_{i-1}(0)f(0, 0) - \int_0^1\! P_{i-1}'(x)f(x, 0) \,\mathrm{d}x \right] &\mbox{ for edge 1} \\ \frac{1}{2i-1} \left[ P_{i-1}(1)f(0, 1) - P_{i-1}(0)f(1, 0) - \int_0^1\! P_{i-1}'(x)f(1-x, x) \,\mathrm{d}{x} \right] &\mbox{ for edge 2} \\ \frac{1}{2i-1} \left[ P_{i-1}(1)f(0, 0) - P_{i-1}(0)f(0, 1) - \int_0^1\! P_{i-1}'(x)f(0, 1-x) \,\mathrm{d}x \right] &\mbox{ for edge 2} \end{cases} \]
The dual basis for the interior shape functions is quite a bit more involved. It is given by
\[ \lambda_{ij}^{\triangle}[f] = \frac{1}{(2i-1)(2i+2j-1)} \left[ \int_0^1\! \int_0^{1-y}\! f(x, y) \left( (1-y)^{i-1}{L_j^{2i}}''(y) - 2i(1-y)^{i-2}{L_j^{2i}}'(y) \right) L_{j+1}''\left(\frac{x}{1-y}\right) \,\mathrm{d}x \,\mathrm{d}y \right] \]
and must additionally be orthogonalized with respect to the dual basis on the vertices and edges of the triangle by Gram Schmidt.
lf::mesh::Orientation
of the according edge, the local coordinate may be flipped to ensure continuity of the function space over the cell interfaces of the mesh. The basis functions and the dual basis must be adjusted accordingly in this case.A complete description of the basis functions and dual basis can be found here.
Definition at line 474 of file hierarchic_fe.h.
|
default |
|
defaultnoexcept |
|
overridedefault |
|
inline |
Definition at line 482 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, and lf::base::RefEl::kSegment().
|
inlineoverridevirtual |
Request the maximal polynomial degree of the basis functions in this finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 507 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, and lf::fe::FeHierarchicTria< SCALAR >::interior_degree_.
|
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 581 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, lf::fe::IJacobi(), lf::fe::ILegendre(), lf::fe::FeHierarchicTria< SCALAR >::interior_degree_, lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions(), lf::mesh::positive, and lf::fe::FeHierarchicTria< SCALAR >::rel_orient_.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToDofs().
|
inlineoverridevirtual |
Evaluation nodes are the vertices, the points of a Gauss quadrature rule for each edge and the points of a quadrature rule on the interior of the triangle.
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 828 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, lf::quad::QuadRule::NumPoints(), lf::quad::QuadRule::Points(), lf::mesh::positive, lf::fe::FeHierarchicTria< SCALAR >::qr_dual_edge_, lf::fe::FeHierarchicTria< SCALAR >::qr_dual_tria_, and lf::fe::FeHierarchicTria< SCALAR >::rel_orient_.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToDofs().
|
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 690 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, lf::fe::IJacobi(), lf::fe::IJacobiDx(), lf::fe::ILegendre(), lf::fe::ILegendreDt(), lf::fe::ILegendreDx(), lf::fe::FeHierarchicTria< SCALAR >::interior_degree_, lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions(), lf::mesh::positive, and lf::fe::FeHierarchicTria< SCALAR >::rel_orient_.
|
inlineoverridevirtual |
Computes the linear combination of reference shape functions matching function values at evaluation nodes.
nodvals | row vector of function values at evaluation nodes The length of this vector must agree with NumEvaluationNodes(). |
If the evaluation nodes are interpolation nodes, that is, if the set of reference shape functions forms a cardinal basis with respect to these nodes, then we have NumEvaluationNodes() == NumRefShapeFunctions() and the linear mapping realized by NodalValuesToDofs() is the identity mapping.
If the vector of values at the evaluation nodes agrees with a vector of function values of a linear combination of reference shape functions at the evaluation nodes, then this method must return the very coefficients of the linear combination.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 891 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, lf::fe::FeHierarchicTria< SCALAR >::EvalReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::EvaluationNodes(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToEdgeDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToVertexDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions().
|
inlineprivate |
Definition at line 945 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, lf::fe::ILegendreDx(), lf::fe::LegendreDx(), lf::quad::QuadRule::NumPoints(), lf::mesh::positive, lf::fe::FeHierarchicTria< SCALAR >::qr_dual_edge_, lf::fe::FeHierarchicTria< SCALAR >::qr_dual_tria_, and lf::fe::FeHierarchicTria< SCALAR >::rel_orient_.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToDofs().
|
inlineprivate |
Definition at line 1026 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, lf::fe::IJacobiDx(), lf::fe::FeHierarchicTria< SCALAR >::interior_degree_, lf::fe::JacobiDx(), lf::fe::LegendreDx(), lf::quad::QuadRule::NumPoints(), lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions(), lf::quad::QuadRule::Points(), lf::mesh::positive, lf::fe::FeHierarchicTria< SCALAR >::qr_dual_edge_, lf::fe::FeHierarchicTria< SCALAR >::qr_dual_tria_, lf::fe::FeHierarchicTria< SCALAR >::rel_orient_, and lf::quad::QuadRule::Weights().
Referenced by lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToDofs().
|
inlineprivate |
Definition at line 932 of file hierarchic_fe.h.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToDofs().
|
inlineoverridevirtual |
Tell the number of evaluation (interpolation) nodes.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 883 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, lf::quad::QuadRule::NumPoints(), lf::fe::FeHierarchicTria< SCALAR >::qr_dual_edge_, and lf::fe::FeHierarchicTria< SCALAR >::qr_dual_tria_.
|
inlineoverridevirtual |
The local shape functions.
Total number of reference shape functions associated with the reference cell.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 516 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions().
Referenced by lf::fe::FeHierarchicTria< SCALAR >::EvalReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs(), lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions(), and lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions().
|
inlineoverridevirtual |
One shape function for each vertex, p-1 shape functions on the edges and max(0, (p-2)*(p-1)/2) shape functions on the triangle.
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 526 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, and lf::fe::FeHierarchicTria< SCALAR >::interior_degree_.
|
inlineoverridevirtual |
One shape function for each vertex, p-1 shape functions on the edges and max(0, (p-2)*(p-1)/2) shape functions on the triangle.
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 558 of file hierarchic_fe.h.
References lf::fe::FeHierarchicTria< SCALAR >::edge_degrees_, and lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions().
|
default |
|
defaultnoexcept |
|
inlineoverridevirtual |
Tells the type of reference cell underlying the parametric finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 503 of file hierarchic_fe.h.
References lf::base::RefEl::kTria().
|
private |
Definition at line 921 of file hierarchic_fe.h.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::Degree(), lf::fe::FeHierarchicTria< SCALAR >::EvalReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::EvaluationNodes(), lf::fe::FeHierarchicTria< SCALAR >::FeHierarchicTria(), lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToEdgeDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs(), lf::fe::FeHierarchicTria< SCALAR >::NumEvaluationNodes(), lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions(), and lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions().
|
private |
Definition at line 920 of file hierarchic_fe.h.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::Degree(), lf::fe::FeHierarchicTria< SCALAR >::EvalReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NumRefShapeFunctions().
|
private |
Definition at line 923 of file hierarchic_fe.h.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::EvaluationNodes(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToEdgeDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NumEvaluationNodes().
|
private |
Definition at line 924 of file hierarchic_fe.h.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::EvaluationNodes(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToEdgeDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NumEvaluationNodes().
|
private |
Definition at line 925 of file hierarchic_fe.h.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::EvalReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::EvaluationNodes(), lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToEdgeDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs().