![]() |
v0.15.0 |
Calculation of base functions at integration points. More...
Classes | |
| struct | MoFEM::BaseFunctionCtx |
| Base class used to exchange data between element data structures and class calculating base functions. More... | |
| struct | MoFEM::BaseFunction |
| Base class if inherited used to calculate base functions. More... | |
| struct | MoFEM::BaseFunction::DofsSideMapData |
| struct | MoFEM::EdgePolynomialBase |
| Calculate base functions on tetrahedral. More... | |
| struct | MoFEM::EntPolynomialBaseCtx |
| Class used to pass element data to calculate base functions on tet,triangle,edge. More... | |
| struct | MoFEM::FatPrismPolynomialBaseCtx |
| Class used to pass element data to calculate base functions on fat prism. More... | |
| struct | MoFEM::FatPrismPolynomialBase |
| Calculate base functions on tetrahedralFIXME: Need moab and mofem finite element structure to work (that not perfect) More... | |
| struct | MoFEM::FlatPrismPolynomialBaseCtx |
| Class used to pass element data to calculate base functions on flat prism. More... | |
| struct | MoFEM::FlatPrismPolynomialBase |
| Calculate base functions on tetrahedralFIXME: Need moab and mofem finite element structure to work (that not perfect) More... | |
| struct | MoFEM::HexPolynomialBase |
| Calculate base functions on tetrahedral. More... | |
| struct | MoFEM::JacobiPolynomialCtx |
| Class used to give arguments to Legendre base functions. More... | |
| struct | MoFEM::JacobiPolynomial |
| Calculating Legendre base functions. More... | |
| struct | MoFEM::LegendrePolynomialCtx |
| Context class for Legendre basis function arguments. More... | |
| struct | MoFEM::LegendrePolynomial |
| Class for calculating Legendre basis functions. More... | |
| struct | MoFEM::LobattoPolynomialCtx |
| Context class for Lobatto basis function arguments. More... | |
| struct | MoFEM::LobattoPolynomial |
| Class for calculating Lobatto basis functions. More... | |
| struct | MoFEM::KernelLobattoPolynomialCtx |
| Context class for Kernel Lobatto basis function arguments. More... | |
| struct | MoFEM::KernelLobattoPolynomial |
| Class for calculating Kernel Lobatto basis functions. More... | |
| struct | MoFEM::QuadPolynomialBase |
| Calculate base functions on triangle. More... | |
| struct | MoFEM::TetPolynomialBase |
| Calculate base functions on tetrahedral. More... | |
| struct | MoFEM::TriPolynomialBase |
| Calculate base functions on triangle. More... | |
Functions | |
| PetscErrorCode | LobattoKernel_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim) |
| Calculate Kernel Lobatto base functions. | |
| template<typename Scalar , typename DsScalar > | |
| MoFEMErrorCode | MoFEM::Legendre_polynomials_T (int p, Scalar &s, DsScalar *diff_s, Scalar *L, Scalar *diffL, int dim) |
| Calculate Legendre approximation basis. | |
| PetscErrorCode | MoFEM::Legendre_polynomials (int p, double s, double *diff_s, double *L, double *diffL, int dim) |
| Calculate Legendre approximation basis. | |
| PetscErrorCode | MoFEM::Lobatto_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim) |
| Calculate Lobatto basis functions [24]. | |
Calculation of base functions at integration points.
| PetscErrorCode MoFEM::Legendre_polynomials | ( | int | p, |
| double | s, | ||
| double * | diff_s, | ||
| double * | L, | ||
| double * | diffL, | ||
| int | dim | ||
| ) |
#include <src/approximation/impl/LegendrePolynomial.cpp>
Calculate Legendre approximation basis.
This function computes Legendre polynomials and their derivatives for use in finite element approximations. Legendre polynomials form a complete orthogonal system on the interval [-1,1] and are widely used in high-order finite element methods due to their excellent numerical properties.
The Legendre polynomials are defined by the recurrence relation:
\[ P_0(s)=1;\quad P_1(s) = s \]
and the following terms are generated inductively using Bonnet's recursion formula:
\[ P_{n+1}(s)=\frac{2n+1}{n+1}sP_n(s)-\frac{n}{n+1}P_{n-1}(s) \]
The parameter \(s\) is the normalized coordinate:
\[ s\in[-1,1] \quad \textrm{and}\; s=s(\xi_0,\xi_1,\xi_2) \]
where \(\xi_i\) are barycentric coordinates of the element.
For more information on Legendre polynomials, see: https://en.wikipedia.org/wiki/Legendre_polynomials
| p | approximation order (degree of highest polynomial) |
| s | position parameter \(s\in[-1,1]\) at which to evaluate polynomials |
| diff_s | derivatives of coordinate transformation, i.e. \(\frac{\partial s}{\partial \xi_i}\) |
| L | [out] array of Legendre polynomial values \(P_0(s), P_1(s), \ldots, P_p(s)\) |
| diffL | [out] array of polynomial derivatives, i.e. \(\frac{\partial P_k}{\partial \xi_i}\) |
| dim | spatial dimension for derivative calculations |
Definition at line 132 of file LegendrePolynomial.cpp.
|
inline |
#include <src/approximation/impl/LegendrePolynomial.cpp>
Calculate Legendre approximation basis.
Lagrange polynomial is given by
\[ L_0(s)=1;\quad L_1(s) = s \]
and following terms are generated inductively
\[ L_{l+1}=\frac{2l+1}{l+1}sL_l(s)-\frac{l}{l+1}L_{l-1}(s) \]
Note that:
\[ s\in[-1,1] \quad \textrm{and}\; s=s(\xi_0,\xi_1,\xi_2) \]
where \(\xi_i\) are barycentric coordinates of element.
| p | is approximation order |
| s | is position \(s\in[-1,1]\) |
| diff_s | derivatives of shape functions, i.e. \(\frac{\partial s}{\partial \xi_i}\) |
| L | approximation functions |
| diffL | derivatives, i.e. \(\frac{\partial L}{\partial \xi_i}\) |
| dim | dimension |
Definition at line 71 of file LegendrePolynomial.cpp.
| PetscErrorCode MoFEM::Lobatto_polynomials | ( | int | p, |
| double | s, | ||
| double * | diff_s, | ||
| double * | L, | ||
| double * | diffL, | ||
| const int | dim | ||
| ) |
#include <src/approximation/impl/LobattoPolynomial.cpp>
Calculate Lobatto basis functions [24].
Lobatto polynomials are hierarchical basis functions that include vertex modes and edge/internal modes. They are constructed from Legendre polynomials and provide excellent conditioning properties for high-order finite element methods.
The first function corresponds to order 2 and the sequence continues up to order p. These polynomials are particularly useful for spectral element methods due to their orthogonality properties on the reference interval [-1,1].
For more information on Lobatto polynomials, see: https://en.wikipedia.org/wiki/Gauss%E2%80%93Lobatto_quadrature
| p | approximation order (maximum polynomial degree) |
| s | coordinate mapping from edge to \([-1, 1]\), i.e., \(s(\xi_1,\cdot,\xi_{dim})\in[-1,1]\) |
| diff_s | [in] Jacobian of the coordinate transformation, i.e. \(\frac{\partial s}{\partial \xi_i}\) |
| L | [out] values of basis functions evaluated at coordinate s |
| diffL | [out] derivatives of basis functions at s, i.e. \(\frac{\partial L}{\partial \xi_i}\) |
| dim | spatial dimension |
Definition at line 80 of file LobattoPolynomial.cpp.
| PetscErrorCode LobattoKernel_polynomials | ( | int | p, |
| double | s, | ||
| double * | diff_s, | ||
| double * | L, | ||
| double * | diffL, | ||
| const int | dim | ||
| ) |
#include <src/approximation/c/base_functions.h>
Calculate Kernel Lobatto base functions.
This is implemented using definitions from Hermes2d https://github.com/hpfem/hermes following book by Pavel Solin et al [solin2003higher].
| p | is approximation order |
| s | is position \(s\in[-1,1]\) |
| diff_s | derivatives of shape functions, i.e. \(\frac{\partial s}{\partial \xi_i}\) |
| L | approximation functions |
| diffL | derivatives, i.e. \(\frac{\partial L}{\partial \xi_i}\) |
| dim | dimension |
Definition at line 229 of file base_functions.c.