v0.14.0
Loading...
Searching...
No Matches
Macros | Functions
base_functions.h File Reference

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

Macro Definition Documentation

◆ LOBATTO_PHI0

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

◆ LOBATTO_PHI0X

#define LOBATTO_PHI0X (   x)    (0)

Derivatives of kernel functions for Lobbatto base.

Definition at line 157 of file base_functions.h.

◆ LOBATTO_PHI1

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

Definition at line 127 of file base_functions.h.

◆ LOBATTO_PHI1X

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

Definition at line 158 of file base_functions.h.

◆ LOBATTO_PHI2

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

Definition at line 128 of file base_functions.h.

◆ LOBATTO_PHI2X

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

Definition at line 159 of file base_functions.h.

◆ LOBATTO_PHI3

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

Definition at line 130 of file base_functions.h.

◆ LOBATTO_PHI3X

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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 151 of file base_functions.c.

154 {
156#ifndef NDEBUG
157 if (dim < 1)
158 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim < 1");
159 if (dim > 3)
160 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim > 3");
161 if (p < 1)
162 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "p < 1");
163#endif // NDEBUG
164 L[0] = x;
165 if (diffL != NULL) {
166 int d = 0;
167 for (; d != dim; ++d) {
168 diffL[d * p + 0] = diff_x[d];
169 }
170 }
171 if (p == 0)
173 double jacobi[(p + 1)];
174 double diff_jacobi[(p + 1) * dim];
175 ierr = Jacobi_polynomials(p, alpha, x, t, diff_x, diff_t, jacobi, diff_jacobi,
176 dim);
177 CHKERRQ(ierr);
178 int l = 1;
179 for (; l < p; l++) {
180 int i = l + 1;
181 const double a = (i + alpha) / ((2 * i + alpha - 1) * (2 * i + alpha));
182 const double b = alpha / ((2 * i + alpha - 2) * (2 * i + alpha));
183 const double c = (i - 1) / ((2 * i + alpha - 2) * (2 * i + alpha - 1));
184 L[l] = a * jacobi[i] + b * t * jacobi[i - 1] - c * t * t * jacobi[i - 2];
185 if (diffL != NULL) {
186 int dd = 0;
187 for (; dd != dim; ++dd) {
188 diffL[dd * p + l] = a * diff_jacobi[dd * (p + 1) + i] +
189 b * (t * diff_jacobi[dd * (p + 1) + i - 1] +
190 diff_t[dd] * jacobi[i - 1]) -
191 c * (t * t * diff_jacobi[dd * (p + 1) + i - 2] +
192 2 * t * diff_t[dd] * jacobi[i - 2]);
193 }
194 }
195 }
197}
static Index< 'L', 3 > L
static Index< 'p', 3 > p
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()
Definition: definitions.h:447
@ 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 ...
Definition: definitions.h:440
const int dim
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
FTensor::Index< 'l', 3 > l
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
constexpr double t
plate stiffness
Definition: plate.cpp:59

◆ 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 [26]

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#endif // NDEBUG
79 L[0] = 1;
80 if (diffL != NULL) {
81 diffL[0 * (p + 1) + 0] = 0;
82 if (dim >= 2) {
83 diffL[1 * (p + 1) + 0] = 0;
84 if (dim == 3) {
85 diffL[2 * (p + 1) + 0] = 0;
86 }
87 }
88 }
89 if (p == 0)
91 L[1] = 2 * x - t + alpha * x;
92 if (diffL != NULL) {
93#ifndef NDEBUG
94 if (diff_x == NULL) {
95 SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "diff_s == NULL");
96 }
97#endif // NDEBUG
98 double d_t = (diff_t) ? diff_t[0] : 0;
99 diffL[0 * (p + 1) + 1] = (2 + alpha) * diff_x[0] - d_t;
100 if (dim >= 2) {
101 double d_t = (diff_t) ? diff_t[1] : 0;
102 diffL[1 * (p + 1) + 1] = (2 + alpha) * diff_x[1] - d_t;
103 if (dim == 3) {
104 double d_t = (diff_t) ? diff_t[2] : 0;
105 diffL[2 * (p + 1) + 1] = (2 + alpha) * diff_x[2] - d_t;
106 }
107 }
108 }
109 if (p == 1)
111 int l = 1;
112 for (; l < p; l++) {
113 int lp1 = l + 1;
114 double a = 2 * lp1 * (lp1 + alpha) * (2 * lp1 + alpha - 2);
115 double b = 2 * lp1 + alpha - 1;
116 double c = (2 * lp1 + alpha) * (2 * lp1 + alpha - 2);
117 double d = 2 * (lp1 + alpha - 1) * (lp1 - 1) * (2 * lp1 + alpha);
118 double A = b * (c * (2 * x - t) + alpha * alpha * t) / a;
119 double B = d * t * t / a;
120 L[lp1] = A * L[l] - B * L[l - 1];
121 if (diffL != NULL) {
122 double d_t = (diff_t) ? diff_t[0] : 0;
123 double diffA = b * (c * (2 * diff_x[0] - d_t) + alpha * alpha * d_t) / a;
124 double diffB = d * 2 * t * d_t / a;
125 diffL[0 * (p + 1) + lp1] = A * diffL[0 * (p + 1) + l] -
126 B * diffL[0 * (p + 1) + l - 1] + diffA * L[l] -
127 diffB * L[l - 1];
128 if (dim >= 2) {
129 double d_t = (diff_t) ? diff_t[1] : 0;
130 double diffA =
131 b * (c * (2 * diff_x[1] - d_t) + alpha * alpha * d_t) / a;
132 double diffB = d * 2 * t * d_t / a;
133 diffL[1 * (p + 1) + lp1] = A * diffL[1 * (p + 1) + l] -
134 B * diffL[1 * (p + 1) + l - 1] +
135 diffA * L[l] - diffB * L[l - 1];
136 if (dim == 3) {
137 double d_t = (diff_t) ? diff_t[2] : 0;
138 double diffA =
139 b * (c * (2 * diff_x[2] - d_t) + alpha * alpha * d_t) / a;
140 double diffB = d * 2 * t * d_t / a;
141 diffL[2 * (p + 1) + lp1] = A * diffL[2 * (p + 1) + l] -
142 B * diffL[2 * (p + 1) + l - 1] +
143 diffA * L[l] - diffB * L[l - 1];
144 }
145 }
146 }
147 }
149}
constexpr AssemblyType A