v0.14.0
Loading...
Searching...
No Matches
h1_hdiv_hcurl_l2.h
Go to the documentation of this file.
1/** \file h1_hdiv_hcurl_l2.h
2\brief Functions to approximate hierarchical spaces
3
4\FIXME: Name Shape Functions is used, in that context is more appropriate
5to use base functions. Need to be changed.
6
7*/
8
9#ifndef __H1_HDIV_HCURL_L2_H__
10#define __H1_HDIV_HCURL_L2_H__
11
12#ifndef DEPRECATED
13#define DEPRECATED
14#endif
15
16#ifdef __cplusplus
17extern "C" {
18#endif
19
20// L2
21
22// Number of dofs for L2 space
23
24/**
25 * @brief Number of base functions on tetrahedron for L2 space
26 */
27#define NBVOLUMETET_L2(P) ((P + 1) * (P + 2) * (P + 3) / 6)
28
29/**
30 * @brief Number of base functions on hexahedron for L2 space
31 */
32#define NBVOLUMEHEX_L2_GENERAL(P, Q, R) ((P + 1) * (Q + 1) * (R + 1))
33
34/**
35 * @brief Number of base functions on hexahedron for L2 space
36 */
37#define NBVOLUMEHEX_L2(P) (NBVOLUMEHEX_L2_GENERAL(P, P, P))
38
39/**
40 * @brief Number of base functions on triangle for L2 space
41 */
42#define NBFACETRI_L2(P) ((P + 1) * (P + 2) / 2)
43
44/**
45 * @brief Number of base functions on edge from L2 space
46 *
47 */
48#define NBEDGE_L2(P) ((P) + 1)
49
50// H1
51
52/**
53 * @brief Number of base function on edge for H1 space
54 */
55#define NBEDGE_H1(P) (((P) > 1) ? (P - 1) : 0)
56
57/**
58 * @brief Number of base function on triangle for H1 space
59 */
60#define NBFACETRI_H1(P) (((P) > 2) ? ((P - 2) * (P - 1) / 2) : 0)
61
62/**
63 * @brief Number of base functions on quad for H1 space
64 */
65#define NBFACEQUAD_H1(P) (((P) > 1) ? ((P - 1) * (P - 1)) : 0)
66
67/**
68 * @brief Number of base functions on quad for L2 space
69 */
70#define NBFACEQUAD_L2(P) (((P) >= 0) ? (P + 1) * (P + 1) : 0)
71
72/**
73 * @brief Number of base functions on tetrahedron for H1 space
74 */
75#define NBVOLUMETET_H1(P) (((P) > 3) ? ((P - 3) * (P - 2) * (P - 1) / 6) : 0)
76
77/**
78 * @brief Number of base functions on prism for H1 space
79 */
80#define NBVOLUMEPRISM_H1(P) ((P > 3) ? ((P - 2) * (P - 2) * (P - 2)) : 0)
81
82/**
83 * @brief Number of base functions on hex for H1 space
84 *
85 */
86#define NBVOLUMEHEX_H1_GENERAL(P, Q, R) \
87 ((((P) > 1) && ((Q) > 1) && ((R) > 1)) ? (((P)-1) * ((Q)-1) * ((R)-1)) : 0)
88
89/**
90 * @brief Number of base functions on hex for H1 space
91 *
92 */
93#define NBVOLUMEHEX_H1(P) (NBVOLUMEHEX_H1_GENERAL(P, P, P))
94
95// H curl
96
97#define NBEDGE_AINSWORTH_HCURL(P) (((P) > 0) ? (P + 1) : 0)
98#define NBFACETRI_AINSWORTH_EDGE_HCURL(P) (((P) > 1) ? P - 1 : 0)
99#define NBFACETRI_AINSWORTH_FACE_HCURL(P) (((P) > 2) ? (P - 1) * (P - 2) : 0)
100#define NBFACETRI_AINSWORTH_HCURL(P) ((P > 1) ? ((P)-1) * (P + 1) : 0)
101#define NBVOLUMETET_AINSWORTH_FACE_HCURL(P) \
102 (((P) > 2) ? (2 * (P - 1) * (P - 2)) : 0)
103#define NBVOLUMETET_AINSWORTH_TET_HCURL(P) \
104 (((P) > 3) ? ((P - 3) * (P - 2) * (P - 1) / 2) : 0)
105#define NBVOLUMETET_AINSWORTH_HCURL(P) \
106 (((P) > 2) ? (P - 2) * (P - 1) * (P + 1) / 2 : 0)
107
108#define NBEDGE_DEMKOWICZ_HCURL(P) (((P) > 0) ? (P) : 0)
109#define NBFACETRI_DEMKOWICZ_HCURL(P) (((P) > 1) ? (P) * ((P)-1) : 0)
110#define NBVOLUMETET_DEMKOWICZ_HCURL(P) \
111 (((P) > 2) ? ((P) * ((P)-1) * ((P)-2) / 2) : 0)
112
113/**
114 * @brief Number of base functions on quad for Hcurl space
115 */
116#define NBFACEQUAD_DEMKOWICZ_FAMILY_HCURL(P, Q) \
117 (((P) > 0 && (Q) > 1) ? P * (Q - 1) : 0)
118#define NBFACEQUAD_DEMKOWICZ_HCURL(P) \
119 (2 * NBFACEQUAD_DEMKOWICZ_FAMILY_HCURL(P, P))
120
121#define NBVOLUMEHEX_DEMKOWICZ_FAMILY_HCURL(P, Q, R) \
122 ((P > 0) && (Q > 1) && (R > 1) ? ((P) * (Q - 1) * (R - 1)) : 0)
123
124#define NBVOLUMEHEX_DEMKOWICZ_HCURL(P) \
125 (3 * NBVOLUMEHEX_DEMKOWICZ_FAMILY_HCURL(P, P, P))
126
127// H div
128
129#define NBEDGE_HDIV(P) (0)
130#define NBFACETRI_AINSWORTH_EDGE_HDIV(P) (((P) > 0) ? (P) : 0)
131#define NBFACETRI_AINSWORTH_FACE_HDIV(P) (((P) > 2) ? (P - 1) * (P - 2) / 2 : 0)
132#define NBFACETRI_AINSWORTH_HDIV(P) (((P) > 0) ? (P + 1) * (P + 2) / 2 : 0)
133#define NBVOLUMETET_AINSWORTH_EDGE_HDIV(P) (((P) > 1) ? (P - 1) : 0)
134#define NBVOLUMETET_AINSWORTH_FACE_HDIV(P) (((P) > 2) ? (P - 1) * (P - 2) : 0)
135#define NBVOLUMETET_AINSWORTH_VOLUME_HDIV(P) \
136 (((P) > 3) ? (P - 3) * (P - 2) * (P - 1) / 2 : 0)
137#define NBVOLUMETET_AINSWORTH_HDIV(P) \
138 (((P) > 1) ? (P - 1) * (P + 1) * (P + 2) / 2 : 0)
139#define NBFACETRI_DEMKOWICZ_HDIV(P) ((P > 0) ? (P) * (P + 1) / 2 : 0)
140#define NBVOLUMETET_DEMKOWICZ_HDIV(P) \
141 (((P) > 1) ? (P) * (P - 1) * (P + 1) / 2 : 0)
142
143#define NBFACEQUAD_DEMKOWICZ_QUAD_HDIV_GENERAL(P, Q) \
144 (((P) > 0 && (Q) > 0) ? ((P) * (Q)) : 0)
145#define NBFACEQUAD_DEMKOWICZ_HDIV(P) \
146 (NBFACEQUAD_DEMKOWICZ_QUAD_HDIV_GENERAL(P, P))
147#define NBVOLUMEHEX_DEMKOWICZ_FAMILY_HDIV(P, Q, R) \
148 ((((P) > 0) && ((Q) > 0) && ((R) > 0)) ? ((P - 1) * Q * R) : 0)
149#define NBVOLUMEHEX_DEMKOWICZ_HDIV(P) \
150 (3 * NBVOLUMEHEX_DEMKOWICZ_FAMILY_HDIV(P, P, P))
151
152/**
153 * @brief Get base functions on triangle for L2 space
154 *
155 * @param p polynomial order
156 * @param N barycentric coordinates (shape functions) at integration points
157 * @param diffN derivatives of barycentric coordinates, i.e. derivatives of
158 * shape functions
159 * @param L2N values of L2 base at integration points
160 * @param diff_L2N dirvatives of base functions at integration points
161 * @param GDIM number of integration points
162 * @param base_polynomials polynomial base used to construct L2 base on
163 * element
164 * @return PetscErrorCode
165 */
167 int p, double *N, double *diffN, double *L2N, double *diff_L2N, int GDIM,
168 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
169 double *L, double *diffL,
170 const int dim));
171/**
172 * @brief Get base functions on tetrahedron for L2 space
173 *
174 * @param p polynomial order
175 * @param N barycentric coordinates (shape functions) at integration points
176 * @param diffN derivatives of barycentric coordinates, i.e. derivatives of
177 * shape functions
178 * @param L2N values of L2 base at integration points
179 * @param diff_L2N dirvatives of base functions at integration points
180 * @param GDIM number of integration points
181 * @param base_polynomials polynomial base used to construct L2 base on element
182 * @return PetscErrorCode
183 */
185 int p, double *N, double *diffN, double *L2N, double *diff_L2N, int GDIM,
186 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
187 double *L, double *diffL,
188 const int dim));
189
190/**
191 * \brief H1_EdgeShapeFunctions_MBTRI
192 *
193 * \param sense of edges, it is array of integers dim 3 (3-edges of triangle)
194 * \param p of edges
195 */
196PetscErrorCode H1_EdgeShapeFunctions_MBTRI(
197 int *sense, int *p, double *N, double *diffN, double *edgeN[3],
198 double *diff_edgeN[3], int GDIM,
199 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
200 double *L, double *diffL,
201 const int dim));
202PetscErrorCode H1_FaceShapeFunctions_MBTRI(
203 const int *face_nodes, int p, double *N, double *diffN, double *faceN,
204 double *diff_faceN, int GDIM,
205 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
206 double *L, double *diffL,
207 const int dim));
208PetscErrorCode H1_EdgeShapeFunctions_MBTET(
209 int *sense, int *p, double *N, double *diffN, double *edgeN[],
210 double *diff_edgeN[], int GDIM,
211 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
212 double *L, double *diffL,
213 const int dim));
214PetscErrorCode H1_FaceShapeFunctions_MBTET(
215 int *faces_nodes, int *p, double *N, double *diffN, double *faceN[],
216 double *diff_faceN[], int GDIM,
217 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
218 double *L, double *diffL,
219 const int dim));
220PetscErrorCode H1_VolumeShapeFunctions_MBTET(
221 int p, double *N, double *diffN, double *volumeN, double *diff_volumeN,
222 int GDIM,
223 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
224 double *L, double *diffL,
225 const int dim));
226PetscErrorCode H1_EdgeShapeDiffMBTETinvJ(int *base_p, int *p,
227 double *edge_diffN[], double *invJac,
228 double *edge_diffNinvJac[], int GDIM);
229PetscErrorCode H1_FaceShapeDiffMBTETinvJ(int *base_p, int *p,
230 double *face_diffN[], double *invJac,
231 double *face_diffNinvJac[], int GDIM);
232PetscErrorCode H1_VolumeShapeDiffMBTETinvJ(int base_p, int p,
233 double *volume_diffN, double *invJac,
234 double *volume_diffNinvJac,
235 int GDIM);
236
237PetscErrorCode H1_EdgeGradientOfDeformation_hierarchical(int p, double *diffN,
238 double *dofs,
239 double *F);
240PetscErrorCode H1_FaceGradientOfDeformation_hierarchical(int p, double *diffN,
241 double *dofs,
242 double *F);
243PetscErrorCode H1_VolumeGradientOfDeformation_hierarchical(int p, double *diffN,
244 double *dofs,
245 double *F);
246
247PetscErrorCode H1_QuadShapeFunctions_MBPRISM(
248 int *faces_nodes, int *p, double *N, double *diffN, double *faceN[],
249 double *diff_faceN[], int GDIM,
250 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
251 double *L, double *diffL,
252 const int dim));
254 int p, double *N, double *diffN, double *volumeN, double *diff_volumeN,
255 int GDIM,
256 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
257 double *L, double *diffL,
258 const int dim));
259PetscErrorCode H1_QuadShapeFunctions_MBQUAD(
260 int *faces_nodes, int p, double *N, double *diffN, double *faceN,
261 double *diff_faceN, int GDIM,
262 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
263 double *L, double *diffL,
264 const int dim));
265PetscErrorCode H1_EdgeShapeFunctions_MBQUAD(
266 int *sense, int *p, double *N, double *diffN, double *edgeN[4],
267 double *diff_edgeN[4], int GDIM,
268 PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
269 double *L, double *diffL,
270 const int dim));
271
272// Hdiv and Hcurl are implemented and declared in other files
273
274#ifdef __cplusplus
275}
276#endif
277
278#endif // __H1_HDIV_HCURL_L2_H__
MatrixDouble invJac
@ F
PetscErrorCode H1_QuadShapeFunctions_MBPRISM(int *faces_nodes, int *p, double *N, double *diffN, double *faceN[], double *diff_faceN[], int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition h1.c:642
PetscErrorCode H1_VolumeShapeDiffMBTETinvJ(int base_p, int p, double *volume_diffN, double *invJac, double *volume_diffNinvJac, int GDIM)
Definition h1.c:592
PetscErrorCode H1_EdgeShapeFunctions_MBTET(int *sense, int *p, double *N, double *diffN, double *edgeN[], double *diff_edgeN[], int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition h1.c:274
PetscErrorCode H1_FaceGradientOfDeformation_hierarchical(int p, double *diffN, double *dofs, double *F)
Definition h1.c:620
PetscErrorCode H1_QuadShapeFunctions_MBQUAD(int *faces_nodes, int p, double *N, double *diffN, double *faceN, double *diff_faceN, int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition h1.c:959
PetscErrorCode H1_VolumeGradientOfDeformation_hierarchical(int p, double *diffN, double *dofs, double *F)
Definition h1.c:631
PetscErrorCode L2_Ainsworth_ShapeFunctions_MBTET(int p, double *N, double *diffN, double *L2N, double *diff_L2N, int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Get base functions on tetrahedron for L2 space.
Definition l2.c:74
PetscErrorCode L2_Ainsworth_ShapeFunctions_MBTRI(int p, double *N, double *diffN, double *L2N, double *diff_L2N, int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Get base functions on triangle for L2 space.
Definition l2.c:19
PetscErrorCode H1_FaceShapeDiffMBTETinvJ(int *base_p, int *p, double *face_diffN[], double *invJac, double *face_diffNinvJac[], int GDIM)
Definition h1.c:574
PetscErrorCode H1_EdgeShapeDiffMBTETinvJ(int *base_p, int *p, double *edge_diffN[], double *invJac, double *edge_diffNinvJac[], int GDIM)
Definition h1.c:556
PetscErrorCode H1_EdgeShapeFunctions_MBQUAD(int *sense, int *p, double *N, double *diffN, double *edgeN[4], double *diff_edgeN[4], int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition h1.c:1091
PetscErrorCode H1_EdgeShapeFunctions_MBTRI(int *sense, int *p, double *N, double *diffN, double *edgeN[3], double *diff_edgeN[3], int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
H1_EdgeShapeFunctions_MBTRI.
Definition h1.c:17
PetscErrorCode H1_FaceShapeFunctions_MBTRI(const int *face_nodes, int p, double *N, double *diffN, double *faceN, double *diff_faceN, int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition h1.c:191
PetscErrorCode H1_EdgeGradientOfDeformation_hierarchical(int p, double *diffN, double *dofs, double *F)
Definition h1.c:609
PetscErrorCode H1_VolumeShapeFunctions_MBTET(int p, double *N, double *diffN, double *volumeN, double *diff_volumeN, int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition h1.c:475
PetscErrorCode H1_FaceShapeFunctions_MBTET(int *faces_nodes, int *p, double *N, double *diffN, double *faceN[], double *diff_faceN[], int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition h1.c:373
PetscErrorCode H1_VolumeShapeFunctions_MBPRISM(int p, double *N, double *diffN, double *volumeN, double *diff_volumeN, int GDIM, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition h1.c:790
const int N
Definition speed_test.cpp:3