v0.14.0
fem_tools.h
Go to the documentation of this file.
1 /** \file fem_tools.h
2  * \brief Loose implementation of some useful functions
3  *
4  * FIXME: Implementation here is very unstructured, need cleaning and pruning
5  *
6  */
7 
8 #ifndef __FEM_H__
9 #define __FEM_H__
10 
11 #include <petscsys.h>
12 
13 #ifdef __cplusplus
14 extern "C" {
15 #endif
16 #include <cblas.h>
17 #include <lapack_wrap.h>
18 #ifdef __cplusplus
19 }
20 #endif
21 
22 #define LAMBDA(E, NU) (E * NU / ((1. + NU) * (1. - 2. * NU)))
23 #define MU(E, NU) (0.5 * E / (1. + NU))
24 #define DELTA(NU_P, NU_PZ, E_P, E_Z) \
25  (((1 + NU_P) * (1 - NU_P - 2 * NU_PZ * (NU_PZ * E_Z / E_P))) / \
26  (E_P * E_P * E_Z))
27 
28 #define N_MBTET0(x, y, z) (1. - x - y - z) ///< tetrahedral shape function
29 #define N_MBTET1(x, y, z) (x) ///< tetrahedral shape function
30 #define N_MBTET2(x, y, z) (y) ///< tetrahedral shape function
31 #define N_MBTET3(x, y, z) (z) ///< tetrahedral shape function
32 #define diffN_MBTET0x (-1.) ///< derivative of tetrahedral shape function
33 #define diffN_MBTET0y (-1.) ///< derivative of tetrahedral shape function
34 #define diffN_MBTET0z (-1.) ///< derivative of tetrahedral shape function
35 #define diffN_MBTET1x (1.) ///< derivative of tetrahedral shape function
36 #define diffN_MBTET1y (0.) ///< derivative of tetrahedral shape function
37 #define diffN_MBTET1z (0.) ///< derivative of tetrahedral shape function
38 #define diffN_MBTET2x (0.) ///< derivative of tetrahedral shape function
39 #define diffN_MBTET2y (1.) ///< derivative of tetrahedral shape function
40 #define diffN_MBTET2z (0.) ///< derivative of tetrahedral shape function
41 #define diffN_MBTET3x (0.) ///< derivative of tetrahedral shape function
42 #define diffN_MBTET3y (0.) ///< derivative of tetrahedral shape function
43 #define diffN_MBTET3z (1.) ///< derivative of tetrahedral shape function
44 
45 // MBTRI
46 #define N_MBTRI0(x, y) (1. - x - y) ///< triangle shape function
47 #define N_MBTRI1(x, y) (x) ///< triangle shape function
48 #define N_MBTRI2(x, y) (y) ///< triangle shape function
49 #define diffN_MBTRI0x (-1.) ///< derivative of triangle shape function
50 #define diffN_MBTRI0y (-1.) ///< derivative of triangle shape function
51 #define diffN_MBTRI1x (1.) ///< derivative of triangle shape function
52 #define diffN_MBTRI1y (0.) ///< derivative of triangle shape function
53 #define diffN_MBTRI2x (0.) ///< derivative of triangle shape function
54 #define diffN_MBTRI2y (1.) ///< derivative of triangle shape function
55 
56 // MBQUAD
57 #define N_MBQUAD0(x, y) ((1. - x) * (1. - y)) ///< quad shape function
58 #define N_MBQUAD1(x, y) ((x) * (1. - y)) ///< quad shape function
59 #define N_MBQUAD2(x, y) ((x) * (y)) ///< quad shape function
60 #define N_MBQUAD3(x, y) ((1. - x) * (y)) ///< quad shape function
61 #define diffN_MBQUAD0x(y) (-(1. - y))
62 #define diffN_MBQUAD0y(x) (-(1. - x))
63 #define diffN_MBQUAD1x(y) ((1. - y))
64 #define diffN_MBQUAD1y(x) (-x)
65 #define diffN_MBQUAD2x(y) (y)
66 #define diffN_MBQUAD2y(x) (x)
67 #define diffN_MBQUAD3x(y) (-y)
68 #define diffN_MBQUAD3y(x) ((1. - x))
69 
70 // MBHEX
71 #define N_MBHEX0(x, y, z) (N_MBQUAD0(x, y) * (1 - z))
72 #define N_MBHEX1(x, y, z) (N_MBQUAD1(x, y) * (1 - z))
73 #define N_MBHEX2(x, y, z) (N_MBQUAD2(x, y) * (1 - z))
74 #define N_MBHEX3(x, y, z) (N_MBQUAD3(x, y) * (1 - z))
75 #define N_MBHEX4(x, y, z) (N_MBQUAD0(x, y) * (z))
76 #define N_MBHEX5(x, y, z) (N_MBQUAD1(x, y) * (z))
77 #define N_MBHEX6(x, y, z) (N_MBQUAD2(x, y) * (z))
78 #define N_MBHEX7(x, y, z) (N_MBQUAD3(x, y) * (z))
79 #define diffN_MBHEX0x(y, z) (diffN_MBQUAD0x(y) * (1 - z))
80 #define diffN_MBHEX1x(y, z) (diffN_MBQUAD1x(y) * (1 - z))
81 #define diffN_MBHEX2x(y, z) (diffN_MBQUAD2x(y) * (1 - z))
82 #define diffN_MBHEX3x(y, z) (diffN_MBQUAD3x(y) * (1 - z))
83 #define diffN_MBHEX4x(y, z) (diffN_MBQUAD0x(y) * (z))
84 #define diffN_MBHEX5x(y, z) (diffN_MBQUAD1x(y) * (z))
85 #define diffN_MBHEX6x(y, z) (diffN_MBQUAD2x(y) * (z))
86 #define diffN_MBHEX7x(y, z) (diffN_MBQUAD3x(y) * (z))
87 #define diffN_MBHEX0y(x, z) (diffN_MBQUAD0y(x) * (1 - z))
88 #define diffN_MBHEX1y(x, z) (diffN_MBQUAD1y(x) * (1 - z))
89 #define diffN_MBHEX2y(x, z) (diffN_MBQUAD2y(x) * (1 - z))
90 #define diffN_MBHEX3y(x, z) (diffN_MBQUAD3y(x) * (1 - z))
91 #define diffN_MBHEX4y(x, z) (diffN_MBQUAD0y(x) * (z))
92 #define diffN_MBHEX5y(x, z) (diffN_MBQUAD1y(x) * (z))
93 #define diffN_MBHEX6y(x, z) (diffN_MBQUAD2y(x) * (z))
94 #define diffN_MBHEX7y(x, z) (diffN_MBQUAD3y(x) * (z))
95 #define diffN_MBHEX0z(x, y) (-N_MBQUAD0(x, y))
96 #define diffN_MBHEX1z(x, y) (-N_MBQUAD1(x, y))
97 #define diffN_MBHEX2z(x, y) (-N_MBQUAD2(x, y))
98 #define diffN_MBHEX3z(x, y) (-N_MBQUAD3(x, y))
99 #define diffN_MBHEX4z(x, y) (N_MBQUAD0(x, y))
100 #define diffN_MBHEX5z(x, y) (N_MBQUAD1(x, y))
101 #define diffN_MBHEX6z(x, y) (N_MBQUAD2(x, y))
102 #define diffN_MBHEX7z(x, y) (N_MBQUAD3(x, y))
103 
104 // MBEDGE
105 #define N_MBEDGE0(x) (1. - (x)) ///< edge shape function
106 #define N_MBEDGE1(x) (x) ///< edge shape function
107 #define diffN_MBEDGE0 (-1.) ///< derivative of edge shape function
108 #define diffN_MBEDGE1 (1.) ///< derivative of edge shape function
109 
110 // MBTRIQ
111 #define N_MBTRIQ0(x, y) ((1. - x - y) * (2 * (1. - x - y) - 1.))
112 #define N_MBTRIQ1(x, y) (x * (2. * x - 1.))
113 #define N_MBTRIQ2(x, y) (y * (2. * y - 1.))
114 #define N_MBTRIQ3(x, y) (4. * (1. - x - y) * x)
115 #define N_MBTRIQ4(x, y) (4. * x * y)
116 #define N_MBTRIQ5(x, y) (4. * (1. - x - y) * y)
117 #define diffN_MBTRIQ0x(x, y) (x + y - 3. * (1. - x - y))
118 #define diffN_MBTRIQ0y(x, y) (x + y - 3. * (1. - x - y))
119 #define diffN_MBTRIQ1x(x, y) (-1. + 4. * x)
120 #define diffN_MBTRIQ1y(x, y) (0.)
121 #define diffN_MBTRIQ2x(x, y) (0.)
122 #define diffN_MBTRIQ2y(x, y) (-1. + 4. * y)
123 #define diffN_MBTRIQ3x(x, y) (4. * ((1. - x - y) - x))
124 #define diffN_MBTRIQ3y(x, y) (-4. * x)
125 #define diffN_MBTRIQ4x(x, y) (4. * y)
126 #define diffN_MBTRIQ4y(x, y) (4. * x)
127 #define diffN_MBTRIQ5x(x, y) (-4. * y)
128 #define diffN_MBTRIQ5y(x, y) (4. * ((1. - x - y) - y))
129 
130 #ifdef __cplusplus
131 extern "C" {
132 #endif
133 
134 /// print matric M
135 void print_mat(double *M, int m, int n);
136 /// print upper part of the symmetric matrix
137 void print_mat_sym_upper(double *M, int m, int n);
138 /// priint complex matrix
139 void print_mat_complex(__CLPK_doublecomplex *M, int m, int n);
140 
141 /// \brief calculate shape functions on triangle
142 /// \param N shape function array
143 /// \param X array of Gauss X coordinates
144 /// \param Y array of Gauss Y coordinates
145 /// \param G_DIM number of Gauss points
146 PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y,
147  const int G_DIM);
148 /// calculate derivatives of shape functions
149 PetscErrorCode ShapeDiffMBTRI(double *diffN);
150 
151 /// calculate face normal
152 /// \param diffN derivatives of shape functions
153 /// \param coords is position of the nodes
154 /// \param normal vector
155 PetscErrorCode ShapeFaceNormalMBTRI(double *diffN, const double *coords,
156  double *normal);
157 PetscErrorCode ShapeFaceBaseMBTRI(double *diffN, const double *coords,
158  double *normal, double *s1, double *s2);
159 
160 /// calculate derivative of normal in respect to nodal positions
161 PetscErrorCode ShapeFaceDiffNormalMBTRI(double *diffN, const double *coords,
162  double *diff_normal);
163 /// calculate jacobioan
164 void ShapeJacMBTRI(double *diffN, const double *coords, double *Jac);
165 /// calculate derivatives of shape functions in space
166 void ShapeDiffMBTRIinvJ(double *diffN, double *invJac, double *diffNinvJac);
167 /// calculate shape functions
168 PetscErrorCode ShapeMBTET(double *N, const double *G_X, const double *G_Y,
169  const double *G_Z, int DIM);
170 /// calculate derivatives of shape functions
171 PetscErrorCode ShapeDiffMBTET(double *diffN);
172 /// determined of jacobian
173 double ShapeDetJacVolume(double *jac);
174 /// calculate jacobian
175 PetscErrorCode ShapeJacMBTET(double *diffN, const double *coords, double *jac);
176 // calculate inverse of jacobian
177 PetscErrorCode ShapeInvJacVolume(double *jac);
178 /// calculate TET volume
179 double ShapeVolumeMBTET(double *diffN, const double *coords);
180 /// calculate shape functions derivatives in space
181 PetscErrorCode ShapeDiffMBTETinvJ(double *diffN, double *invJac,
182  double *diffNinvJac);
183 
184 /// calculate spin matrix from vector
185 // \param spinOmega is a spin matrix
186 // \param vecOmega is a spin vector
187 PetscErrorCode Spin(double *spinOmega, double *vecOmega);
188 
189 /// Compose complex matrix (3x3) from two real matrices
190 PetscErrorCode make_complex_matrix(double *reA, double *imA,
192 /// Complex normal
193 PetscErrorCode
194 Normal_hierarchical(int order_approx, int *order_edge_approx, int order,
195  int *order_edge, double *diffN, double *diffN_face,
196  double *diffN_edge[], double *dofs, double *dofs_edge[],
197  double *dofs_face, double *idofs, double *idofs_edge[],
198  double *idofs_face, __CLPK_doublecomplex *xnormal,
199  __CLPK_doublecomplex *s1, __CLPK_doublecomplex *s2, int gg);
200 PetscErrorCode Base_scale(__CLPK_doublecomplex *xnormal,
202 
203 /**
204  * \brief calculate local coordinates for given global coordinates
205  *
206  * new version for multiple points need to be implemented
207  */
208 PetscErrorCode ShapeMBTET_inverse(double *N, double *diffN,
209  const double *elem_coords,
210  const double *glob_coords,
211  double *loc_coords);
212 
213 /**
214  * \brief calculate local coordinates of triangle element for given global
215 coordinates in 2D (Assume e.g. z=0) \f[ \left[\begin{array}{cc}
216 \frac{\partial N_{1}}{\partial\xi}x_{N_{1}}+\frac{\partial
217 N_{2}}{\partial\xi}x_{N_{2}}+\frac{\partial N_{3}}{\partial\xi}x_{N_{3}} &
218 \frac{\partial N_{1}}{\partial\eta}x_{N_{1}}+\frac{\partial
219 N_{2}}{\partial\eta}x_{N_{2}}+\frac{\partial N_{3}}{\partial\eta}x_{N_{3}}\\
220 \frac{\partial N_{1}}{\partial\xi}y_{N_{1}}+\frac{\partial
221 N_{2}}{\partial\xi}y_{N_{2}}+\frac{\partial N_{3}}{\partial\xi}y_{N_{3}} &
222 \frac{\partial N_{1}}{\partial\eta}y_{N_{1}}+\frac{\partial
223 N_{2}}{\partial\eta}y_{N_{2}}+\frac{\partial N_{3}}{\partial\eta}y_{N_{3}}
224 \end{array}\right]\left\{ \begin{array}{c}
225 \xi\\
226 \eta
227 \end{array}\right\} =\left\{ \begin{array}{c}
228 x_{gp}-\left(N_{1}x_{N_{1}}+N_{2}x_{N_{2}}+N_{3}x_{N_{3}}\right)\\
229 y_{gp}-\left(N_{1}y_{N_{1}}+N_{2}y_{N_{2}}+N_{3}y_{N_{3}}\right)
230 \end{array}\right\}
231  \f]
232 
233  /// \param N shape function array
234  /// \param diffN array of shape function derivative w.r.t to local coordinates
235  /// \param elem_coords global coordinates of element
236  /// \param glob_coords global coordinates of required point
237  /// \param loc_coords local coordinates of required point
238  */
239 PetscErrorCode ShapeMBTRI_inverse(double *N, double *diffN,
240  const double *elem_coords,
241  const double *glob_coords,
242  double *loc_coords);
243 
244 /// calculate gradient of deformation
245 PetscErrorCode GradientOfDeformation(double *diffN, double *dofs, double *F);
246 
247 // 2 Node edge
248 PetscErrorCode ShapeMBEDGE(double *N, const double *G_X, int DIM);
249 PetscErrorCode ShapeDiffMBEDGE(double *diffN);
250 
251 // 10 Node Tet
252 PetscErrorCode ShapeMBTRIQ(double *N, const double *X, const double *Y,
253  const int G_DIM);
254 PetscErrorCode ShapeDiffMBTRIQ(double *diffN, const double *X, const double *Y,
255  const int G_DIM);
256 PetscErrorCode ShapeMBTETQ(double *N, const double x, const double y,
257  const double z);
258 PetscErrorCode ShapeDiffMBTETQ(double *diffN, const double x, const double y,
259  const double z);
260 PetscErrorCode ShapeMBTETQ_GAUSS(double *N, const double *X, const double *Y,
261  const double *Z, const int G_DIM);
262 PetscErrorCode ShapeDiffMBTETQ_GAUSS(double *diffN, const double *X,
263  const double *Y, const double *Z,
264  const int G_DIM);
265 PetscErrorCode ShapeJacMBTETQ(const double *diffN, const double *coords,
266  double *Jac);
267 PetscErrorCode
268 ShapeMBTETQ_detJac_at_Gauss_Points(double *detJac_at_Gauss_Points,
269  const double *diffN, const double *coords,
270  int G_DIM);
271 double ShapeVolumeMBTETQ(const double *diffN, const double *coords, int G_DIM,
272  double *G_TET_W);
273 PetscErrorCode ShapeMBTETQ_inverse(double *N, double *diffN,
274  const double *elem_coords,
275  const double *glob_coords,
276  double *loc_coords, const double eps);
277 
278 // complex part
280  __CLPK_doublecomplex *diffNinvJac,
281  enum CBLAS_TRANSPOSE Trans);
282 PetscErrorCode ShapeFaceNormalMBTRI_complex(double *diffN,
283  __CLPK_doublecomplex *xcoords,
284  __CLPK_doublecomplex *xnormal);
285 PetscErrorCode MakeComplexTensor(double *reA, double *imA,
287 PetscErrorCode InvertComplexGradient(__CLPK_doublecomplex *xF);
289  __CLPK_doublecomplex *det_xF);
290 
291 // integration
292 /// Compute weights and integration points for edge using Grundmann_Moeller rule
293 PetscErrorCode Grundmann_Moeller_integration_points_1D_EDGE(int rule,
294  double *G_TRI_X,
295  double *G_TRI_W);
296 /// Compute weights and integration points for 2D Triangle using
297 /// Grundmann_Moeller rule
298 PetscErrorCode Grundmann_Moeller_integration_points_2D_TRI(int rule,
299  double *G_TRI_X,
300  double *G_TRI_Y,
301  double *G_TRI_W);
302 /// Compute weights and integration points for 3D Tet using Grundmann_Moeller
303 /// rule
304 PetscErrorCode Grundmann_Moeller_integration_points_3D_TET(int rule,
305  double *G_TET_X,
306  double *G_TET_Y,
307  double *G_TET_Z,
308  double *G_TET_W);
309 
310 // http://people.sc.fsu.edu/~jburkardt/datasets/quadrature_rules_tet/quadrature_rules_tet.html
311 // TRI
312 static const double G_TRI_X1[] = {3.3333333333333331e-01};
313 static const double G_TRI_Y1[] = {3.3333333333333331e-01};
314 static const double G_TRI_W1[] = {1.};
315 static const double G_TRI_X3[] = {0.5, 0., 0.5};
316 static const double G_TRI_Y3[] = {0., 0.5, 0.5};
317 static const double G_TRI_W3[] = {
318  3.3333333333333331e-01, 3.3333333333333331e-01, 3.3333333333333331e-01};
319 static const double G_TRI_X4[] = {
320  7.503111022260811058e-02, 1.785587282636164064e-01,
321  2.800199154990741235e-01, 6.663902460147014262e-01};
322 static const double G_TRI_Y4[] = {
323  2.800199154990741235e-01, 6.663902460147014262e-01,
324  7.503111022260811058e-02, 1.785587282636164064e-01};
325 static const double G_TRI_W4[] = {
326  1.8195861825602258066e-01, 3.1804138174397683647e-01,
327  1.8195861825602258066e-01, 3.1804138174397683647e-01};
328 static const double G_TRI_X7[] = {
329  0.333333333333333, 0.736712498968435, 0.736712498968435, 0.237932366472434,
330  0.237932366472434, 0.025355134551932, 0.025355134551932};
331 static const double G_TRI_Y7[] = {
332  0.333333333333333, 0.237932366472434, 0.025355134551932, 0.736712498968435,
333  0.025355134551932, 0.736712498968435, 0.237932366472434};
334 static const double G_TRI_W7[] = {
335  0.375000000000000, 0.104166666666667, 0.104166666666667, 0.104166666666667,
336  0.104166666666667, 0.104166666666667, 0.104166666666667};
337 static const double G_TRI_X13[] = {
338  0.333333333333333, 0.479308067841923, 0.260345966079038, 0.260345966079038,
339  0.869739794195568, 0.065130102902216, 0.065130102902216, 0.638444188569809,
340  0.638444188569809, 0.312865496004875, 0.312865496004875, 0.048690315425316,
341  0.048690315425316};
342 static const double G_TRI_Y13[] = {
343  0.333333333333333, 0.260345966079038, 0.479308067841923, 0.260345966079038,
344  0.065130102902216, 0.869739794195568, 0.065130102902216, 0.312865496004875,
345  0.048690315425316, 0.638444188569809, 0.048690315425316, 0.638444188569809,
346  0.312865496004875};
347 static const double G_TRI_W13[] = {
348  -0.149570044467670, 0.175615257433204, 0.175615257433204, 0.175615257433204,
349  0.053347235608839, 0.053347235608839, 0.053347235608839, 0.077113760890257,
350  0.077113760890257, 0.077113760890257, 0.077113760890257, 0.077113760890257,
351  0.077113760890257};
352 
353 static const double G_TRI_X19[] = {
354  0.333333333333333, 0.797426985353087, 0.101286507323456, 0.101286507323456,
355  0.059715871789770, 0.470142064105115, 0.470142064105115, 0.535795346449899,
356  0.232102326775050, 0.232102326775050, 0.941038278231121, 0.029480860884440,
357  0.029480860884440, 0.738416812340510, 0.738416812340510, 0.232102326775050,
358  0.232102326775050, 0.029480860884440, 0.029480860884440};
359 static const double G_TRI_Y19[] = {
360  0.333333333333333, 0.101286507323456, 0.797426985353087, 0.101286507323456,
361  0.470142064105115, 0.059715871789770, 0.470142064105115, 0.232102326775050,
362  0.535795346449899, 0.232102326775050, 0.029480860884440, 0.941038278231121,
363  0.029480860884440, 0.232102326775050, 0.029480860884440, 0.738416812340510,
364  0.029480860884440, 0.738416812340510, 0.232102326775050};
365 static const double G_TRI_W19[] = {
366  9.71357962827961025E-002, 3.13347002271398278E-002,
367  3.13347002271398278E-002, 3.13347002271398278E-002,
368  7.78275410047754301E-002, 7.78275410047754301E-002,
369  7.78275410047754301E-002, 7.96477389272090969E-002,
370  7.96477389272090969E-002, 7.96477389272090969E-002,
371  2.55776756586981006E-002, 2.55776756586981006E-002,
372  2.55776756586981006E-002, 4.32835393772893970E-002,
373  4.32835393772893970E-002, 4.32835393772893970E-002,
374  4.32835393772893970E-002, 4.32835393772893970E-002,
375  4.32835393772893970E-002};
376 
377 static const double G_TRI_X28[] = {
378  0.333333333333333, 0.948021718143423, 0.025989140928288,
379  0.025989140928288, 0.811424994704155, 0.094287502647923,
380  0.094287502647923, 0.010726449965571, 0.494636775017215,
381  0.494636775017215, 0.585313234770972, 0.207343382614514,
382  0.207343382614514, 0.122184388599019, 0.438907805700491,
383  0.438907805700491, 0.677937654882590, 0.677937654882590,
384  0.044841677589131, 0.044841677589131, 0.277220667528279,
385  0.277220667528279, 0.858870281282636, 0.858870281282636,
386  0.0000000000000000, 0.0000000000000000, 0.141129718717364,
387  0.141129718717364};
388 static const double G_TRI_Y28[] = {
389  0.333333333333333, 0.025989140928288, 0.948021718143423, 0.025989140928288,
390  0.094287502647923, 0.811424994704155, 0.094287502647923, 0.494636775017215,
391  0.010726449965571, 0.494636775017215, 0.207343382614514, 0.585313234770972,
392  0.207343382614514, 0.438907805700491, 0.122184388599019, 0.438907805700491,
393  0.044841677589131, 0.277220667528279, 0.677937654882590, 0.277220667528279,
394  0.677937654882590, 0.044841677589131, 0.000000000000000, 0.141129718717364,
395  0.858870281282636, 0.141129718717364, 0.858870281282636, 0.000000000000000};
396 static const double G_TRI_W28[] = {
397  0.08797730116222190, 0.008744311553736190, 0.008744311553736190,
398  0.008744311553736190, 0.03808157199393533, 0.03808157199393533,
399  0.03808157199393533, 0.01885544805613125, 0.01885544805613125,
400  0.01885544805613125, 0.07215969754474100, 0.07215969754474100,
401  0.07215969754474100, 0.06932913870553720, 0.06932913870553720,
402  0.06932913870553720, 0.04105631542928860, 0.04105631542928860,
403  0.04105631542928860, 0.04105631542928860, 0.04105631542928860,
404  0.04105631542928860, 0.007362383783300573, 0.007362383783300573,
405  0.007362383783300573, 0.007362383783300573, 0.007362383783300573,
406  0.007362383783300573};
407 
408 static const double G_TRI_X37[] = {
409  0.333333333333333, 0.950275662924106, 0.024862168537947, 0.024862168537947,
410  0.171614914923835, 0.414192542538082, 0.414192542538082, 0.539412243677190,
411  0.230293878161405, 0.230293878161405, 0.772160036676533, 0.113919981661734,
412  0.113919981661734, 0.009085399949835, 0.495457300025082, 0.495457300025082,
413  0.062277290305887, 0.468861354847056, 0.468861354847056, 0.022076289653624,
414  0.022076289653624, 0.851306504174348, 0.851306504174348, 0.126617206172027,
415  0.126617206172027, 0.018620522802521, 0.018620522802521, 0.689441970728591,
416  0.689441970728591, 0.291937506468888, 0.291937506468888, 0.096506481292159,
417  0.096506481292159, 0.635867859433873, 0.635867859433873, 0.267625659273968,
418  0.267625659273968};
419 static const double G_TRI_Y37[] = {
420  0.333333333333333, 0.024862168537947, 0.950275662924106, 0.024862168537947,
421  0.414192542538082, 0.171614914923835, 0.414192542538082, 0.230293878161405,
422  0.539412243677190, 0.230293878161405, 0.113919981661734, 0.772160036676533,
423  0.113919981661734, 0.495457300025082, 0.009085399949835, 0.495457300025082,
424  0.468861354847056, 0.062277290305887, 0.468861354847056, 0.851306504174348,
425  0.126617206172027, 0.022076289653624, 0.126617206172027, 0.022076289653624,
426  0.851306504174348, 0.689441970728591, 0.291937506468888, 0.018620522802521,
427  0.291937506468888, 0.018620522802521, 0.689441970728591, 0.635867859433873,
428  0.267625659273968, 0.096506481292159, 0.267625659273968, 0.096506481292159,
429  0.635867859433873};
430 static const double G_TRI_W37[] = {
431  0.051739766065744, 0.008007799555565, 0.008007799555565, 0.008007799555565,
432  0.046868898981822, 0.046868898981822, 0.046868898981822, 0.046590940183976,
433  0.046590940183976, 0.046590940183976, 0.031016943313796, 0.031016943313796,
434  0.031016943313796, 0.010791612736631, 0.010791612736631, 0.010791612736631,
435  0.032195534242432, 0.032195534242432, 0.032195534242432, 0.015445834210702,
436  0.015445834210702, 0.015445834210702, 0.015445834210702, 0.015445834210702,
437  0.015445834210702, 0.017822989923179, 0.017822989923179, 0.017822989923179,
438  0.017822989923179, 0.017822989923179, 0.017822989923179, 0.037038683681385,
439  0.037038683681385, 0.037038683681385, 0.037038683681385, 0.037038683681385,
440  0.037038683681385};
441 static const double G_TRI_X286[] = {0.04347826086956522,
442  0.1304347826086956,
443  0.2173913043478261,
444  0.3043478260869565,
445  0.391304347826087,
446  0.4782608695652174,
447  0.5652173913043478,
448  0.6521739130434783,
449  0.7391304347826086,
450  0.8260869565217391,
451  0.9130434782608695,
452  0.04347826086956522,
453  0.1304347826086956,
454  0.2173913043478261,
455  0.3043478260869565,
456  0.391304347826087,
457  0.4782608695652174,
458  0.5652173913043478,
459  0.6521739130434783,
460  0.7391304347826086,
461  0.8260869565217391,
462  0.04347826086956522,
463  0.1304347826086956,
464  0.2173913043478261,
465  0.3043478260869565,
466  0.391304347826087,
467  0.4782608695652174,
468  0.5652173913043478,
469  0.6521739130434783,
470  0.7391304347826086,
471  0.04347826086956522,
472  0.1304347826086956,
473  0.2173913043478261,
474  0.3043478260869565,
475  0.391304347826087,
476  0.4782608695652174,
477  0.5652173913043478,
478  0.6521739130434783,
479  0.04347826086956522,
480  0.1304347826086956,
481  0.2173913043478261,
482  0.3043478260869565,
483  0.391304347826087,
484  0.4782608695652174,
485  0.5652173913043478,
486  0.04347826086956522,
487  0.1304347826086956,
488  0.2173913043478261,
489  0.3043478260869565,
490  0.391304347826087,
491  0.4782608695652174,
492  0.04347826086956522,
493  0.1304347826086956,
494  0.2173913043478261,
495  0.3043478260869565,
496  0.391304347826087,
497  0.04347826086956522,
498  0.1304347826086956,
499  0.2173913043478261,
500  0.3043478260869565,
501  0.04347826086956522,
502  0.1304347826086956,
503  0.2173913043478261,
504  0.04347826086956522,
505  0.1304347826086956,
506  0.04347826086956522,
507  0.04761904761904762,
508  0.1428571428571428,
509  0.2380952380952381,
510  0.3333333333333333,
511  0.4285714285714285,
512  0.5238095238095238,
513  0.6190476190476191,
514  0.7142857142857143,
515  0.8095238095238095,
516  0.9047619047619048,
517  0.04761904761904762,
518  0.1428571428571428,
519  0.2380952380952381,
520  0.3333333333333333,
521  0.4285714285714285,
522  0.5238095238095238,
523  0.6190476190476191,
524  0.7142857142857143,
525  0.8095238095238095,
526  0.04761904761904762,
527  0.1428571428571428,
528  0.2380952380952381,
529  0.3333333333333333,
530  0.4285714285714285,
531  0.5238095238095238,
532  0.6190476190476191,
533  0.7142857142857143,
534  0.04761904761904762,
535  0.1428571428571428,
536  0.2380952380952381,
537  0.3333333333333333,
538  0.4285714285714285,
539  0.5238095238095238,
540  0.6190476190476191,
541  0.04761904761904762,
542  0.1428571428571428,
543  0.2380952380952381,
544  0.3333333333333333,
545  0.4285714285714285,
546  0.5238095238095238,
547  0.04761904761904762,
548  0.1428571428571428,
549  0.2380952380952381,
550  0.3333333333333333,
551  0.4285714285714285,
552  0.04761904761904762,
553  0.1428571428571428,
554  0.2380952380952381,
555  0.3333333333333333,
556  0.04761904761904762,
557  0.1428571428571428,
558  0.2380952380952381,
559  0.04761904761904762,
560  0.1428571428571428,
561  0.04761904761904762,
562  0.05263157894736842,
563  0.1578947368421053,
564  0.2631578947368421,
565  0.3684210526315789,
566  0.4736842105263158,
567  0.5789473684210527,
568  0.6842105263157895,
569  0.7894736842105263,
570  0.8947368421052632,
571  0.05263157894736842,
572  0.1578947368421053,
573  0.2631578947368421,
574  0.3684210526315789,
575  0.4736842105263158,
576  0.5789473684210527,
577  0.6842105263157895,
578  0.7894736842105263,
579  0.05263157894736842,
580  0.1578947368421053,
581  0.2631578947368421,
582  0.3684210526315789,
583  0.4736842105263158,
584  0.5789473684210527,
585  0.6842105263157895,
586  0.05263157894736842,
587  0.1578947368421053,
588  0.2631578947368421,
589  0.3684210526315789,
590  0.4736842105263158,
591  0.5789473684210527,
592  0.05263157894736842,
593  0.1578947368421053,
594  0.2631578947368421,
595  0.3684210526315789,
596  0.4736842105263158,
597  0.05263157894736842,
598  0.1578947368421053,
599  0.2631578947368421,
600  0.3684210526315789,
601  0.05263157894736842,
602  0.1578947368421053,
603  0.2631578947368421,
604  0.05263157894736842,
605  0.1578947368421053,
606  0.05263157894736842,
607  0.05882352941176471,
608  0.1764705882352941,
609  0.2941176470588235,
610  0.4117647058823529,
611  0.5294117647058824,
612  0.6470588235294118,
613  0.7647058823529411,
614  0.8823529411764706,
615  0.05882352941176471,
616  0.1764705882352941,
617  0.2941176470588235,
618  0.4117647058823529,
619  0.5294117647058824,
620  0.6470588235294118,
621  0.7647058823529411,
622  0.05882352941176471,
623  0.1764705882352941,
624  0.2941176470588235,
625  0.4117647058823529,
626  0.5294117647058824,
627  0.6470588235294118,
628  0.05882352941176471,
629  0.1764705882352941,
630  0.2941176470588235,
631  0.4117647058823529,
632  0.5294117647058824,
633  0.05882352941176471,
634  0.1764705882352941,
635  0.2941176470588235,
636  0.4117647058823529,
637  0.05882352941176471,
638  0.1764705882352941,
639  0.2941176470588235,
640  0.05882352941176471,
641  0.1764705882352941,
642  0.05882352941176471,
643  0.06666666666666667,
644  0.2,
645  0.3333333333333333,
646  0.4666666666666667,
647  0.6,
648  0.7333333333333333,
649  0.8666666666666667,
650  0.06666666666666667,
651  0.2,
652  0.3333333333333333,
653  0.4666666666666667,
654  0.6,
655  0.7333333333333333,
656  0.06666666666666667,
657  0.2,
658  0.3333333333333333,
659  0.4666666666666667,
660  0.6,
661  0.06666666666666667,
662  0.2,
663  0.3333333333333333,
664  0.4666666666666667,
665  0.06666666666666667,
666  0.2,
667  0.3333333333333333,
668  0.06666666666666667,
669  0.2,
670  0.06666666666666667,
671  0.07692307692307693,
672  0.2307692307692308,
673  0.3846153846153846,
674  0.5384615384615384,
675  0.6923076923076923,
676  0.8461538461538461,
677  0.07692307692307693,
678  0.2307692307692308,
679  0.3846153846153846,
680  0.5384615384615384,
681  0.6923076923076923,
682  0.07692307692307693,
683  0.2307692307692308,
684  0.3846153846153846,
685  0.5384615384615384,
686  0.07692307692307693,
687  0.2307692307692308,
688  0.3846153846153846,
689  0.07692307692307693,
690  0.2307692307692308,
691  0.07692307692307693,
692  0.09090909090909091,
693  0.2727272727272727,
694  0.4545454545454545,
695  0.6363636363636364,
696  0.8181818181818182,
697  0.09090909090909091,
698  0.2727272727272727,
699  0.4545454545454545,
700  0.6363636363636364,
701  0.09090909090909091,
702  0.2727272727272727,
703  0.4545454545454545,
704  0.09090909090909091,
705  0.2727272727272727,
706  0.09090909090909091,
707  0.1111111111111111,
708  0.3333333333333333,
709  0.5555555555555556,
710  0.7777777777777778,
711  0.1111111111111111,
712  0.3333333333333333,
713  0.5555555555555556,
714  0.1111111111111111,
715  0.3333333333333333,
716  0.1111111111111111,
717  0.1428571428571428,
718  0.4285714285714285,
719  0.7142857142857143,
720  0.1428571428571428,
721  0.4285714285714285,
722  0.1428571428571428,
723  0.2,
724  0.6,
725  0.2,
726  0.3333333333333333};
727 static const double G_TRI_Y286[] = {0.04347826086956522,
728  0.04347826086956522,
729  0.04347826086956522,
730  0.04347826086956522,
731  0.04347826086956522,
732  0.04347826086956522,
733  0.04347826086956522,
734  0.04347826086956522,
735  0.04347826086956522,
736  0.04347826086956522,
737  0.04347826086956522,
738  0.1304347826086956,
739  0.1304347826086956,
740  0.1304347826086956,
741  0.1304347826086956,
742  0.1304347826086956,
743  0.1304347826086956,
744  0.1304347826086956,
745  0.1304347826086956,
746  0.1304347826086956,
747  0.1304347826086956,
748  0.2173913043478261,
749  0.2173913043478261,
750  0.2173913043478261,
751  0.2173913043478261,
752  0.2173913043478261,
753  0.2173913043478261,
754  0.2173913043478261,
755  0.2173913043478261,
756  0.2173913043478261,
757  0.3043478260869565,
758  0.3043478260869565,
759  0.3043478260869565,
760  0.3043478260869565,
761  0.3043478260869565,
762  0.3043478260869565,
763  0.3043478260869565,
764  0.3043478260869565,
765  0.391304347826087,
766  0.391304347826087,
767  0.391304347826087,
768  0.391304347826087,
769  0.391304347826087,
770  0.391304347826087,
771  0.391304347826087,
772  0.4782608695652174,
773  0.4782608695652174,
774  0.4782608695652174,
775  0.4782608695652174,
776  0.4782608695652174,
777  0.4782608695652174,
778  0.5652173913043478,
779  0.5652173913043478,
780  0.5652173913043478,
781  0.5652173913043478,
782  0.5652173913043478,
783  0.6521739130434783,
784  0.6521739130434783,
785  0.6521739130434783,
786  0.6521739130434783,
787  0.7391304347826086,
788  0.7391304347826086,
789  0.7391304347826086,
790  0.8260869565217391,
791  0.8260869565217391,
792  0.9130434782608695,
793  0.04761904761904762,
794  0.04761904761904762,
795  0.04761904761904762,
796  0.04761904761904762,
797  0.04761904761904762,
798  0.04761904761904762,
799  0.04761904761904762,
800  0.04761904761904762,
801  0.04761904761904762,
802  0.04761904761904762,
803  0.1428571428571428,
804  0.1428571428571428,
805  0.1428571428571428,
806  0.1428571428571428,
807  0.1428571428571428,
808  0.1428571428571428,
809  0.1428571428571428,
810  0.1428571428571428,
811  0.1428571428571428,
812  0.2380952380952381,
813  0.2380952380952381,
814  0.2380952380952381,
815  0.2380952380952381,
816  0.2380952380952381,
817  0.2380952380952381,
818  0.2380952380952381,
819  0.2380952380952381,
820  0.3333333333333333,
821  0.3333333333333333,
822  0.3333333333333333,
823  0.3333333333333333,
824  0.3333333333333333,
825  0.3333333333333333,
826  0.3333333333333333,
827  0.4285714285714285,
828  0.4285714285714285,
829  0.4285714285714285,
830  0.4285714285714285,
831  0.4285714285714285,
832  0.4285714285714285,
833  0.5238095238095238,
834  0.5238095238095238,
835  0.5238095238095238,
836  0.5238095238095238,
837  0.5238095238095238,
838  0.6190476190476191,
839  0.6190476190476191,
840  0.6190476190476191,
841  0.6190476190476191,
842  0.7142857142857143,
843  0.7142857142857143,
844  0.7142857142857143,
845  0.8095238095238095,
846  0.8095238095238095,
847  0.9047619047619048,
848  0.05263157894736842,
849  0.05263157894736842,
850  0.05263157894736842,
851  0.05263157894736842,
852  0.05263157894736842,
853  0.05263157894736842,
854  0.05263157894736842,
855  0.05263157894736842,
856  0.05263157894736842,
857  0.1578947368421053,
858  0.1578947368421053,
859  0.1578947368421053,
860  0.1578947368421053,
861  0.1578947368421053,
862  0.1578947368421053,
863  0.1578947368421053,
864  0.1578947368421053,
865  0.2631578947368421,
866  0.2631578947368421,
867  0.2631578947368421,
868  0.2631578947368421,
869  0.2631578947368421,
870  0.2631578947368421,
871  0.2631578947368421,
872  0.3684210526315789,
873  0.3684210526315789,
874  0.3684210526315789,
875  0.3684210526315789,
876  0.3684210526315789,
877  0.3684210526315789,
878  0.4736842105263158,
879  0.4736842105263158,
880  0.4736842105263158,
881  0.4736842105263158,
882  0.4736842105263158,
883  0.5789473684210527,
884  0.5789473684210527,
885  0.5789473684210527,
886  0.5789473684210527,
887  0.6842105263157895,
888  0.6842105263157895,
889  0.6842105263157895,
890  0.7894736842105263,
891  0.7894736842105263,
892  0.8947368421052632,
893  0.05882352941176471,
894  0.05882352941176471,
895  0.05882352941176471,
896  0.05882352941176471,
897  0.05882352941176471,
898  0.05882352941176471,
899  0.05882352941176471,
900  0.05882352941176471,
901  0.1764705882352941,
902  0.1764705882352941,
903  0.1764705882352941,
904  0.1764705882352941,
905  0.1764705882352941,
906  0.1764705882352941,
907  0.1764705882352941,
908  0.2941176470588235,
909  0.2941176470588235,
910  0.2941176470588235,
911  0.2941176470588235,
912  0.2941176470588235,
913  0.2941176470588235,
914  0.4117647058823529,
915  0.4117647058823529,
916  0.4117647058823529,
917  0.4117647058823529,
918  0.4117647058823529,
919  0.5294117647058824,
920  0.5294117647058824,
921  0.5294117647058824,
922  0.5294117647058824,
923  0.6470588235294118,
924  0.6470588235294118,
925  0.6470588235294118,
926  0.7647058823529411,
927  0.7647058823529411,
928  0.8823529411764706,
929  0.06666666666666667,
930  0.06666666666666667,
931  0.06666666666666667,
932  0.06666666666666667,
933  0.06666666666666667,
934  0.06666666666666667,
935  0.06666666666666667,
936  0.2,
937  0.2,
938  0.2,
939  0.2,
940  0.2,
941  0.2,
942  0.3333333333333333,
943  0.3333333333333333,
944  0.3333333333333333,
945  0.3333333333333333,
946  0.3333333333333333,
947  0.4666666666666667,
948  0.4666666666666667,
949  0.4666666666666667,
950  0.4666666666666667,
951  0.6,
952  0.6,
953  0.6,
954  0.7333333333333333,
955  0.7333333333333333,
956  0.8666666666666667,
957  0.07692307692307693,
958  0.07692307692307693,
959  0.07692307692307693,
960  0.07692307692307693,
961  0.07692307692307693,
962  0.07692307692307693,
963  0.2307692307692308,
964  0.2307692307692308,
965  0.2307692307692308,
966  0.2307692307692308,
967  0.2307692307692308,
968  0.3846153846153846,
969  0.3846153846153846,
970  0.3846153846153846,
971  0.3846153846153846,
972  0.5384615384615384,
973  0.5384615384615384,
974  0.5384615384615384,
975  0.6923076923076923,
976  0.6923076923076923,
977  0.8461538461538461,
978  0.09090909090909091,
979  0.09090909090909091,
980  0.09090909090909091,
981  0.09090909090909091,
982  0.09090909090909091,
983  0.2727272727272727,
984  0.2727272727272727,
985  0.2727272727272727,
986  0.2727272727272727,
987  0.4545454545454545,
988  0.4545454545454545,
989  0.4545454545454545,
990  0.6363636363636364,
991  0.6363636363636364,
992  0.8181818181818182,
993  0.1111111111111111,
994  0.1111111111111111,
995  0.1111111111111111,
996  0.1111111111111111,
997  0.3333333333333333,
998  0.3333333333333333,
999  0.3333333333333333,
1000  0.5555555555555556,
1001  0.5555555555555556,
1002  0.7777777777777778,
1003  0.1428571428571428,
1004  0.1428571428571428,
1005  0.1428571428571428,
1006  0.4285714285714285,
1007  0.4285714285714285,
1008  0.7142857142857143,
1009  0.2,
1010  0.2,
1011  0.6,
1012  0.3333333333333333};
1013 static const double G_TRI_W286[] = {
1014  2.912193380035668, 2.912193380035668, 2.912193380035668,
1015  2.912193380035668, 2.912193380035668, 2.912193380035668,
1016  2.912193380035668, 2.912193380035668, 2.912193380035668,
1017  2.912193380035668, 2.912193380035668, 2.912193380035668,
1018  2.912193380035668, 2.912193380035668, 2.912193380035668,
1019  2.912193380035668, 2.912193380035668, 2.912193380035668,
1020  2.912193380035668, 2.912193380035668, 2.912193380035668,
1021  2.912193380035668, 2.912193380035668, 2.912193380035668,
1022  2.912193380035668, 2.912193380035668, 2.912193380035668,
1023  2.912193380035668, 2.912193380035668, 2.912193380035668,
1024  2.912193380035668, 2.912193380035668, 2.912193380035668,
1025  2.912193380035668, 2.912193380035668, 2.912193380035668,
1026  2.912193380035668, 2.912193380035668, 2.912193380035668,
1027  2.912193380035668, 2.912193380035668, 2.912193380035668,
1028  2.912193380035668, 2.912193380035668, 2.912193380035668,
1029  2.912193380035668, 2.912193380035668, 2.912193380035668,
1030  2.912193380035668, 2.912193380035668, 2.912193380035668,
1031  2.912193380035668, 2.912193380035668, 2.912193380035668,
1032  2.912193380035668, 2.912193380035668, 2.912193380035668,
1033  2.912193380035668, 2.912193380035668, 2.912193380035668,
1034  2.912193380035668, 2.912193380035668, 2.912193380035668,
1035  2.912193380035668, 2.912193380035668, 2.912193380035668,
1036  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1037  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1038  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1039  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1040  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1041  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1042  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1043  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1044  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1045  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1046  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1047  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1048  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1049  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1050  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1051  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1052  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1053  -9.914451197589852, -9.914451197589852, -9.914451197589852,
1054  -9.914451197589852, 13.33158527957992, 13.33158527957992,
1055  13.33158527957992, 13.33158527957992, 13.33158527957992,
1056  13.33158527957992, 13.33158527957992, 13.33158527957992,
1057  13.33158527957992, 13.33158527957992, 13.33158527957992,
1058  13.33158527957992, 13.33158527957992, 13.33158527957992,
1059  13.33158527957992, 13.33158527957992, 13.33158527957992,
1060  13.33158527957992, 13.33158527957992, 13.33158527957992,
1061  13.33158527957992, 13.33158527957992, 13.33158527957992,
1062  13.33158527957992, 13.33158527957992, 13.33158527957992,
1063  13.33158527957992, 13.33158527957992, 13.33158527957992,
1064  13.33158527957992, 13.33158527957992, 13.33158527957992,
1065  13.33158527957992, 13.33158527957992, 13.33158527957992,
1066  13.33158527957992, 13.33158527957992, 13.33158527957992,
1067  13.33158527957992, 13.33158527957992, 13.33158527957992,
1068  13.33158527957992, 13.33158527957992, 13.33158527957992,
1069  13.33158527957992, -9.027792408986382, -9.027792408986382,
1070  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1071  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1072  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1073  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1074  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1075  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1076  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1077  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1078  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1079  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1080  -9.027792408986382, -9.027792408986382, -9.027792408986382,
1081  -9.027792408986382, 3.258672079964582, 3.258672079964582,
1082  3.258672079964582, 3.258672079964582, 3.258672079964582,
1083  3.258672079964582, 3.258672079964582, 3.258672079964582,
1084  3.258672079964582, 3.258672079964582, 3.258672079964582,
1085  3.258672079964582, 3.258672079964582, 3.258672079964582,
1086  3.258672079964582, 3.258672079964582, 3.258672079964582,
1087  3.258672079964582, 3.258672079964582, 3.258672079964582,
1088  3.258672079964582, 3.258672079964582, 3.258672079964582,
1089  3.258672079964582, 3.258672079964582, 3.258672079964582,
1090  3.258672079964582, 3.258672079964582, -0.6133639040302452,
1091  -0.6133639040302452, -0.6133639040302452, -0.6133639040302452,
1092  -0.6133639040302452, -0.6133639040302452, -0.6133639040302452,
1093  -0.6133639040302452, -0.6133639040302452, -0.6133639040302452,
1094  -0.6133639040302452, -0.6133639040302452, -0.6133639040302452,
1095  -0.6133639040302452, -0.6133639040302452, -0.6133639040302452,
1096  -0.6133639040302452, -0.6133639040302452, -0.6133639040302452,
1097  -0.6133639040302452, -0.6133639040302452, 0.05511571669513555,
1098  0.05511571669513555, 0.05511571669513555, 0.05511571669513555,
1099  0.05511571669513555, 0.05511571669513555, 0.05511571669513555,
1100  0.05511571669513555, 0.05511571669513555, 0.05511571669513555,
1101  0.05511571669513555, 0.05511571669513555, 0.05511571669513555,
1102  0.05511571669513555, 0.05511571669513555, -0.001979122382447095,
1103  -0.001979122382447095, -0.001979122382447095, -0.001979122382447095,
1104  -0.001979122382447095, -0.001979122382447095, -0.001979122382447095,
1105  -0.001979122382447095, -0.001979122382447095, -0.001979122382447095,
1106  2.02054621415273e-05, 2.02054621415273e-05, 2.02054621415273e-05,
1107  2.02054621415273e-05, 2.02054621415273e-05, 2.02054621415273e-05,
1108  -2.874940020535803e-08, -2.874940020535803e-08, -2.874940020535803e-08,
1109  8.829438425435718e-13};
1110 
1111 // TET
1112 static const double G_TET_X1[] = {0.25};
1113 static const double G_TET_Y1[] = {0.25};
1114 static const double G_TET_Z1[] = {0.25};
1115 static const double G_TET_W1[] = {1.};
1116 static const double G_TET_X4[] = {0.1757281246520584, 0.2445310270213291,
1117  0.5556470949048655, 0.0240937534217468};
1118 static const double G_TET_Y4[] = {0.5656137776620919, 0.0501800797762026,
1119  0.1487681308666864, 0.2354380116950194};
1120 static const double G_TET_Z4[] = {0.2180665126782654, 0.5635595064952189,
1121  0.0350112499848832, 0.1833627308416330};
1122 static const double G_TET_W4[] = {0.25, 0.25, 0.25, 0.25};
1123 static const double G_TET_X5[] = {0.25000000000000000, 0.50000000000000000,
1124  0.16666666666666667, 0.16666666666666667,
1125  0.16666666666666667};
1126 static const double G_TET_Y5[] = {0.25000000000000000, 0.16666666666666667,
1127  0.50000000000000000, 0.16666666666666667,
1128  0.16666666666666667};
1129 static const double G_TET_Z5[] = {0.25000000000000000, 0.16666666666666667,
1130  0.16666666666666667, 0.50000000000000000,
1131  0.16666666666666667};
1132 static const double G_TET_W5[] = {-0.80000000000000000, 0.45000000000000000,
1133  0.45000000000000000, 0.45000000000000000,
1134  0.45000000000000000};
1135 static const double G_TET_X10[] = {0.5684305841968444, 0.1438564719343852,
1136  0.1438564719343852, 0.1438564719343852,
1137  0.0000000000000000, 0.5000000000000000,
1138  0.5000000000000000, 0.5000000000000000,
1139  0.0000000000000000, 0.0000000000000000};
1140 static const double G_TET_Y10[] = {0.1438564719343852, 0.1438564719343852,
1141  0.1438564719343852, 0.5684305841968444,
1142  0.5000000000000000, 0.0000000000000000,
1143  0.5000000000000000, 0.0000000000000000,
1144  0.5000000000000000, 0.0000000000000000};
1145 static const double G_TET_Z10[] = {0.1438564719343852, 0.1438564719343852,
1146  0.5684305841968444, 0.1438564719343852,
1147  0.5000000000000000, 0.5000000000000000,
1148  0.0000000000000000, 0.0000000000000000,
1149  0.0000000000000000, 0.5000000000000000};
1150 static const double G_TET_W10[] = {0.2177650698804054, 0.2177650698804054,
1151  0.2177650698804054, 0.2177650698804054,
1152  0.0214899534130631, 0.0214899534130631,
1153  0.0214899534130631, 0.0214899534130631,
1154  0.0214899534130631, 0.0214899534130631};
1155 
1156 static const double G_TET_X45[] = {
1157  0.2500000000000000, 0.6175871903000830, 0.1274709365666390,
1158  0.1274709365666390, 0.1274709365666390, 0.9037635088221031,
1159  0.0320788303926323, 0.0320788303926323, 0.0320788303926323,
1160  0.4502229043567190, 0.0497770956432810, 0.0497770956432810,
1161  0.0497770956432810, 0.4502229043567190, 0.4502229043567190,
1162  0.3162695526014501, 0.1837304473985499, 0.1837304473985499,
1163  0.1837304473985499, 0.3162695526014501, 0.3162695526014501,
1164  0.0229177878448171, 0.2319010893971509, 0.2319010893971509,
1165  0.5132800333608811, 0.2319010893971509, 0.2319010893971509,
1166  0.2319010893971509, 0.0229177878448171, 0.5132800333608811,
1167  0.2319010893971509, 0.0229177878448171, 0.5132800333608811,
1168  0.7303134278075384, 0.0379700484718286, 0.0379700484718286,
1169  0.1937464752488044, 0.0379700484718286, 0.0379700484718286,
1170  0.0379700484718286, 0.7303134278075384, 0.1937464752488044,
1171  0.0379700484718286, 0.7303134278075384, 0.1937464752488044};
1172 static const double G_TET_Y45[] = {
1173  0.2500000000000000, 0.1274709365666390, 0.1274709365666390,
1174  0.1274709365666390, 0.6175871903000830, 0.0320788303926323,
1175  0.0320788303926323, 0.0320788303926323, 0.9037635088221031,
1176  0.0497770956432810, 0.4502229043567190, 0.0497770956432810,
1177  0.4502229043567190, 0.0497770956432810, 0.4502229043567190,
1178  0.1837304473985499, 0.3162695526014501, 0.1837304473985499,
1179  0.3162695526014501, 0.1837304473985499, 0.3162695526014501,
1180  0.2319010893971509, 0.0229177878448171, 0.2319010893971509,
1181  0.2319010893971509, 0.5132800333608811, 0.2319010893971509,
1182  0.0229177878448171, 0.5132800333608811, 0.2319010893971509,
1183  0.5132800333608811, 0.2319010893971509, 0.0229177878448171,
1184  0.0379700484718286, 0.7303134278075384, 0.0379700484718286,
1185  0.0379700484718286, 0.1937464752488044, 0.0379700484718286,
1186  0.7303134278075384, 0.1937464752488044, 0.0379700484718286,
1187  0.1937464752488044, 0.0379700484718286, 0.7303134278075384};
1188 static const double G_TET_Z45[] = {
1189  0.2500000000000000, 0.1274709365666390, 0.1274709365666390,
1190  0.6175871903000830, 0.1274709365666390, 0.0320788303926323,
1191  0.0320788303926323, 0.9037635088221031, 0.0320788303926323,
1192  0.0497770956432810, 0.0497770956432810, 0.4502229043567190,
1193  0.4502229043567190, 0.4502229043567190, 0.0497770956432810,
1194  0.1837304473985499, 0.1837304473985499, 0.3162695526014501,
1195  0.3162695526014501, 0.3162695526014501, 0.1837304473985499,
1196  0.2319010893971509, 0.2319010893971509, 0.0229177878448171,
1197  0.2319010893971509, 0.2319010893971509, 0.5132800333608811,
1198  0.5132800333608811, 0.2319010893971509, 0.0229177878448171,
1199  0.0229177878448171, 0.5132800333608811, 0.2319010893971509,
1200  0.0379700484718286, 0.0379700484718286, 0.7303134278075384,
1201  0.0379700484718286, 0.0379700484718286, 0.1937464752488044,
1202  0.1937464752488044, 0.0379700484718286, 0.7303134278075384,
1203  0.7303134278075384, 0.1937464752488044, 0.0379700484718286};
1204 static const double G_TET_W45[] = {
1205  -0.2359620398477557, 0.0244878963560562, 0.0244878963560562,
1206  0.0244878963560562, 0.0244878963560562, 0.0039485206398261,
1207  0.0039485206398261, 0.0039485206398261, 0.0039485206398261,
1208  0.0263055529507371, 0.0263055529507371, 0.0263055529507371,
1209  0.0263055529507371, 0.0263055529507371, 0.0263055529507371,
1210  0.0829803830550589, 0.0829803830550589, 0.0829803830550589,
1211  0.0829803830550589, 0.0829803830550589, 0.0829803830550589,
1212  0.0254426245481023, 0.0254426245481023, 0.0254426245481023,
1213  0.0254426245481023, 0.0254426245481023, 0.0254426245481023,
1214  0.0254426245481023, 0.0254426245481023, 0.0254426245481023,
1215  0.0254426245481023, 0.0254426245481023, 0.0254426245481023,
1216  0.0134324384376852, 0.0134324384376852, 0.0134324384376852,
1217  0.0134324384376852, 0.0134324384376852, 0.0134324384376852,
1218  0.0134324384376852, 0.0134324384376852, 0.0134324384376852,
1219  0.0134324384376852, 0.0134324384376852, 0.0134324384376852};
1220 static const double NC_TET_X84[] = {
1221  0.1000000000000000, 0.1000000000000000, 0.1000000000000000,
1222  0.7000000000000000, 0.1000000000000000, 0.1000000000000000,
1223  0.1000000000000000, 0.1000000000000000, 0.1000000000000000,
1224  0.1000000000000000, 0.6000000000000000, 0.6000000000000000,
1225  0.6000000000000000, 0.2000000000000000, 0.2000000000000000,
1226  0.2000000000000000, 0.1000000000000000, 0.1000000000000000,
1227  0.1000000000000000, 0.1000000000000000, 0.1000000000000000,
1228  0.1000000000000000, 0.5000000000000000, 0.5000000000000000,
1229  0.5000000000000000, 0.3000000000000000, 0.3000000000000000,
1230  0.3000000000000000, 0.2000000000000000, 0.2000000000000000,
1231  0.2000000000000000, 0.2000000000000000, 0.2000000000000000,
1232  0.2000000000000000, 0.1000000000000000, 0.1000000000000000,
1233  0.1000000000000000, 0.5000000000000000, 0.5000000000000000,
1234  0.5000000000000000, 0.4000000000000000, 0.4000000000000000,
1235  0.4000000000000000, 0.1000000000000000, 0.1000000000000000,
1236  0.1000000000000000, 0.4000000000000000, 0.4000000000000000,
1237  0.4000000000000000, 0.4000000000000000, 0.4000000000000000,
1238  0.4000000000000000, 0.3000000000000000, 0.3000000000000000,
1239  0.3000000000000000, 0.3000000000000000, 0.3000000000000000,
1240  0.3000000000000000, 0.2000000000000000, 0.2000000000000000,
1241  0.2000000000000000, 0.2000000000000000, 0.2000000000000000,
1242  0.2000000000000000, 0.1000000000000000, 0.1000000000000000,
1243  0.1000000000000000, 0.1000000000000000, 0.1000000000000000,
1244  0.1000000000000000, 0.2000000000000000, 0.2000000000000000,
1245  0.2000000000000000, 0.4000000000000000, 0.3000000000000000,
1246  0.3000000000000000, 0.3000000000000000, 0.1000000000000000,
1247  0.3000000000000000, 0.3000000000000000, 0.3000000000000000,
1248  0.2000000000000000, 0.2000000000000000, 0.2000000000000000};
1249 static const double NC_TET_Y84[] = {
1250  0.1000000000000000, 0.1000000000000000, 0.7000000000000000,
1251  0.1000000000000000, 0.1000000000000000, 0.1000000000000000,
1252  0.6000000000000000, 0.6000000000000000, 0.2000000000000000,
1253  0.2000000000000000, 0.1000000000000000, 0.1000000000000000,
1254  0.2000000000000000, 0.1000000000000000, 0.1000000000000000,
1255  0.6000000000000000, 0.1000000000000000, 0.1000000000000000,
1256  0.5000000000000000, 0.5000000000000000, 0.3000000000000000,
1257  0.3000000000000000, 0.1000000000000000, 0.1000000000000000,
1258  0.3000000000000000, 0.1000000000000000, 0.1000000000000000,
1259  0.5000000000000000, 0.2000000000000000, 0.2000000000000000,
1260  0.1000000000000000, 0.1000000000000000, 0.5000000000000000,
1261  0.5000000000000000, 0.2000000000000000, 0.2000000000000000,
1262  0.5000000000000000, 0.2000000000000000, 0.2000000000000000,
1263  0.1000000000000000, 0.4000000000000000, 0.1000000000000000,
1264  0.1000000000000000, 0.4000000000000000, 0.4000000000000000,
1265  0.1000000000000000, 0.3000000000000000, 0.3000000000000000,
1266  0.2000000000000000, 0.2000000000000000, 0.1000000000000000,
1267  0.1000000000000000, 0.4000000000000000, 0.4000000000000000,
1268  0.2000000000000000, 0.2000000000000000, 0.1000000000000000,
1269  0.1000000000000000, 0.4000000000000000, 0.4000000000000000,
1270  0.3000000000000000, 0.3000000000000000, 0.1000000000000000,
1271  0.1000000000000000, 0.4000000000000000, 0.4000000000000000,
1272  0.3000000000000000, 0.3000000000000000, 0.2000000000000000,
1273  0.2000000000000000, 0.2000000000000000, 0.2000000000000000,
1274  0.4000000000000000, 0.2000000000000000, 0.3000000000000000,
1275  0.3000000000000000, 0.1000000000000000, 0.3000000000000000,
1276  0.3000000000000000, 0.2000000000000000, 0.2000000000000000,
1277  0.3000000000000000, 0.3000000000000000, 0.2000000000000000};
1278 static const double NC_TET_Z84[] = {
1279  0.1000000000000000, 0.7000000000000000, 0.1000000000000000,
1280  0.1000000000000000, 0.6000000000000000, 0.2000000000000000,
1281  0.1000000000000000, 0.2000000000000000, 0.1000000000000000,
1282  0.6000000000000000, 0.1000000000000000, 0.2000000000000000,
1283  0.1000000000000000, 0.1000000000000000, 0.6000000000000000,
1284  0.1000000000000000, 0.5000000000000000, 0.3000000000000000,
1285  0.1000000000000000, 0.3000000000000000, 0.1000000000000000,
1286  0.5000000000000000, 0.1000000000000000, 0.3000000000000000,
1287  0.1000000000000000, 0.1000000000000000, 0.5000000000000000,
1288  0.1000000000000000, 0.1000000000000000, 0.5000000000000000,
1289  0.2000000000000000, 0.5000000000000000, 0.2000000000000000,
1290  0.1000000000000000, 0.2000000000000000, 0.5000000000000000,
1291  0.2000000000000000, 0.2000000000000000, 0.1000000000000000,
1292  0.2000000000000000, 0.1000000000000000, 0.4000000000000000,
1293  0.1000000000000000, 0.4000000000000000, 0.1000000000000000,
1294  0.4000000000000000, 0.2000000000000000, 0.1000000000000000,
1295  0.3000000000000000, 0.1000000000000000, 0.3000000000000000,
1296  0.2000000000000000, 0.2000000000000000, 0.1000000000000000,
1297  0.4000000000000000, 0.1000000000000000, 0.4000000000000000,
1298  0.2000000000000000, 0.3000000000000000, 0.1000000000000000,
1299  0.4000000000000000, 0.1000000000000000, 0.4000000000000000,
1300  0.3000000000000000, 0.3000000000000000, 0.2000000000000000,
1301  0.4000000000000000, 0.2000000000000000, 0.4000000000000000,
1302  0.3000000000000000, 0.2000000000000000, 0.4000000000000000,
1303  0.2000000000000000, 0.2000000000000000, 0.3000000000000000,
1304  0.1000000000000000, 0.3000000000000000, 0.3000000000000000,
1305  0.2000000000000000, 0.3000000000000000, 0.2000000000000000,
1306  0.3000000000000000, 0.2000000000000000, 0.3000000000000000};
1307 static const double NC_TET_W84[] = {
1308  0.2843915343915344, 0.2843915343915344, 0.2843915343915344,
1309  0.2843915343915344, -0.3882275132275133, -0.3882275132275133,
1310  -0.3882275132275133, -0.3882275132275133, -0.3882275132275133,
1311  -0.3882275132275133, -0.3882275132275133, -0.3882275132275133,
1312  -0.3882275132275133, -0.3882275132275133, -0.3882275132275133,
1313  -0.3882275132275133, 0.8776455026455027, 0.8776455026455027,
1314  0.8776455026455027, 0.8776455026455027, 0.8776455026455027,
1315  0.8776455026455027, 0.8776455026455027, 0.8776455026455027,
1316  0.8776455026455027, 0.8776455026455027, 0.8776455026455027,
1317  0.8776455026455027, 0.1236772486772487, 0.1236772486772487,
1318  0.1236772486772487, 0.1236772486772487, 0.1236772486772487,
1319  0.1236772486772487, 0.1236772486772487, 0.1236772486772487,
1320  0.1236772486772487, 0.1236772486772487, 0.1236772486772487,
1321  0.1236772486772487, -0.8584656084656085, -0.8584656084656085,
1322  -0.8584656084656085, -0.8584656084656085, -0.8584656084656085,
1323  -0.8584656084656085, -0.2632275132275133, -0.2632275132275133,
1324  -0.2632275132275133, -0.2632275132275133, -0.2632275132275133,
1325  -0.2632275132275133, -0.2632275132275133, -0.2632275132275133,
1326  -0.2632275132275133, -0.2632275132275133, -0.2632275132275133,
1327  -0.2632275132275133, -0.2632275132275133, -0.2632275132275133,
1328  -0.2632275132275133, -0.2632275132275133, -0.2632275132275133,
1329  -0.2632275132275133, -0.2632275132275133, -0.2632275132275133,
1330  -0.2632275132275133, -0.2632275132275133, -0.2632275132275133,
1331  -0.2632275132275133, 0.0145502645502645, 0.0145502645502645,
1332  0.0145502645502645, 0.0145502645502645, 1.0165343915343916,
1333  1.0165343915343916, 1.0165343915343916, 1.0165343915343916,
1334  -0.0251322751322751, -0.0251322751322751, -0.0251322751322751,
1335  -0.0251322751322751, -0.0251322751322751, -0.0251322751322751};
1336 
1337 #ifdef __cplusplus
1338 }
1339 #endif
1340 
1341 #endif
G_TET_X5
static const double G_TET_X5[]
Definition: fem_tools.h:1123
NC_TET_X84
static const double NC_TET_X84[]
Definition: fem_tools.h:1220
G_TRI_X37
static const double G_TRI_X37[]
Definition: fem_tools.h:408
ShapeMBTETQ_detJac_at_Gauss_Points
PetscErrorCode ShapeMBTETQ_detJac_at_Gauss_Points(double *detJac_at_Gauss_Points, const double *diffN, const double *coords, int G_DIM)
Definition: fem_tools.c:979
G_TRI_X1
static const double G_TRI_X1[]
Definition: fem_tools.h:312
ShapeJacMBTRI
void ShapeJacMBTRI(double *diffN, const double *coords, double *Jac)
calculate jacobioan
G_TRI_W3
static const double G_TRI_W3[]
Definition: fem_tools.h:317
G_TRI_W19
static const double G_TRI_W19[]
Definition: fem_tools.h:365
G_TET_Z5
static const double G_TET_Z5[]
Definition: fem_tools.h:1129
G_TRI_W13
static const double G_TRI_W13[]
Definition: fem_tools.h:347
ShapeMBEDGE
PetscErrorCode ShapeMBEDGE(double *N, const double *G_X, int DIM)
Definition: fem_tools.c:761
print_mat
void print_mat(double *M, int m, int n)
print matric M
G_TET_Z4
static const double G_TET_Z4[]
Definition: fem_tools.h:1120
ShapeMBTET_inverse
PetscErrorCode ShapeMBTET_inverse(double *N, double *diffN, const double *elem_coords, const double *glob_coords, double *loc_coords)
calculate local coordinates for given global coordinates
Definition: fem_tools.c:335
ShapeMBTETQ_GAUSS
PetscErrorCode ShapeMBTETQ_GAUSS(double *N, const double *X, const double *Y, const double *Z, const int G_DIM)
Definition: fem_tools.c:908
print_mat_sym_upper
void print_mat_sym_upper(double *M, int m, int n)
print upper part of the symmetric matrix
Grundmann_Moeller_integration_points_2D_TRI
PetscErrorCode Grundmann_Moeller_integration_points_2D_TRI(int rule, double *G_TRI_X, double *G_TRI_Y, double *G_TRI_W)
Definition: fem_tools.c:98
Base_scale
PetscErrorCode Base_scale(__CLPK_doublecomplex *xnormal, __CLPK_doublecomplex *xs1, __CLPK_doublecomplex *xs2)
Definition: fem_tools.c:741
G_TET_Y5
static const double G_TET_Y5[]
Definition: fem_tools.h:1126
G_TET_Y4
static const double G_TET_Y4[]
Definition: fem_tools.h:1118
GradientOfDeformation
PetscErrorCode GradientOfDeformation(double *diffN, double *dofs, double *F)
calculate gradient of deformation
Definition: fem_tools.c:425
NC_TET_Y84
static const double NC_TET_Y84[]
Definition: fem_tools.h:1249
G_TET_X10
static const double G_TET_X10[]
Definition: fem_tools.h:1135
G_TET_W10
static const double G_TET_W10[]
Definition: fem_tools.h:1150
G_TET_W1
static const double G_TET_W1[]
Definition: fem_tools.h:1115
G_TRI_W4
static const double G_TRI_W4[]
Definition: fem_tools.h:325
order
constexpr int order
Definition: dg_projection.cpp:18
G_TRI_X7
static const double G_TRI_X7[]
Definition: fem_tools.h:328
G_TET_X45
static const double G_TET_X45[]
Definition: fem_tools.h:1156
G_TRI_Y286
static const double G_TRI_Y286[]
Definition: fem_tools.h:727
ShapeDiffMBTRI
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:194
PlasticOps::M
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
ShapeDiffMBTETQ_GAUSS
PetscErrorCode ShapeDiffMBTETQ_GAUSS(double *diffN, const double *X, const double *Y, const double *Z, const int G_DIM)
Definition: fem_tools.c:927
G_TET_Z1
static const double G_TET_Z1[]
Definition: fem_tools.h:1114
ShapeDiffMBTETinvJ_complex
void ShapeDiffMBTETinvJ_complex(double *diffN, __CLPK_doublecomplex *invJac, __CLPK_doublecomplex *diffNinvJac, enum CBLAS_TRANSPOSE Trans)
Definition: fem_tools.c:437
G_TRI_Y37
static const double G_TRI_Y37[]
Definition: fem_tools.h:419
G_TRI_X3
static const double G_TRI_X3[]
Definition: fem_tools.h:315
Spin
PetscErrorCode Spin(double *spinOmega, double *vecOmega)
calculate spin matrix from vector
Definition: fem_tools.c:546
G_TRI_W286
static const double G_TRI_W286[]
Definition: fem_tools.h:1013
NC_TET_W84
static const double NC_TET_W84[]
Definition: fem_tools.h:1307
ShapeFaceDiffNormalMBTRI
PetscErrorCode ShapeFaceDiffNormalMBTRI(double *diffN, const double *coords, double *diff_normal)
calculate derivative of normal in respect to nodal positions
Definition: fem_tools.c:237
ShapeDiffMBTETQ
PetscErrorCode ShapeDiffMBTETQ(double *diffN, const double x, const double y, const double z)
Definition: fem_tools.c:873
G_TRI_W7
static const double G_TRI_W7[]
Definition: fem_tools.h:334
G_TET_Y1
static const double G_TET_Y1[]
Definition: fem_tools.h:1113
ShapeMBTRI_inverse
PetscErrorCode ShapeMBTRI_inverse(double *N, double *diffN, const double *elem_coords, const double *glob_coords, double *loc_coords)
calculate local coordinates of triangle element for given global coordinates in 2D (Assume e....
Definition: fem_tools.c:380
ShapeMBTRIQ
PetscErrorCode ShapeMBTRIQ(double *N, const double *X, const double *Y, const int G_DIM)
Definition: fem_tools.c:780
ShapeFaceNormalMBTRI
PetscErrorCode ShapeFaceNormalMBTRI(double *diffN, const double *coords, double *normal)
Definition: fem_tools.c:229
ShapeDiffMBTETinvJ
PetscErrorCode ShapeDiffMBTETinvJ(double *diffN, double *invJac, double *diffNinvJac)
calculate shape functions derivatives in space
Definition: fem_tools.c:415
NC_TET_Z84
static const double NC_TET_Z84[]
Definition: fem_tools.h:1278
G_TET_X1
static const double G_TET_X1[]
Definition: fem_tools.h:1112
lapack_wrap.h
G_TRI_Y1
static const double G_TRI_Y1[]
Definition: fem_tools.h:313
G_TRI_W28
static const double G_TRI_W28[]
Definition: fem_tools.h:396
make_complex_matrix
PetscErrorCode make_complex_matrix(double *reA, double *imA, __CLPK_doublecomplex *xA)
Compose complex matrix (3x3) from two real matrices.
Definition: fem_tools.c:557
G_TRI_X28
static const double G_TRI_X28[]
Definition: fem_tools.h:377
print_mat_complex
void print_mat_complex(__CLPK_doublecomplex *M, int m, int n)
priint complex matrix
ShapeVolumeMBTET
double ShapeVolumeMBTET(double *diffN, const double *coords)
calculate TET volume
Definition: fem_tools.c:298
Grundmann_Moeller_integration_points_1D_EDGE
PetscErrorCode Grundmann_Moeller_integration_points_1D_EDGE(int rule, double *G_TRI_X, double *G_TRI_W)
Compute weights and integration points for edge using Grundmann_Moeller rule.
Definition: fem_tools.c:55
G_TRI_Y19
static const double G_TRI_Y19[]
Definition: fem_tools.h:359
G_TRI_Y13
static const double G_TRI_Y13[]
Definition: fem_tools.h:342
G_TRI_X19
static const double G_TRI_X19[]
Definition: fem_tools.h:353
InvertComplexGradient
PetscErrorCode InvertComplexGradient(__CLPK_doublecomplex *xF)
Definition: fem_tools.c:499
G_TRI_X4
static const double G_TRI_X4[]
Definition: fem_tools.h:319
G_TET_W45
static const double G_TET_W45[]
Definition: fem_tools.h:1204
convert.n
n
Definition: convert.py:82
G_TRI_Y4
static const double G_TRI_Y4[]
Definition: fem_tools.h:322
N
const int N
Definition: speed_test.cpp:3
G_TET_X4
static const double G_TET_X4[]
Definition: fem_tools.h:1116
G_TET_W4
static const double G_TET_W4[]
Definition: fem_tools.h:1122
G_TET_Z45
static const double G_TET_Z45[]
Definition: fem_tools.h:1188
ShapeDetJacVolume
double ShapeDetJacVolume(double *jac)
determined of jacobian
Definition: fem_tools.c:22
ShapeJacMBTET
PetscErrorCode ShapeJacMBTET(double *diffN, const double *coords, double *jac)
calculate jacobian
Definition: fem_tools.c:288
G_TRI_X13
static const double G_TRI_X13[]
Definition: fem_tools.h:337
G_TRI_W1
static const double G_TRI_W1[]
Definition: fem_tools.h:314
G_TET_Y10
static const double G_TET_Y10[]
Definition: fem_tools.h:1140
G_TRI_Y7
static const double G_TRI_Y7[]
Definition: fem_tools.h:331
ShapeMBTETQ_inverse
PetscErrorCode ShapeMBTETQ_inverse(double *N, double *diffN, const double *elem_coords, const double *glob_coords, double *loc_coords, const double eps)
Definition: fem_tools.c:1007
G_TET_W5
static const double G_TET_W5[]
Definition: fem_tools.h:1132
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
ShapeDiffMBTET
PetscErrorCode ShapeDiffMBTET(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:319
G_TRI_W37
static const double G_TRI_W37[]
Definition: fem_tools.h:430
ShapeFaceBaseMBTRI
PetscErrorCode ShapeFaceBaseMBTRI(double *diffN, const double *coords, double *normal, double *s1, double *s2)
Definition: fem_tools.c:204
ShapeJacMBTETQ
PetscErrorCode ShapeJacMBTETQ(const double *diffN, const double *coords, double *Jac)
Definition: fem_tools.c:967
Grundmann_Moeller_integration_points_3D_TET
PetscErrorCode Grundmann_Moeller_integration_points_3D_TET(int rule, double *G_TET_X, double *G_TET_Y, double *G_TET_Z, double *G_TET_W)
Definition: fem_tools.c:143
DeterminantComplexGradient
PetscErrorCode DeterminantComplexGradient(__CLPK_doublecomplex *xF, __CLPK_doublecomplex *det_xF)
Definition: fem_tools.c:526
Normal_hierarchical
PetscErrorCode Normal_hierarchical(int order_approx, int *order_edge_approx, int order, int *order_edge, double *diffN, double *diffN_face, double *diffN_edge[], double *dofs, double *dofs_edge[], double *dofs_face, double *idofs, double *idofs_edge[], double *idofs_face, __CLPK_doublecomplex *xnormal, __CLPK_doublecomplex *s1, __CLPK_doublecomplex *s2, int gg)
Complex normal.
Definition: fem_tools.c:569
ShapeFaceNormalMBTRI_complex
PetscErrorCode ShapeFaceNormalMBTRI_complex(double *diffN, __CLPK_doublecomplex *xcoords, __CLPK_doublecomplex *xnormal)
Definition: fem_tools.c:452
ShapeVolumeMBTETQ
double ShapeVolumeMBTETQ(const double *diffN, const double *coords, int G_DIM, double *G_TET_W)
Definition: fem_tools.c:993
ShapeDiffMBTRIinvJ
void ShapeDiffMBTRIinvJ(double *diffN, double *invJac, double *diffNinvJac)
calculate derivatives of shape functions in space
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
G_TRI_X286
static const double G_TRI_X286[]
Definition: fem_tools.h:441
ShapeDiffMBTRIQ
PetscErrorCode ShapeDiffMBTRIQ(double *diffN, const double *X, const double *Y, const int G_DIM)
Definition: fem_tools.c:795
G_TRI_Y28
static const double G_TRI_Y28[]
Definition: fem_tools.h:388
ShapeInvJacVolume
PetscErrorCode ShapeInvJacVolume(double *jac)
Definition: fem_tools.c:39
invJac
MatrixDouble invJac
Definition: HookeElement.hpp:683
ShapeDiffMBEDGE
PetscErrorCode ShapeDiffMBEDGE(double *diffN)
Definition: fem_tools.c:771
MakeComplexTensor
PetscErrorCode MakeComplexTensor(double *reA, double *imA, __CLPK_doublecomplex *xA)
Definition: fem_tools.c:487
G_TRI_Y3
static const double G_TRI_Y3[]
Definition: fem_tools.h:316
__CLPK_doublecomplex
Definition: lapack_wrap.h:34
ShapeMBTETQ
PetscErrorCode ShapeMBTETQ(double *N, const double x, const double y, const double z)
Definition: fem_tools.c:858
G_TET_Y45
static const double G_TET_Y45[]
Definition: fem_tools.h:1172
G_TET_Z10
static const double G_TET_Z10[]
Definition: fem_tools.h:1145
F
@ F
Definition: free_surface.cpp:394
ShapeMBTET
PetscErrorCode ShapeMBTET(double *N, const double *G_X, const double *G_Y, const double *G_Z, int DIM)
calculate shape functions
Definition: fem_tools.c:306
ShapeMBTRI
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:182