v0.14.0
Hdiv.hpp
Go to the documentation of this file.
1 /** \file Hdiv.hpp
2 
3  \brief Implementation of H-curl base function.
4 
5  Based on Hierarchic Finite Element Bases on Unstructured Tetrahedral
6  Meshes, by Mark Ainsworth and Joe Coyle
7  Shape functions for MBTRI/MBTET and HCurl space
8 
9 */
10 
11 #ifndef __HDIV_HPP__
12 #define __HDIV_HPP__
13 
14 namespace MoFEM {
15 
16 /** \brief Broken base Ainsworth subentries order change hooks.
17  *
18  * Hooks enabling change of DOFs numbers/polynomial order on subentries, when
19  * base on element is constructed.
20  *
21  * Note that this functionality is global, that all functions are static.
22  *
23  **/
25  static boost::function<int(int)> broken_nbfacetri_edge_hdiv;
26  static boost::function<int(int)> broken_nbfacetri_face_hdiv;
27  static boost::function<int(int)> broken_nbvolumetet_edge_hdiv;
28  static boost::function<int(int)> broken_nbvolumetet_face_hdiv;
29  static boost::function<int(int)> broken_nbvolumetet_volume_hdiv;
30 };
31 
32 /**
33  * \brief Hdiv base functions, Edge-based face functions by Ainsworth \cite
34  * NME:NME847
35  * @param faces_nodes Face nodes on tetrahedral
36  * @param p Approximation order on faces
37  * @param N Shape functions
38  * @param diffN Derivatives of shape functions
39  * @param phi_f_e Base functions (returned)
40  * @param diff_phi_f_e Derivatives of base functions (returned)
41  * @param gdim Number of integration pts
42  * @param base_polynomials base function (Legendre/Lobbatto-Gauss)
43  * @return error code
44  */
46  int *faces_nodes, int *p, double *N, double *diffN, double *phi_f_e[4][3],
47  double *diff_phi_f_e[4][3], int gdim,
48  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
49  double *L, double *diffL,
50  const int dim));
51 
52 /**
53  * \brief Hdiv base functions, Edge-based face functions by Ainsworth \cite
54  * NME:NME847
55  * @param faces_nodes Face nodes on face
56  * @param p Approximation order on faces
57  * @param N Shape functions
58  * @param diffN Derivatives of shape functions
59  * @param phi_f_e Base functions (returned)
60  * @param diff_phi_f_e Derivatives of base functions (returned)
61  * @param gdim Number of integration pts
62  * @param base_polynomials base function (Legendre/Lobbatto-Gauss)
63  * @param nb Number of nodes on entity (4 if tet, 3 if triangle)
64  * @return error code
65  */
67  int *faces_nodes, int p, double *N, double *diffN, double *phi_f_e[3],
68  double *diff_phi_f_e[3], int gdim, int nb,
69  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
70  double *L, double *diffL,
71  const int dim));
72 
73 /**
74  * \brief Face bubble functions by Ainsworth \cite NME:NME847
75  * @param faces_nodes Face nodes on tetrahedral
76  * @param p Approx. order on faces
77  * @param N Shape function
78  * @param diffN Derivatives of shape functions
79  * @param phi_f Base functions
80  * @param diff_phi_f Derivatives of base functions
81  * @param gdim Number of integration pts
82  * @param base_polynomials Base function (Legendre/Lobbatto-Gauss)
83  * @return error code
84  */
86  int *faces_nodes, int *p, double *N, double *diffN, double *phi_f[],
87  double *diff_phi_f[], int gdim,
88  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
89  double *L, double *diffL,
90  const int dim));
91 
92 /**
93  * \brief Face bubble functions by Ainsworth \cite NME:NME847
94  * @param faces_nodes Face nodes on face
95  * @param p Approx. order on face
96  * @param N Shape function
97  * @param diffN Derivatives of shape functions
98  * @param phi_f Base functions
99  * @param diff_phi_f Derivatives of base functions
100  * @param gdim Number of integration pts
101  * @param nb Number of nodes on entity (4 if tet, 3 if triangle)
102  * @param base_polynomials Base function (Legendre/Lobbatto-Gauss)
103  * @return error code
104  */
106  int *faces_nodes, int p, double *N, double *diffN, double *phi_f,
107  double *diff_phi_f, int gdim, int nb,
108  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
109  double *L, double *diffL,
110  const int dim));
111 
112 /**
113  * \brief Hdiv base function, Edge-based interior (volume) functions by
114  * Ainsworth \cite NME:NME847
115  * @param p volume approx. order
116  * @param N shape functions
117  * @param diffN derivatives of shape functions
118  * @param phi_v_e base functions (return)
119  * @param diff_phi_v_e derivatives of base functions (returned)
120  * @param gdim number of integration points
121  * @param base_polynomials base function (Legendre/Lobbatto-Gauss)
122  * @return error code
123  */
125  int p, double *N, double *diffN, double *phi_v_e[6],
126  double *diff_phi_v_e[6], int gdim,
127  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
128  double *L, double *diffL,
129  const int dim));
130 
131 /** Hdiv Face-based interior functions by Ainsworth \cite NME:NME847
132  * @param p Approximation order on face
133  * @param N Shape functions on face
134  * @param diffN Derivatives of shape functions of face
135  * @param phi_v_f Base functions (returned)
136  * @param diff_phi_v_f Derivatives of base functions (returned)
137  * @param gdim Number of integration points
138  * @param base_polynomials Base function (Legendre/Lobbatto-Gauss)
139  * @return Error code
140  */
142  int p, double *N, double *diffN, double *phi_v_f[], double *diff_phi_v_f[],
143  int gdim,
144  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
145  double *L, double *diffL,
146  const int dim));
147 
148 /**
149  * \brief Interior bubble functions by Ainsworth \cite NME:NME847
150  * @param p Volume order
151  * @param N Shape functions
152  * @param diffN Derivatives of shape functions
153  * @param phi_v Base functions
154  * @param diff_phi_v Derivatives of shape functions
155  * @param gdim Number of integration points
156  * @param base_polynomials Base function (Legendre/Lobbatto-Gauss)
157  * @return Error code
158  */
160  int p, double *N, double *diffN, double *phi_v, double *diff_phi_v,
161  int gdim,
162  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
163  double *L, double *diffL,
164  const int dim));
165 
166 /**
167  * \brirf HDiv base finuctions on triangle by Demkowicz \cite
168  * fuentes2015orientation
169  *
170  * @param faces_nodes Nodes on face
171  * @param p Approx. order
172  * @param N Shape functions
173  * @param diffN Derivatives of shape functions
174  * @param phi_f Returned base functions on face
175  * @param diff_phi_f Returned derivatives of base functions on face
176  * @param gdim Number of integration points
177  * @param nb nb is 4 for face on tetrahedral and 3 for face
178  * @return error code
179  */
180 MoFEMErrorCode Hdiv_Demkowicz_Face_MBTET_ON_FACE(int *faces_nodes, int p,
181  double *N, double *diffN,
182  double *phi_f,
183  double *diff_phi_f, int gdim,
184  int nb);
185 
186 /**
187  * \brirf HDiv base finuctions in tetrahedral interior by Demkowicz \cite
188  * fuentes2015orientation
189  *
190  * @param p Approximation order
191  * @param N Shape functions
192  * @param diffN Derivatives of base functions
193  * @param p_face Max order on faces
194  * @param phi_f Precalculated face base functions
195  * @param diff_phi_f Precalculated derivatives of face base functions
196  * @param phi_v Returned base functions in volume
197  * @param diff_phi_v Returned derivatives of base functions in volume
198  * @param gdim Number of integration points
199  * @return error code
200  */
201 MoFEMErrorCode Hdiv_Demkowicz_Interior_MBTET(int p, double *N, double *diffN,
202  int p_face[], double *phi_f[4],
203  double *diff_phi_f[4],
204  double *phi_v, double *diff_phi_v,
205  int gdim);
206 
207 } // namespace MoFEM
208 
209 #endif // __HDIV_HPP__
MoFEM::Hdiv_Ainsworth_FaceBasedVolumeShapeFunctions_MBTET
MoFEMErrorCode Hdiv_Ainsworth_FaceBasedVolumeShapeFunctions_MBTET(int p, double *N, double *diffN, double *phi_v_f[], double *diff_phi_v_f[], int gdim, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Definition: Hdiv.cpp:401
MoFEM::Hdiv_Demkowicz_Face_MBTET_ON_FACE
MoFEMErrorCode Hdiv_Demkowicz_Face_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int gdim, int nb)
Definition: Hdiv.cpp:634
MoFEM::AinsworthOrderHooks::broken_nbvolumetet_volume_hdiv
static boost::function< int(int)> broken_nbvolumetet_volume_hdiv
Definition: Hdiv.hpp:29
MoFEM::Hdiv_Ainsworth_EdgeFaceShapeFunctions_MBTET
MoFEMErrorCode Hdiv_Ainsworth_EdgeFaceShapeFunctions_MBTET(int *faces_nodes, int *p, double *N, double *diffN, double *phi_f_e[4][3], double *diff_phi_f_e[4][3], int gdim, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Hdiv base functions, Edge-based face functions by Ainsworth .
Definition: Hdiv.cpp:24
MoFEM::Hdiv_Ainsworth_EdgeFaceShapeFunctions_MBTET_ON_FACE
MoFEMErrorCode Hdiv_Ainsworth_EdgeFaceShapeFunctions_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f_e[3], double *diff_phi_f_e[3], int gdim, int nb, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Hdiv base functions, Edge-based face functions by Ainsworth .
Definition: Hdiv.cpp:47
MoFEM::AinsworthOrderHooks::broken_nbvolumetet_face_hdiv
static boost::function< int(int)> broken_nbvolumetet_face_hdiv
Definition: Hdiv.hpp:28
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Hdiv_Demkowicz_Interior_MBTET
MoFEMErrorCode Hdiv_Demkowicz_Interior_MBTET(int p, double *N, double *diffN, int p_face[], double *phi_f[4], double *diff_phi_f[4], double *phi_v, double *diff_phi_v, int gdim)
Definition: Hdiv.cpp:780
MoFEM::AinsworthOrderHooks::broken_nbfacetri_edge_hdiv
static boost::function< int(int)> broken_nbfacetri_edge_hdiv
Definition: Hdiv.hpp:25
MoFEM::AinsworthOrderHooks::broken_nbvolumetet_edge_hdiv
static boost::function< int(int)> broken_nbvolumetet_edge_hdiv
Definition: Hdiv.hpp:27
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
MoFEM::Hdiv_Ainsworth_FaceBubbleShapeFunctions_ON_FACE
MoFEMErrorCode Hdiv_Ainsworth_FaceBubbleShapeFunctions_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int gdim, int nb, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face bubble functions by Ainsworth .
Definition: Hdiv.cpp:174
MoFEM::AinsworthOrderHooks::broken_nbfacetri_face_hdiv
static boost::function< int(int)> broken_nbfacetri_face_hdiv
Definition: Hdiv.hpp:26
N
const int N
Definition: speed_test.cpp:3
MoFEM::AinsworthOrderHooks
Broken base Ainsworth subentries order change hooks.
Definition: Hdiv.hpp:24
MoFEM::Hdiv_Ainsworth_VolumeBubbleShapeFunctions_MBTET
MoFEMErrorCode Hdiv_Ainsworth_VolumeBubbleShapeFunctions_MBTET(int p, double *N, double *diffN, double *phi_v, double *diff_phi_v, int gdim, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Interior bubble functions by Ainsworth .
Definition: Hdiv.cpp:507
MoFEM::Hdiv_Ainsworth_FaceBubbleShapeFunctions
MoFEMErrorCode Hdiv_Ainsworth_FaceBubbleShapeFunctions(int *faces_nodes, int *p, double *N, double *diffN, double *phi_f[], double *diff_phi_f[], int gdim, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Face bubble functions by Ainsworth .
Definition: Hdiv.cpp:151
MoFEM::Hdiv_Ainsworth_EdgeBasedVolumeShapeFunctions_MBTET
MoFEMErrorCode Hdiv_Ainsworth_EdgeBasedVolumeShapeFunctions_MBTET(int p, double *N, double *diffN, double *phi_v_e[6], double *diff_phi_v_e[6], int gdim, PetscErrorCode(*base_polynomials)(int p, double s, double *diff_s, double *L, double *diffL, const int dim))
Hdiv base function, Edge-based interior (volume) functions by Ainsworth .
Definition: Hdiv.cpp:307
convert.int
int
Definition: convert.py:64