v0.14.0
EdgeQuadHexPolynomials.hpp
Go to the documentation of this file.
1 /** \file EdgeQuadHexPolynomials.hpp
2 
3  \brief Implementation of base Demkowicz functions on quad and hex
4 
5  Based on Hierarchic Finite Element Bases on Unstructured QUADS (2D)
6  and HEXES (3D) Meshes
7 
8 */
9 
10 #ifndef __EDGE_QUAD_HEX_HPP__
11 #define __EDGE_QUAD_HEX_HPP__
12 
13 namespace MoFEM {
15 
16 MoFEMErrorCode Legendre_polynomials01(int p, double s, double *L);
17 
18 MoFEMErrorCode Integrated_Legendre01(int p, double s, double *L, double *diffL);
19 
20 MoFEMErrorCode Face_orientMat(int *face_nodes, double orientMat[2][2]);
21 
22 /*
23 Segment (1d) basis functions
24 
25  0-----------1
26 */
28  double *diffN, double *bubbleN,
29  double *diff_bubbleN,
30  int nb_integration_pts);
31 
32 MoFEMErrorCode L2_ShapeFunctions_ONSEGMENT(int p, double *N, double *diffN,
33  double *funN, double *funDiffN,
34  int nb_integration_pts);
35 
36 /*
37  Quads
38  3-------2------2
39  | | y
40  | | ^
41  3 1 |
42  | | |
43  | | 0----- > x
44  0-------0------1
45 */
46 
47 /**
48 * \brief H1 Edge base functions on Quad
49 
50 Function generates hierarchical base of H1 comforting functions on a 2D
51 Quad.
52 
53 * @param sense array of orientation of edges (take 1 or -1)
54 * @param p array of orders (in each direction) of base
55 functions
56 * @param N array vertex shape functions evaluated at each
57 integration point
58 * @param diffN derivatives of vertex shape functions
59 * @return edgeN base functions on edges at qd pts
60 * @return diff_edgeN derivatives of edge shape functions at qd pts
61 * @param nb_integration_pts number of integration points
62 * @param base_polynomials polynomial base function (f.e. Legendre of
63 Integrated Legendre)
64 * @return error code
65 */
66 MoFEMErrorCode H1_EdgeShapeFunctions_ONQUAD(int *sense, int *p, double *N,
67  double *diffN, double *edgeN[4],
68  double *edgeDiffN[4],
69  int nb_integration_pts);
70 
71 /**
72 * \brief H1 Face bubble functions on Quad.
73 
74 Function generates hierarchical base of H1 comforting functions on a 2D quad.
75 
76 * @param face_nodes face nodes order
77 * @param p array of orders (in each direction) of base functions
78 * @param N array vertex shape functions evaluated at each
79 integration point
80 * @param diffN derivatives of vertex shape functions
81 * @return edgeN base functions on edges at qd pts
82 * @return diff_edgeN derivatives of edge shape functions at qd pts
83 * @param nb_integration_pts number of integration points
84 * @param base_polynomials polynomial base function (f.e. Legendre of Integrated
85 Legendre)
86 * @return error code
87 */
88 MoFEMErrorCode H1_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N,
89  double *diffN, double *faceN,
90  double *diff_faceN,
91  int nb_integration_pts);
92 
93 /**
94 * \brief L2 Face base functions on Quad
95 
96 Function generates hierarchical base of H1 comforting functions on a 2D quad.
97 
98 * @param p array of orders (in each direction) of base functions
99 * @param N array vertex shape functions evaluated at each
100 integration point
101 * @param diffN derivatives of vertex shape functions
102 * @return edgeN base functions on edges at qd pts
103 * @return diff_edgeN derivatives of edge shape functions at qd pts
104 * @param nb_integration_pts number of integration points
105 * @param base_polynomials polynomial base function (f.e. Legendre of Integrated
106 Legendre)
107 * @return error code
108 */
109 MoFEMErrorCode L2_FaceShapeFunctions_ONQUAD(int *p, double *N, double *diffN,
110  double *face_buble,
111  double *diff_face_bubble,
112  int nb_integration_pts);
113 
114 MoFEMErrorCode Hcurl_EdgeShapeFunctions_ONQUAD(int *sense, int *p, double *N,
115  double *diffN, double *edgeN[4],
116  double *curl_edgeN[4],
117  int nb_integration_pts);
118 
119 MoFEMErrorCode Hcurl_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p,
120  double *N, double *diffN,
121  double *faceN[],
122  double *diff_faceN[],
123  int nb_integration_pts);
124 
125 MoFEMErrorCode Hdiv_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p,
126  double *N, double *diffN,
127  double *faceN, double *diff_faceN,
128  int nb_integration_pts);
129 
130 /* Reference Hex and its canonical vertex and edge numbering
131  7---------10----------6
132  /| /|
133  / | / |
134  11 | 9 | x3
135  / 7 / | |
136  / | / 6 | x2
137  4----------8----------5 | | /
138  | | | | | /
139  | 3----------2----------2 | /
140  | / | / o -------- x1
141  4 / 5 /
142  | 3 | 1
143  | / | /
144  |/ |/
145  0---------0-----------1
146 
147  Hex Face Canonical numbering
148 
149  1. 1 2 3 4
150  2. 1 2 6 5
151  3. 2 3 7 6
152  4. 4 3 7 8
153  5. 1 4 8 5
154  6. 5 6 7 8
155 */
156 
157 MoFEMErrorCode H1_EdgeShapeFunctions_ONHEX(int *sense, int *p, double *N,
158  double *N_diff, double *edgeN[12],
159  double *diff_edgeN[12],
160  int nb_integration_pts);
161 
163 H1_FaceShapeFunctions_ONHEX(int *face_nodes, int *face_nodes_order, int *p,
164  double *N, double *N_diff, double *faceN[6],
165  double *diff_faceN[6], int nb_integration_pts);
166 
167 MoFEMErrorCode H1_InteriorShapeFunctions_ONHEX(const int *p, double *N,
168  double *N_diff, double *faceN,
169  double *diff_faceN,
170  int nb_integration_pts);
171 
172 MoFEMErrorCode L2_InteriorShapeFunctions_ONHEX(const int *p, double *N,
173  double *N_diff, double *volN,
174  double *diff_volN,
175  int nb_integration_pts);
176 
177 MoFEMErrorCode Hcurl_EdgeShapeFunctions_ONHEX(int *sense, int *p, double *N,
178  double *N_diff, double *edgeN[12],
179  double *diff_edgeN[12],
180  int nb_integration_pts);
181 
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);
185 
187  double *N_diff,
188  double *volN[3],
189  double *diff_volN[3],
190  int nb_integration_pts);
191 
193 Hdiv_FaceShapeFunctions_ONHEX(int *face_nodes, int *face_nodes_order, int *p,
194  double *N, double *diffN, double *faceN[6],
195  double *div_faceN[6], int nb_integration_pts);
196 
198  double *N_diff,
199  double *bubleN[3],
200  double *div_bubleN[3],
201  int nb_integration_pts);
202 
203 } // namespace DemkowiczHexAndQuad
204 
205 } // namespace MoFEM
206 
207 #endif // __EDGE_QUAD_HEX_HPP__
MoFEM::DemkowiczHexAndQuad::Legendre_polynomials01
MoFEMErrorCode Legendre_polynomials01(int p, double s, double *L)
MoFEM::DemkowiczHexAndQuad::L2_FaceShapeFunctions_ONQUAD
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.
Definition: EdgeQuadHexPolynomials.cpp:335
MoFEM::DemkowiczHexAndQuad::Face_orientMat
MoFEMErrorCode Face_orientMat(int *face_nodes, double orientMat[2][2])
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::DemkowiczHexAndQuad::L2_ShapeFunctions_ONSEGMENT
MoFEMErrorCode L2_ShapeFunctions_ONSEGMENT(int p, double *N, double *diffN, double *funN, double *funDiffN, int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:166
MoFEM::DemkowiczHexAndQuad::H1_EdgeShapeFunctions_ONQUAD
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.
Definition: EdgeQuadHexPolynomials.cpp:202
MoFEM::DemkowiczHexAndQuad::Hcurl_EdgeShapeFunctions_ONHEX
MoFEMErrorCode Hcurl_EdgeShapeFunctions_ONHEX(int *sense, int *p, double *N, double *N_diff, double *edgeN[12], double *diff_edgeN[12], int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:1019
MoFEM::DemkowiczHexAndQuad::Hdiv_InteriorShapeFunctions_ONHEX
MoFEMErrorCode Hdiv_InteriorShapeFunctions_ONHEX(int *p, double *N, double *N_diff, double *bubleN[3], double *div_bubleN[3], int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:1544
MoFEM::DemkowiczHexAndQuad::Integrated_Legendre01
MoFEMErrorCode Integrated_Legendre01(int p, double s, double *L, double *diffL)
MoFEM::DemkowiczHexAndQuad::Hcurl_FaceShapeFunctions_ONHEX
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)
Definition: EdgeQuadHexPolynomials.cpp:1139
MoFEM::DemkowiczHexAndQuad::L2_InteriorShapeFunctions_ONHEX
MoFEMErrorCode L2_InteriorShapeFunctions_ONHEX(const int *p, double *N, double *N_diff, double *volN, double *diff_volN, int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:966
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::DemkowiczHexAndQuad::Hdiv_FaceShapeFunctions_ONQUAD
MoFEMErrorCode Hdiv_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN, double *diff_faceN, int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:590
MoFEM::DemkowiczHexAndQuad::Hdiv_FaceShapeFunctions_ONHEX
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)
Definition: EdgeQuadHexPolynomials.cpp:1404
MoFEM::DemkowiczHexAndQuad::Hcurl_InteriorShapeFunctions_ONHEX
MoFEMErrorCode Hcurl_InteriorShapeFunctions_ONHEX(int *p, double *N, double *N_diff, double *volN[3], double *diff_volN[3], int nb_integration_pts)
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
DemkowiczHexAndQuad
Definition: EdgeQuadHexPolynomials.cpp:9
MoFEM::DemkowiczHexAndQuad::H1_FaceShapeFunctions_ONQUAD
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.
Definition: EdgeQuadHexPolynomials.cpp:275
MoFEM::DemkowiczHexAndQuad::H1_InteriorShapeFunctions_ONHEX
MoFEMErrorCode H1_InteriorShapeFunctions_ONHEX(const int *p, double *N, double *N_diff, double *faceN, double *diff_faceN, int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:905
MoFEM::DemkowiczHexAndQuad::H1_EdgeShapeFunctions_ONHEX
MoFEMErrorCode H1_EdgeShapeFunctions_ONHEX(int *sense, int *p, double *N, double *N_diff, double *edgeN[12], double *diff_edgeN[12], int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:707
N
const int N
Definition: speed_test.cpp:3
MoFEM::DemkowiczHexAndQuad::Hcurl_FaceShapeFunctions_ONQUAD
MoFEMErrorCode Hcurl_FaceShapeFunctions_ONQUAD(int *face_nodes, int *p, double *N, double *diffN, double *faceN[], double *diff_faceN[], int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:484
MoFEM::DemkowiczHexAndQuad::Hcurl_EdgeShapeFunctions_ONQUAD
MoFEMErrorCode Hcurl_EdgeShapeFunctions_ONQUAD(int *sense, int *p, double *N, double *diffN, double *edgeN[4], double *curl_edgeN[4], int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:393
MoFEM::DemkowiczHexAndQuad::H1_FaceShapeFunctions_ONHEX
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)
Definition: EdgeQuadHexPolynomials.cpp:795
MoFEM::DemkowiczHexAndQuad::H1_BubbleShapeFunctions_ONSEGMENT
MoFEMErrorCode H1_BubbleShapeFunctions_ONSEGMENT(int p, double *N, double *diffN, double *bubbleN, double *diff_bubbleN, int nb_integration_pts)
Definition: EdgeQuadHexPolynomials.cpp:141