v0.14.0
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
5 to 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
17 extern "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 fro L2 space
46  *
47  */
48 #define NBEDGE_L2(P) ((P) + 1)
49 
50 // H1
51 
52 /**
53  * @brief Numer 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_GEMERAL(P, Q) \
144  (((P) > 0 && (Q) > 0) ? ((P) * (Q)) : 0)
145 #define NBFACEQUAD_DEMKOWICZ_HDIV(P) \
146  (NBFACEQUAD_DEMKOWICZ_QUAD_HDIV_GEMERAL(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 // Bubbles for H div space
153 
154 /**
155  * @brief Get base functions on triangle for L2 space
156  *
157  * @param p polynomial order
158  * @param N barycentric coordinates (shape functions) at integration points
159  * @param diffN derivatives of barycentric coordinates, i.e. derivatives of
160  * shape functions
161  * @param L2N values of L2 base at integration points
162  * @param diff_L2N dirvatives of base functions at integration points
163  * @param GDIM number of integration points
164  * @param base_polynomials polynomial base used to construct L2 base on
165  * element
166  * @return PetscErrorCode
167  */
168 PetscErrorCode L2_Ainsworth_ShapeFunctions_MBTRI(
169  int p, double *N, double *diffN, double *L2N, double *diff_L2N, int GDIM,
170  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
171  double *L, double *diffL,
172  const int dim));
173 /**
174  * @brief Get base functions on tetrahedron for L2 space
175  *
176  * @param p polynomial order
177  * @param N barycentric coordinates (shape functions) at integration points
178  * @param diffN derivatives of barycentric coordinates, i.e. derivatives of
179  * shape functions
180  * @param L2N values of L2 base at integration points
181  * @param diff_L2N dirvatives of base functions at integration points
182  * @param GDIM number of integration points
183  * @param base_polynomials polynomial base used to construct L2 base on element
184  * @return PetscErrorCode
185  */
186 PetscErrorCode L2_Ainsworth_ShapeFunctions_MBTET(
187  int p, double *N, double *diffN, double *L2N, double *diff_L2N, int GDIM,
188  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
189  double *L, double *diffL,
190  const int dim));
191 
192 /**
193  * \brief H1_EdgeShapeFunctions_MBTRI
194  *
195  * \param sense of edges, it is array of integers dim 3 (3-edges of triangle)
196  * \param p of edges
197  */
198 PetscErrorCode H1_EdgeShapeFunctions_MBTRI(
199  int *sense, int *p, double *N, double *diffN, double *edgeN[3],
200  double *diff_edgeN[3], int GDIM,
201  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
202  double *L, double *diffL,
203  const int dim));
204 PetscErrorCode H1_FaceShapeFunctions_MBTRI(
205  const int *face_nodes, int p, double *N, double *diffN, double *faceN,
206  double *diff_faceN, int GDIM,
207  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
208  double *L, double *diffL,
209  const int dim));
210 PetscErrorCode H1_EdgeShapeFunctions_MBTET(
211  int *sense, int *p, double *N, double *diffN, double *edgeN[],
212  double *diff_edgeN[], int GDIM,
213  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
214  double *L, double *diffL,
215  const int dim));
216 PetscErrorCode H1_FaceShapeFunctions_MBTET(
217  int *faces_nodes, int *p, double *N, double *diffN, double *faceN[],
218  double *diff_faceN[], int GDIM,
219  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
220  double *L, double *diffL,
221  const int dim));
222 PetscErrorCode H1_VolumeShapeFunctions_MBTET(
223  int p, double *N, double *diffN, double *volumeN, double *diff_volumeN,
224  int GDIM,
225  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
226  double *L, double *diffL,
227  const int dim));
228 PetscErrorCode H1_EdgeShapeDiffMBTETinvJ(int *base_p, int *p,
229  double *edge_diffN[], double *invJac,
230  double *edge_diffNinvJac[], int GDIM);
231 PetscErrorCode H1_FaceShapeDiffMBTETinvJ(int *base_p, int *p,
232  double *face_diffN[], double *invJac,
233  double *face_diffNinvJac[], int GDIM);
234 PetscErrorCode H1_VolumeShapeDiffMBTETinvJ(int base_p, int p,
235  double *volume_diffN, double *invJac,
236  double *volume_diffNinvJac,
237  int GDIM);
238 
239 PetscErrorCode H1_EdgeGradientOfDeformation_hierarchical(int p, double *diffN,
240  double *dofs,
241  double *F);
242 PetscErrorCode H1_FaceGradientOfDeformation_hierarchical(int p, double *diffN,
243  double *dofs,
244  double *F);
245 PetscErrorCode H1_VolumeGradientOfDeformation_hierarchical(int p, double *diffN,
246  double *dofs,
247  double *F);
248 
249 PetscErrorCode H1_QuadShapeFunctions_MBPRISM(
250  int *faces_nodes, int *p, double *N, double *diffN, double *faceN[],
251  double *diff_faceN[], int GDIM,
252  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
253  double *L, double *diffL,
254  const int dim));
255 PetscErrorCode H1_VolumeShapeFunctions_MBPRISM(
256  int p, double *N, double *diffN, double *volumeN, double *diff_volumeN,
257  int GDIM,
258  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
259  double *L, double *diffL,
260  const int dim));
261 PetscErrorCode H1_QuadShapeFunctions_MBQUAD(
262  int *faces_nodes, int p, double *N, double *diffN, double *faceN,
263  double *diff_faceN, int GDIM,
264  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
265  double *L, double *diffL,
266  const int dim));
267 PetscErrorCode H1_EdgeShapeFunctions_MBQUAD(
268  int *sense, int *p, double *N, double *diffN, double *edgeN[4],
269  double *diff_edgeN[4], int GDIM,
270  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
271  double *L, double *diffL,
272  const int dim));
273 
274 // Hdiv and Hcurl are implemented and declared in other files
275 
276 #ifdef __cplusplus
277 }
278 #endif
279 
280 #endif // __H1_HDIV_HCURL_L2_H__
L2_Ainsworth_ShapeFunctions_MBTET
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
L2_Ainsworth_ShapeFunctions_MBTRI
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
H1_EdgeShapeFunctions_MBQUAD
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
H1_EdgeShapeFunctions_MBTRI
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
H1_QuadShapeFunctions_MBQUAD
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
H1_VolumeShapeDiffMBTETinvJ
PetscErrorCode H1_VolumeShapeDiffMBTETinvJ(int base_p, int p, double *volume_diffN, double *invJac, double *volume_diffNinvJac, int GDIM)
Definition: h1.c:592
H1_FaceShapeDiffMBTETinvJ
PetscErrorCode H1_FaceShapeDiffMBTETinvJ(int *base_p, int *p, double *face_diffN[], double *invJac, double *face_diffNinvJac[], int GDIM)
Definition: h1.c:574
H1_VolumeGradientOfDeformation_hierarchical
PetscErrorCode H1_VolumeGradientOfDeformation_hierarchical(int p, double *diffN, double *dofs, double *F)
Definition: h1.c:631
H1_EdgeShapeFunctions_MBTET
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
H1_QuadShapeFunctions_MBPRISM
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
H1_VolumeShapeFunctions_MBPRISM
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
H1_FaceShapeFunctions_MBTRI
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
H1_FaceGradientOfDeformation_hierarchical
PetscErrorCode H1_FaceGradientOfDeformation_hierarchical(int p, double *diffN, double *dofs, double *F)
Definition: h1.c:620
H1_EdgeGradientOfDeformation_hierarchical
PetscErrorCode H1_EdgeGradientOfDeformation_hierarchical(int p, double *diffN, double *dofs, double *F)
Definition: h1.c:609
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
N
const int N
Definition: speed_test.cpp:3
H1_FaceShapeFunctions_MBTET
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
H1_VolumeShapeFunctions_MBTET
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
invJac
MatrixDouble invJac
Definition: HookeElement.hpp:683
F
@ F
Definition: free_surface.cpp:394
H1_EdgeShapeDiffMBTETinvJ
PetscErrorCode H1_EdgeShapeDiffMBTETinvJ(int *base_p, int *p, double *edge_diffN[], double *invJac, double *edge_diffNinvJac[], int GDIM)
Definition: h1.c:556