Dune Core Modules (unstable)
Lagrange finite element for simplices with arbitrary compile-time dimension and polynomial order. More...
#include <dune/localfunctions/lagrange/lagrangesimplex.hh>
Public Types | |
using | Traits = LocalFiniteElementTraits< Impl::LagrangeSimplexLocalBasis< D, R, d, k >, Impl::LagrangeSimplexLocalCoefficients< d, k >, Impl::LagrangeSimplexLocalInterpolation< Impl::LagrangeSimplexLocalBasis< D, R, d, k > > > |
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. | |
Detailed Description
class Dune::LagrangeSimplexLocalFiniteElement< D, R, d, k >
Lagrange finite element for simplices with arbitrary compile-time dimension and polynomial order.
- Template Parameters
-
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 \(\phi_i\) of order \(k>1\) on the unit simplex \(G = \{ x \in [0,1]^{d} \,|\, \sum_{j=1}^d x_j \leq 1\}\) are implemented as \(\phi_i(x) = \hat{\phi}_i(kx)\) where \(\hat{\phi}_i\) is the \(i\)-th basis function on the scaled simplex \(\hat{G} = kG = \{ x \in [0,k]^{d} \,|\, \sum_{j=1}^d x_i \leq k\}\). On \(\hat{G}\) the uniform Lagrange points of order \(k\) are exactly the points \((i_1,\dots,i_d) \in \mathbb{N}_0^d \cap \hat{G}\) with non-negative integer coordinates in \(\hat{G}\). Using the lexicographic enumeration of those points we can identify each \(i=(i_1,...,i_d)\) with the flat index \(i\) of the Lagrange node and associated basis function.
Since we can map any point \(\hat{x} \in \hat{G}\) to its barycentric coordinates \((\hat{x}_1, \dots, \hat{x}_d, k-\sum_{j=1}^d \hat{x}_j)\) we define for any such point its auxiliary \((d+1)\)-th coordinate \(\hat{x}_{d+1} = k-\sum_{j=1}^d \hat{x}_j\) and use this in particular for Lagrange points \(i=(i_1,...,i_d)\). Then the Lagrange basis function on \(\hat{G}\) associated to the Lagrange point \(i=(i_1,\dots,i_d)\) is given by
\[ \hat{\phi}_i(\hat{x}) = \prod_{j=1}^{d+1} L_{i_j}(\hat{x}_j) \]
where we used the barycentric coordinates of \(\hat{x}\) and \(i\) and the univariate Lagrange polynomials
\[ L_n(t) = \prod_{m=0}^{n-1} \frac{t-m}{n-m} = \frac{1}{n!}\prod_{m=0}^{n-1}(t-m). \]
This factorization can be interpreted geometrically. To this end note that any component \(1,\dots,d+1\) of the barycentric coordinates is associated to a vertex of the simplex and thus to the opposite facet. Furthermore the Lagrange nodes are all located on hyperplanes parallel to those facets. In the factorized form the term \(L_{i_j}(\hat{x}_j)\) guarantees that \(\hat{\phi}_i\) is zero on any of these hyperplane that is parallel to the facet associated with \(j\) and lying between the Lagrange node \(i\) and this facet (excluding \(i\) and including the facet).
Using the factorized form, evaluating all basis functions \(\hat{\phi}_i\) at a point \(\hat{x}\) requires to first evaluate \(\alpha_{m,j}=L_m(\hat{x}_j)\) for all \(j\in \{1,\dots,d+1\}\) and \(m \in \{0,\dots,k\}\) and then computing the products \(\hat{\phi}_i(\hat{x})=\prod_{j=1}^{d+1} \alpha_{i_j,j}\) for all admissible multi indices (or, equivalently, Lagrange nodes in barycentric coordinates) \((i_1,\dots,i_{d+1})\) with \(\sum_{j=1}^{d+1} i_j = k\). The evaluation of the univariate Lagrange polynomials can be done efficiently using the recursion
\[ L_0(t) = 1, \qquad L_{n+1}(t) = L_n(t)\frac{t-n}{n+1} \qquad n\geq 0. \]
Constructor & Destructor Documentation
◆ LagrangeSimplexLocalFiniteElement()
|
inline |
Default-construct the finite element
The documentation for this class was generated from the following file:
- dune/localfunctions/lagrange/lagrangesimplex.hh
