|
| v0.14.0
|
Go to the documentation of this file.
9 #ifndef __H1_HDIV_HCURL_L2_H__
10 #define __H1_HDIV_HCURL_L2_H__
27 #define NBVOLUMETET_L2(P) ((P + 1) * (P + 2) * (P + 3) / 6)
32 #define NBVOLUMEHEX_L2_GENERAL(P, Q, R) ((P + 1) * (Q + 1) * (R + 1))
37 #define NBVOLUMEHEX_L2(P) (NBVOLUMEHEX_L2_GENERAL(P, P, P))
42 #define NBFACETRI_L2(P) ((P + 1) * (P + 2) / 2)
48 #define NBEDGE_L2(P) ((P) + 1)
55 #define NBEDGE_H1(P) (((P) > 1) ? (P - 1) : 0)
60 #define NBFACETRI_H1(P) (((P) > 2) ? ((P - 2) * (P - 1) / 2) : 0)
65 #define NBFACEQUAD_H1(P) (((P) > 1) ? ((P - 1) * (P - 1)) : 0)
70 #define NBFACEQUAD_L2(P) (((P) >= 0) ? (P + 1) * (P + 1) : 0)
75 #define NBVOLUMETET_H1(P) (((P) > 3) ? ((P - 3) * (P - 2) * (P - 1) / 6) : 0)
80 #define NBVOLUMEPRISM_H1(P) ((P > 3) ? ((P - 2) * (P - 2) * (P - 2)) : 0)
86 #define NBVOLUMEHEX_H1_GENERAL(P, Q, R) \
87 ((((P) > 1) && ((Q) > 1) && ((R) > 1)) ? (((P)-1) * ((Q)-1) * ((R)-1)) : 0)
93 #define NBVOLUMEHEX_H1(P) (NBVOLUMEHEX_H1_GENERAL(P, P, P))
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)
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)
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))
121 #define NBVOLUMEHEX_DEMKOWICZ_FAMILY_HCURL(P, Q, R) \
122 ((P > 0) && (Q > 1) && (R > 1) ? ((P) * (Q - 1) * (R - 1)) : 0)
124 #define NBVOLUMEHEX_DEMKOWICZ_HCURL(P) \
125 (3 * NBVOLUMEHEX_DEMKOWICZ_FAMILY_HCURL(P, P, P))
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)
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))
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,
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,
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,
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,
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,
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,
221 int p,
double *
N,
double *diffN,
double *volumeN,
double *diff_volumeN,
223 PetscErrorCode (*base_polynomials)(
int p,
double s,
double *diff_s,
224 double *
L,
double *diffL,
227 double *edge_diffN[],
double *
invJac,
228 double *edge_diffNinvJac[],
int GDIM);
230 double *face_diffN[],
double *
invJac,
231 double *face_diffNinvJac[],
int GDIM);
233 double *volume_diffN,
double *
invJac,
234 double *volume_diffNinvJac,
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,
254 int p,
double *
N,
double *diffN,
double *volumeN,
double *diff_volumeN,
256 PetscErrorCode (*base_polynomials)(
int p,
double s,
double *diff_s,
257 double *
L,
double *diffL,
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,
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,
278 #endif // __H1_HDIV_HCURL_L2_H__
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.
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.
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))
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.
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))
PetscErrorCode H1_VolumeShapeDiffMBTETinvJ(int base_p, int p, double *volume_diffN, double *invJac, double *volume_diffNinvJac, int GDIM)
PetscErrorCode H1_FaceShapeDiffMBTETinvJ(int *base_p, int *p, double *face_diffN[], double *invJac, double *face_diffNinvJac[], int GDIM)
PetscErrorCode H1_VolumeGradientOfDeformation_hierarchical(int p, double *diffN, double *dofs, double *F)
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))
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))
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))
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))
PetscErrorCode H1_FaceGradientOfDeformation_hierarchical(int p, double *diffN, double *dofs, double *F)
PetscErrorCode H1_EdgeGradientOfDeformation_hierarchical(int p, double *diffN, double *dofs, double *F)
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))
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))
PetscErrorCode H1_EdgeShapeDiffMBTETinvJ(int *base_p, int *p, double *edge_diffN[], double *invJac, double *edge_diffNinvJac[], int GDIM)