v0.14.0
base_functions.c
Go to the documentation of this file.
1 /** \file base_functions.c
2 
3 */
4 
5 #include <cblas.h>
6 #include <petscsys.h>
7 #include <phg-quadrule/quad.h>
8 
9 #include <definitions.h>
10 
11 #include <base_functions.h>
12 
13 static PetscErrorCode ierr;
14 
15 PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L,
16  double *diffL, const int dim) {
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 }
66 
67 PetscErrorCode Jacobi_polynomials(int p, double alpha, double x, double t,
68  double *diff_x, double *diff_t, double *L,
69  double *diffL, const int dim) {
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 }
150 
151 PetscErrorCode IntegratedJacobi_polynomials(int p, double alpha, double x,
152  double t, double *diff_x,
153  double *diff_t, double *L,
154  double *diffL, const int dim) {
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 }
198 
199 PetscErrorCode Lobatto_polynomials(int p, double s, double *diff_s, double *L,
200  double *diffL, const int dim) {
201 
203 #ifndef NDEBUG
204  if (dim < 1)
205  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim < 1");
206  if (dim > 3)
207  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim > 3");
208  if (p < 2)
209  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "p < 2");
210 #endif // NDEBUG
211  double l[p + 1];
212  ierr = Legendre_polynomials(p, s, NULL, l, NULL, 1);
213  CHKERRQ(ierr);
214 
215  L[0] = 1;
216  if (diffL != NULL) {
217  for (int d = 0; d != dim; ++d) {
218  diffL[d * (p + 1) + 0] = 0;
219  }
220  }
221  L[1] = s;
222  if (diffL != NULL) {
223 #ifndef NDEBUG
224  if (diff_s == NULL) {
225  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "diff_s == NULL");
226  }
227 #endif // NDEBUG
228  for (int d = 0; d != dim; ++d) {
229  diffL[d * (p + 1) + 1] = diff_s[d];
230  }
231  }
232 
233  // Integrated Legendre
234  for (int k = 2; k <= p; k++) {
235  const double factor = 2 * (2 * k - 1);
236  L[k] = 1.0 / factor * (l[k] - l[k - 2]);
237  }
238 
239  if (diffL != NULL) {
240  for (int k = 2; k <= p; k++) {
241  double a = l[k - 1] / 2.;
242  for (int d = 0; d != dim; ++d) {
243  diffL[d * (p + 1) + k] = a * diff_s[d];
244  }
245  }
246  }
248 }
249 
250 static double f_phi0(double x) { return LOBATTO_PHI0(x); }
251 static double f_phi1(double x) { return LOBATTO_PHI1(x); }
252 static double f_phi2(double x) { return LOBATTO_PHI2(x); }
253 static double f_phi3(double x) { return LOBATTO_PHI3(x); }
254 static double f_phi4(double x) { return LOBATTO_PHI4(x); }
255 static double f_phi5(double x) { return LOBATTO_PHI5(x); }
256 static double f_phi6(double x) { return LOBATTO_PHI6(x); }
257 static double f_phi7(double x) { return LOBATTO_PHI7(x); }
258 static double f_phi8(double x) { return LOBATTO_PHI8(x); }
259 static double f_phi9(double x) { return LOBATTO_PHI9(x); }
260 
261 static double (*f_phi[])(double x) = {f_phi0, f_phi1, f_phi2, f_phi3, f_phi4,
263 
264 static double f_phi0x(double x) { return LOBATTO_PHI0X(x); }
265 static double f_phi1x(double x) { return LOBATTO_PHI1X(x); }
266 static double f_phi2x(double x) { return LOBATTO_PHI2X(x); }
267 static double f_phi3x(double x) { return LOBATTO_PHI3X(x); }
268 static double f_phi4x(double x) { return LOBATTO_PHI4X(x); }
269 static double f_phi5x(double x) { return LOBATTO_PHI5X(x); }
270 static double f_phi6x(double x) { return LOBATTO_PHI6X(x); }
271 static double f_phi7x(double x) { return LOBATTO_PHI7X(x); }
272 static double f_phi8x(double x) { return LOBATTO_PHI8X(x); }
273 static double f_phi9x(double x) { return LOBATTO_PHI9X(x); }
274 
275 static double (*f_phix[])(double x) = {f_phi0x, f_phi1x, f_phi2x, f_phi3x,
277  f_phi8x, f_phi9x};
278 
279 // /// Legendre polynomials
280 // #define Legendre0(x) (1.0)
281 // #define Legendre1(x) (x)
282 // #define Legendre2(x) (1.0 / 2.0 * (3 * (x) * (x) - 1))
283 // #define Legendre3(x) (1.0 / 2.0 * (5 * (x) * (x) - 3) * (x))
284 // #define Legendre4(x) (1.0 / 8.0 * ((35 * (x) * (x) - 30) * (x) * (x) + 3))
285 // #define Legendre5(x) (1.0 / 8.0 * ((63 * (x) * (x) - 70) * (x) * (x) + 15) *
286 // (x)) #define Legendre6(x) (1.0 / 16.0 * (((231 * (x) * (x) - 315) * (x) * (x)
287 // + 105) * (x) * (x) - 5)) #define Legendre7(x) (1.0 / 16.0 * (((429 * (x) *
288 // (x) - 693) * (x) * (x) + 315) * (x) * (x) - 35) * (x)) #define Legendre8(x)
289 // (1.0 / 128.0 * ((((6435 * (x) * (x) - 12012) * (x) * (x) + 6930) * (x) * (x)
290 // - 1260) * (x) * (x) + 35)) #define Legendre9(x) (1.0 / 128.0 * ((((12155 *
291 // (x) * (x) - 25740) * (x) * (x) + 18018) * (x) * (x) - 4620) * (x) * (x) +
292 // 315) * (x)) #define Legendre10(x) (1.0 / 256.0 * (((((46189 * (x) * (x) -
293 // 109395) * (x) * (x) + 90090) * (x) * (x) - 30030) * (x) * (x) + 3465) * (x) *
294 // (x) - 63))
295 //
296 // /// derivatives of Legendre polynomials
297 // #define Legendre0x(x) (0.0)
298 // #define Legendre1x(x) (1.0)
299 // #define Legendre2x(x) (3.0 * (x))
300 // #define Legendre3x(x) (15.0 / 2.0 * (x) * (x) - 3.0 / 2.0)
301 // #define Legendre4x(x) (5.0 / 2.0 * (x) * (7.0 * (x) * (x) - 3.0))
302 // #define Legendre5x(x) ((315.0 / 8.0 * (x) * (x) - 105.0 / 4.0) * (x) * (x)
303 // + 15.0 / 8.0) #define Legendre6x(x) (21.0 / 8.0 * (x) * ((33.0 * (x) * (x)
304 // - 30.0) * (x) * (x) + 5.0)) #define Legendre7x(x) (((3003.0 / 16.0 * (x) *
305 // (x) - 3465.0 / 16.0) * (x) * (x) + 945.0 / 16.0) * (x) * (x) - 35.0 / 16.0)
306 // #define Legendre8x(x) (9.0 / 16.0 * (x) * (((715.0 * (x) * (x) - 1001.0) *
307 // (x) * (x) + 385.0) * (x) * (x) - 35.0)) #define Legendre9x(x) ((((109395.0 /
308 // 128.0 * (x) * (x) - 45045.0 / 32.0) * (x) * (x) + 45045.0 / 64.0) * (x) * (x)
309 // - 3465.0 / 32.0) * (x) * (x) + 315.0 / 128.0) #define Legendre10x(x) (2.0 /
310 // 256.0 * (x) * ((((230945.0 * (x) * (x) - 437580.0) * (x) * (x) + 270270.0) *
311 // (x) * (x) - 60060.0) * (x) * (x) + 3465.0))
312 //
313 // /// first two Lobatto shape functions
314 // #define l0(x) ((1.0 - (x)) * 0.5)
315 // #define l1(x) ((1.0 + (x)) * 0.5)
316 //
317 // #define l0l1(x) ((1.0 - (x)*(x)) * 0.25)
318 //
319 // /// other Lobatto shape functions
320 // #define l2(x) (phi0(x) * l0l1(x))
321 // #define l3(x) (phi1(x) * l0l1(x))
322 // #define l4(x) (phi2(x) * l0l1(x))
323 // #define l5(x) (phi3(x) * l0l1(x))
324 // #define l6(x) (phi4(x) * l0l1(x))
325 // #define l7(x) (phi5(x) * l0l1(x))
326 // #define l8(x) (phi6(x) * l0l1(x))
327 // #define l9(x) (phi7(x) * l0l1(x))
328 // #define l10(x) (phi8(x) * l0l1(x))
329 // #define l11(x) (phi9(x) * l0l1(x))
330 //
331 // /// derivatives of Lobatto functions
332 // #define dl0(x) (-0.5)
333 // #define dl1(x) (0.5)
334 // #define dl2(x) (sqrt(3.0/2.0) * Legendre1(x))
335 // #define dl3(x) (sqrt(5.0/2.0) * Legendre2(x))
336 // #define dl4(x) (sqrt(7.0/2.0) * Legendre3(x))
337 // #define dl5(x) (sqrt(9.0/2.0) * Legendre4(x))
338 // #define dl6(x) (sqrt(11.0/2.0) * Legendre5(x))
339 // #define dl7(x) (sqrt(13.0/2.0) * Legendre6(x))
340 // #define dl8(x) (sqrt(15.0/2.0) * Legendre7(x))
341 // #define dl9(x) (sqrt(17.0/2.0) * Legendre8(x))
342 // #define dl10(x) (sqrt(19.0/2.0) * Legendre9(x))
343 // #define dl11(x) (sqrt(21.0/2.0) * Legendre10(x))
344 
345 PetscErrorCode LobattoKernel_polynomials(int p, double s, double *diff_s,
346  double *L, double *diffL,
347  const int dim) {
349 #ifndef NDEBUG
350  if (dim < 1)
351  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim < 1");
352  if (dim > 3)
353  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "dim > 3");
354  if (p < 0)
355  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "p < 0");
356  if (p > 9)
357  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
358  "Polynomial beyond order 9 is not implemented");
359 #endif // NDEBUG
360  if (L) {
361  int l = 0;
362  for (; l != p + 1; l++) {
363  L[l] = f_phi[l](s);
364  }
365  }
366  if (diffL != NULL) {
367 #ifndef NDEBUG
368  if (diff_s == NULL) {
369  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA, "diff_s == NULL");
370  }
371 #endif // NDEBUG
372  int l = 0;
373  for (; l != p + 1; l++) {
374  double a = f_phix[l](s);
375  diffL[0 * (p + 1) + l] = diff_s[0] * a;
376  if (dim >= 2) {
377  diffL[1 * (p + 1) + l] = diff_s[1] * a;
378  if (dim == 3) {
379  diffL[2 * (p + 1) + l] = diff_s[2] * a;
380  }
381  }
382  }
383  }
385 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
LOBATTO_PHI4
#define LOBATTO_PHI4(x)
Definition: base_functions.h:132
f_phi5x
static double f_phi5x(double x)
Definition: base_functions.c:269
sdf_hertz.d
float d
Definition: sdf_hertz.py:5
LOBATTO_PHI1X
#define LOBATTO_PHI1X(x)
Definition: base_functions.h:158
f_phi7
static double f_phi7(double x)
Definition: base_functions.c:257
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.
Definition: base_functions.c:151
LOBATTO_PHI7
#define LOBATTO_PHI7(x)
Definition: base_functions.h:141
LOBATTO_PHI9X
#define LOBATTO_PHI9X(x)
Definition: base_functions.h:179
f_phi7x
static double f_phi7x(double x)
Definition: base_functions.c:271
ierr
static PetscErrorCode ierr
Definition: base_functions.c:13
LOBATTO_PHI8
#define LOBATTO_PHI8(x)
Definition: base_functions.h:144
f_phi1x
static double f_phi1x(double x)
Definition: base_functions.c:265
f_phi4x
static double f_phi4x(double x)
Definition: base_functions.c:268
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
base_functions.h
LOBATTO_PHI6X
#define LOBATTO_PHI6X(x)
Definition: base_functions.h:169
LOBATTO_PHI1
#define LOBATTO_PHI1(x)
Definition: base_functions.h:127
LOBATTO_PHI2X
#define LOBATTO_PHI2X(x)
Definition: base_functions.h:159
f_phi8x
static double f_phi8x(double x)
Definition: base_functions.c:272
Lobatto_polynomials
PetscErrorCode Lobatto_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Lobatto base functions .
Definition: base_functions.c:199
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
f_phi3
static double f_phi3(double x)
Definition: base_functions.c:253
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.
Definition: base_functions.c:67
f_phi4
static double f_phi4(double x)
Definition: base_functions.c:254
f_phi
static double(* f_phi[])(double x)
Definition: base_functions.c:261
LOBATTO_PHI9
#define LOBATTO_PHI9(x)
Definition: base_functions.h:149
a
constexpr double a
Definition: approx_sphere.cpp:30
double
f_phi8
static double f_phi8(double x)
Definition: base_functions.c:258
f_phi2x
static double f_phi2x(double x)
Definition: base_functions.c:266
LOBATTO_PHI5X
#define LOBATTO_PHI5X(x)
Definition: base_functions.h:166
f_phi9
static double f_phi9(double x)
Definition: base_functions.c:259
LOBATTO_PHI8X
#define LOBATTO_PHI8X(x)
Definition: base_functions.h:175
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
f_phi9x
static double f_phi9x(double x)
Definition: base_functions.c:273
f_phi0
static double f_phi0(double x)
Definition: base_functions.c:250
Legendre_polynomials
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Legendre approximation basis.
Definition: base_functions.c:15
definitions.h
useful compiler directives and definitions
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
f_phi6x
static double f_phi6x(double x)
Definition: base_functions.c:270
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
f_phi3x
static double f_phi3x(double x)
Definition: base_functions.c:267
LobattoKernel_polynomials
PetscErrorCode LobattoKernel_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Kernel Lobatto base functions.
Definition: base_functions.c:345
f_phi6
static double f_phi6(double x)
Definition: base_functions.c:256
LOBATTO_PHI0
#define LOBATTO_PHI0(x)
Definitions taken from Hermes2d code.
Definition: base_functions.h:126
FTensor::dd
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
LOBATTO_PHI3
#define LOBATTO_PHI3(x)
Definition: base_functions.h:130
f_phi1
static double f_phi1(double x)
Definition: base_functions.c:251
LOBATTO_PHI4X
#define LOBATTO_PHI4X(x)
Definition: base_functions.h:163
f_phi0x
static double f_phi0x(double x)
Definition: base_functions.c:264
LOBATTO_PHI2
#define LOBATTO_PHI2(x)
Definition: base_functions.h:128
LOBATTO_PHI7X
#define LOBATTO_PHI7X(x)
Definition: base_functions.h:172
f_phi5
static double f_phi5(double x)
Definition: base_functions.c:255
f_phi2
static double f_phi2(double x)
Definition: base_functions.c:252
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
LOBATTO_PHI5
#define LOBATTO_PHI5(x)
Definition: base_functions.h:135
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
f_phix
static double(* f_phix[])(double x)
Definition: base_functions.c:275
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
LOBATTO_PHI6
#define LOBATTO_PHI6(x)
Definition: base_functions.h:138
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
LOBATTO_PHI3X
#define LOBATTO_PHI3X(x)
Definition: base_functions.h:161
LOBATTO_PHI0X
#define LOBATTO_PHI0X(x)
Derivatives of kernel functions for Lobbatto base.
Definition: base_functions.h:157
quad.h