v0.14.0
Hcurl.hpp
Go to the documentation of this file.
1 /** \file Hcurl.hpp
2 
3  \brief Implementation of H-curl base function.
4 
5  Based on Hierarchic Finite Element Bases on Unstructured Tetrahedral
6  Meshes, by \cite NME:NME847 and \cite fuentes2015orientation
7  Shape functions for MBTRI/MBTET and HCurl space
8 
9 */
10 
11 #ifndef __HCURL_HPP__
12 #define __HCURL_HPP__
13 
14 namespace MoFEM {
15 
16 /**
17  * \brief Edge based H-curl base functions on tetrahedral
18 
19  Function generates hierarchical base of h-curl comforting functions on
20  tetrahedral edge. For more details see \cite ainsworth2011bernstein.
21 
22  On each tetrahedral's edge we have P+1 functions. See NBEDGE_AINSWORTH_HCURL
23 
24  * @param sense sense fo edge (i.e. unique orientation)
25  * @param p array of oder for each edge
26  * @param N array shape functions evaluated at each integration
27  point
28  * @param diffN derivatives of shape functions
29  * @param edgeN base functions on edges
30  * @param diff_edgeN derivatives of edge shape functions
31  * @param nb_integration_pts number of integration points
32  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
33  * @return error code
34  */
36  int *sense, int *p, double *N, double *diffN, double *edgeN[],
37  double *diff_edgeN[], int nb_integration_pts,
38  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
39  double *L, double *diffL,
40  const int dim));
41 
42 /**
43  * \brief Edge based H-curl base functions on edge
44 
45  Function generates hierarchical base of h-curl comforting functions on
46  tetrahedral edge. For more details see \cite ainsworth2011bernstein.
47 
48  On each edge we have P+1 functions. See NBEDGE_AINSWORTH_HCURL
49 
50  * @param sense sense fo edge (i.e. unique orientation)
51  * @param p array of oder for each edge
52  * @param N array shape functions evaluated at each integration
53  point
54  * @param diffN derivatives of shape functions
55  * @param edgeN base functions on edges
56  * @param diff_edgeN derivatives of edge shape functions
57  * @param nb_integration_pts number of integration points
58  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
59  * @return error code
60  */
62  int sense, int p, double *N, double *diffN, double *edgeN,
63  double *diff_edgeN, int nb_integration_pts,
64  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
65  double *L, double *diffL,
66  const int dim));
67 
68 /**
69  * \brief Edge based H-curl base functions on face
70 
71  Function generates hierarchical base of h-curl comforting functions on
72  tetrahedral edge. For more details see \cite ainsworth2011bernstein.
73 
74  On each face's edge we have P+1 functions. See NBEDGE_AINSWORTH_HCURL
75 
76  * @param sense sense fo edge (i.e. unique orientation)
77  * @param p array of oder for each edge
78  * @param N array shape functions evaluated at each integration
79  point
80  * @param diffN derivatives of shape functions
81  * @param edgeN base functions on edges
82  * @param diff_edgeN derivatives of edge shape functions
83  * @param nb_integration_pts number of integration points
84  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
85  * @return error code
86  */
88  int *sense, int *p, double *N, double *diffN, double *edgeN[],
89  double *diff_edgeN[], int nb_integration_pts,
90  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
91  double *L, double *diffL,
92  const int dim));
93 
94 /** \brief Face edge base functions of Hcurl space on tetrahedral.
95 
96  On each edge we have (P-1) base functions, and each face has 3 edges and are 4
97  faces on tetrahedral.
98 
99  See NBFACETRI_AINSWORTH_EDGE_HCURL
100 
101  * @param face_nodes array [4*3] of local indices of face nodes
102  * @param p approximation order
103  * @param N array shape functions evaluated at each integration
104  point
105  * @param diffN derivatives of nodal shape functions
106  * @param phi_f[4] calculated shape functions for each face
107  * @param diff_phi_v[4] derivatives of shape functions for each face
108  * @param nb_integration_pts number of shape functions
109  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
110  * @return error code
111 
112 */
114  int *faces_nodes, int *p, double *N, double *diffN, double *phi_f_e[4][3],
115  double *diff_phi_f_e[4][3], int nb_integration_pts,
116  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
117  double *L, double *diffL,
118  const int dim));
119 
120 /** \brief Face edge base functions of Hcurl space.
121 
122  On each edge we have (P-1) base functions, and each face has 3 edges and are 4
123  faces on tetrahedral.
124 
125  See NBFACETRI_AINSWORTH_EDGE_HCURL
126 
127  * @param face_nodes array [4*3] of local indices of face nodes
128  * @param p approximation order
129  * @param N array shape functions evaluated at each integration
130  point
131  * @param diffN derivatives of nodal shape functions
132  * @param phi_f[4] calculated shape functions for each face
133  * @param diff_phi_v[4] derivatives of shape functions for each face
134  * @param nb_integration_pts number of shape functions
135  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
136  * @return error code
137 
138 */
140  int *faces_nodes, int p, double *N, double *diffN, double *phi_f_e[3],
141  double *diff_phi_f_e[3], int nb_integration_pts,
142  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
143  double *L, double *diffL,
144  const int dim));
145 
146 /** \brief Face edge base functions of Hcurl space on face on tetrahedral.
147 
148  On each face we have P*(P-1) base functions and are 4 faces.
149 
150  See NBFACETRI_AINSWORTH_EDGE_HCURL
151 
152  * @param face_nodes array [4*3] of local indices of face nodes
153  * @param p approximation order
154  * @param N array shape functions evaluated at each integration
155  point
156  * @param diffN derivatives of nodal shape functions
157  * @param phi_f[4] calculated shape functions for each face
158  * @param diff_phi_v[4] derivatives of shape functions for each face
159  * @param nb_integration_pts number of shape functions
160  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
161  * @return error code
162 
163 */
165  int *faces_nodes, int *p, double *N, double *diffN, double *phi_f[4],
166  double *diff_phi_f[4], int nb_integration_pts,
167  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
168  double *L, double *diffL,
169  const int dim));
170 
171 /** \brief Face edge base functions of Hcurl space on face.
172 
173  On each face we have P*(P-1) base functions and are 4 faces.
174 
175  See NBFACETRI_AINSWORTH_EDGE_HCURL
176 
177  * @param face_nodes array [4*3] of local indices of face nodes
178  * @param p approximation order
179  * @param N array shape functions evaluated at each integration
180  point
181  * @param diffN derivatives of nodal shape functions
182  * @param phi_f[4] calculated shape functions for each face
183  * @param diff_phi_v[4] derivatives of shape functions for each face
184  * @param nb_integration_pts number of shape functions
185  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
186  * @return error code
187 
188 */
190  int *faces_nodes, int p, double *N, double *diffN, double *phi_f,
191  double *diff_phi_f, int nb_integration_pts,
192  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
193  double *L, double *diffL,
194  const int dim));
195 
196 /** \brief Face base interior function
197 
198 On each face we have P*(P-1)/2 and are 4 faces on Tetrahedral.
199 
200 See NBFACETRI_AINSWORTH_FACE_HCURL
201 
202 * @param face_nodes array [4*3] of local indices of face nodes
203 * @param p approximation order
204 * @param N nodal shape functions
205 * @param diffN derivatives of nodal shape functions
206 * @param phi_v calculated shape functions
207 * @param diff_phi_v derivatives of shape functions
208 * @param nb_integration_pts number of shape functions
209 * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
210 * @return error code
211 
212 
213 */
215  int *faces_nodes, int p, double *N, double *diffN, double *phi_v,
216  double *diff_phi_v, int nb_integration_pts,
217  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
218  double *L, double *diffL,
219  const int dim));
220 
221 /** \brief Volume interior function
222 
223 On volume have (P-3)*(P-2)*(P-1)/2.
224 
225 See NBVOLUMETET_AINSWORTH_TET_HCURL
226 
227 * @param p approximation order
228 * @param N nodal shape functions
229 * @param diffN derivatives of nodal shape functions
230 * @param phi_v calculated shape functions
231 * @param diff_phi_v derivatives of shape functions
232 * @param nb_integration_pts number of shape functions
233 * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
234 * @return error code
235 
236 */
238  int p, double *N, double *diffN, double *phi_v, double *diff_phi_v,
239  int nb_integration_pts,
240  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
241  double *L, double *diffL,
242  const int dim));
243 
244 /** \brief Face H-curl functions
245 
246  Face H-curl functions are set of Eddge-Based Face functions
247  (Hcurl_Ainsworth_EdgeBasedFaceFunctions_MBTET) and Bubble-Face functions
248  (Hcurl_Ainsworth_BubbleFaceFunctions_MBTET).
249 
250  See NBVOLUMETET_AINSWORTH_FACE_HCURL
251 
252  * @param face_nodes array [4*3] of local indices of face nodes
253  * @param p approximation order
254  * @param N nodal shape functions
255  * @param diffN derivatives of nodal shape functions
256  * @param phi_f[4] calculated shape functions for each face
257  * @param diff_phi_v[4] derivatives of shape functions for each face
258  * @param nb_integration_pts number of shape functions
259  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
260  * @return error code
261 
262 */
264  int *face_nodes, int *p, double *N, double *diffN, double *phi_f[4],
265  double *diff_phi_f[4], int nb_integration_pts,
266  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
267  double *L, double *diffL,
268  const int dim));
269 
270 /** \brief Face H-curl functions
271 
272  Face H-curl functions are set of Eddge-Based Face functions
273  (Hcurl_Ainsworth_EdgeBasedFaceFunctions_MBTET) and Bubble-Face functions
274  (Hcurl_Ainsworth_BubbleFaceFunctions_MBTET).
275 
276  See NBVOLUMETET_AINSWORTH_FACE_HCURL
277 
278  * @param face_nodes array [4*3] of local indices of face nodes
279  * @param p approximation order
280  * @param N nodal shape functions
281  * @param diffN derivatives of nodal shape functions
282  * @param phi_f[4] calculated shape functions for each face
283  * @param diff_phi_v[4] derivatives of shape functions for each face
284  * @param nb_integration_pts number of shape functions
285  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
286  * @return error code
287 
288 */
290  int *faces_nodes, int p, double *N, double *diffN, double *phi_f,
291  double *diff_phi_f, int nb_integration_pts,
292  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
293  double *L, double *diffL,
294  const int dim));
295 
296 /**
297  \brief H-curl volume base functions
298 
299  Volume base functions are collection of iace interior base functions and
300  volume interior base functions.
301 
302  * @param p approximation order
303  * @param N nodal shape functions
304  * @param diffN derivatives of nodal shape functions
305  * @param phi_v calculated shape functions
306  * @param diff_phi_v derivatives of shape functions
307  * @param nb_integration_pts number of shape functions
308  * @param base_polynomials polynomial base function (f.e. Legendre of Lobatto)
309  * @return error code
310 
311  */
313  int p, double *N, double *diffN, double *phi_v, double *diff_phi_v,
314  int nb_integration_pts,
315  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
316  double *L, double *diffL,
317  const int dim));
318 
319 /**
320  * \brief Edge based H-curl base functions on tetrahedral
321 
322  Function generates hierarchical base of h-curl comforting functions on
323  tetrahedral edge. For more details see \cite ainsworth2011bernstein.
324 
325  * @param sense sense fo edge (i.e. unique orientation)
326  * @param p array of oder for each edge
327  * @param n array shape functions evaluated at each integration
328  point
329  * @param diff_n derivatives of shape functions
330  * @param phi base functions on edges
331  * @param diff_phi derivatives of edge shape functions
332  * @param nb_integration_pts number of integration points
333  * @return error code
334  */
336  int *sense, int *p, double *n, double *diff_n, double *phi[],
337  double *diff_phi[], int nb_integration_pts);
338 
339 /**
340  * \brief Edge based H-curl base functions on teriangle
341 
342  Function generates hierarchical base of h-curl comforting functions on
343  tetrahedral edge. For more details see \cite ainsworth2011bernstein.
344 
345  * @param sense sense fo edge (i.e. unique orientation)
346  * @param p array of oder for each edge
347  * @param n array shape functions evaluated at each integration
348  point
349  * @param diff_n derivatives of shape functions
350  * @param phi base functions on edges
351  * @param diff_phi derivatives of edge shape functions
352  * @param nb_integration_pts number of integration points
353  * @return error code
354  */
356  int *sense, int *p, double *n, double *diff_n, double *phi[],
357  double *diff_phi[], int nb_integration_pts);
358 
359 /**
360  * \brief Edge based H-curl base functions on edge
361 
362  Function generates hierarchical base of h-curl comforting functions on
363  tetrahedral edge. For more details see \cite ainsworth2011bernstein.
364 
365  * @param sense sense fo edge (i.e. unique orientation)
366  * @param p array of oder for each edge
367  * @param n array shape functions evaluated at each integration
368  point
369  * @param diff_n derivatives of shape functions
370  * @param phi base functions on edges
371  * @param diff_phi derivatives of edge shape functions
372  * @param nb_integration_pts number of integration points
373  * @return error code
374  */
376  int sense, int p, double *n, double *diff_n, double *phi,
377  double *diff_phi, int nb_integration_pts);
378 
379 /** \brief Face base interior function
380 
381 * @param face_nodes array [4*3] of local indices of face nodes
382 * @param p approximation order
383 * @param n nodal shape functions
384 * @param diff_n derivatives of nodal shape functions
385 * @param phi calculated shape functions
386 * @param diff_phi derivatives of shape functions
387 * @param nb_integration_pts number of shape functions
388 * @return error code
389 
390 */
392  int *faces_nodes, int *p, double *n, double *diff_n, double *phi[],
393  double *diff_phi[], int nb_integration_pts);
394 
395 /** \brief Face base interior function
396 
397 * @param face_nodes array [4*3] of local indices of face nodes
398 * @param p approximation order
399 * @param n nodal shape functions
400 * @param diff_n derivatives of nodal shape functions
401 * @param phi calculated shape functions
402 * @param diff_phi derivatives of shape functions
403 * @param nb_integration_pts number of shape functions
404 * @return error code
405 
406 */
408  int *faces_nodes, int p, double *n, double *diff_n, double *phi,
409  double *diff_phi, int nb_integration_pts);
410 
411 /** \brief Volume base interior function
412 
413 * @param p approximation order
414 * @param n nodal shape functions
415 * @param diff_n derivatives of nodal shape functions
416 * @param phi calculated shape functions
417 * @param diff_phi derivatives of shape functions
418 * @param nb_integration_pts number of shape functions
419 * @return error code
420 
421 
422 */
424  int p, double *n, double *diff_n, double *phi,
425  double *diff_phi, int nb_integration_pts);
426 
427 
428 } // namespace MoFEM
429 
430 #endif // __HCURL_HPP__
MoFEM::Hcurl_Ainsworth_FaceFunctions_MBTET
MoFEMErrorCode Hcurl_Ainsworth_FaceFunctions_MBTET(int *face_nodes, int *p, double *N, double *diffN, double *phi_f[4], double *diff_phi_f[4], int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face H-curl functions.
Definition: Hcurl.cpp:1052
MoFEM::Hcurl_Ainsworth_EdgeBaseFunctions_MBTET_ON_EDGE
MoFEMErrorCode Hcurl_Ainsworth_EdgeBaseFunctions_MBTET_ON_EDGE(int sense, int p, double *N, double *diffN, double *edgeN, double *diff_edgeN, int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Edge based H-curl base functions on edge.
Definition: Hcurl.cpp:175
phi
static double phi
Definition: poisson_2d_dis_galerkin.cpp:29
MoFEM::Hcurl_Demkowicz_EdgeBaseFunctions_MBTET
MoFEMErrorCode Hcurl_Demkowicz_EdgeBaseFunctions_MBTET(int *sense, int *p, double *n, double *diff_n, double *phi[], double *diff_phi[], int nb_integration_pts)
Edge based H-curl base functions on tetrahedral.
Definition: Hcurl.cpp:2068
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Hcurl_Demkowicz_EdgeBaseFunctions_MBEDGE
MoFEMErrorCode Hcurl_Demkowicz_EdgeBaseFunctions_MBEDGE(int sense, int p, double *n, double *diff_n, double *phi, double *diff_phi, int nb_integration_pts)
Edge based H-curl base functions on edge.
Definition: Hcurl.cpp:2144
MoFEM::Hcurl_Ainsworth_FaceInteriorFunctions_MBTET
MoFEMErrorCode Hcurl_Ainsworth_FaceInteriorFunctions_MBTET(int *faces_nodes, int p, double *N, double *diffN, double *phi_v, double *diff_phi_v, int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face base interior function.
Definition: Hcurl.cpp:775
MoFEM::Hcurl_Ainsworth_FaceFunctions_MBTET_ON_FACE
MoFEMErrorCode Hcurl_Ainsworth_FaceFunctions_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face H-curl functions.
Definition: Hcurl.cpp:1249
MoFEM::Hcurl_Ainsworth_EdgeBasedFaceFunctions_MBTET
MoFEMErrorCode Hcurl_Ainsworth_EdgeBasedFaceFunctions_MBTET(int *faces_nodes, int *p, double *N, double *diffN, double *phi_f_e[4][3], double *diff_phi_f_e[4][3], int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face edge base functions of Hcurl space on tetrahedral.
Definition: Hcurl.cpp:363
MoFEM::Hcurl_Demkowicz_EdgeBaseFunctions_MBTRI
MoFEMErrorCode Hcurl_Demkowicz_EdgeBaseFunctions_MBTRI(int *sense, int *p, double *n, double *diff_n, double *phi[], double *diff_phi[], int nb_integration_pts)
Edge based H-curl base functions on teriangle.
Definition: Hcurl.cpp:2105
MoFEM::Hcurl_Ainsworth_EdgeBaseFunctions_MBTET_ON_FACE
MoFEMErrorCode Hcurl_Ainsworth_EdgeBaseFunctions_MBTET_ON_FACE(int *sense, int *p, double *N, double *diffN, double *edgeN[], double *diff_edgeN[], int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Edge based H-curl base functions on face.
Definition: Hcurl.cpp:237
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::Hcurl_Ainsworth_EdgeBasedFaceFunctions_MBTET_ON_FACE
MoFEMErrorCode Hcurl_Ainsworth_EdgeBasedFaceFunctions_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f_e[3], double *diff_phi_f_e[3], int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face edge base functions of Hcurl space.
Definition: Hcurl.cpp:458
MoFEM::Hcurl_Demkowicz_VolumeBaseFunctions_MBTET
MoFEMErrorCode Hcurl_Demkowicz_VolumeBaseFunctions_MBTET(int p, double *n, double *diff_n, double *phi, double *diff_phi, int nb_integration_pts)
Volume base interior function.
Definition: Hcurl.cpp:2468
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
convert.n
n
Definition: convert.py:82
N
const int N
Definition: speed_test.cpp:3
MoFEM::Hcurl_Ainsworth_EdgeBaseFunctions_MBTET
MoFEMErrorCode Hcurl_Ainsworth_EdgeBaseFunctions_MBTET(int *sense, int *p, double *N, double *diffN, double *edgeN[], double *diff_edgeN[], int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Edge based H-curl base functions on tetrahedral.
Definition: Hcurl.cpp:16
MoFEM::Hcurl_Ainsworth_VolumeInteriorFunctions_MBTET
MoFEMErrorCode Hcurl_Ainsworth_VolumeInteriorFunctions_MBTET(int p, double *N, double *diffN, double *phi_v, double *diff_phi_v, int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Volume interior function.
Definition: Hcurl.cpp:909
MoFEM::Hcurl_Demkowicz_FaceBaseFunctions_MBTET
MoFEMErrorCode Hcurl_Demkowicz_FaceBaseFunctions_MBTET(int *faces_nodes, int *p, double *n, double *diff_n, double *phi[], double *diff_phi[], int nb_integration_pts)
Face base interior function.
Definition: Hcurl.cpp:2395
MoFEM::Hcurl_Ainsworth_BubbleFaceFunctions_MBTET
MoFEMErrorCode Hcurl_Ainsworth_BubbleFaceFunctions_MBTET(int *faces_nodes, int *p, double *N, double *diffN, double *phi_f[4], double *diff_phi_f[4], int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face edge base functions of Hcurl space on face on tetrahedral.
Definition: Hcurl.cpp:545
MoFEM::Hcurl_Demkowicz_FaceBaseFunctions_MBTRI
MoFEMErrorCode Hcurl_Demkowicz_FaceBaseFunctions_MBTRI(int *faces_nodes, int p, double *n, double *diff_n, double *phi, double *diff_phi, int nb_integration_pts)
Face base interior function.
Definition: Hcurl.cpp:2433
MoFEM::Hcurl_Ainsworth_VolumeFunctions_MBTET
MoFEMErrorCode Hcurl_Ainsworth_VolumeFunctions_MBTET(int p, double *N, double *diffN, double *phi_v, double *diff_phi_v, int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
H-curl volume base functions.
Definition: Hcurl.cpp:1403
MoFEM::Hcurl_Ainsworth_BubbleFaceFunctions_MBTET_ON_FACE
MoFEMErrorCode Hcurl_Ainsworth_BubbleFaceFunctions_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int nb_integration_pts, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face edge base functions of Hcurl space on face.
Definition: Hcurl.cpp:663