![]() |
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 |
| Class used to give arguments to Legendre base functions. More... | |
| struct | MoFEM::LegendrePolynomial |
| Calculating Legendre base functions. More... | |
| struct | MoFEM::LobattoPolynomialCtx |
| Class used to give arguments to Lobatto base functions. More... | |
| struct | MoFEM::LobattoPolynomial |
| Calculating Lobatto base functions. More... | |
| struct | MoFEM::KernelLobattoPolynomialCtx |
| Class used to give arguments to Kernel Lobatto base functions. More... | |
| struct | MoFEM::KernelLobattoPolynomial |
| Calculating Lobatto base 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 | Legendre_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim) |
| Calculate Legendre approximation basis. | |
| PetscErrorCode | Lobatto_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim) |
| Calculate Lobatto base functions [FUENTES2015353]. | |
| PetscErrorCode | LobattoKernel_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim) |
| Calculate Kernel Lobatto base functions. | |
Calculation of base functions at integration points.
| PetscErrorCode Legendre_polynomials | ( | int | p, |
| double | s, | ||
| double * | diff_s, | ||
| double * | L, | ||
| double * | diffL, | ||
| const int | dim | ||
| ) |
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 15 of file base_functions.c.
| PetscErrorCode Lobatto_polynomials | ( | int | p, |
| double | s, | ||
| double * | diff_s, | ||
| double * | L, | ||
| double * | diffL, | ||
| const int | dim | ||
| ) |
Calculate Lobatto base functions [FUENTES2015353].
Order of first function is 2 and goes to p.
| p | is approximation order |
| s | is a mapping of coordinates of edge to \([-1, 1]\), i.e., \(s(\xi_1,\cdot,\xi_{dim})\in[-1,1]\) |
| diff_s | jacobian of the transformation, i.e. \(\frac{\partial
s}{\partial \xi_i}\)
|
| L | values basis functions at s |
| diffL | derivatives of basis functions at s, i.e. \(\frac{\partial L}{\partial \xi_i}\) |
| dim | dimension |
Definition at line 182 of file base_functions.c.
| PetscErrorCode LobattoKernel_polynomials | ( | int | p, |
| double | s, | ||
| double * | diff_s, | ||
| double * | L, | ||
| double * | diffL, | ||
| const int | dim | ||
| ) |
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 328 of file base_functions.c.