ParMooN
 All Classes Functions Variables Friends Pages
BF_C_T_P5_2D.h
1 // ***********************************************************************
2 // P5 element, conforming, 2D
3 // ***********************************************************************
4 
5 // base function values
6 static void C_T_P5_2D_Funct(double xi, double eta, double *values)
7 {
8  double t3, t5, t7, t9, t11, t13, t15, t18, t20, t22, t24, t26, t28, t29;
9  double t30, t31, t32, t33, t34, t35, t36, t37, t38, t39, t40, t42, t44;
10  double t46, t50, t52, t53, t54, t55, t56, t58, t63, t64, t65, t66, t67;
11  double t70, t73, t76, t77, t80, t82, t84, t90, t93, t98, t99, t100, t101;
12  double t102, t109, t110, t112, t116, t117, t119, t123, t129, t132, t133;
13  double t134, t143, t148, t153, t162;
14 
15  t3 = xi*xi;
16  t5 = xi*eta;
17  t7 = eta*eta;
18  t9 = t3*xi;
19  t11 = t3*eta;
20  t13 = xi*t7;
21  t15 = t7*eta;
22  t18 = t3*t3;
23  t20 = t9*eta;
24  t22 = t3*t7;
25  t24 = xi*t15;
26  t26 = t7*t7;
27  t28 = t18*xi;
28  t29 = 625.0/24.0*t28;
29  t30 = t18*eta;
30  t31 = 3125.0/24.0*t30;
31  t32 = t9*t7;
32  t33 = 3125.0/12.0*t32;
33  t34 = t3*t15;
34  t35 = 3125.0/12.0*t34;
35  t36 = xi*t26;
36  t37 = 3125.0/24.0*t36;
37  t38 = t26*eta;
38  t39 = 625.0/24.0*t38;
39  t40 = 625.0/8.0*t18+625.0/2.0*t20+1875.0/4.0*t22+625.0/2.0*t24+625.0/8.0*t26-t29-t31-t33-t35-t37-t39;
40  t42 = 25.0*xi;
41  t44 = 1925.0/12.0*t5;
42  t46 = 8875.0/12.0*t11;
43  t50 = 4375.0/4.0*t22;
44  t52 = 3125.0/24.0*t28;
45  t53 = 3125.0/6.0*t30;
46  t54 = 3125.0/4.0*t32;
47  t55 = 3125.0/6.0*t34;
48  t56 = t42-1925.0/12.0*t3-t44+8875.0/24.0*t9+t46+8875.0/24.0*t13-4375.0/12.0*t18-4375.0/4.0*t20-t50-4375.0/12.0*t24+t52+t53+t54+t55+t37;
49  t58 = 1175.0/12.0*t5;
50  t63 = 3125.0/4.0*t22;
51  t64 = 625.0/12.0*t24;
52  t65 = 3125.0/12.0*t28;
53  t66 = 3125.0/4.0*t30;
54  t67 = -t42+2675.0/12.0*t3+t58-7375.0/12.0*t9-t46-125.0*t13+8125.0/12.0*t18+5625.0/4.0*t20+t63+t64-t65-t66-t54-t35;
55  t70 = 75.0/2.0*t5;
56  t73 = 125.0/6.0*t13;
57  t76 = 625.0/4.0*t22;
58  t77 = 50.0/3.0*xi-325.0/2.0*t3-t70+6125.0/12.0*t9+3875.0/12.0*t11+t73-625.0*t18-3125.0/4.0*t20-t76+t65+t53+t33;
59  t80 = 25.0/4.0*t5;
60  t82 = 1375.0/24.0*t11;
61  t84 = 625.0/4.0*t20;
62  t90 = 25.0*eta;
63  t93 = 8875.0/12.0*t13;
64  t98 = 3125.0/6.0*t32;
65  t99 = 3125.0/4.0*t34;
66  t100 = 3125.0/6.0*t36;
67  t101 = 3125.0/24.0*t38;
68  t102 = t90-t44-1925.0/12.0*t7+8875.0/24.0*t11+t93+8875.0/24.0*t15-4375.0/12.0*t20-t50-4375.0/4.0*t24-4375.0/12.0*t26+t31+t98+t99+t100+t101;
69  t109 = 3125.0/2.0*t32;
70  t110 = 3125.0/2.0*t34;
71  t112 = 125.0*t5;
72  t116 = 6875.0/4.0*t22;
73  t117 = 625.0/4.0*t24;
74  t119 = 125.0/3.0*t5;
75  t123 = 625.0/2.0*t22;
76  t129 = 625.0/12.0*t20;
77  t132 = 3125.0/4.0*t36;
78  t133 = 3125.0/12.0*t38;
79  t134 = -t90+t58+2675.0/12.0*t7-125.0*t11-t93-7375.0/12.0*t15+t129+t63+5625.0/4.0*t24+8125.0/12.0*t26-t33-t99-t132-t133;
80  t143 = 25.0/6.0*t5;
81  t148 = 125.0/6.0*t11;
82  t153 = 50.0/3.0*eta-t70-325.0/2.0*t7+t148+3875.0/12.0*t13+6125.0/12.0*t15-t76-3125.0/4.0*t24-625.0*t26+t35+t100+t133;
83  t162 = 1375.0/24.0*t13;
84 
85  values[0] = 1.0-137.0/12.0*xi-137.0/12.0*eta+375.0/8.0*t3+375.0/4.0*t5+375.0/8.0*t7-2125.0/24.0*t9-2125.0/8.0*t11-2125.0/8.0*t13-2125.0/24.0*t15+t40;
86  values[1] = t56;
87  values[2] = t67;
88  values[3] = t77;
89  values[4] = -25.0/4.0*xi+1525.0/24.0*t3+t80-5125.0/24.0*t9-t82+6875.0/24.0*t18+t84-t52-t31;
90  values[5] = xi-125.0/12.0*t3+875.0/24.0*t9-625.0/12.0*t18+t29;
91  values[6] = t102;
92  values[7] = 250.0*t5-5875.0/6.0*t11-5875.0/6.0*t13+1250.0*t20+2500.0*t22+1250.0*t24-t53-t109-t110-t100;
93  values[8] = -t112+3625.0/4.0*t11+1125.0/4.0*t13-3125.0/2.0*t20-t116-t117+t66+t109+t99;
94  values[9] = t119-2125.0/6.0*t11-125.0/3.0*t13+2500.0/3.0*t20+t123-t53-t98;
95  values[10] = -t80+t82-t84+t31;
96  values[11] = t134;
97  values[12] = -t112+1125.0/4.0*t11+3625.0/4.0*t13-t84-t116-3125.0/2.0*t24+t54+t110+t132;
98  values[13] = 125.0/4.0*t5-375.0/2.0*t11-375.0/2.0*t13+t84+t50+t117-t54-t99;
99  values[14] = -t143+125.0/4.0*t11+t73-t129-t76+t33;
100  values[15] = t153;
101  values[16] = t119-125.0/3.0*t11-2125.0/6.0*t13+t123+2500.0/3.0*t24-t55-t100;
102  values[17] = -t143+t148+125.0/4.0*t13-t76-t64+t35;
103  values[18] = -25.0/4.0*eta+t80+1525.0/24.0*t7-t162-5125.0/24.0*t15+t117+6875.0/24.0*t26-t37-t101;
104  values[19] = -t80+t162-t117+t37;
105  values[20] = eta-125.0/12.0*t7+875.0/24.0*t15-625.0/12.0*t26+t39;
106 }
107 
108 // values of the derivatives in xi direction
109 static void C_T_P5_2D_DeriveXi(double xi, double eta, double *values)
110 {
111  double t3, t5, t7, t9, t11, t13, t15, t17, t18, t19, t20, t21, t22, t23;
112  double t24, t25, t26, t27, t29, t31, t35, t37, t38, t39, t40, t41, t43;
113  double t48, t49, t50, t51, t52, t54, t57, t60, t61, t63, t65, t67, t74;
114  double t77, t78, t79, t87, t88, t90, t94, t95, t97, t101, t105, t107;
115  double t117, t120, t131;
116 
117  t3 = xi*xi;
118  t5 = xi*eta;
119  t7 = eta*eta;
120  t9 = t3*xi;
121  t11 = t3*eta;
122  t13 = xi*t7;
123  t15 = t7*eta;
124  t17 = t3*t3;
125  t18 = 3125.0/24.0*t17;
126  t19 = t9*eta;
127  t20 = 3125.0/6.0*t19;
128  t21 = t3*t7;
129  t22 = 3125.0/4.0*t21;
130  t23 = xi*t15;
131  t24 = 3125.0/6.0*t23;
132  t25 = t7*t7;
133  t26 = 3125.0/24.0*t25;
134  t27 = -137.0/12.0+375.0/4.0*xi+375.0/4.0*eta-2125.0/8.0*t3-2125.0/4.0*t5-2125.0/8.0*t7+625.0/2.0*t9+1875.0/2.0*t11+1875.0/2.0*t13+625.0/2.0*t15-t18-t20-t22-t24-t26;
135  t29 = 1925.0/12.0*eta;
136  t31 = 8875.0/6.0*t5;
137  t35 = 4375.0/2.0*t13;
138  t37 = 15625.0/24.0*t17;
139  t38 = 6250.0/3.0*t19;
140  t39 = 9375.0/4.0*t21;
141  t40 = 3125.0/3.0*t23;
142  t41 = 25.0-1925.0/6.0*xi-t29+8875.0/8.0*t3+t31+8875.0/24.0*t7-4375.0/3.0*t9-13125.0/4.0*t11-t35-4375.0/12.0*t15+t37+t38+t39+t40+t26;
143  t43 = 1175.0/12.0*eta;
144  t48 = 3125.0/2.0*t13;
145  t49 = 625.0/12.0*t15;
146  t50 = 15625.0/12.0*t17;
147  t51 = 3125.0*t19;
148  t52 = -25.0+2675.0/6.0*xi+t43-7375.0/4.0*t3-t31-125.0*t7+8125.0/3.0*t9+16875.0/4.0*t11+t48+t49-t50-t51-t39-t24;
149  t54 = 75.0/2.0*eta;
150  t57 = 125.0/6.0*t7;
151  t60 = 625.0/2.0*t13;
152  t61 = 50.0/3.0-325.0*xi-t54+6125.0/4.0*t3+3875.0/6.0*t5+t57-2500.0*t9-9375.0/4.0*t11-t60+t50+t38+t22;
153  t63 = 25.0/4.0*eta;
154  t65 = 1375.0/12.0*t5;
155  t67 = 1875.0/4.0*t11;
156  t74 = 8875.0/12.0*t7;
157  t77 = 3125.0/2.0*t21;
158  t78 = 3125.0/2.0*t23;
159  t79 = 3125.0/6.0*t25;
160  t87 = 9375.0/2.0*t21;
161  t88 = 3125.0*t23;
162  t90 = 125.0*eta;
163  t94 = 6875.0/2.0*t13;
164  t95 = 625.0/4.0*t15;
165  t97 = 125.0/3.0*eta;
166  t101 = 625.0*t13;
167  t105 = 625.0/4.0*t11;
168  t107 = 3125.0/4.0*t25;
169  t117 = 25.0/6.0*eta;
170  t120 = 125.0/3.0*t5;
171  t131 = t63-1375.0/24.0*t7+t95-t26;
172 
173  values[0] = t27;
174  values[1] = t41;
175  values[2] = t52;
176  values[3] = t61;
177  values[4] = -25.0/4.0+1525.0/12.0*xi+t63-5125.0/8.0*t3-t65+6875.0/6.0*t9+t67-t37-t20;
178  values[5] = 1.0-125.0/6.0*xi+875.0/8.0*t3-625.0/3.0*t9+t18;
179  values[6] = -t29+8875.0/12.0*t5+t74-4375.0/4.0*t11-t35-4375.0/4.0*t15+t20+t77+t78+t79;
180  values[7] = 250.0*eta-5875.0/3.0*t5-5875.0/6.0*t7+3750.0*t11+5000.0*t13+1250.0*t15-t38-t87-t88-t79;
181  values[8] = -t90+3625.0/2.0*t5+1125.0/4.0*t7-9375.0/2.0*t11-t94-t95+t51+t87+t78;
182  values[9] = t97-2125.0/3.0*t5-125.0/3.0*t7+2500.0*t11+t101-t38-t77;
183  values[10] = -t63+t65-t67+t20;
184  values[11] = t43-250.0*t5-t74+t105+t48+5625.0/4.0*t15-t22-t78-t107;
185  values[12] = -t90+1125.0/2.0*t5+3625.0/4.0*t7-t67-t94-3125.0/2.0*t15+t39+t88+t107;
186  values[13] = 125.0/4.0*eta-375.0*t5-375.0/2.0*t7+t67+t35+t95-t39-t78;
187  values[14] = -t117+125.0/2.0*t5+t57-t105-t60+t22;
188  values[15] = -t54+t120+3875.0/12.0*t7-t60-3125.0/4.0*t15+t24+t79;
189  values[16] = t97-250.0/3.0*t5-2125.0/6.0*t7+t101+2500.0/3.0*t15-t40-t79;
190  values[17] = -t117+t120+125.0/4.0*t7-t60-t49+t24;
191  values[18] = t131;
192  values[19] = -t131;
193  values[20] = 0.0;
194 }
195 
196 // values of the derivatives in eta direction
197 static void C_T_P5_2D_DeriveEta(double xi, double eta, double *values)
198 {
199  double t3, t5, t7, t9, t11, t13, t15, t17, t18, t19, t20, t21, t22, t23;
200  double t24, t25, t26, t27, t28, t29, t32, t34, t35, t36, t38, t41, t42;
201  double t43, t45, t47, t49, t51, t53, t54, t57, t62, t63, t64, t65, t66;
202  double t73, t74, t76, t80, t81, t83, t87, t92, t95, t96, t97, t106, t110;
203  double t115, t123;
204 
205  t3 = xi*xi;
206  t5 = xi*eta;
207  t7 = eta*eta;
208  t9 = t3*xi;
209  t11 = t3*eta;
210  t13 = xi*t7;
211  t15 = t7*eta;
212  t17 = t3*t3;
213  t18 = 3125.0/24.0*t17;
214  t19 = t9*eta;
215  t20 = 3125.0/6.0*t19;
216  t21 = t3*t7;
217  t22 = 3125.0/4.0*t21;
218  t23 = xi*t15;
219  t24 = 3125.0/6.0*t23;
220  t25 = t7*t7;
221  t26 = 3125.0/24.0*t25;
222  t27 = -137.0/12.0+375.0/4.0*xi+375.0/4.0*eta-2125.0/8.0*t3-2125.0/4.0*t5-2125.0/8.0*t7+625.0/2.0*t9+1875.0/2.0*t11+1875.0/2.0*t13+625.0/2.0*t15-t18-t20-t22-t24-t26;
223  t28 = 1925.0/12.0*xi;
224  t29 = 8875.0/12.0*t3;
225  t32 = 4375.0/2.0*t11;
226  t34 = 3125.0/6.0*t17;
227  t35 = 3125.0/2.0*t19;
228  t36 = 3125.0/2.0*t21;
229  t38 = 1175.0/12.0*xi;
230  t41 = 3125.0/2.0*t11;
231  t42 = 625.0/4.0*t13;
232  t43 = 3125.0/4.0*t17;
233  t45 = 75.0/2.0*xi;
234  t47 = 125.0/3.0*t5;
235  t49 = 625.0/2.0*t11;
236  t51 = 25.0/4.0*xi;
237  t53 = 625.0/4.0*t9;
238  t54 = t51-1375.0/24.0*t3+t53-t18;
239  t57 = 8875.0/6.0*t5;
240  t62 = 3125.0/3.0*t19;
241  t63 = 9375.0/4.0*t21;
242  t64 = 6250.0/3.0*t23;
243  t65 = 15625.0/24.0*t25;
244  t66 = 25.0-t28-1925.0/6.0*eta+8875.0/24.0*t3+t57+8875.0/8.0*t7-4375.0/12.0*t9-t32-13125.0/4.0*t13-4375.0/3.0*t15+t18+t62+t63+t64+t65;
245  t73 = 3125.0*t19;
246  t74 = 9375.0/2.0*t21;
247  t76 = 125.0*xi;
248  t80 = 6875.0/2.0*t11;
249  t81 = 1875.0/4.0*t13;
250  t83 = 125.0/3.0*xi;
251  t87 = 625.0*t11;
252  t92 = 625.0/12.0*t9;
253  t95 = 3125.0*t23;
254  t96 = 15625.0/12.0*t25;
255  t97 = -25.0+t38+2675.0/6.0*eta-125.0*t3-t57-7375.0/4.0*t7+t92+t41+16875.0/4.0*t13+8125.0/3.0*t15-t20-t63-t95-t96;
256  t106 = 25.0/6.0*xi;
257  t110 = 125.0/6.0*t3;
258  t115 = 50.0/3.0-t45-325.0*eta+t110+3875.0/6.0*t5+6125.0/4.0*t7-t49-9375.0/4.0*t13-2500.0*t15+t22+t64+t96;
259  t123 = 1375.0/12.0*t5;
260 
261  values[0] = t27;
262  values[1] = -t28+t29+8875.0/12.0*t5-4375.0/4.0*t9-t32-4375.0/4.0*t13+t34+t35+t36+t24;
263  values[2] = t38-t29-250.0*t5+5625.0/4.0*t9+t41+t42-t43-t35-t22;
264  values[3] = -t45+3875.0/12.0*t3+t47-3125.0/4.0*t9-t49+t34+t20;
265  values[4] = t54;
266  values[5] = 0.0;
267  values[6] = t66;
268  values[7] = 250.0*xi-5875.0/6.0*t3-5875.0/3.0*t5+1250.0*t9+5000.0*t11+3750.0*t13-t34-t73-t74-t64;
269  values[8] = -t76+3625.0/4.0*t3+1125.0/2.0*t5-3125.0/2.0*t9-t80-t81+t43+t73+t63;
270  values[9] = t83-2125.0/6.0*t3-250.0/3.0*t5+2500.0/3.0*t9+t87-t34-t62;
271  values[10] = -t54;
272  values[11] = t97;
273  values[12] = -t76+1125.0/4.0*t3+3625.0/2.0*t5-t53-t80-9375.0/2.0*t13+t35+t74+ t95;
274  values[13] = 125.0/4.0*xi-375.0/2.0*t3-375.0*t5+t53+t32+t81-t35-t63;
275  values[14] = -t106+125.0/4.0*t3+t47-t92-t49+t20;
276  values[15] = t115;
277  values[16] = t83-125.0/3.0*t3-2125.0/3.0*t5+t87+2500.0*t13-t36-t64;
278  values[17] = -t106+t110+125.0/2.0*t5-t49-t42+t22;
279  values[18] = -25.0/4.0+t51+1525.0/12.0*eta-t123-5125.0/8.0*t7+t81+6875.0/6.0*t15-t24-t65;
280  values[19] = -t51+t123-t81+t24;
281  values[20] = 1.0-125.0/6.0*eta+875.0/8.0*t7-625.0/3.0*t15+t26;
282 }
283 
284 // values of the derivatives in xi-xi direction
285 static void C_T_P5_2D_DeriveXiXi(double xi, double eta, double *values)
286 {
287  double t3, t5, t7, t9, t10, t11, t12, t13, t14, t15, t16, t19, t22, t23;
288  double t24, t25, t26, t31, t32, t33, t39, t42, t44, t51, t52, t57, t58;
289  double t62, t66, t70, t79;
290 
291  t3 = xi*xi;
292  t5 = xi*eta;
293  t7 = eta*eta;
294  t9 = t3*xi;
295  t10 = 3125.0/6.0*t9;
296  t11 = t3*eta;
297  t12 = 3125.0/2.0*t11;
298  t13 = xi*t7;
299  t14 = 3125.0/2.0*t13;
300  t15 = t7*eta;
301  t16 = 3125.0/6.0*t15;
302  t19 = 8875.0/6.0*eta;
303  t22 = 4375.0/2.0*t7;
304  t23 = 15625.0/6.0*t9;
305  t24 = 6250.0*t11;
306  t25 = 9375.0/2.0*t13;
307  t26 = 3125.0/3.0*t15;
308  t31 = 3125.0/2.0*t7;
309  t32 = 15625.0/3.0*t9;
310  t33 = 9375.0*t11;
311  t39 = 625.0/2.0*t7;
312  t42 = 1375.0/12.0*eta;
313  t44 = 1875.0/2.0*t5;
314  t51 = 3125.0*t13;
315  t52 = 3125.0/2.0*t15;
316  t57 = 9375.0*t13;
317  t58 = 3125.0*t15;
318  t62 = 6875.0/2.0*t7;
319  t66 = 625.0*t7;
320  t70 = 625.0/2.0*t5;
321  t79 = 125.0/3.0*eta-t39+t16;
322 
323  values[0] = 375.0/4.0-2125.0/4.0*xi-2125.0/4.0*eta+1875.0/2.0*t3+1875.0*t5+1875.0/2.0*t7-t10-t12-t14-t16;
324  values[1] = -1925.0/6.0+8875.0/4.0*xi+t19-4375.0*t3-13125.0/2.0*t5-t22+t23+t24+t25+t26;
325  values[2] = 2675.0/6.0-7375.0/2.0*xi-t19+8125.0*t3+16875.0/2.0*t5+t31-t32-t33-t25-t16;
326  values[3] = -325.0+6125.0/2.0*xi+3875.0/6.0*eta-7500.0*t3-9375.0/2.0*t5-t39+t32+t24+t14;
327  values[4] = 1525.0/12.0-5125.0/4.0*xi-t42+6875.0/2.0*t3+t44-t23-t12;
328  values[5] = -125.0/6.0+875.0/4.0*xi-625.0*t3+t10;
329  values[6] = 8875.0/12.0*eta-4375.0/2.0*t5-t22+t12+t51+t52;
330  values[7] = -5875.0/3.0*eta+7500.0*t5+5000.0*t7-t24-t57-t58;
331  values[8] = 3625.0/2.0*eta-9375.0*t5-t62+t33+t57+t52;
332  values[9] = -2125.0/3.0*eta+5000.0*t5+t66-t24-t51;
333  values[10] = t42-t44+t12;
334  values[11] = -250.0*eta+t70+t31-t14-t52;
335  values[12] = 1125.0/2.0*eta-t44-t62+t25+t58;
336  values[13] = -375.0*eta+t44+t22-t25-t52;
337  values[14] = 125.0/2.0*eta-t70-t39+t14;
338  values[15] = t79;
339  values[16] = -250.0/3.0*eta+t66-t26;
340  values[17] = t79;
341  values[18] = 0.0;
342  values[19] = 0.0;
343  values[20] = 0.0;
344 }
345 
346 // values of the derivatives in xi-eta direction
347 static void C_T_P5_2D_DeriveXiEta(double xi, double eta, double *values)
348 {
349  double t3, t5, t7, t9, t10, t11, t12, t13, t14, t15, t16, t18, t21, t23;
350  double t24, t25, t29, t30, t31, t34, t36, t39, t40, t42, t45, t46, t47;
351  double t54, t55, t60, t61, t66, t69, t71, t82, t93;
352 
353  t3 = xi*xi;
354  t5 = xi*eta;
355  t7 = eta*eta;
356  t9 = t3*xi;
357  t10 = 3125.0/6.0*t9;
358  t11 = t3*eta;
359  t12 = 3125.0/2.0*t11;
360  t13 = xi*t7;
361  t14 = 3125.0/2.0*t13;
362  t15 = t7*eta;
363  t16 = 3125.0/6.0*t15;
364  t18 = 8875.0/6.0*xi;
365  t21 = 4375.0*t5;
366  t23 = 6250.0/3.0*t9;
367  t24 = 9375.0/2.0*t11;
368  t25 = 3125.0*t13;
369  t29 = 3125.0*t5;
370  t30 = 625.0/4.0*t7;
371  t31 = 3125.0*t9;
372  t34 = 125.0/3.0*eta;
373  t36 = 625.0*t5;
374  t39 = 1875.0/4.0*t3;
375  t40 = 25.0/4.0-1375.0/12.0*xi+t39-t10;
376  t42 = 8875.0/6.0*eta;
377  t45 = 3125.0*t11;
378  t46 = 9375.0/2.0*t13;
379  t47 = 6250.0/3.0*t15;
380  t54 = 9375.0*t11;
381  t55 = 9375.0*t13;
382  t60 = 6875.0*t5;
383  t61 = 1875.0/4.0*t7;
384  t66 = 1250.0*t5;
385  t69 = 625.0/4.0*t3;
386  t71 = 3125.0*t15;
387  t82 = 125.0/3.0*xi;
388  t93 = 25.0/4.0-1375.0/12.0*eta+t61-t16;
389 
390  values[0] = 375.0/4.0-2125.0/4.0*xi-2125.0/4.0*eta+1875.0/2.0*t3+1875.0*t5+1875.0/2.0*t7-t10-t12-t14-t16;
391  values[1] = -1925.0/12.0+t18+8875.0/12.0*eta-13125.0/4.0*t3-t21-4375.0/4.0*t7+t23+t24+t25+t16;
392  values[2] = 1175.0/12.0-t18-250.0*eta+16875.0/4.0*t3+t29+t30-t31-t24-t14;
393  values[3] = -75.0/2.0+3875.0/6.0*xi+t34-9375.0/4.0*t3-t36+t23+t12;
394  values[4] = t40;
395  values[5] = 0.0;
396  values[6] = -1925.0/12.0+8875.0/12.0*xi+t42-4375.0/4.0*t3-t21-13125.0/4.0*t7+t10+t45+t46+t47;
397  values[7] = 250.0-5875.0/3.0*xi-5875.0/3.0*eta+3750.0*t3+10000.0*t5+3750.0*t7-t23-t54-t55-t47;
398  values[8] = -125.0+3625.0/2.0*xi+1125.0/2.0*eta-9375.0/2.0*t3-t60-t61+t31+t54+t46;
399  values[9] = 125.0/3.0-2125.0/3.0*xi-250.0/3.0*eta+2500.0*t3+t66-t23-t45;
400  values[10] = -t40;
401  values[11] = 1175.0/12.0-250.0*xi-t42+t69+t29+16875.0/4.0*t7-t12-t46-t71;
402  values[12] = -125.0+1125.0/2.0*xi+3625.0/2.0*eta-t39-t60-9375.0/2.0*t7+t24+t55+t71;
403  values[13] = 125.0/4.0-375.0*xi-375.0*eta+t39+t21+t61-t24-t46;
404  values[14] = -25.0/6.0+125.0/2.0*xi+t34-t69-t36+t12;
405  values[15] = -75.0/2.0+t82+3875.0/6.0*eta-t36-9375.0/4.0*t7+t14+t47;
406  values[16] = 125.0/3.0-250.0/3.0*xi-2125.0/3.0*eta+t66+2500.0*t7-t25-t47;
407  values[17] = -25.0/6.0+t82+125.0/2.0*eta-t36-t30+t14;
408  values[18] = t93;
409  values[19] = -t93;
410  values[20] = 0.0;
411 }
412 
413 // values of the derivatives in eta-eta direction
414 static void C_T_P5_2D_DeriveEtaEta(double xi, double eta, double *values)
415 {
416  double t3, t5, t7, t9, t10, t11, t12, t13, t14, t15, t16, t19, t21, t22;
417  double t25, t26, t29, t30, t31, t35, t36, t37, t38, t43, t44, t47, t48;
418  double t51, t56, t57, t74;
419 
420  t3 = xi*xi;
421  t5 = xi*eta;
422  t7 = eta*eta;
423  t9 = t3*xi;
424  t10 = 3125.0/6.0*t9;
425  t11 = t3*eta;
426  t12 = 3125.0/2.0*t11;
427  t13 = xi*t7;
428  t14 = 3125.0/2.0*t13;
429  t15 = t7*eta;
430  t16 = 3125.0/6.0*t15;
431  t19 = 4375.0/2.0*t3;
432  t21 = 3125.0/2.0*t9;
433  t22 = 3125.0*t11;
434  t25 = 3125.0/2.0*t3;
435  t26 = 625.0/2.0*t5;
436  t29 = 625.0/2.0*t3;
437  t30 = 125.0/3.0*xi-t29+t10;
438  t31 = 8875.0/6.0*xi;
439  t35 = 3125.0/3.0*t9;
440  t36 = 9375.0/2.0*t11;
441  t37 = 6250.0*t13;
442  t38 = 15625.0/6.0*t15;
443  t43 = 3125.0*t9;
444  t44 = 9375.0*t11;
445  t47 = 6875.0/2.0*t3;
446  t48 = 1875.0/2.0*t5;
447  t51 = 625.0*t3;
448  t56 = 9375.0*t13;
449  t57 = 15625.0/3.0*t15;
450  t74 = 1375.0/12.0*xi;
451 
452  values[0] = 375.0/4.0-2125.0/4.0*xi-2125.0/4.0*eta+1875.0/2.0*t3+1875.0*t5+1875.0/2.0*t7-t10-t12-t14-t16;
453  values[1] = 8875.0/12.0*xi-t19-4375.0/2.0*t5+t21+t22+t14;
454  values[2] = -250.0*xi+t25+t26-t21-t12;
455  values[3] = t30;
456  values[4] = 0.0;
457  values[5] = 0.0;
458  values[6] = -1925.0/6.0+t31+8875.0/4.0*eta-t19-13125.0/2.0*t5-4375.0*t7+t35+t36+t37+t38;
459  values[7] = -5875.0/3.0*xi+5000.0*t3+7500.0*t5-t43-t44-t37;
460  values[8] = 1125.0/2.0*xi-t47-t48+t43+t36;
461  values[9] = -250.0/3.0*xi+t51-t35;
462  values[10] = 0.0;
463  values[11] = 2675.0/6.0-t31-7375.0/2.0*eta+t25+16875.0/2.0*t5+8125.0*t7-t10-t36-t56-t57;
464  values[12] = 3625.0/2.0*xi-t47-9375.0*t5+t21+t44+t56;
465  values[13] = -375.0*xi+t19+t48-t21-t36;
466  values[14] = t30;
467  values[15] = -325.0+3875.0/6.0*xi+6125.0/2.0*eta-t29-9375.0/2.0*t5-7500.0*t7+t12+t37+t57;
468  values[16] = -2125.0/3.0*xi+t51+5000.0*t5-t22-t37;
469  values[17] = 125.0/2.0*xi-t29-t26+t12;
470  values[18] = 1525.0/12.0-t74-5125.0/4.0*eta+t48+6875.0/2.0*t7-t14-t38;
471  values[19] = t74-t48+t14;
472  values[20] = -125.0/6.0+875.0/4.0*eta-625.0*t7+t16;
473 }
474 
475 // ***********************************************************************
476 
477 TBaseFunct2D *BF_C_T_P5_2D_Obj = new TBaseFunct2D
478  (21, BF_C_T_P5_2D, BFUnitTriangle,
479  C_T_P5_2D_Funct, C_T_P5_2D_DeriveXi,
480  C_T_P5_2D_DeriveEta, C_T_P5_2D_DeriveXiXi,
481  C_T_P5_2D_DeriveXiEta, C_T_P5_2D_DeriveEtaEta, 5, 5,
482  0, NULL);
Definition: BaseFunct2D.h:27