v0.15.0
Loading...
Searching...
No Matches
base_functions.h
Go to the documentation of this file.
1/** \file base_functions.c
2
3*/
4
5#ifndef __BASE_FUNCTIONS_H__
6#define __BASE_FUNCTIONS_H__
7
8#ifdef __cplusplus
9extern "C" {
10#endif
11
12// Note: Implementation of Legendre_polynomials is moved to LegendrePolynomial.hpp
13
14/**
15\brief Calculate Jacobi approximation basis
16
17For more details see \cite fuentes2015orientation
18
19\param p is approximation order
20\param alpha polynomial parameter
21\param x is position \f$s\in[0,t]\f$
22\param t range of polynomial
23\param diff_x derivatives of shape functions, i.e. \f$\frac{\partial x}{\partial
24\xi_i}\f$ \param diff_t derivatives of shape functions, i.e. \f$\frac{\partial
25t}{\partial \xi_i}\f$ \retval L approximation functions \retval diffL
26derivatives, i.e. \f$\frac{\partial L}{\partial \xi_i}\f$ \param dim dimension
27\return error code
28
29*/
30PetscErrorCode Jacobi_polynomials(int p, double alpha, double x, double t,
31 double *diff_x, double *diff_t, double *L,
32 double *diffL, const int dim);
33
34/**
35\brief Calculate integrated Jacobi approximation basis
36
37For more details see \cite fuentes2015orientation
38
39\param p is approximation order
40\param alpha polynomial parameter
41\param x is position \f$s\in[0,t]\f$
42\param t range of polynomial
43\param diff_x derivatives of shape functions, i.e. \f$\frac{\partial x}{\partial
44\xi_i}\f$ \param diff_t derivatives of shape functions, i.e. \f$\frac{\partial
45t}{\partial \xi_i}\f$ \retval L approximation functions \retval diffL
46derivatives, i.e. \f$\frac{\partial L}{\partial \xi_i}\f$ \param dim dimension
47\return error code
48
49*/
50PetscErrorCode IntegratedJacobi_polynomials(int p, double alpha, double x,
51 double t, double *diff_x,
52 double *diff_t, double *L,
53 double *diffL, const int dim);
54
55// Note: Implementation of Lobatto_polynomials is moved to LobattoPolynomial.hpp
56
57/**
58 \brief Calculate Kernel Lobatto base functions.
59
60 \ingroup mofem_base_functions
61
62 This is implemented using definitions from Hermes2d
63 <https://github.com/hpfem/hermes> following book by Pavel Solin et al \cite
64 solin2003higher.
65
66 \param p is approximation order
67 \param s is position \f$s\in[-1,1]\f$
68 \param diff_s derivatives of shape functions, i.e. \f$\frac{\partial
69 s}{\partial \xi_i}\f$ \retval L approximation functions \retval diffL
70 derivatives, i.e. \f$\frac{\partial L}{\partial \xi_i}\f$ \param dim dimension
71 \return error code
72
73*/
74PetscErrorCode LobattoKernel_polynomials(int p, double s, double *diff_s,
75 double *L, double *diffL,
76 const int dim);
77
78/// Definitions taken from Hermes2d code
79
80/// kernel functions for Lobatto base
81#define LOBATTO_PHI0(x) (-2.0 * 1.22474487139158904909864203735)
82#define LOBATTO_PHI1(x) (-2.0 * 1.58113883008418966599944677222 * (x))
83#define LOBATTO_PHI2(x) \
84 (-1.0 / 2.0 * 1.87082869338697069279187436616 * (5 * (x) * (x)-1))
85#define LOBATTO_PHI3(x) \
86 (-1.0 / 2.0 * 2.12132034355964257320253308631 * (7 * (x) * (x)-3) * (x))
87#define LOBATTO_PHI4(x) \
88 (-1.0 / 4.0 * 2.34520787991171477728281505677 * \
89 (21 * (x) * (x) * (x) * (x)-14 * (x) * (x) + 1))
90#define LOBATTO_PHI5(x) \
91 (-1.0 / 4.0 * 2.54950975679639241501411205451 * \
92 ((33 * (x) * (x)-30) * (x) * (x) + 5) * (x))
93#define LOBATTO_PHI6(x) \
94 (-1.0 / 32.0 * 2.73861278752583056728484891400 * \
95 (((429 * (x) * (x)-495) * (x) * (x) + 135) * (x) * (x)-5))
96#define LOBATTO_PHI7(x) \
97 (-1.0 / 32.0 * 2.91547594742265023543707643877 * \
98 (((715 * (x) * (x)-1001) * (x) * (x) + 385) * (x) * (x)-35) * (x))
99#define LOBATTO_PHI8(x) \
100 (-1.0 / 64.0 * 3.08220700148448822512509619073 * \
101 ((((2431 * (x) * (x)-4004) * (x) * (x) + 2002) * (x) * (x)-308) * (x) * \
102 (x) + \
103 7))
104#define LOBATTO_PHI9(x) \
105 (-1.0 / 128.0 * 6.4807406984078603784382721642 * \
106 ((((4199 * (x) * (x)-7956) * (x) * (x) + 4914) * (x) * (x)-1092) * (x) * \
107 (x) + \
108 63) * \
109 (x))
110
111/// Derivatives of kernel functions for Lobbatto base
112#define LOBATTO_PHI0X(x) (0)
113#define LOBATTO_PHI1X(x) (-2.0 * 1.58113883008418966599944677222)
114#define LOBATTO_PHI2X(x) \
115 (-1.0 / 2.0 * 1.87082869338697069279187436616 * (10 * (x)))
116#define LOBATTO_PHI3X(x) \
117 (-1.0 / 2.0 * 2.12132034355964257320253308631 * (21.0 * (x) * (x)-3.0))
118#define LOBATTO_PHI4X(x) \
119 (-1.0 / 4.0 * 2.34520787991171477728281505677 * \
120 ((84.0 * (x) * (x)-28.0) * (x)))
121#define LOBATTO_PHI5X(x) \
122 (-1.0 / 4.0 * 2.54950975679639241501411205451 * \
123 ((165.0 * (x) * (x)-90.0) * (x) * (x) + 5.0))
124#define LOBATTO_PHI6X(x) \
125 (-1.0 / 32.0 * 2.73861278752583056728484891400 * \
126 (((2574.0 * (x) * (x)-1980.0) * (x) * (x) + 270.0) * (x)))
127#define LOBATTO_PHI7X(x) \
128 (-1.0 / 32.0 * 2.91547594742265023543707643877 * \
129 (((5005.0 * (x) * (x)-5005.0) * (x) * (x) + 1155.0) * (x) * (x)-35.0))
130#define LOBATTO_PHI8X(x) \
131 (-1.0 / 64.0 * 3.08220700148448822512509619073 * \
132 ((((19448.0 * (x) * (x)-24024.0) * (x) * (x) + 8008.0) * (x) * (x)-616.0) * \
133 (x)))
134#define LOBATTO_PHI9X(x) \
135 (-1.0 / 128.0 * 6.4807406984078603784382721642 * \
136 ((((37791.0 * (x) * (x)-55692.0) * (x) * (x) + 24570.0) * (x) * \
137 (x)-3276.0) * \
138 (x) * (x)-63.0))
139
140#ifdef __cplusplus
141}
142#endif
143
144#endif //__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.
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 LobattoKernel_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Kernel Lobatto base functions.
constexpr double t
plate stiffness
Definition plate.cpp:58