LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
lf::fe::FeHierarchicQuad< SCALAR > Class Template Referencefinal

Hierarchic Finite Elements of arbitrary degree on quadrilaterals. More...

#include <lf/fe/fe.h>

Inheritance diagram for lf::fe::FeHierarchicQuad< SCALAR >:
lf::fe::ScalarReferenceFiniteElement< SCALAR >

Public Member Functions

 FeHierarchicQuad (const FeHierarchicQuad &)=default
 
 FeHierarchicQuad (FeHierarchicQuad &&) noexcept=default
 
FeHierarchicQuadoperator= (const FeHierarchicQuad &)=default
 
FeHierarchicQuadoperator= (FeHierarchicQuad &&) noexcept=default
 
 ~FeHierarchicQuad () override=default
 
 FeHierarchicQuad (unsigned interior_degree, std::array< unsigned, 4 > 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 (p-1)^2 shape functions on the quadrilateral.
 
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 (p-1)^2 shape functions on the quadrilateral.
 
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 quadrature rule on each edge and the points of a quadrature rule on the interior of the quadrilateral.
 
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.
 
- 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 interior_degree_
 
std::array< unsigned, 4 > edge_degrees_
 
std::array< const quad::QuadRule *, 4 > qr_dual_edge_
 
FeHierarchicSegment< SCALAR > fe1d_
 
std::span< const lf::mesh::Orientationrel_orient_
 

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
 
ScalarReferenceFiniteElementoperator= (const ScalarReferenceFiniteElement &)=default
 
ScalarReferenceFiniteElementoperator= (ScalarReferenceFiniteElement &&) noexcept=default
 

Detailed Description

template<typename SCALAR>
class lf::fe::FeHierarchicQuad< SCALAR >

Hierarchic Finite Elements of arbitrary degree on quadrilaterals.

The basis functions on the quadrilateral has a tensor product structure and can thus be represented by products of basis functions on segments. The vertex basis functions on the reference quad are therefore given by

\[ \begin{align*} \widehat{b^{\cdot}}^0(x, y) &:= (1 - x)(1 - y) \\ \widehat{b^{\cdot}}^1(x, y) &:= x(1 - y) \\ \widehat{b^{\cdot}}^2(x, y) &:= xy \\ \widehat{b^{\cdot}}^3(x, y) &:= (1 - x)y. \end{align*} \]

The edge basis functions can be written as

\[ \begin{align*} \widehat{b^{-}}^{0,n} &:= (1-y)L_n(x) \\ \widehat{b^{-}}^{1,n} &:= xL_n(y) \\ \widehat{b^{-}}^{2,n} &:= yL_n(1-x) \\ \widehat{b^{-}}^{3,n} &:= (1-x)L_n(1-y) \end{align*} \]

where \( n \geq 2 \) is the degree of the basis function. Finally, the face bubbles are given by

\[ \widehat{b^{\square}}^{n,m}(x, y) := L_n(x)L_m(y) \]

where \( n \geq 2 \), \( m \geq 2 \).

The dual basis is therefore also quite simple, as we can recycle the one from the segments by first applying the dual basis along the \(x\)-axis and then apply the dual basis to the resulting 1d function.

Attention
Note that for the basis functions associated with the edges, depending on the 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.

See also
ScalarReferenceFiniteElement

Definition at line 1141 of file hierarchic_fe.h.

Constructor & Destructor Documentation

◆ FeHierarchicQuad() [1/3]

template<typename SCALAR >
lf::fe::FeHierarchicQuad< SCALAR >::FeHierarchicQuad ( const FeHierarchicQuad< SCALAR > & )
default

◆ FeHierarchicQuad() [2/3]

template<typename SCALAR >
lf::fe::FeHierarchicQuad< SCALAR >::FeHierarchicQuad ( FeHierarchicQuad< SCALAR > && )
defaultnoexcept

◆ ~FeHierarchicQuad()

template<typename SCALAR >
lf::fe::FeHierarchicQuad< SCALAR >::~FeHierarchicQuad ( )
overridedefault

◆ FeHierarchicQuad() [3/3]

template<typename SCALAR >
lf::fe::FeHierarchicQuad< SCALAR >::FeHierarchicQuad ( unsigned interior_degree,
std::array< unsigned, 4 > edge_degrees,
const quad::QuadRuleCache & qr_cache,
std::span< const lf::mesh::Orientation > rel_orient )
inline

Member Function Documentation

◆ Degree()

template<typename SCALAR >
unsigned lf::fe::FeHierarchicQuad< SCALAR >::Degree ( ) const
inlineoverridevirtual

Request the maximal polynomial degree of the basis functions in this finite element.

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 1179 of file hierarchic_fe.h.

References lf::fe::FeHierarchicQuad< SCALAR >::edge_degrees_, and lf::fe::FeHierarchicQuad< SCALAR >::interior_degree_.

Referenced by lf::fe::FeHierarchicQuad< SCALAR >::NodalValuesToDofs().

◆ EvalReferenceShapeFunctions()

template<typename SCALAR >
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > lf::fe::FeHierarchicQuad< SCALAR >::EvalReferenceShapeFunctions ( const Eigen::MatrixXd & refcoords) const
inlineoverridevirtual

Evaluation of all reference shape functions in a number of points.

Parameters
refcoordscoordinates 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
Returns
An Eigen Matrix of size 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.

Note
shape functions are assumed to be real-valued.

Example: Triangular Linear Lagrangian finite elements

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 1248 of file hierarchic_fe.h.

References lf::fe::FeHierarchicQuad< SCALAR >::edge_degrees_, lf::fe::FeHierarchicQuad< SCALAR >::fe1d_, lf::fe::FeHierarchicQuad< SCALAR >::interior_degree_, lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions(), lf::mesh::positive, and lf::fe::FeHierarchicQuad< SCALAR >::rel_orient_.

◆ EvaluationNodes()

template<typename SCALAR >
Eigen::MatrixXd lf::fe::FeHierarchicQuad< SCALAR >::EvaluationNodes ( ) const
inlineoverridevirtual

Evaluation nodes are the vertices, the points of a quadrature rule on each edge and the points of a quadrature rule on the interior of the quadrilateral.

Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom, nodal interpolation in the simplest case.

Returns
A d x N matrix, where d is the dimension of the reference cell, and N the number of interpolation nodes. The columns of this matrix contain their reference coordinates.

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.

Unisolvence

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.

Note
It is not required that any vector a values at evaluation nodes can be produced by a suitable linear combination of reference shape functions. For instance, this will not be possible, if there are more evaluation points than reference shape functions. If both sets have the same size, however, the interpolation problem has a solution for any vector of values at the evluation points.

Example: Principal lattice

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\).

Moment-based degrees of freedom

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 1452 of file hierarchic_fe.h.

References lf::fe::FeHierarchicQuad< SCALAR >::edge_degrees_, lf::fe::FeHierarchicQuad< SCALAR >::fe1d_, lf::fe::FeHierarchicQuad< SCALAR >::interior_degree_, lf::mesh::positive, lf::fe::FeHierarchicQuad< SCALAR >::qr_dual_edge_, and lf::fe::FeHierarchicQuad< SCALAR >::rel_orient_.

◆ GradientsReferenceShapeFunctions()

template<typename SCALAR >
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > lf::fe::FeHierarchicQuad< SCALAR >::GradientsReferenceShapeFunctions ( const Eigen::MatrixXd & refcoords) const
inlineoverridevirtual

Computation of the gradients of all reference shape functions in a number of points.

Parameters
refcoordscoordinates of N points in the reference cell passed as columns of a matrix of size dim x N.
Returns
An Eigen Matrix of size 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.

Example: Triangular Linear Lagrangian finite elements

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 1335 of file hierarchic_fe.h.

References lf::fe::FeHierarchicQuad< SCALAR >::edge_degrees_, lf::fe::FeHierarchicQuad< SCALAR >::fe1d_, lf::fe::FeHierarchicQuad< SCALAR >::interior_degree_, lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions(), lf::mesh::positive, and lf::fe::FeHierarchicQuad< SCALAR >::rel_orient_.

◆ NodalValuesToDofs()

template<typename SCALAR >
Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > lf::fe::FeHierarchicQuad< SCALAR >::NodalValuesToDofs ( const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > & nodvals) const
inlineoverridevirtual

Computes the linear combination of reference shape functions matching function values at evaluation nodes.

Parameters
nodvalsrow vector of function values at evaluation nodes The length of this vector must agree with NumEvaluationNodes().
Returns
The coefficients of the local "nodal interpolant" with respect to the reference shape functions. This is a row vector of length NumRefShapeFunctions().

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.

Note
default implementation is the identity mapping

Requirement: reproduction of local finite element functions

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 1545 of file hierarchic_fe.h.

References lf::fe::FeHierarchicQuad< SCALAR >::Degree(), lf::fe::FeHierarchicQuad< SCALAR >::edge_degrees_, lf::fe::FeHierarchicQuad< SCALAR >::fe1d_, lf::fe::ILegendreDx(), lf::fe::FeHierarchicQuad< SCALAR >::interior_degree_, lf::fe::LegendreDx(), lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions(), lf::mesh::positive, lf::fe::FeHierarchicQuad< SCALAR >::qr_dual_edge_, and lf::fe::FeHierarchicQuad< SCALAR >::rel_orient_.

◆ NumEvaluationNodes()

template<typename SCALAR >
lf::base::size_type lf::fe::FeHierarchicQuad< SCALAR >::NumEvaluationNodes ( ) const
inlineoverridevirtual

◆ NumRefShapeFunctions() [1/3]

template<typename SCALAR >
lf::base::size_type lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions ( ) const
inlineoverridevirtual

The local shape functions.

Total number of reference shape functions associated with the reference cell.

Note
the total number of shape functions is the sum of the number of interior shape functions for all sub-entities and the entity itself.

Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 1188 of file hierarchic_fe.h.

References lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions().

Referenced by lf::fe::FeHierarchicQuad< SCALAR >::EvalReferenceShapeFunctions(), lf::fe::FeHierarchicQuad< SCALAR >::GradientsReferenceShapeFunctions(), lf::fe::FeHierarchicQuad< SCALAR >::NodalValuesToDofs(), lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions(), and lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions().

◆ NumRefShapeFunctions() [2/3]

template<typename SCALAR >
lf::base::size_type lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions ( dim_t codim) const
inlineoverridevirtual

One shape function for each vertex, p-1 shape functions on the edges and (p-1)^2 shape functions on the quadrilateral.

The number of interior reference shape functions for sub-entities of a particular co-dimension.

Parameters
codimco-dimension of the subentity
Returns
number of interior reference shape function belonging to the specified sub-entity.
Note
this method will throw an exception when different numbers of reference shape functions are assigned to different sub-entities of the same co-dimension

Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 1199 of file hierarchic_fe.h.

References lf::fe::FeHierarchicQuad< SCALAR >::edge_degrees_, and lf::fe::FeHierarchicQuad< SCALAR >::interior_degree_.

◆ NumRefShapeFunctions() [3/3]

template<typename SCALAR >
lf::base::size_type lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions ( dim_t codim,
sub_idx_t subidx ) const
inlineoverridevirtual

One shape function for each vertex, p-1 shape functions on the edges and (p-1)^2 shape functions on the quadrilateral.

The number of interior reference shape functions for every sub-entity.

Parameters
codimdo-dimension of the subentity
subidxlocal index of the sub-entity
Returns
number of interior reference shape function belonging to the specified sub-entity.

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 1228 of file hierarchic_fe.h.

References lf::fe::FeHierarchicQuad< SCALAR >::edge_degrees_, and lf::fe::FeHierarchicQuad< SCALAR >::NumRefShapeFunctions().

◆ operator=() [1/2]

template<typename SCALAR >
FeHierarchicQuad & lf::fe::FeHierarchicQuad< SCALAR >::operator= ( const FeHierarchicQuad< SCALAR > & )
default

◆ operator=() [2/2]

template<typename SCALAR >
FeHierarchicQuad & lf::fe::FeHierarchicQuad< SCALAR >::operator= ( FeHierarchicQuad< SCALAR > && )
defaultnoexcept

◆ RefEl()

template<typename SCALAR >
lf::base::RefEl lf::fe::FeHierarchicQuad< SCALAR >::RefEl ( ) const
inlineoverridevirtual

Tells the type of reference cell underlying the parametric finite element.

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 1175 of file hierarchic_fe.h.

References lf::base::RefEl::kQuad().

Member Data Documentation

◆ edge_degrees_

template<typename SCALAR >
std::array<unsigned, 4> lf::fe::FeHierarchicQuad< SCALAR >::edge_degrees_
private

◆ fe1d_

template<typename SCALAR >
FeHierarchicSegment<SCALAR> lf::fe::FeHierarchicQuad< SCALAR >::fe1d_
private

◆ interior_degree_

template<typename SCALAR >
unsigned lf::fe::FeHierarchicQuad< SCALAR >::interior_degree_
private

◆ qr_dual_edge_

template<typename SCALAR >
std::array<const quad::QuadRule *, 4> lf::fe::FeHierarchicQuad< SCALAR >::qr_dual_edge_
private

◆ rel_orient_

template<typename SCALAR >
std::span<const lf::mesh::Orientation> lf::fe::FeHierarchicQuad< SCALAR >::rel_orient_
private

The documentation for this class was generated from the following file: