|
dune-localfunctions 2.11
|
Lagrange finite element for simplices with arbitrary compile-time dimension and polynomial order. More...
#include <dune/localfunctions/lagrange/lagrangesimplex.hh>

Public Types | |
| using | Traits |
| Export number types, dimensions, etc. | |
Public Member Functions | |
| LagrangeSimplexLocalFiniteElement () | |
| template<typename VertexMap> | |
| LagrangeSimplexLocalFiniteElement (const VertexMap &vertexmap) | |
| Constructs a finite element given a vertex reordering. | |
| const Traits::LocalBasisType & | localBasis () const |
| Returns the local basis, i.e., the set of shape functions. | |
| const Traits::LocalCoefficientsType & | localCoefficients () const |
| Returns the assignment of the degrees of freedom to the element subentities. | |
| const Traits::LocalInterpolationType & | localInterpolation () const |
| Returns object that evaluates degrees of freedom. | |
Static Public Member Functions | |
| static constexpr std::size_t | size () |
| The number of shape functions. | |
| static constexpr GeometryType | type () |
| The reference element that the local finite element is defined on. | |
Lagrange finite element for simplices with arbitrary compile-time dimension and polynomial order.
| D | type used for domain coordinates |
| R | type used for function values |
| d | dimension of the reference element |
| k | polynomial order |
The Lagrange basis functions 

![$G = \{ x \in [0,1]^{d} \,|\, \sum_{j=1}^d x_j \leq 1\}$](form_32.png)



![$\hat{G} = kG = \{ x \in [0,k]^{d} \,|\, \sum_{j=1}^d x_i \leq k\}$](form_36.png)






Since we can map any point 






![\[ \hat{\phi}_i(\hat{x}) = \prod_{j=1}^{d+1} L_{i_j}(\hat{x}_j)
\]](form_46.png)
where we used the barycentric coordinates of 

![\[ L_n(t) = \prod_{m=0}^{n-1} \frac{t-m}{n-m} = \frac{1}{n!}\prod_{m=0}^{n-1}(t-m).
\]](form_48.png)
This factorization can be interpreted geometrically. To this end note that any component 





Using the factorized form, evaluating all basis functions 







![\[ L_0(t) = 1, \qquad L_{n+1}(t) = L_n(t)\frac{t-n}{n+1} \qquad n\geq 0.
\]](form_58.png)
| using Dune::LagrangeSimplexLocalFiniteElement< D, R, d, k >::Traits |
Export number types, dimensions, etc.
|
inline |
Default-construct the finite element
|
inline |
Constructs a finite element given a vertex reordering.
|
inline |
Returns the local basis, i.e., the set of shape functions.
|
inline |
Returns the assignment of the degrees of freedom to the element subentities.
|
inline |
Returns object that evaluates degrees of freedom.
|
inlinestaticconstexpr |
The number of shape functions.
|
inlinestaticconstexpr |
The reference element that the local finite element is defined on.