v0.14.0
Loading...
Searching...
No Matches
Base functions

Calculation of base functions at integration points. More...

Collaboration diagram for Base functions:

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 [28].
 
PetscErrorCode LobattoKernel_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim)
 Calculate Kernel Lobatto base functions.
 

Detailed Description

Calculation of base functions at integration points.

Function Documentation

◆ Legendre_polynomials()

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.

Parameters
pis approximation order
sis position \(s\in[-1,1]\)
diff_sderivatives of shape functions, i.e. \(\frac{\partial s}{\partial \xi_i}\)
Return values
Lapproximation functions
diffLderivatives, i.e. \(\frac{\partial L}{\partial \xi_i}\)
Parameters
dimdimension
Returns
error code
Examples
EshelbianPlasticity.cpp, and edge_and_bubble_shape_functions_on_quad.cpp.

Definition at line 15 of file base_functions.c.

16 {
18#ifndef NDEBUG
19 if (dim < 1)
20 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim < 1");
21 if (dim > 3)
22 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim > 3");
23 if (p < 0)
24 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "p < 0");
25#endif // NDEBUG
26
27 L[0] = 1;
28 if (diffL != NULL) {
29 for (int d = 0; d != dim; ++d) {
30 diffL[d * (p + 1) + 0] = 0;
31 }
32 }
33 if (p == 0)
35
36 L[1] = s;
37 if (diffL != NULL) {
38#ifndef NDEBUG
39 if (diff_s == NULL) {
40 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "diff_s == NULL");
41 }
42#endif // NDEBUG
43 for (int d = 0; d != dim; ++d) {
44 diffL[d * (p + 1) + 1] = diff_s[d];
45 }
46 }
47 if (p == 1)
49
50 int l = 1;
51 for (; l < p; l++) {
52 double A = ((2 * (double)l + 1) / ((double)l + 1));
53 double B = ((double)l / ((double)l + 1));
54 L[l + 1] = A * s * L[l] - B * L[l - 1];
55 if (diffL != NULL) {
56 for (int d = 0; d != dim; ++d) {
57 diffL[d * (p + 1) + l + 1] =
58 A * (diff_s[d] * L[l] + s * diffL[d * (p + 1) + l]) -
59 B * diffL[d * (p + 1) + l - 1];
60 }
61 }
62 }
63
65}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ MOFEM_INVALID_DATA
Definition definitions.h:36
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
FTensor::Index< 'l', 3 > l

◆ Lobatto_polynomials()

PetscErrorCode Lobatto_polynomials ( int p,
double s,
double * diff_s,
double * L,
double * diffL,
const int dim )

Calculate Lobatto base functions [28].

Order of first function is 2 and goes to p.

Parameters
pis approximation order
sis a mapping of coordinates of edge to \([-1, 1]\), i.e., \(s(\xi_1,\cdot,\xi_{dim})\in[-1,1]\)
diff_sjacobian of the transformation, i.e. \(\frac{\partial s}{\partial \xi_i}\)
  • output
Return values
Lvalues basis functions at s
diffLderivatives of basis functions at s, i.e. \(\frac{\partial L}{\partial \xi_i}\)
Parameters
dimdimension
Returns
error code

Definition at line 182 of file base_functions.c.

183 {
184
186#ifndef NDEBUG
187 if (dim < 1)
188 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim < 1");
189 if (dim > 3)
190 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim > 3");
191 if (p < 2)
192 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "p < 2");
193#endif // NDEBUG
194 double l[p + 1];
195 ierr = Legendre_polynomials(p, s, NULL, l, NULL, 1);
196 CHKERRQ(ierr);
197
198 L[0] = 1;
199 if (diffL != NULL) {
200 for (int d = 0; d != dim; ++d) {
201 diffL[d * (p + 1) + 0] = 0;
202 }
203 }
204 L[1] = s;
205 if (diffL != NULL) {
206#ifndef NDEBUG
207 if (diff_s == NULL) {
208 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "diff_s == NULL");
209 }
210#endif // NDEBUG
211 for (int d = 0; d != dim; ++d) {
212 diffL[d * (p + 1) + 1] = diff_s[d];
213 }
214 }
215
216 // Integrated Legendre
217 for (int k = 2; k <= p; k++) {
218 const double factor = 2 * (2 * k - 1);
219 L[k] = 1.0 / factor * (l[k] - l[k - 2]);
220 }
221
222 if (diffL != NULL) {
223 for (int k = 2; k <= p; k++) {
224 double a = l[k - 1] / 2.;
225 for (int d = 0; d != dim; ++d) {
226 diffL[d * (p + 1) + k] = a * diff_s[d];
227 }
228 }
229 }
231}
constexpr double a
static PetscErrorCode ierr
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Legendre approximation basis.
FTensor::Index< 'k', 3 > k

◆ LobattoKernel_polynomials()

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].

Parameters
pis approximation order
sis position \(s\in[-1,1]\)
diff_sderivatives of shape functions, i.e. \(\frac{\partial s}{\partial \xi_i}\)
Return values
Lapproximation functions
diffLderivatives, i.e. \(\frac{\partial L}{\partial \xi_i}\)
Parameters
dimdimension
Returns
error code

Definition at line 328 of file base_functions.c.

330 {
332#ifndef NDEBUG
333 if (dim < 1)
334 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim < 1");
335 if (dim > 3)
336 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim > 3");
337 if (p < 0)
338 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "p < 0");
339 if (p > 9)
340 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
341 "Polynomial beyond order 9 is not implemented");
342 if (diffL != NULL) {
343 if (diff_s == NULL) {
344 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "diff_s == NULL");
345 }
346 }
347
348#endif // NDEBUG
349 if (L) {
350 int l = 0;
351 for (; l != p + 1; l++) {
352 L[l] = f_phi[l](s);
353 }
354 }
355 if (diffL != NULL) {
356 int l = 0;
357 for (; l != p + 1; l++) {
358 double a = f_phix[l](s);
359 int d = 0;
360 for (; d < dim; ++d) {
361 diffL[d * (p + 1) + l] = diff_s[d] * a;
362 }
363 }
364 }
366}
static double(* f_phi[])(double x)
static double(* f_phix[])(double x)
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32