|
| v0.14.0
|
Go to the documentation of this file.
10 #ifndef __EDGE_QUAD_HEX_HPP__
11 #define __EDGE_QUAD_HEX_HPP__
28 double *diffN,
double *bubbleN,
30 int nb_integration_pts);
33 double *funN,
double *funDiffN,
34 int nb_integration_pts);
67 double *diffN,
double *edgeN[4],
69 int nb_integration_pts);
89 double *diffN,
double *faceN,
91 int nb_integration_pts);
111 double *diff_face_bubble,
112 int nb_integration_pts);
115 double *diffN,
double *edgeN[4],
116 double *curl_edgeN[4],
117 int nb_integration_pts);
120 double *
N,
double *diffN,
122 double *diff_faceN[],
123 int nb_integration_pts);
126 double *
N,
double *diffN,
127 double *faceN,
double *diff_faceN,
128 int nb_integration_pts);
158 double *N_diff,
double *edgeN[12],
159 double *diff_edgeN[12],
160 int nb_integration_pts);
164 double *
N,
double *N_diff,
double *faceN[6],
165 double *diff_faceN[6],
int nb_integration_pts);
168 double *N_diff,
double *faceN,
170 int nb_integration_pts);
173 double *N_diff,
double *volN,
175 int nb_integration_pts);
178 double *N_diff,
double *edgeN[12],
179 double *diff_edgeN[12],
180 int nb_integration_pts);
183 int *face_nodes,
int *face_nodes_order,
int *p,
double *
N,
double *N_diff,
184 double *faceN[6][2],
double *diff_faceN[6][2],
int nb_integration_pts);
189 double *diff_volN[3],
190 int nb_integration_pts);
194 double *
N,
double *diffN,
double *faceN[6],
195 double *div_faceN[6],
int nb_integration_pts);
200 double *div_bubleN[3],
201 int nb_integration_pts);
207 #endif // __EDGE_QUAD_HEX_HPP__
MoFEMErrorCode Legendre_polynomials01(int p, double s, double *L)
MoFEMErrorCode L2_FaceShapeFunctions_ONQUAD(int *p, double *N, double *diffN, double *face_buble, double *diff_face_bubble, int nb_integration_pts)
L2 Face base functions on Quad.
MoFEMErrorCode Face_orientMat(int *face_nodes, double orientMat[2][2])
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
MoFEMErrorCode L2_ShapeFunctions_ONSEGMENT(int p, double *N, double *diffN, double *funN, double *funDiffN, int nb_integration_pts)
MoFEMErrorCode H1_EdgeShapeFunctions_ONQUAD(int *sense, int *p, double *N, double *diffN, double *edgeN[4], double *edgeDiffN[4], int nb_integration_pts)
H1 Edge base functions on Quad.
MoFEMErrorCode Hcurl_EdgeShapeFunctions_ONHEX(int *sense, int *p, double *N, double *N_diff, double *edgeN[12], double *diff_edgeN[12], int nb_integration_pts)
MoFEMErrorCode Hdiv_InteriorShapeFunctions_ONHEX(int *p, double *N, double *N_diff, double *bubleN[3], double *div_bubleN[3], int nb_integration_pts)
MoFEMErrorCode Integrated_Legendre01(int p, double s, double *L, double *diffL)
MoFEMErrorCode Hcurl_FaceShapeFunctions_ONHEX(int *face_nodes, int *face_nodes_order, int *p, double *N, double *N_diff, double *faceN[6][2], double *diff_faceN[6][2], int nb_integration_pts)
MoFEMErrorCode L2_InteriorShapeFunctions_ONHEX(const int *p, double *N, double *N_diff, double *volN, double *diff_volN, int nb_integration_pts)
implementation of Data Operators for Forces and Sources
MoFEMErrorCode Hdiv_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN, double *diff_faceN, int nb_integration_pts)
MoFEMErrorCode Hdiv_FaceShapeFunctions_ONHEX(int *face_nodes, int *face_nodes_order, int *p, double *N, double *diffN, double *faceN[6], double *div_faceN[6], int nb_integration_pts)
MoFEMErrorCode Hcurl_InteriorShapeFunctions_ONHEX(int *p, double *N, double *N_diff, double *volN[3], double *diff_volN[3], int nb_integration_pts)
MoFEMErrorCode H1_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN, double *diff_faceN, int nb_integration_pts)
H1 Face bubble functions on Quad.
MoFEMErrorCode H1_InteriorShapeFunctions_ONHEX(const int *p, double *N, double *N_diff, double *faceN, double *diff_faceN, int nb_integration_pts)
MoFEMErrorCode H1_EdgeShapeFunctions_ONHEX(int *sense, int *p, double *N, double *N_diff, double *edgeN[12], double *diff_edgeN[12], int nb_integration_pts)
MoFEMErrorCode Hcurl_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN[], double *diff_faceN[], int nb_integration_pts)
MoFEMErrorCode Hcurl_EdgeShapeFunctions_ONQUAD(int *sense, int *p, double *N, double *diffN, double *edgeN[4], double *curl_edgeN[4], int nb_integration_pts)
MoFEMErrorCode H1_FaceShapeFunctions_ONHEX(int *face_nodes, int *face_nodes_order, int *p, double *N, double *N_diff, double *faceN[6], double *diff_faceN[6], int nb_integration_pts)
MoFEMErrorCode H1_BubbleShapeFunctions_ONSEGMENT(int p, double *N, double *diffN, double *bubbleN, double *diff_bubbleN, int nb_integration_pts)