v0.14.0
Loading...
Searching...
No Matches
base_functions.h File Reference

Go to the source code of this file.

Macros

#define LOBATTO_PHI0(x)
 Definitions taken from Hermes2d code.
 
#define LOBATTO_PHI1(x)
 
#define LOBATTO_PHI2(x)
 
#define LOBATTO_PHI3(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)
 Derivatives of kernel functions for Lobbatto base.
 
#define LOBATTO_PHI1X(x)
 
#define LOBATTO_PHI2X(x)
 
#define LOBATTO_PHI3X(x)
 
#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.
 
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.
 
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.
 
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.
 

Macro Definition Documentation

◆ LOBATTO_PHI0

#define LOBATTO_PHI0 ( x)
Value:
(-2.0 * 1.22474487139158904909864203735)

Definitions taken from Hermes2d code.

kernel functions for Lobatto base

Definition at line 126 of file base_functions.h.

◆ LOBATTO_PHI0X

#define LOBATTO_PHI0X ( x)
Value:
(0)

Derivatives of kernel functions for Lobbatto base.

Definition at line 157 of file base_functions.h.

◆ LOBATTO_PHI1

#define LOBATTO_PHI1 ( x)
Value:
(-2.0 * 1.58113883008418966599944677222 * (x))

Definition at line 127 of file base_functions.h.

◆ LOBATTO_PHI1X

#define LOBATTO_PHI1X ( x)
Value:
(-2.0 * 1.58113883008418966599944677222)

Definition at line 158 of file base_functions.h.

◆ LOBATTO_PHI2

#define LOBATTO_PHI2 ( x)
Value:
(-1.0 / 2.0 * 1.87082869338697069279187436616 * (5 * (x) * (x)-1))

Definition at line 128 of file base_functions.h.

128#define LOBATTO_PHI2(x) \
129 (-1.0 / 2.0 * 1.87082869338697069279187436616 * (5 * (x) * (x)-1))

◆ LOBATTO_PHI2X

#define LOBATTO_PHI2X ( x)
Value:
(-1.0 / 2.0 * 1.87082869338697069279187436616 * (10 * (x)))

Definition at line 159 of file base_functions.h.

159#define LOBATTO_PHI2X(x) \
160 (-1.0 / 2.0 * 1.87082869338697069279187436616 * (10 * (x)))

◆ LOBATTO_PHI3

#define LOBATTO_PHI3 ( x)
Value:
(-1.0 / 2.0 * 2.12132034355964257320253308631 * (7 * (x) * (x)-3) * (x))

Definition at line 130 of file base_functions.h.

130#define LOBATTO_PHI3(x) \
131 (-1.0 / 2.0 * 2.12132034355964257320253308631 * (7 * (x) * (x)-3) * (x))

◆ LOBATTO_PHI3X

#define LOBATTO_PHI3X ( x)
Value:
(-1.0 / 2.0 * 2.12132034355964257320253308631 * (21.0 * (x) * (x)-3.0))

Definition at line 161 of file base_functions.h.

161#define LOBATTO_PHI3X(x) \
162 (-1.0 / 2.0 * 2.12132034355964257320253308631 * (21.0 * (x) * (x)-3.0))

◆ LOBATTO_PHI4

#define LOBATTO_PHI4 ( x)
Value:
(-1.0 / 4.0 * 2.34520787991171477728281505677 * \
(21 * (x) * (x) * (x) * (x)-14 * (x) * (x) + 1))

Definition at line 132 of file base_functions.h.

132#define LOBATTO_PHI4(x) \
133 (-1.0 / 4.0 * 2.34520787991171477728281505677 * \
134 (21 * (x) * (x) * (x) * (x)-14 * (x) * (x) + 1))

◆ LOBATTO_PHI4X

#define LOBATTO_PHI4X ( x)
Value:
(-1.0 / 4.0 * 2.34520787991171477728281505677 * \
((84.0 * (x) * (x)-28.0) * (x)))

Definition at line 163 of file base_functions.h.

163#define LOBATTO_PHI4X(x) \
164 (-1.0 / 4.0 * 2.34520787991171477728281505677 * \
165 ((84.0 * (x) * (x)-28.0) * (x)))

◆ LOBATTO_PHI5

#define LOBATTO_PHI5 ( x)
Value:
(-1.0 / 4.0 * 2.54950975679639241501411205451 * \
((33 * (x) * (x)-30) * (x) * (x) + 5) * (x))

Definition at line 135 of file base_functions.h.

135#define LOBATTO_PHI5(x) \
136 (-1.0 / 4.0 * 2.54950975679639241501411205451 * \
137 ((33 * (x) * (x)-30) * (x) * (x) + 5) * (x))

◆ LOBATTO_PHI5X

#define LOBATTO_PHI5X ( x)
Value:
(-1.0 / 4.0 * 2.54950975679639241501411205451 * \
((165.0 * (x) * (x)-90.0) * (x) * (x) + 5.0))

Definition at line 166 of file base_functions.h.

166#define LOBATTO_PHI5X(x) \
167 (-1.0 / 4.0 * 2.54950975679639241501411205451 * \
168 ((165.0 * (x) * (x)-90.0) * (x) * (x) + 5.0))

◆ LOBATTO_PHI6

#define LOBATTO_PHI6 ( x)
Value:
(-1.0 / 32.0 * 2.73861278752583056728484891400 * \
(((429 * (x) * (x)-495) * (x) * (x) + 135) * (x) * (x)-5))

Definition at line 138 of file base_functions.h.

138#define LOBATTO_PHI6(x) \
139 (-1.0 / 32.0 * 2.73861278752583056728484891400 * \
140 (((429 * (x) * (x)-495) * (x) * (x) + 135) * (x) * (x)-5))

◆ LOBATTO_PHI6X

#define LOBATTO_PHI6X ( x)
Value:
(-1.0 / 32.0 * 2.73861278752583056728484891400 * \
(((2574.0 * (x) * (x)-1980.0) * (x) * (x) + 270.0) * (x)))

Definition at line 169 of file base_functions.h.

169#define LOBATTO_PHI6X(x) \
170 (-1.0 / 32.0 * 2.73861278752583056728484891400 * \
171 (((2574.0 * (x) * (x)-1980.0) * (x) * (x) + 270.0) * (x)))

◆ LOBATTO_PHI7

#define LOBATTO_PHI7 ( x)
Value:
(-1.0 / 32.0 * 2.91547594742265023543707643877 * \
(((715 * (x) * (x)-1001) * (x) * (x) + 385) * (x) * (x)-35) * (x))

Definition at line 141 of file base_functions.h.

141#define LOBATTO_PHI7(x) \
142 (-1.0 / 32.0 * 2.91547594742265023543707643877 * \
143 (((715 * (x) * (x)-1001) * (x) * (x) + 385) * (x) * (x)-35) * (x))

◆ LOBATTO_PHI7X

#define LOBATTO_PHI7X ( x)
Value:
(-1.0 / 32.0 * 2.91547594742265023543707643877 * \
(((5005.0 * (x) * (x)-5005.0) * (x) * (x) + 1155.0) * (x) * (x)-35.0))

Definition at line 172 of file base_functions.h.

172#define LOBATTO_PHI7X(x) \
173 (-1.0 / 32.0 * 2.91547594742265023543707643877 * \
174 (((5005.0 * (x) * (x)-5005.0) * (x) * (x) + 1155.0) * (x) * (x)-35.0))

◆ LOBATTO_PHI8

#define LOBATTO_PHI8 ( x)
Value:
(-1.0 / 64.0 * 3.08220700148448822512509619073 * \
((((2431 * (x) * (x)-4004) * (x) * (x) + 2002) * (x) * (x)-308) * (x) * \
(x) + \
7))

Definition at line 144 of file base_functions.h.

144#define LOBATTO_PHI8(x) \
145 (-1.0 / 64.0 * 3.08220700148448822512509619073 * \
146 ((((2431 * (x) * (x)-4004) * (x) * (x) + 2002) * (x) * (x)-308) * (x) * \
147 (x) + \
148 7))

◆ LOBATTO_PHI8X

#define LOBATTO_PHI8X ( x)
Value:
(-1.0 / 64.0 * 3.08220700148448822512509619073 * \
((((19448.0 * (x) * (x)-24024.0) * (x) * (x) + 8008.0) * (x) * (x)-616.0) * \
(x)))

Definition at line 175 of file base_functions.h.

175#define LOBATTO_PHI8X(x) \
176 (-1.0 / 64.0 * 3.08220700148448822512509619073 * \
177 ((((19448.0 * (x) * (x)-24024.0) * (x) * (x) + 8008.0) * (x) * (x)-616.0) * \
178 (x)))

◆ LOBATTO_PHI9

#define LOBATTO_PHI9 ( x)
Value:
(-1.0 / 128.0 * 6.4807406984078603784382721642 * \
((((4199 * (x) * (x)-7956) * (x) * (x) + 4914) * (x) * (x)-1092) * (x) * \
(x) + \
63) * \
(x))

Definition at line 149 of file base_functions.h.

149#define LOBATTO_PHI9(x) \
150 (-1.0 / 128.0 * 6.4807406984078603784382721642 * \
151 ((((4199 * (x) * (x)-7956) * (x) * (x) + 4914) * (x) * (x)-1092) * (x) * \
152 (x) + \
153 63) * \
154 (x))

◆ LOBATTO_PHI9X

#define LOBATTO_PHI9X ( x)
Value:
(-1.0 / 128.0 * 6.4807406984078603784382721642 * \
((((37791.0 * (x) * (x)-55692.0) * (x) * (x) + 24570.0) * (x) * \
(x)-3276.0) * \
(x) * (x)-63.0))

Definition at line 179 of file base_functions.h.

179#define LOBATTO_PHI9X(x) \
180 (-1.0 / 128.0 * 6.4807406984078603784382721642 * \
181 ((((37791.0 * (x) * (x)-55692.0) * (x) * (x) + 24570.0) * (x) * \
182 (x)-3276.0) * \
183 (x) * (x)-63.0))

Function Documentation

◆ IntegratedJacobi_polynomials()

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 [29]

Parameters
pis approximation order
alphapolynomial parameter
xis position \(s\in[0,t]\)
trange of polynomial
diff_xderivatives of shape functions, i.e. \(\frac{\partial x}{\partial \xi_i}\)
diff_tderivatives of shape functions, i.e. \(\frac{\partial t}{\partial \xi_i}\)
Return values
Lapproximation functions
diffLderivatives, i.e. \(\frac{\partial L}{\partial \xi_i}\)
Parameters
dimdimension
Returns
error code

Definition at line 134 of file base_functions.c.

137 {
139#ifndef NDEBUG
140 if (dim < 1)
141 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim < 1");
142 if (dim > 3)
143 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim > 3");
144 if (p < 1)
145 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "p < 1");
146#endif // NDEBUG
147 L[0] = x;
148 if (diffL != NULL) {
149 int d = 0;
150 for (; d != dim; ++d) {
151 diffL[d * p + 0] = diff_x[d];
152 }
153 }
154 if (p == 0)
156 double jacobi[(p + 1)];
157 double diff_jacobi[(p + 1) * dim];
158 ierr = Jacobi_polynomials(p, alpha, x, t, diff_x, diff_t, jacobi, diff_jacobi,
159 dim);
160 CHKERRQ(ierr);
161 int l = 1;
162 for (; l < p; l++) {
163 int i = l + 1;
164 const double a = (i + alpha) / ((2 * i + alpha - 1) * (2 * i + alpha));
165 const double b = alpha / ((2 * i + alpha - 2) * (2 * i + alpha));
166 const double c = (i - 1) / ((2 * i + alpha - 2) * (2 * i + alpha - 1));
167 L[l] = a * jacobi[i] + b * t * jacobi[i - 1] - c * t * t * jacobi[i - 2];
168 if (diffL != NULL) {
169 int d = 0;
170 for (; d != dim; ++d) {
171 diffL[d * p + l] = a * diff_jacobi[d * (p + 1) + i] +
172 b * (t * diff_jacobi[d * (p + 1) + i - 1] +
173 diff_t[d] * jacobi[i - 1]) -
174 c * (t * t * diff_jacobi[d * (p + 1) + i - 2] +
175 2 * t * diff_t[d] * jacobi[i - 2]);
176 }
177 }
178 }
180}
constexpr double a
static PetscErrorCode ierr
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.
#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< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
FTensor::Index< 'l', 3 > l
constexpr double t
plate stiffness
Definition plate.cpp:58

◆ Jacobi_polynomials()

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 [29]

Parameters
pis approximation order
alphapolynomial parameter
xis position \(s\in[0,t]\)
trange of polynomial
diff_xderivatives of shape functions, i.e. \(\frac{\partial x}{\partial \xi_i}\)
diff_tderivatives of shape functions, i.e. \(\frac{\partial t}{\partial \xi_i}\)
Return values
Lapproximation functions
diffLderivatives, i.e. \(\frac{\partial L}{\partial \xi_i}\)
Parameters
dimdimension
Returns
error code

Definition at line 67 of file base_functions.c.

69 {
71#ifndef NDEBUG
72 if (dim < 1)
73 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim < 1");
74 if (dim > 3)
75 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim > 3");
76 if (p < 0)
77 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "p < 0");
78
79 if (diffL != NULL) {
80 if (diff_x == NULL) {
81 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "diff_s == NULL");
82 }
83 }
84
85#endif // NDEBUG
86 L[0] = 1;
87 if (diffL != NULL) {
88 diffL[0 * (p + 1) + 0] = 0;
89 if (dim >= 2) {
90 diffL[1 * (p + 1) + 0] = 0;
91 if (dim == 3) {
92 diffL[2 * (p + 1) + 0] = 0;
93 }
94 }
95 }
96 if (p == 0)
98 L[1] = 2 * x - t + alpha * x;
99 if (diffL != NULL) {
100 int d = 0;
101 for (; d < dim; ++d) {
102 double d_t = (diff_t) ? diff_t[d] : 0;
103 diffL[d * (p + 1) + 1] = (2 + alpha) * diff_x[d] - d_t;
104 }
105 }
106 if (p == 1)
108 int l = 1;
109 for (; l < p; l++) {
110 int lp1 = l + 1;
111 double a = 2 * lp1 * (lp1 + alpha) * (2 * lp1 + alpha - 2);
112 double b = 2 * lp1 + alpha - 1;
113 double c = (2 * lp1 + alpha) * (2 * lp1 + alpha - 2);
114 double d = 2 * (lp1 + alpha - 1) * (lp1 - 1) * (2 * lp1 + alpha);
115 double A = b * (c * (2 * x - t) + alpha * alpha * t) / a;
116 double B = d * t * t / a;
117 L[lp1] = A * L[l] - B * L[l - 1];
118 if (diffL != NULL) {
119 int z = 0;
120 for (; z < dim; ++z) {
121 double d_t = (diff_t) ? diff_t[z] : 0;
122 double diffA =
123 b * (c * (2 * diff_x[z] - d_t) + alpha * alpha * d_t) / a;
124 double diffB = d * 2 * t * d_t / a;
125 diffL[z * (p + 1) + lp1] = A * diffL[z * (p + 1) + l] -
126 B * diffL[z * (p + 1) + l - 1] +
127 diffA * L[l] - diffB * L[l - 1];
128 }
129 }
130 }
132}