![]() |
v0.14.0 |
Go to the source code of this file.
Macros | |
#define | LOBATTO_PHI0(x) (-2.0 * 1.22474487139158904909864203735) |
Definitions taken from Hermes2d code. More... | |
#define | LOBATTO_PHI1(x) (-2.0 * 1.58113883008418966599944677222 * (x)) |
#define | LOBATTO_PHI2(x) (-1.0 / 2.0 * 1.87082869338697069279187436616 * (5 * (x) * (x)-1)) |
#define | LOBATTO_PHI3(x) (-1.0 / 2.0 * 2.12132034355964257320253308631 * (7 * (x) * (x)-3) * (x)) |
#define | LOBATTO_PHI4(x) |
#define | LOBATTO_PHI5(x) |
#define | LOBATTO_PHI6(x) |
#define | LOBATTO_PHI7(x) |
#define | LOBATTO_PHI8(x) |
#define | LOBATTO_PHI9(x) |
#define | LOBATTO_PHI0X(x) (0) |
Derivatives of kernel functions for Lobbatto base. More... | |
#define | LOBATTO_PHI1X(x) (-2.0 * 1.58113883008418966599944677222) |
#define | LOBATTO_PHI2X(x) (-1.0 / 2.0 * 1.87082869338697069279187436616 * (10 * (x))) |
#define | LOBATTO_PHI3X(x) (-1.0 / 2.0 * 2.12132034355964257320253308631 * (21.0 * (x) * (x)-3.0)) |
#define | LOBATTO_PHI4X(x) |
#define | LOBATTO_PHI5X(x) |
#define | LOBATTO_PHI6X(x) |
#define | LOBATTO_PHI7X(x) |
#define | LOBATTO_PHI8X(x) |
#define | LOBATTO_PHI9X(x) |
Functions | |
PetscErrorCode | Legendre_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim) |
Calculate Legendre approximation basis. More... | |
PetscErrorCode | Jacobi_polynomials (int p, double alpha, double x, double t, double *diff_x, double *diff_t, double *L, double *diffL, const int dim) |
Calculate Jacobi approximation basis. More... | |
PetscErrorCode | IntegratedJacobi_polynomials (int p, double alpha, double x, double t, double *diff_x, double *diff_t, double *L, double *diffL, const int dim) |
Calculate integrated Jacobi approximation basis. More... | |
PetscErrorCode | Lobatto_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim) |
Calculate Lobatto base functions [25]. More... | |
PetscErrorCode | LobattoKernel_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim) |
Calculate Kernel Lobatto base functions. More... | |
#define LOBATTO_PHI0 | ( | x | ) | (-2.0 * 1.22474487139158904909864203735) |
Definitions taken from Hermes2d code.
kernel functions for Lobatto base
Definition at line 126 of file base_functions.h.
#define LOBATTO_PHI0X | ( | x | ) | (0) |
Derivatives of kernel functions for Lobbatto base.
Definition at line 157 of file base_functions.h.
#define LOBATTO_PHI1 | ( | x | ) | (-2.0 * 1.58113883008418966599944677222 * (x)) |
Definition at line 127 of file base_functions.h.
#define LOBATTO_PHI1X | ( | x | ) | (-2.0 * 1.58113883008418966599944677222) |
Definition at line 158 of file base_functions.h.
#define LOBATTO_PHI2 | ( | x | ) | (-1.0 / 2.0 * 1.87082869338697069279187436616 * (5 * (x) * (x)-1)) |
Definition at line 128 of file base_functions.h.
#define LOBATTO_PHI2X | ( | x | ) | (-1.0 / 2.0 * 1.87082869338697069279187436616 * (10 * (x))) |
Definition at line 159 of file base_functions.h.
#define LOBATTO_PHI3 | ( | x | ) | (-1.0 / 2.0 * 2.12132034355964257320253308631 * (7 * (x) * (x)-3) * (x)) |
Definition at line 130 of file base_functions.h.
#define LOBATTO_PHI3X | ( | x | ) | (-1.0 / 2.0 * 2.12132034355964257320253308631 * (21.0 * (x) * (x)-3.0)) |
Definition at line 161 of file base_functions.h.
#define LOBATTO_PHI4 | ( | x | ) |
Definition at line 132 of file base_functions.h.
#define LOBATTO_PHI4X | ( | x | ) |
Definition at line 163 of file base_functions.h.
#define LOBATTO_PHI5 | ( | x | ) |
Definition at line 135 of file base_functions.h.
#define LOBATTO_PHI5X | ( | x | ) |
Definition at line 166 of file base_functions.h.
#define LOBATTO_PHI6 | ( | x | ) |
Definition at line 138 of file base_functions.h.
#define LOBATTO_PHI6X | ( | x | ) |
Definition at line 169 of file base_functions.h.
#define LOBATTO_PHI7 | ( | x | ) |
Definition at line 141 of file base_functions.h.
#define LOBATTO_PHI7X | ( | x | ) |
Definition at line 172 of file base_functions.h.
#define LOBATTO_PHI8 | ( | x | ) |
Definition at line 144 of file base_functions.h.
#define LOBATTO_PHI8X | ( | x | ) |
Definition at line 175 of file base_functions.h.
#define LOBATTO_PHI9 | ( | x | ) |
Definition at line 149 of file base_functions.h.
#define LOBATTO_PHI9X | ( | x | ) |
Definition at line 179 of file base_functions.h.
PetscErrorCode IntegratedJacobi_polynomials | ( | int | p, |
double | alpha, | ||
double | x, | ||
double | t, | ||
double * | diff_x, | ||
double * | diff_t, | ||
double * | L, | ||
double * | diffL, | ||
const int | dim | ||
) |
Calculate integrated Jacobi approximation basis.
For more details see [26]
p | is approximation order |
alpha | polynomial parameter |
x | is position \(s\in[0,t]\) |
t | range of polynomial |
diff_x | derivatives of shape functions, i.e. \(\frac{\partial x}{\partial \xi_i}\) |
diff_t | derivatives of shape functions, i.e. \(\frac{\partial t}{\partial \xi_i}\) |
L | approximation functions |
diffL | derivatives, i.e. \(\frac{\partial L}{\partial \xi_i}\) |
dim | dimension |
Definition at line 151 of file base_functions.c.
PetscErrorCode Jacobi_polynomials | ( | int | p, |
double | alpha, | ||
double | x, | ||
double | t, | ||
double * | diff_x, | ||
double * | diff_t, | ||
double * | L, | ||
double * | diffL, | ||
const int | dim | ||
) |
Calculate Jacobi approximation basis.
For more details see [26]
p | is approximation order |
alpha | polynomial parameter |
x | is position \(s\in[0,t]\) |
t | range of polynomial |
diff_x | derivatives of shape functions, i.e. \(\frac{\partial x}{\partial \xi_i}\) |
diff_t | derivatives of shape functions, i.e. \(\frac{\partial t}{\partial \xi_i}\) |
L | approximation functions |
diffL | derivatives, i.e. \(\frac{\partial L}{\partial \xi_i}\) |
dim | dimension |
Definition at line 67 of file base_functions.c.