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 /**
17  * \brief Hdiv base functions, Edge-based face functions by Ainsworth \cite
18  * NME:NME847
19  * @param faces_nodes Face nodes on tetrahedral
20  * @param p Approximation order on faces
21  * @param N Shape functions
22  * @param diffN Derivatives of shape functions
23  * @param phi_f_e Base functions (returned)
24  * @param diff_phi_f_e Derivatives of base functions (returned)
25  * @param gdim Number of integration pts
26  * @param base_polynomials base function (Legendre/Lobbatto-Gauss)
27  * @return error code
28  */
30  int *faces_nodes, int *p, double *N, double *diffN, double *phi_f_e[4][3],
31  double *diff_phi_f_e[4][3], int gdim,
32  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
33  double *L, double *diffL,
34  const int dim));
35 
36 /**
37  * \brief Hdiv base functions, Edge-based face functions by Ainsworth \cite
38  * NME:NME847
39  * @param faces_nodes Face nodes on face
40  * @param p Approximation order on faces
41  * @param N Shape functions
42  * @param diffN Derivatives of shape functions
43  * @param phi_f_e Base functions (returned)
44  * @param diff_phi_f_e Derivatives of base functions (returned)
45  * @param gdim Number of integration pts
46  * @param base_polynomials base function (Legendre/Lobbatto-Gauss)
47  * @param nb Number of nodes on entity (4 if tet, 3 if triangle)
48  * @return error code
49  */
51  int *faces_nodes, int p, double *N, double *diffN, double *phi_f_e[3],
52  double *diff_phi_f_e[3], int gdim, int nb,
53  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
54  double *L, double *diffL,
55  const int dim));
56 
57 /**
58  * \brief Face bubble functions by Ainsworth \cite NME:NME847
59  * @param faces_nodes Face nodes on tetrahedral
60  * @param p Approx. order on faces
61  * @param N Shape function
62  * @param diffN Derivatives of shape functions
63  * @param phi_f Base functions
64  * @param diff_phi_f Derivatives of base functions
65  * @param gdim Number of integration pts
66  * @param base_polynomials Base function (Legendre/Lobbatto-Gauss)
67  * @return error code
68  */
70  int *faces_nodes, int *p, double *N, double *diffN, double *phi_f[],
71  double *diff_phi_f[], int gdim,
72  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
73  double *L, double *diffL,
74  const int dim));
75 
76 /**
77  * \brief Face bubble functions by Ainsworth \cite NME:NME847
78  * @param faces_nodes Face nodes on face
79  * @param p Approx. order on face
80  * @param N Shape function
81  * @param diffN Derivatives of shape functions
82  * @param phi_f Base functions
83  * @param diff_phi_f Derivatives of base functions
84  * @param gdim Number of integration pts
85  * @param nb Number of nodes on entity (4 if tet, 3 if triangle)
86  * @param base_polynomials Base function (Legendre/Lobbatto-Gauss)
87  * @return error code
88  */
90  int *faces_nodes, int p, double *N, double *diffN, double *phi_f,
91  double *diff_phi_f, int gdim, int nb,
92  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
93  double *L, double *diffL,
94  const int dim));
95 
96 /**
97  * \brief Hdiv base function, Edge-based interior (volume) functions by
98  * Ainsworth \cite NME:NME847
99  * @param p volume approx. order
100  * @param N shape functions
101  * @param diffN derivatives of shape functions
102  * @param phi_v_e base functions (return)
103  * @param diff_phi_v_e derivatives of base functions (returned)
104  * @param gdim number of integration points
105  * @param base_polynomials base function (Legendre/Lobbatto-Gauss)
106  * @return error code
107  */
109  int p, double *N, double *diffN, double *phi_v_e[6],
110  double *diff_phi_v_e[6], int gdim,
111  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
112  double *L, double *diffL,
113  const int dim));
114 
115 /** Hdiv Face-based interior functions by Ainsworth \cite NME:NME847
116  * @param p Approximation order on face
117  * @param N Shape functions on face
118  * @param diffN Derivatives of shape functions of face
119  * @param phi_v_f Base functions (returned)
120  * @param diff_phi_v_f Derivatives of base functions (returned)
121  * @param gdim Number of integration points
122  * @param base_polynomials Base function (Legendre/Lobbatto-Gauss)
123  * @return Error code
124  */
126  int p, double *N, double *diffN, double *phi_v_f[], double *diff_phi_v_f[],
127  int gdim,
128  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
129  double *L, double *diffL,
130  const int dim));
131 
132 /**
133  * \brief Interior bubble functions by Ainsworth \cite NME:NME847
134  * @param p Volume order
135  * @param N Shape functions
136  * @param diffN Derivatives of shape functions
137  * @param phi_v Base functions
138  * @param diff_phi_v Derivatives of shape functions
139  * @param gdim Number of integration points
140  * @param base_polynomials Base function (Legendre/Lobbatto-Gauss)
141  * @return Error code
142  */
144  int p, double *N, double *diffN, double *phi_v, double *diff_phi_v,
145  int gdim,
146  PetscErrorCode (*base_polynomials)(int p, double s, double *diff_s,
147  double *L, double *diffL,
148  const int dim));
149 
150 /**
151  * \brirf HDiv base finuctions on triangle by Demkowicz \cite
152  * fuentes2015orientation
153  *
154  * @param faces_nodes Nodes on face
155  * @param p Approx. order
156  * @param N Shape functions
157  * @param diffN Derivatives of shape functions
158  * @param phi_f Returned base functions on face
159  * @param diff_phi_f Returned derivatives of base functions on face
160  * @param gdim Number of integration points
161  * @param nb nb is 4 for face on tetrahedral and 3 for face
162  * @return error code
163  */
164 MoFEMErrorCode Hdiv_Demkowicz_Face_MBTET_ON_FACE(int *faces_nodes, int p,
165  double *N, double *diffN,
166  double *phi_f,
167  double *diff_phi_f, int gdim,
168  int nb);
169 
170 /**
171  * \brirf HDiv base finuctions in tetrahedral interior by Demkowicz \cite
172  * fuentes2015orientation
173  *
174  * @param p Approximation order
175  * @param N Shape functions
176  * @param diffN Derivatives of base functions
177  * @param p_face Max order on faces
178  * @param phi_f Precalculated face base functions
179  * @param diff_phi_f Precalculated derivatives of face base functions
180  * @param phi_v Returned base functions in volume
181  * @param diff_phi_v Returned derivatives of base functions in volume
182  * @param gdim Number of integration points
183  * @return error code
184  */
185 MoFEMErrorCode Hdiv_Demkowicz_Interior_MBTET(int p, double *N, double *diffN,
186  int p_face[], double *phi_f[4],
187  double *diff_phi_f[4],
188  double *phi_v, double *diff_phi_v,
189  int gdim);
190 
191 } // namespace MoFEM
192 
193 #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:368
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:617
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:13
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:37
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:762
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:160
N
const int N
Definition: speed_test.cpp:3
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:484
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:137
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:280