v0.14.0
Functions | Variables
base_functions.c File Reference
#include <cblas.h>
#include <petscsys.h>
#include <phg-quadrule/quad.h>
#include <definitions.h>
#include <base_functions.h>

Go to the source code of this file.

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 [29]. More...
 
static double f_phi0 (double x)
 
static double f_phi1 (double x)
 
static double f_phi2 (double x)
 
static double f_phi3 (double x)
 
static double f_phi4 (double x)
 
static double f_phi5 (double x)
 
static double f_phi6 (double x)
 
static double f_phi7 (double x)
 
static double f_phi8 (double x)
 
static double f_phi9 (double x)
 
static double f_phi0x (double x)
 
static double f_phi1x (double x)
 
static double f_phi2x (double x)
 
static double f_phi3x (double x)
 
static double f_phi4x (double x)
 
static double f_phi5x (double x)
 
static double f_phi6x (double x)
 
static double f_phi7x (double x)
 
static double f_phi8x (double x)
 
static double f_phi9x (double x)
 
PetscErrorCode LobattoKernel_polynomials (int p, double s, double *diff_s, double *L, double *diffL, const int dim)
 Calculate Kernel Lobatto base functions. More...
 

Variables

static PetscErrorCode ierr
 
static double(* f_phi [])(double x)
 
static double(* f_phix [])(double x)
 

Function Documentation

◆ f_phi0()

static double f_phi0 ( double  x)
static

Definition at line 250 of file base_functions.c.

250 { return LOBATTO_PHI0(x); }

◆ f_phi0x()

static double f_phi0x ( double  x)
static

Definition at line 264 of file base_functions.c.

264 { return LOBATTO_PHI0X(x); }

◆ f_phi1()

static double f_phi1 ( double  x)
static

Definition at line 251 of file base_functions.c.

251 { return LOBATTO_PHI1(x); }

◆ f_phi1x()

static double f_phi1x ( double  x)
static

Definition at line 265 of file base_functions.c.

265 { return LOBATTO_PHI1X(x); }

◆ f_phi2()

static double f_phi2 ( double  x)
static

Definition at line 252 of file base_functions.c.

252 { return LOBATTO_PHI2(x); }

◆ f_phi2x()

static double f_phi2x ( double  x)
static

Definition at line 266 of file base_functions.c.

266 { return LOBATTO_PHI2X(x); }

◆ f_phi3()

static double f_phi3 ( double  x)
static

Definition at line 253 of file base_functions.c.

253 { return LOBATTO_PHI3(x); }

◆ f_phi3x()

static double f_phi3x ( double  x)
static

Definition at line 267 of file base_functions.c.

267 { return LOBATTO_PHI3X(x); }

◆ f_phi4()

static double f_phi4 ( double  x)
static

Definition at line 254 of file base_functions.c.

254 { return LOBATTO_PHI4(x); }

◆ f_phi4x()

static double f_phi4x ( double  x)
static

Definition at line 268 of file base_functions.c.

268 { return LOBATTO_PHI4X(x); }

◆ f_phi5()

static double f_phi5 ( double  x)
static

Definition at line 255 of file base_functions.c.

255 { return LOBATTO_PHI5(x); }

◆ f_phi5x()

static double f_phi5x ( double  x)
static

Definition at line 269 of file base_functions.c.

269 { return LOBATTO_PHI5X(x); }

◆ f_phi6()

static double f_phi6 ( double  x)
static

Definition at line 256 of file base_functions.c.

256 { return LOBATTO_PHI6(x); }

◆ f_phi6x()

static double f_phi6x ( double  x)
static

Definition at line 270 of file base_functions.c.

270 { return LOBATTO_PHI6X(x); }

◆ f_phi7()

static double f_phi7 ( double  x)
static

Definition at line 257 of file base_functions.c.

257 { return LOBATTO_PHI7(x); }

◆ f_phi7x()

static double f_phi7x ( double  x)
static

Definition at line 271 of file base_functions.c.

271 { return LOBATTO_PHI7X(x); }

◆ f_phi8()

static double f_phi8 ( double  x)
static

Definition at line 258 of file base_functions.c.

258 { return LOBATTO_PHI8(x); }

◆ f_phi8x()

static double f_phi8x ( double  x)
static

Definition at line 272 of file base_functions.c.

272 { return LOBATTO_PHI8X(x); }

◆ f_phi9()

static double f_phi9 ( double  x)
static

Definition at line 259 of file base_functions.c.

259 { return LOBATTO_PHI9(x); }

◆ f_phi9x()

static double f_phi9x ( double  x)
static

Definition at line 273 of file base_functions.c.

273 { return LOBATTO_PHI9X(x); }

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

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 }

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

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 }

Variable Documentation

◆ f_phi

double(* f_phi[])(double x)
static
Initial value:

Definition at line 261 of file base_functions.c.

◆ f_phix

double(* f_phix[])(double x)
static
Initial value:

Definition at line 275 of file base_functions.c.

◆ ierr

PetscErrorCode ierr
static

Definition at line 13 of file base_functions.c.

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
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
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
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
LOBATTO_PHI9
#define LOBATTO_PHI9(x)
Definition: base_functions.h:149
a
constexpr double a
Definition: approx_sphere.cpp:30
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
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
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
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