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

Linear Lagrange finite element on triangular reference element. More...

#include <lf/uscalfe/uscalfe.h>

Inheritance diagram for lf::uscalfe::FeLagrangeO1Tria< SCALAR >:
lf::fe::ScalarReferenceFiniteElement< SCALAR >

Public Member Functions

 FeLagrangeO1Tria (const FeLagrangeO1Tria &)=default
 
 FeLagrangeO1Tria (FeLagrangeO1Tria &&) noexcept=default
 
FeLagrangeO1Triaoperator= (const FeLagrangeO1Tria &)=default
 
FeLagrangeO1Triaoperator= (FeLagrangeO1Tria &&) noexcept=default
 
 FeLagrangeO1Tria ()=default
 
 ~FeLagrangeO1Tria () 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
 The local shape functions: barycentric coordinate functions.
 
size_type NumRefShapeFunctions (dim_t codim) const override
 Exactly one shape function attached to each node, none to other sub-entities.
 
size_type NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const override
 Exactly one shape function attached to each node, none to other sub-entities.
 
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
 Evalutation nodes are just the vertices of the triangle.
 
size_type NumEvaluationNodes () const override
 Three 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.
 
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

- 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::uscalfe::FeLagrangeO1Tria< SCALAR >

Linear Lagrange finite element on triangular reference element.

This is a specialization of fe::ScalarReferenceFiniteElement. Refer to its documentation.

The reference shape functions are

\[ \widehat{b}^1(\widehat{\mathbf{x}}) = 1-\widehat{x}_1-\widehat{x}_2\quad,\quad \widehat{b}^2(\widehat{\mathbf{x}}) = \widehat{x}_1\quad,\quad \widehat{b}^3(\widehat{\mathbf{x}}) = \widehat{x}_2 \;. \]

The basis function \(\widehat{b}^i\) is associated with vertex \(i\) of the triangle.

See also
fe::ScalarReferenceFiniteElement

Definition at line 57 of file lagr_fe.h.

Constructor & Destructor Documentation

◆ FeLagrangeO1Tria() [1/3]

template<typename SCALAR >
lf::uscalfe::FeLagrangeO1Tria< SCALAR >::FeLagrangeO1Tria ( const FeLagrangeO1Tria< SCALAR > & )
default

◆ FeLagrangeO1Tria() [2/3]

template<typename SCALAR >
lf::uscalfe::FeLagrangeO1Tria< SCALAR >::FeLagrangeO1Tria ( FeLagrangeO1Tria< SCALAR > && )
defaultnoexcept

◆ FeLagrangeO1Tria() [3/3]

template<typename SCALAR >
lf::uscalfe::FeLagrangeO1Tria< SCALAR >::FeLagrangeO1Tria ( )
default

◆ ~FeLagrangeO1Tria()

template<typename SCALAR >
lf::uscalfe::FeLagrangeO1Tria< SCALAR >::~FeLagrangeO1Tria ( )
overridedefault

Member Function Documentation

◆ Degree()

template<typename SCALAR >
unsigned lf::uscalfe::FeLagrangeO1Tria< 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 71 of file lagr_fe.h.

◆ EvalReferenceShapeFunctions()

template<typename SCALAR >
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > lf::uscalfe::FeLagrangeO1Tria< 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 104 of file lagr_fe.h.

◆ EvaluationNodes()

template<typename SCALAR >
Eigen::MatrixXd lf::uscalfe::FeLagrangeO1Tria< SCALAR >::EvaluationNodes ( ) const
inlineoverridevirtual

Evalutation nodes are just the vertices of the triangle.

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 135 of file lagr_fe.h.

References lf::base::RefEl::NodeCoords(), and lf::uscalfe::FeLagrangeO1Tria< SCALAR >::RefEl().

◆ GradientsReferenceShapeFunctions()

template<typename SCALAR >
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > lf::uscalfe::FeLagrangeO1Tria< 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 117 of file lagr_fe.h.

◆ NumEvaluationNodes()

template<typename SCALAR >
size_type lf::uscalfe::FeLagrangeO1Tria< SCALAR >::NumEvaluationNodes ( ) const
inlineoverridevirtual

Three evaluation nodes.

Tell the number of evaluation (interpolation) nodes.

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

Definition at line 142 of file lagr_fe.h.

References lf::base::RefEl::NumNodes(), and lf::uscalfe::FeLagrangeO1Tria< SCALAR >::RefEl().

◆ NumRefShapeFunctions() [1/3]

template<typename SCALAR >
size_type lf::uscalfe::FeLagrangeO1Tria< SCALAR >::NumRefShapeFunctions ( ) const
inlineoverridevirtual

The local shape functions: barycentric coordinate 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 76 of file lagr_fe.h.

◆ NumRefShapeFunctions() [2/3]

template<typename SCALAR >
size_type lf::uscalfe::FeLagrangeO1Tria< SCALAR >::NumRefShapeFunctions ( dim_t codim) const
inlineoverridevirtual

Exactly one shape function attached to each node, none to other sub-entities.

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 84 of file lagr_fe.h.

◆ NumRefShapeFunctions() [3/3]

template<typename SCALAR >
size_type lf::uscalfe::FeLagrangeO1Tria< SCALAR >::NumRefShapeFunctions ( dim_t codim,
sub_idx_t subidx ) const
inlineoverridevirtual

Exactly one shape function attached to each node, none to other sub-entities.

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 94 of file lagr_fe.h.

◆ operator=() [1/2]

template<typename SCALAR >
FeLagrangeO1Tria & lf::uscalfe::FeLagrangeO1Tria< SCALAR >::operator= ( const FeLagrangeO1Tria< SCALAR > & )
default

◆ operator=() [2/2]

template<typename SCALAR >
FeLagrangeO1Tria & lf::uscalfe::FeLagrangeO1Tria< SCALAR >::operator= ( FeLagrangeO1Tria< SCALAR > && )
defaultnoexcept

◆ RefEl()

template<typename SCALAR >
base::RefEl lf::uscalfe::FeLagrangeO1Tria< SCALAR >::RefEl ( ) const
inlineoverridevirtual

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

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

Definition at line 67 of file lagr_fe.h.

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

Referenced by lf::uscalfe::FeLagrangeO1Tria< SCALAR >::EvaluationNodes(), and lf::uscalfe::FeLagrangeO1Tria< SCALAR >::NumEvaluationNodes().


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