xref: /petsc/src/ksp/ksp/tutorials/ex43-solcx.h (revision d8e47b638cf8f604a99e9678e1df24f82d959cd7) !
1 #pragma once
2 
evaluate_solCx(PetscReal pos[],PetscReal _eta_A,PetscReal _eta_B,PetscReal _x_c,PetscInt _n,PetscReal vel[],PetscReal * pressure,PetscReal total_stress[],PetscReal strain_rate[])3 static void evaluate_solCx(PetscReal pos[], PetscReal _eta_A, PetscReal _eta_B, /* Input parameters: density, viscosity A, viscosity B */
4                            PetscReal _x_c, PetscInt _n,                         /* Input parameters: viscosity jump location, wavenumber in z */
5                            PetscReal vel[], PetscReal *pressure, PetscReal total_stress[], PetscReal strain_rate[])
6 {
7   PetscReal Z, u1, u2, u3, u4, u5, u6, ZA, ZB;
8   PetscReal sum1, sum2, sum3, sum4, sum5, sum6, x, z, xc;
9   PetscReal _PC1A, _PC2A, _PC3A, _PC4A, _PC1B, _PC2B, _PC3B, _PC4B, _PC1, _PC2, _PC3, _PC4;
10   PetscInt  n, nx;
11 
12   PetscReal t1, t2, t3, t4, t5, t6, t7, t8, t9, t10, t11, t12, t13, t14, t15, t16, t17, t18, t19, t20, t21, t22, t23, t24, t25, t26, t27, t28, t29, t30, t31, t32, t33, t34, t35, t36, t37, t38, t39, t40;
13   PetscReal t41, t42, t43, t44, t45, t46, t47, t48, t49, t50, t51, t52, t53, t54, t55, t56, t57, t58, t59, t60, t61, t62, t63, t64, t65, t66, t67, t68, t69, t70, t71, t72, t73, t74, t75, t76, t77, t78, t79, t80;
14   PetscReal t81, t82, t83, t84, t85, t86, t87, t88, t89, t90, t91, t92, t93, t94, t95, t96, t97, t98, t99, t100, t101, t102, t103, t104, t105, t106, t107, t108, t109, t110, t111, t112, t113, t115, t116, t117, t118, t119, t120;
15   PetscReal t121, t122, t123, t124, t125, t126, t127, t128, t129, t130, t131, t132, t133, t134, t135, t136, t137, t138, t139, t140, t141, t142, t143, t144, t145, t146, t147, t148, t149, t150, t151, t152, t153, t154, t155, t156, t157, t158, t159, t160;
16   PetscReal t161, t162, t163, t164, t165, t166, t167, t168, t169, t170, t171, t172, t173, t174, t175, t176, t177, t178, t179, t180, t181, t182, t183, t184, t186, t187, t188, t189, t190, t191, t192, t193, t194, t195, t196, t197, t198, t199;
17   PetscReal t201, t202, t203, t204, t206, t207, t208, t209, t210, t211, t212, t213, t215, t216, t217, t218, t219, t220, t221, t222, t223, t224, t225, t226, t227, t228, t229, t230, t231, t232, t233, t234, t235, t236, t237, t238, t239, t240;
18   PetscReal t241, t242, t243, t244, t245, t246, t247, t248, t249, t250, t251, t252, t253, t254, t255, t256, t257, t258, t259, t260, t261, t262, t263, t264, t265, t267, t268, t269, t270, t272, t273, t274, t275, t276, t277, t278, t279, t280;
19   PetscReal t281, t282, t283, t284, t285, t286, t288, t289, t290, t291, t292, t295, t296, t297, t298, t299, t300, t301, t303, t304, t305, t307, t308, t310, t311, t312, t313, t314, t315, t316, t317, t318, t319, t320;
20   PetscReal t321, t322, t323, t324, t325, t326, t327, t328, t329, t330, t331, t332, t334, t335, t336, t337, t338, t339, t340, t341, t342, t344, t345, t346, t347, t348, t349, t350, t351, t352, t353, t354, t355, t356, t358, t359, t360;
21   PetscReal t361, t362, t363, t364, t365, t366, t367, t368, t369, t370, t371, t372, t373, t374, t375, t376, t377, t378, t379, t380, t381, t382, t383, t384, t385, t386, t387, t389, t390, t391, t393, t394, t395, t396, t397, t398;
22   PetscReal t401, t402, t403, t404, t405, t406, t407, t408, t409, t410, t411, t412, t413, t414, t415, t416, t417, t418, t419, t421, t422, t423, t424, t425, t426, t427, t428, t429, t430, t431, t432, t433, t434, t436, t437, t438, t439, t440;
23   PetscReal t441, t442, t443, t444, t445, t446, t447, t448, t450, t451, t453, t454, t455, t456, t457, t458, t459, t461, t462, t463, t464, t465, t466, t468, t469, t470, t471, t474, t475, t478, t480;
24   PetscReal t482, t483, t484, t485, t488, t489, t490, t492, t493, t495, t497, t498, t499, t501, t502, t503, t504, t505, t507, t508, t509, t510, t511, t512, t513, t515, t518, t520;
25   PetscReal t522, t525, t527, t528, t529, t530, t532, t533, t534, t535, t536, t538, t539, t541, t542, t544, t545, t546, t547, t548, t549, t550, t551, t552, t553, t554, t555, t556, t557, t560;
26   PetscReal t561, t562, t563, t564, t567, t568, t571, t573, t575, t576, t578, t579, t583, t590, t591, t594, t595, t596, t597, t598, t600;
27   PetscReal t601, t602, t604, t606, t607, t608, t611, t613, t615, t616, t617, t619, t621, t623, t624, t625, t626, t627, t629, t630, t632, t633, t634, t638, t639, t640;
28   PetscReal t641, t642, t643, t644, t645, t647, t648, t649, t650, t651, t652, t653, t654, t655, t656, t657, t658, t659, t660, t662, t663, t665, t666, t667, t668, t670, t671, t672, t673, t674, t675, t676, t679, t680;
29   PetscReal t682, t683, t684, t685, t686, t688, t689, t690, t691, t693, t694, t695, t696, t697, t698, t699, t700, t701, t702, t704, t705, t708, t709, t711, t712, t713, t714, t717, t718, t719;
30   PetscReal t721, t722, t723, t726, t727, t728, t730, t733, t734, t735, t736, t737, t738, t739, t740, t741, t744, t745, t746, t749, t750, t752, t753, t754, t755, t757, t758, t759, t760;
31   PetscReal t761, t762, t763, t764, t766, t767, t768, t770, t771, t772, t773, t774, t775, t776, t777, t778, t780, t781, t782, t785, t786, t789, t790, t791, t792, t793, t794, t795, t796, t797, t798, t800;
32   PetscReal t801, t806, t807, t808, t809, t811, t812, t817, t818, t819, t821, t822, t824, t827, t828, t830, t834, t835, t837, t840;
33   PetscReal t842, t843, t844, t845, t846, t849, t850, t853, t854, t855, t857, t858, t859, t860, t863, t864, t867, t868, t869, t873, t874, t877, t878, t879, t880;
34   PetscReal t884, t888, t891, t894, t900, t901, t903, t904, t907, t908, t909, t911, t914, t915, t916, t919, t920;
35   PetscReal t923, t924, t925, t926, t927, t929, t932, t935, t937, t939, t942, t943, t944, t945, t947, t948, t949, t950, t952, t953, t954, t955, t956, t957;
36   PetscReal t961, t964, t965, t966, t967, t968, t969, t971, t972, t974, t977, t978, t981, t983, t987, t988, t992, t993, t994, t997, t998;
37   PetscReal t1001, t1003, t1005, t1006, t1009, t1010, t1012, t1013, t1015, t1016, t1017, t1018, t1020, t1021, t1029, t1031, t1032, t1033, t1040;
38   PetscReal t1041, t1042, t1044, t1047, t1050, t1054, t1055, t1057, t1058, t1063, t1068, t1069, t1070, t1079, t1080;
39   PetscReal t1088, t1089, t1091, t1092, t1094, t1096, t1101, t1102, t1103, t1104, t1105, t1108, t1112, t1113, t1118, t1119, t1120;
40   PetscReal t1121, t1122, t1123, t1124, t1125, t1126, t1127, t1128, t1129, t1130, t1132, t1133, t1134, t1135, t1138, t1139, t1140, t1141, t1142, t1145, t1146, t1148, t1149, t1150, t1153, t1154, t1156, t1157, t1158, t1159;
41   PetscReal t1161, t1162, t1165, t1166, t1170, t1171, t1172, t1173, t1175, t1176, t1178, t1180, t1181, t1182, t1185, t1189, t1192, t1193, t1195, t1196, t1199;
42   PetscReal t1201, t1203, t1209, t1210, t1211, t1213, t1214, t1218, t1221, t1224, t1225, t1226, t1228, t1233, t1234, t1235, t1236, t1237, t1240;
43   PetscReal t1241, t1242, t1243, t1244, t1245, t1248, t1251, t1252, t1257, t1258, t1259, t1260, t1263, t1268, t1269, t1272, t1280;
44   PetscReal t1282, t1283, t1284, t1285, t1287, t1288, t1289, t1292, t1293, t1296, t1297, t1300, t1304, t1307, t1310, t1311, t1312, t1316, t1317, t1320;
45   PetscReal t1321, t1323, t1328, t1330, t1331, t1332, t1333, t1336, t1338, t1343, t1344, t1346, t1349, t1350, t1354;
46   PetscReal t1366, t1369, t1370, t1371, t1376, t1378, t1380, t1383, t1386, t1387, t1388, t1391, t1393, t1399;
47   PetscReal t1411, t1412, t1420, t1427;
48   PetscReal t1450, t1456, t1468, t1472, t1474, t1478;
49   PetscReal t1504, t1511;
50   PetscReal t1545;
51   PetscReal t1564, t1583;
52 
53   /* del_rho = sin(n*Pi*z)*cos(nx*Pi*x)  n=nx gives only non-zero terms*/
54   n  = _n; /* only one n in Fourier series because  del_rho has cos term */
55   nx = 1;  /* can set these two to whatever we like */
56 
57   ZA = _eta_A; /* left column viscosity */
58   ZB = _eta_B; /* right column viscosity */
59   xc = _x_c;
60 
61   x = pos[0];
62   z = pos[1];
63 
64   sum1 = 0.0;
65   sum2 = 0.0;
66   sum3 = 0.0;
67   sum4 = 0.0;
68   sum5 = 0.0;
69   sum6 = 0.0;
70 
71   /* Note that there is no Fourier sum here. */
72   _PC1A = 0;
73   t1    = nx * 0.3141592654e1;
74   t2    = PetscSinReal(t1);
75   t3    = nx * t2;
76   t4    = n * n;
77   t5    = t4 * t4;
78   t6    = 0.3141592654e1 * 0.3141592654e1;
79   t8    = t3 * t5 * t6;
80   t9    = ZA * xc;
81   t12   = PetscExpReal(xc * n * 0.3141592654e1);
82   t13   = t12 * t12;
83   t15   = n * 0.3141592654e1;
84   t16   = PetscExpReal(t15);
85   t17   = t16 * t16;
86   t18   = t17 * t16;
87   t19   = ZB * t13 * t18;
88   t20   = t9 * t19;
89   t23   = ZA * ZA;
90   t24   = nx * nx;
91   t25   = t24 * nx;
92   t26   = t23 * t25;
93   t28   = t13 * t13;
94   t29   = t28 * t13;
95   t33   = nx * ZB;
96   t34   = t1 * xc;
97   t35   = PetscSinReal(t34);
98   t36   = t4 * n;
99   t37   = t35 * t36;
100   t38   = t33 * t37;
101   t39   = 0.3141592654e1 * ZA;
102   t40   = t13 * t12;
103   t41   = t17 * t40;
104   t45   = ZB * ZB;
105   t46   = t45 * t24;
106   t47   = t46 * t4;
107   t48   = 0.3141592654e1 * xc;
108   t49   = t13 * t17;
109   t53   = xc * xc;
110   t54   = t36 * t53;
111   t56   = t54 * t6 * t45;
112   t57   = PetscCosReal(t34);
113   t58   = t57 * t24;
114   t59   = t28 * t12;
115   t60   = t17 * t59;
116   t61   = t58 * t60;
117   t64   = t25 * t2;
118   t65   = t64 * t15;
119   t72   = nx * t23;
120   t74   = t72 * t2 * t5;
121   t75   = t6 * t53;
122   t76   = t16 * t29;
123   t80   = t23 * n;
124   t81   = t80 * 0.3141592654e1;
125   t82   = t18 * t28;
126   t86   = nx * t5;
127   t87   = t23 * t6;
128   t89   = xc * t2;
129   t90   = t13 * t18;
130   t91   = t89 * t90;
131   t94   = t28 * t28;
132   t96   = t24 * n;
133   t98   = t4 * t45;
134   t99   = t98 * 0.3141592654e1;
135   t100  = t58 * t41;
136   t104  = 0.3141592654e1 * t25;
137   t105  = ZA * n * t104;
138   t106  = t2 * ZB;
139   t110  = t17 * t17;
140   t111  = ZA * t110;
141   t116  = n * t28;
142   t122  = t64 * t4 * t6;
143   t126  = t23 * t29 * t4;
144   t128  = t24 * xc;
145   t132  = t36 * t23;
146   t133  = t6 * t57;
147   t135  = t128 * t41;
148   t138  = t6 * xc;
149   t142  = t72 * t2;
150   t147 = 0.4e1 * t8 * t20 - 0.2e1 * t26 * t2 * t16 * t29 - 0.8e1 * t38 * t39 * t41 + 0.4e1 * t47 * t48 * t49 - 0.8e1 * t56 * t61 - 0.4e1 * t65 * t20 + 0.2e1 * t26 * t2 * t18 * t28 - 0.4e1 * t74 * t75 * t76 - 0.2e1 * t81 * t64 * t82 - 0.4e1 * t86 * t87 * t91 - t23 * t94 * t96 + 0.8e1 * t99 * t100 - 0.2e1 * t105 * t106 * t82 - 0.4e1 * t38 * t48 * t111 * t12 + 0.2e1 * t116 * ZB * t111 * t24 + 0.4e1 * t122 * t20 + 0.4e1 * t126 * 0.3141592654e1 * t17 * t128 + 0.8e1 * t132 * t133 * t135 + 0.4e1 * t74 * t138 * t76 - 0.2e1 * t142 * t4 * t18 * t28;
151   t149 = ZA * t25 * t2;
152   t150 = ZB * t28;
153   t154 = t35 * t5;
154   t155 = t72 * t154;
155   t156 = t75 * t41;
156   t159 = nx * ZA;
157   t160 = t2 * t36;
158   t161 = t159 * t160;
159   t162 = 0.3141592654e1 * ZB;
160   t163 = t28 * t16;
161   t167 = t23 * t57;
162   t168 = t167 * t24;
163   t169 = n * t110;
164   t170 = t169 * t40;
165   t173 = ZA * ZB;
166   t174 = t173 * t90;
167   t177 = t36 * 0.3141592654e1;
168   t181 = t80 * t104;
169   t184 = n * t17;
170   t188 = t17 * t29;
171   t190 = t4 * 0.3141592654e1;
172   t191 = t190 * t24;
173   t206 = t138 * t60;
174   t209 = t23 * t4;
175   t211 = t209 * t6 * t25;
176   t212 = t89 * t76;
177   t216 = ZB * t16 * t29;
178   t217 = t9 * t216;
179   t220 = ZB * t110;
180   t221 = ZA * t24;
181   t222 = t221 * n;
182   t225 = t132 * t75;
183   t232 = t45 * t28;
184   t233 = t110 * t24;
185   t234 = t233 * n;
186   t236 = t209 * 0.3141592654e1;
187   t237 = t17 * xc;
188   t239 = t237 * t13 * t24;
189   t242 = -0.2e1 * t149 * t150 * t16 - 0.8e1 * t155 * t156 - 0.2e1 * t161 * t162 * t163 + 0.2e1 * t168 * t170 + 0.2e1 * t65 * t174 - 0.2e1 * t142 * t177 * t76 + 0.4e1 * t181 * t91 - 0.4e1 * t168 * t184 * t59 - 0.4e1 * t188 * t23 * t191 + 0.4e1 * t38 * t48 * ZA * t17 * t40 + 0.4e1 * t49 * t23 * t191 + 0.2e1 * t26 * t2 * t13 * t18 - 0.8e1 * t155 * t206 + 0.4e1 * t211 * t212 - 0.4e1 * t8 * t217 + 0.2e1 * t220 * t222 - 0.8e1 * t225 * t100 + 0.2e1 * t142 * t4 * t16 * t29 + t232 * t234 - 0.4e1 * t236 * t239;
190   t244 = nx * t45;
191   t245 = t244 * t37;
192   t246 = t110 * t40;
193   t251 = t237 * t59;
194   t256 = t64 * t90;
195   t260 = t36 * t45 * t133;
196   t263 = t45 * t57;
197   t264 = t263 * t24;
198   t265 = t169 * t12;
199   t269 = t6 * t36;
200   t270 = t17 * t24;
201   t274 = t110 * t13;
202   t276 = t190 * t128;
203   t279 = nx * t36;
204   t281 = t28 * t40;
205   t282 = t281 * t35;
206   t286 = t138 * t41;
207   t289 = t75 * t60;
208   t296 = t190 * t173;
209   t305 = t86 * t45 * t35;
210   t312 = t33 * t154;
211   t313 = t6 * ZA;
212   t324 = t232 * t270;
213   t327 = -0.2e1 * t245 * t48 * t246 + 0.4e1 * t159 * t37 * t162 * t251 + 0.4e1 * t209 * t75 * t256 + 0.8e1 * t260 * t135 + 0.2e1 * t264 * t265 + 0.32e2 * t9 * t150 * t269 * t270 + 0.4e1 * t274 * t23 * t276 + 0.2e1 * t279 * t45 * t282 * t48 + 0.8e1 * t155 * t286 + 0.8e1 * t155 * t289 - 0.8e1 * t150 * ZA * t96 * t17 + 0.8e1 * t296 * t61 - 0.2e1 * t105 * t106 * t163 - 0.2e1 * t81 * t256 - 0.8e1 * t305 * t156 - 0.4e1 * t33 * t282 * t177 * t9 - 0.16e2 * t312 * t313 * t237 * t40 - 0.4e1 * t168 * t184 * t40 + 0.2e1 * t168 * t265 + 0.16e2 * t269 * t53 * t324;
214   t328 = t3 * t4;
215   t331 = t72 * t37;
216   t332 = t48 * t60;
217   t335 = n * t94;
218   t345 = t72 * t35;
219   t349 = t173 * t57;
220   t355 = t53 * t17;
221   t364 = t54 * t6 * ZB;
222   t365 = t28 * t17;
223   t369 = xc * ZB;
224   t370 = t269 * t369;
225   t371 = ZA * t57;
226   t373 = t371 * t270 * t40;
227   t385 = nx * t35;
228   t396 = t4 * xc;
229   t397 = t396 * t162;
230   t415 = t37 * t48;
231   t418 = -0.32e2 * t364 * t365 * t221 - 0.16e2 * t370 * t373 - 0.4e1 * t331 * t48 * t41 + 0.4e1 * t86 * t23 * t53 * t6 * t2 * t90 + 0.2e1 * t385 * t177 * t23 * xc * t246 + 0.16e2 * t132 * t53 * t6 * t28 * t270 - 0.4e1 * t397 * t371 * t233 * t12 - 0.12e2 * t173 * t58 * t190 * t251 + 0.2e1 * t385 * t36 * 0.3141592654e1 * t23 * xc * t59 - 0.8e1 * t99 * t61 - 0.2e1 * t244 * t59 * t415;
232   t427 = t371 * t270 * t59;
233   t439 = t209 * t48;
234   t440 = t110 * t12;
235   t441 = t58 * t440;
236   t447 = t36 * xc;
237   t455 = t48 * t440;
238   t471 = ZB * t17;
239   t492 = 0.12e2 * t397 * t373 - 0.4e1 * t122 * t217 + 0.16e2 * t364 * t427 + 0.16e2 * t312 * t313 * t355 * t40 - 0.8e1 * t279 * t39 * t35 * ZB * t60 + 0.2e1 * t439 * t441 - 0.2e1 * t81 * t64 * t163 + 0.8e1 * t447 * t87 * t61 + 0.2e1 * t23 * t59 * t57 * t276 + 0.2e1 * t245 * t455 - 0.4e1 * t349 * t96 * t440 - 0.16e2 * t370 * t427 + 0.4e1 * t181 * t212 - 0.16e2 * t365 * t23 * t269 * t128 + 0.16e2 * t86 * t138 * ZA * t35 * t471 * t59 + 0.8e1 * t305 * t289 - 0.4e1 * t439 * t100 + 0.2e1 * ZB * t25 * t2 * ZA * t18 * t28 + 0.2e1 * t142 * t4 * t28 * t16 - 0.8e1 * t56 * t100;
240   t499 = ZA * t53 * t19;
241   t505 = t396 * 0.3141592654e1;
242   t518 = t173 * t53 * t16 * t29;
243   t533 = t23 * t28;
244   t535 = t188 * t45;
245   t538 = t24 * t4;
246   t545 = t3 * t177;
247   t546 = t173 * t76;
248   t555 = t45 * t110;
249   t557 = t72 * t160;
250   t561 = -0.8e1 * t225 * t61 - 0.2e1 * t161 * t162 * t82 + t533 * t234 + 0.4e1 * t535 * t191 + 0.4e1 * t167 * t538 * t332 + 0.4e1 * t349 * t96 * t60 + 0.2e1 * t545 * t546 - 0.2e1 * t264 * t170 + 0.4e1 * t397 * t281 * ZA * t58 - t555 * t96 - 0.4e1 * t557 * t48 * t76;
251   t567 = t396 * 0.3141592654e1 * t45;
252   t568 = t58 * t246;
253   t597 = t58 * n;
254   t615 = t13 * t45;
255   t616 = t615 * t233;
256   t619 = t94 * t45;
257   t621 = t45 * t59;
258   t625 = 0.2e1 * t149 * t216 + 0.2e1 * t567 * t568 - 0.16e2 * t269 * xc * t324 - 0.2e1 * t236 * xc * t281 * t58 - 0.2e1 * t142 * t177 * t90 - 0.8e1 * t567 * t100 + 0.2e1 * t65 * t546 - 0.8e1 * t305 * t206 + 0.2e1 * n * t45 * t281 * t57 * t24 - t23 * t110 * t96 - 0.8e1 * t296 * t100 + 0.2e1 * t23 * t281 * t597 + 0.4e1 * t545 * t20 + 0.2e1 * t159 * t2 * t4 * ZB * t163 - 0.4e1 * t557 * t48 * t90 + 0.4e1 * t122 * t518 + 0.8e1 * t263 * t538 * t332 - 0.4e1 * t505 * t616 - t619 * t96 - 0.2e1 * t621 * t57 * t276;
259   t626 = t49 * t45;
260   t660 = t29 * t45;
261   t685 = 0.2e1 * t545 * t174 - 0.4e1 * t126 * 0.3141592654e1 * t24 * xc - 0.4e1 * t47 * t48 * t188 + 0.4e1 * t505 * t660 * t24 - 0.2e1 * t142 * t177 * t163 - 0.2e1 * t142 * t4 * t13 * t18 + 0.8e1 * t260 * t128 * t60 - 0.2e1 * t328 * t546 - 0.2e1 * t26 * t2 * t28 * t16 + 0.4e1 * t545 * t217 - 0.4e1 * t209 * t138 * t256;
262   t690 = t6 * 0.3141592654e1;
263   t691 = ZA * t690;
264   t693 = t24 * t24;
265   t694 = t693 * xc;
266   t695 = t188 * t694;
267   t698 = t23 * ZA;
268   t699 = t698 * t690;
269   t700 = t699 * t5;
270   t704 = t5 * t4;
271   t705 = t691 * t704;
272   t709 = t691 * t5;
273   t713 = t5 * n;
274   t714 = t713 * ZB;
275   t718 = t698 * t6;
276   t719 = t713 * t28;
277   t722 = t699 * t704;
278   t726 = t713 * t94;
279   t733 = t713 * t45;
280   t736 = t87 * t36;
281   t740 = -0.4e1 * t691 * t98 * t695 + 0.8e1 * t700 * t270 * t13 + 0.4e1 * t705 * t660 * xc + 0.8e1 * t709 * t660 * t128 + 0.2e1 * t87 * t714 * t110 + t718 * t719 * t110 - 0.4e1 * t722 * t237 * t13 - t313 * t726 * t45 - 0.4e1 * t699 * t704 * xc * t29 + t313 * t733 * t28 + 0.4e1 * t736 * t150 * t233;
282   t746 = t313 * t36;
283   t752 = t6 * t6;
284   t753 = t23 * t752;
285   t759 = t698 * t752;
286   t760 = t759 * t36;
287   t761 = t17 * t693;
288   t762 = xc * t28;
289   t763 = t761 * t762;
290   t766 = t87 * t713;
291   t773 = t699 * t4;
292   t774 = t110 * t693;
293   t775 = xc * t13;
294   t785 = t704 * t17;
295   t789 = -0.16e2 * t736 * t150 * t270 + t718 * t116 * t693 - 0.2e1 * t746 * t555 * t24 + 0.4e1 * t705 * t535 + 0.64e2 * t753 * t713 * t17 * t150 * t128 - 0.16e2 * t760 * t763 + 0.2e1 * t766 * t150 * t110 + 0.4e1 * t722 * t274 * xc + 0.4e1 * t773 * t774 * t775 - 0.8e1 * t766 * t150 * t17 + 0.8e1 * t700 * t233 * t775 + 0.4e1 * t699 * t785 * t13;
296   t791 = t691 * t4;
297   t792 = t45 * t693;
298   t793 = t49 * t792;
299   t796 = t759 * t713;
300   t797 = t53 * t28;
301   t798 = t270 * t797;
302   t801 = t87 * n;
303   t818 = t5 * t36;
304   t819 = t753 * t818;
305   t827 = t753 * t36 * ZB;
306   t830 = xc * t45;
307   t834 = -0.4e1 * t791 * t793 + 0.32e2 * t796 * t798 + 0.2e1 * t801 * ZB * t693 * t110 + 0.2e1 * t718 * t36 * t28 * t24 - 0.8e1 * t700 * t128 * t29 - 0.8e1 * t700 * t239 - 0.8e1 * t801 * t150 * t761 + 0.32e2 * t819 * t365 * t369 - 0.64e2 * t753 * t714 * t798 + 0.32e2 * t827 * t763 + 0.4e1 * t705 * t830 * t49;
308   t842 = xc * t29;
309   t843 = t270 * t842;
310   t849 = t759 * t818;
311   t853 = t691 * t396;
312   t857 = t691 * t5 * t45;
313   t869 = t313 * n;
314   t874 = -0.2e1 * t718 * t36 * t94 * t24 - 0.4e1 * t773 * t761 * t29 + 0.8e1 * t700 * t843 + 0.2e1 * t87 * t726 * ZB + 0.16e2 * t849 * t797 * t17 + 0.4e1 * t853 * t793 + 0.8e1 * t857 * t239 + 0.2e1 * t801 * t150 * t693 - 0.8e1 * t700 * t270 * t29 - 0.8e1 * t709 * t49 * t46 - t869 * t619 * t693 + t869 * t232 * t693;
315   t877 = ZA * t752;
316   t878 = t877 * t818;
317   t911 = 0.16e2 * t878 * t53 * t45 * t365 - 0.4e1 * t699 * t785 * t29 - 0.4e1 * t705 * t188 * t830 + 0.2e1 * t801 * t94 * t693 * ZB - 0.8e1 * t857 * t843 - t718 * t726 + 0.4e1 * t773 * t761 * t13 - 0.4e1 * t705 * t775 * t555 + 0.2e1 * t746 * t232 * t233 - 0.16e2 * t878 * t830 * t365 - 0.2e1 * t746 * t619 * t24;
318   t916 = t110 * t28;
319   t945 = t28 * t693 * t45 * t17;
320   t948 = 0.32e2 * t877 * t733 * t798 + 0.2e1 * t718 * t36 * t916 * t24 - 0.4e1 * t705 * t626 + t718 * n * t916 * t693 - t869 * t792 * t110 - 0.4e1 * t773 * t761 * t775 + t718 * t719 + 0.2e1 * t746 * t232 * t24 - 0.16e2 * t849 * t365 * xc - t718 * t713 * t110 - 0.4e1 * t773 * t694 * t29 + 0.16e2 * t877 * t54 * t945;
321   t974 = t761 * t797;
322   t987 = 0.4e1 * t773 * t695 + 0.4e1 * t736 * t150 * t24 + 0.4e1 * t722 * t842 * t17 - 0.16e2 * t877 * t447 * t945 + 0.2e1 * t87 * t714 * t28 + t313 * t713 * t916 * t45 - 0.4e1 * t853 * t615 * t774 - 0.32e2 * t877 * t713 * xc * t324 + 0.16e2 * t760 * t974 + 0.4e1 * t736 * t94 * t24 * ZB + t869 * t792 * t916 - 0.8e1 * t691 * t5 * xc * t616;
323   t1021 = -t718 * t169 * t693 - 0.32e2 * t827 * t974 + 0.2e1 * t801 * t150 * t774 + 0.4e1 * t791 * t188 * t792 + 0.4e1 * t736 * t220 * t24 + 0.4e1 * t791 * t842 * t792 + 0.8e1 * t709 * t660 * t270 - t718 * t335 * t693 - 0.2e1 * t718 * t36 * t110 * t24 - 0.32e2 * t819 * t797 * t471 - t313 * t733 * t110 - 0.32e2 * t796 * t270 * t762;
324 
325   _PC2A = (t147 - 0.4e1 * t65 * t217 + t418 + 0.2e1 * t150 * t222 + t327 - 0.2e1 * t149 * t19 + 0.2e1 * t335 * ZB * t24 * ZA - 0.16e2 * t312 * t313 * t355 * t59 - 0.4e1 * t281 * ZB * ZA * t597 - 0.2e1 * t505 * t45 * t281 * t58 - 0.4e1 * t211 * t2 * t53 * t76 + 0.8e1 * t305 * t286 - 0.4e1 * t122 * t499 - 0.4e1 * t331 * t332 + 0.8e1 * t345 * t177 * t60 - 0.2e1 * t142 * t177 * t82 + 0.2e1 * t72 * t281 * t415 + 0.4e1 * t349 * t96 * t41 - 0.2e1 * t81 * t64 * t76 + 0.2e1 * t58 * t80 * t59 + 0.8e1 * t345 * t177 * t41 - 0.4e1 * t8 * t499 + t242 + 0.4e1 * t8 * t518 + t625 + t685 + 0.2e1 * t328 * t174 + 0.2e1 * t331 * t455 - 0.2e1 * t33 * t2 * t4 * ZA * t82 - 0.4e1 * t626 * t191 + 0.16e2 * t364 * t373 - 0.2e1 * t621 * t597 - 0.2e1 * t439 * t568 + t492 + t533 * t96 + t232 * t96 + 0.2e1 * t567 * t441 + t561) / (t740 + t789 + t834 + t874 + t911 + t948 + t987 + t1021);
326 
327   t1   = n * n;
328   t2   = t1 * n;
329   t3   = t2 * 0.3141592654e1;
330   t4   = t3 * xc;
331   t5   = ZB * ZB;
332   t7   = PetscExpReal(n * 0.3141592654e1);
333   t8   = t7 * t7;
334   t9   = t5 * t8;
335   t12  = PetscExpReal(xc * n * 0.3141592654e1);
336   t13  = t12 * t12;
337   t14  = t13 * t13;
338   t15  = t14 * t13;
339   t19  = nx * nx;
340   t21  = nx * 0.3141592654e1;
341   t22  = PetscSinReal(t21);
342   t23  = t19 * nx * t22;
343   t24  = t23 * 0.3141592654e1;
344   t25  = ZA * ZB;
345   t26  = t7 * t15;
346   t27  = t25 * t26;
347   t30  = t21 * xc;
348   t31  = PetscSinReal(t30);
349   t32  = t31 * nx;
350   t33  = t32 * n;
351   t34  = ZA * ZA;
352   t35  = t8 * t8;
353   t36  = t34 * t35;
354   t40  = t2 * t34;
355   t41  = 0.3141592654e1 * t8;
356   t42  = t41 * t15;
357   t45  = t1 * t5;
358   t46  = t14 * t14;
359   t49  = t19 * t5;
360   t51  = t19 * t46;
361   t53  = t19 * t34;
362   t55  = t8 * t7;
363   t56  = t13 * t55;
364   t57  = t25 * t56;
365   t60  = t2 * nx;
366   t61  = 0.3141592654e1 * 0.3141592654e1;
367   t63  = t60 * t31 * t61;
368   t64  = xc * xc;
369   t65  = ZA * t64;
370   t66  = ZB * t8;
371   t67  = t14 * t12;
372   t68  = t66 * t67;
373   t69  = t65 * t68;
374   t72  = -0.4e1 * t4 * t9 * t15 + 0.4e1 * t24 * t27 + 0.4e1 * t33 * t36 * t12 - 0.4e1 * t40 * t42 - t45 * t46 + t45 * t14 - t49 * t14 + t51 * t5 - t53 * t14 + 0.4e1 * t24 * t57 + 0.32e2 * t63 * t69;
375   t73  = t1 * nx;
376   t75  = t73 * t31 * 0.3141592654e1;
377   t76  = t8 * t67;
378   t77  = t25 * t76;
379   t80  = t1 * t1;
380   t81  = t80 * t34;
381   t83  = t61 * t14;
382   t87  = t1 * t19;
383   t88  = PetscCosReal(t30);
384   t90  = t87 * t88 * t61;
385   t91  = t5 * t64;
386   t92  = t13 * t12;
387   t93  = t8 * t92;
388   t94  = t91 * t93;
389   t100 = ZB * t64 * ZA * t8 * t92;
390   t103 = n * t19;
391   t105 = t103 * t88 * 0.3141592654e1;
392   t106 = ZA * xc;
393   t107 = ZB * t35;
394   t109 = t106 * t107 * t12;
395   t112 = t34 * xc;
396   t113 = t112 * t93;
397   t116 = t35 * t14;
398   t118 = t1 * ZA;
399   t119 = ZB * t14;
400   t122 = t1 * t46;
401   t125 = t19 * ZB;
402   t126 = t35 * ZA;
403   t127 = t125 * t126;
404   t129 = t1 * ZB;
405   t132 = -0.16e2 * t75 * t77 + 0.16e2 * t81 * t64 * t83 * t8 + 0.16e2 * t90 * t94 - 0.32e2 * t90 * t100 + 0.8e1 * t105 * t109 - 0.8e1 * t75 * t113 + t45 * t116 + 0.2e1 * t118 * t119 + 0.2e1 * t122 * t25 - 0.2e1 * t127 + 0.2e1 * t129 * t126;
406   t134 = t1 * t34;
407   t136 = t34 * t64;
408   t137 = t136 * t76;
409   t141 = t91 * t76;
410   t145 = t103 * t34;
411   t146 = 0.3141592654e1 * xc;
412   t147 = t8 * t13;
413   t153 = t14 * ZA;
414   t156 = xc * t5;
415   t157 = t156 * t93;
416   t160 = t103 * t5;
417   t162 = t146 * t8 * t15;
418   t166 = t34 * t7 * t15;
419   t169 = t134 * t116 - 0.16e2 * t63 * t137 - t49 * t116 - 0.16e2 * t63 * t141 - t53 * t116 + 0.4e1 * t145 * t146 * t147 - 0.2e1 * t51 * t25 - 0.2e1 * t125 * t153 - 0.16e2 * t75 * t157 + 0.4e1 * t160 * t162 - 0.4e1 * t24 * t166;
420   t170 = t106 * t68;
421   t177 = t35 * t92;
422   t178 = t112 * t177;
423   t181 = t156 * t76;
424   t186 = t35 * t12;
425   t187 = t112 * t186;
426   t193 = t5 * 0.3141592654e1;
427   t206 = t34 * t14;
428   t207 = t206 * t7;
429   t210 = -0.32e2 * t63 * t170 + 0.32e2 * t90 * t170 + 0.8e1 * t75 * t109 + 0.4e1 * t105 * t178 - 0.16e2 * t75 * t181 - 0.16e2 * t90 * t113 - 0.4e1 * t75 * t187 + 0.16e2 * t90 * t141 - 0.4e1 * t103 * t15 * t193 * xc + 0.16e2 * t73 * t22 * t34 * t146 * t26 + 0.4e1 * t32 * n * t34 * t67 + 0.4e1 * t24 * t207;
430   t217 = t106 * t66 * t92;
431   t226 = t88 * t19 * n;
432   t227 = 0.3141592654e1 * t34;
433   t229 = t227 * xc * t67;
434   t232 = t73 * t31;
435   t234 = t146 * t5 * t67;
436   t238 = t61 * ZB;
437   t239 = t14 * t8;
438   t240 = t238 * t239;
439   t243 = t136 * t93;
440   t246 = -0.8e1 * t33 * t25 * t186 + 0.32e2 * t90 * t217 - t45 * t35 + t53 * t35 - t134 * t35 - t134 * t46 + t134 * t14 - 0.4e1 * t226 * t229 + 0.4e1 * t232 * t234 + 0.32e2 * t87 * t65 * t240 + 0.16e2 * t63 * t243;
441   t247 = t14 * t92;
442   t249 = t227 * t247 * xc;
443   t254 = t73 * t22;
444   t259 = t60 * t22 * t61;
445   t260 = t112 * t26;
446   t264 = t146 * t247 * t5;
447   t268 = xc * t14;
448   t274 = t5 * t14;
449   t275 = t274 * t8;
450   t280 = n * nx;
451   t281 = t280 * t22;
452   t282 = t55 * t14;
453   t283 = t25 * t282;
454   t290 = ZA * t247 * xc * ZB;
455   t295 = t22 * nx * t1 * 0.3141592654e1;
456   t298 = -0.4e1 * t232 * t249 + 0.8e1 * t105 * t217 - 0.4e1 * t254 * t227 * t26 - 0.8e1 * t259 * t260 - 0.4e1 * t232 * t264 - 0.16e2 * t81 * t61 * t268 * t8 + 0.16e2 * t80 * t64 * t61 * t275 - 0.4e1 * t232 * t229 + 0.8e1 * t281 * t283 - 0.4e1 * t105 * t187 + 0.8e1 * t75 * t290 + 0.4e1 * t295 * t27;
457   t301 = t61 * t5;
458   t307 = t87 * t34;
459   t312 = t61 * xc;
460   t313 = t312 * t239;
461   t317 = t34 * t55 * t14;
462   t329 = ZB * t13 * t55;
463   t330 = t65 * t329;
464   t337 = -0.16e2 * t87 * t64 * t301 * t239 - 0.32e2 * t90 * t69 - 0.16e2 * t307 * t64 * t61 * t239 + 0.16e2 * t307 * t313 + 0.4e1 * t24 * t317 + t53 * t46 + t49 * t35 - 0.32e2 * t63 * t100 - 0.4e1 * t280 * t31 * t34 * t247 + 0.8e1 * t259 * t330 - 0.4e1 * t280 * t31 * t247 * t5;
465   t340 = t5 * t35;
466   t344 = t25 * t93;
467   t356 = t41 * t13;
468   t360 = t23 * n * t61;
469   t363 = t25 * t64 * t7 * t15;
470   t366 = t156 * t177;
471   t369 = t14 * t7;
472   t370 = t25 * t369;
473   t373 = t156 * t186;
474   t378 = 0.4e1 * t24 * t283 + 0.4e1 * t33 * t340 * t12 - 0.16e2 * t75 * t344 - 0.4e1 * t280 * t31 * t5 * t67 + 0.8e1 * t33 * t25 * t247 + 0.32e2 * t63 * t217 + 0.4e1 * t40 * t356 - 0.8e1 * t360 * t363 + 0.4e1 * t75 * t366 + 0.4e1 * t295 * t370 - 0.4e1 * t75 * t373 - 0.4e1 * t105 * t366;
475   t382 = t112 * t76;
476   t387 = t80 * t61;
477   t391 = t136 * t26;
478   t409 = 0.16e2 * t63 * t382 + 0.4e1 * t226 * t234 - 0.16e2 * t387 * xc * t275 + 0.8e1 * t259 * t391 - 0.16e2 * t105 * t344 + 0.4e1 * t226 * t264 - 0.8e1 * t105 * t170 + 0.16e2 * t232 * t193 * t76 + 0.8e1 * t360 * t330 - 0.8e1 * t105 * t290 + 0.16e2 * t90 * t243;
479   t423 = t153 * t8;
480   t426 = t34 * t13;
481   t427 = t426 * t55;
482   t430 = t34 * t8;
483   t437 = t80 * ZA;
484   t441 = 0.4e1 * t145 * t42 - 0.16e2 * t90 * t157 + 0.24e2 * t75 * t217 + 0.4e1 * t226 * t249 + 0.4e1 * t254 * t227 * t282 + 0.4e1 * t160 * t356 - 0.8e1 * t129 * t423 - 0.8e1 * t281 * t427 - 0.8e1 * t33 * t430 * t67 + 0.8e1 * t33 * t430 * t92 + 0.32e2 * t437 * ZB * t313;
485   t453 = t106 * ZB * t7 * t15;
486   t456 = t2 * t5;
487   t459 = t112 * t56;
488   t462 = t126 * t14;
489   t474 = t40 * 0.3141592654e1;
490   t475 = xc * t8;
491   t480 = t146 * t13 * t35;
492   t483 = -0.4e1 * t103 * xc * t193 * t147 + 0.16e2 * t87 * t61 * t156 * t239 + 0.8e1 * t259 * t453 - 0.4e1 * t456 * t356 + 0.8e1 * t259 * t459 - 0.2e1 * t125 * t462 - 0.8e1 * t281 * t207 + 0.16e2 * t295 * t459 - 0.8e1 * t60 * t22 * ZA * t312 * t329 + 0.4e1 * t474 * t475 * t15 + 0.4e1 * t160 * t480;
493   t497 = t136 * t56;
494   t504 = t9 * t13;
495   t509 = t475 * t13;
496   t512 = -0.8e1 * t105 * t113 - 0.4e1 * t254 * t227 * t56 + 0.8e1 * t281 * t57 + 0.4e1 * t295 * t283 + 0.2e1 * t129 * t462 + 0.4e1 * t24 * t370 - 0.8e1 * t360 * t497 - 0.4e1 * t24 * t427 - 0.4e1 * t145 * t162 + 0.4e1 * t4 * t504 - 0.8e1 * t281 * t370 - 0.4e1 * t474 * t509;
497   t528 = t5 * t13;
498   t529 = t528 * t35;
499   t532 = t106 * t329;
500   t542 = -0.16e2 * t295 * t453 - 0.32e2 * t437 * t64 * t240 + 0.8e1 * t281 * t317 + 0.24e2 * t75 * t170 - 0.4e1 * t75 * t178 + 0.8e1 * t360 * t453 - 0.4e1 * t4 * t529 - 0.16e2 * t295 * t532 - 0.8e1 * t33 * t344 - 0.16e2 * t90 * t181 + 0.4e1 * t33 * t340 * t92;
501   t557 = t146 * t15;
502   t562 = xc * t15;
503   t563 = t562 * t5;
504   t573 = 0.16e2 * t232 * t193 * t93 - 0.8e1 * t259 * t363 - 0.8e1 * t259 * t497 + 0.8e1 * t33 * t77 + 0.8e1 * t360 * t391 + 0.4e1 * t254 * t227 * t369 + 0.4e1 * t145 * t557 + 0.8e1 * t281 * t166 + 0.4e1 * t3 * t563 + 0.8e1 * t105 * t382 - 0.4e1 * t145 * t480 - 0.4e1 * t33 * t36 * t92;
505   t600 = 0.4e1 * t456 * t42 - 0.8e1 * t360 * t260 - 0.4e1 * t40 * t557 - 0.4e1 * t105 * t373 + 0.16e2 * t226 * t227 * t93 - 0.16e2 * t90 * t382 - 0.4e1 * t145 * t356 - 0.16e2 * t63 * t157 - 0.32e2 * t87 * t25 * t313 - 0.16e2 * t226 * t227 * t76 - 0.16e2 * t63 * t113;
506   t623 = xc * t13;
507   t627 = 0.8e1 * t125 * t423 - 0.8e1 * t360 * t532 + 0.16e2 * t90 * t137 - 0.4e1 * t160 * t42 + 0.16e2 * t63 * t94 + 0.16e2 * t63 * t181 - 0.8e1 * t281 * t27 - 0.8e1 * t75 * t382 + 0.8e1 * t360 * t459 + 0.4e1 * t295 * t57 + 0.16e2 * t105 * t77 + 0.4e1 * t474 * t623 * t35;
508   t632 = t61 * 0.3141592654e1;
509   t633 = t632 * t8;
510   t634 = t80 * n;
511   t638 = t632 * t634;
512   t639 = t638 * xc;
513   t642 = t61 * t34;
514   t643 = t122 * t19;
515   t649 = t61 * t61;
516   t650 = t649 * t1;
517   t652 = t19 * t19;
518   t653 = t14 * t652;
519   t654 = t653 * t9;
520   t657 = t14 * t1;
521   t658 = t657 * t19;
522   t665 = t632 * t34;
523   t666 = t665 * t2;
524   t667 = t8 * t19;
525   t668 = t667 * t623;
526   t674 = t665 * n;
527   t675 = t652 * xc;
528   t682 = 0.8e1 * t633 * t426 * t634 - 0.8e1 * t639 * t529 - 0.4e1 * t642 * t643 + 0.2e1 * t642 * t116 * t80 + 0.32e2 * t650 * t64 * t654 + 0.4e1 * t301 * t658 + 0.4e1 * t387 * t46 * ZA * ZB - 0.16e2 * t666 * t668 - 0.16e2 * t666 * t667 * t15 - 0.8e1 * t674 * t675 * t15 + 0.4e1 * t238 * t153 * t80;
529   t683 = t46 * t652;
530   t686 = t633 * t15;
531   t691 = t35 * t80;
532   t698 = t35 * t652;
533   t705 = t14 * t80;
534   t708 = t61 * t35;
535   t717 = -0.2e1 * t642 * t683 - 0.8e1 * t686 * t5 * t634 * xc - 0.2e1 * t301 * t691 + 0.8e1 * t638 * t563 - 0.2e1 * t642 * t691 - 0.2e1 * t642 * t698 - 0.2e1 * t301 * t698 - 0.2e1 * t301 * t683 + 0.2e1 * t642 * t705 + 0.2e1 * t708 * t274 * t80 + 0.2e1 * t301 * t653 - 0.2e1 * t642 * t80 * t46;
536   t727 = t61 * t46;
537   t737 = t649 * t34;
538   t738 = t737 * t1;
539   t739 = t8 * t652;
540   t740 = t739 * t268;
541   t746 = t61 * ZA;
542   t754 = t632 * n * xc;
543   t758 = 0.2e1 * t301 * t705 + 0.2e1 * t642 * t653 - 0.8e1 * t665 * xc * t634 * t15 - 0.2e1 * t727 * t5 * t80 - 0.32e2 * t650 * xc * t654 + 0.2e1 * t301 * t698 * t14 - 0.32e2 * t738 * t740 + 0.8e1 * t674 * t739 * t562 + 0.4e1 * t746 * t119 * t652 + 0.8e1 * t674 * t698 * t623 - 0.8e1 * t754 * t528 * t698;
544   t762 = t633 * t13;
545   t764 = t5 * n * t652;
546   t767 = t80 * t1;
547   t768 = t649 * t767;
548   t772 = t649 * ZA;
549   t773 = t772 * t129;
550   t777 = t35 * t1 * t19;
551   t780 = t632 * t5;
552   t781 = t780 * t15;
553   t786 = t698 * ZA;
554   t790 = t64 * t14;
555   t800 = t649 * t8;
556   t809 = 0.4e1 * t238 * t126 * t80 - 0.8e1 * t762 * t764 - 0.32e2 * t768 * xc * t275 + 0.64e2 * t773 * t740 - 0.4e1 * t301 * t777 - 0.8e1 * t781 * n * t8 * t675 + 0.4e1 * t238 * t786 + 0.32e2 * t768 * t34 * t790 * t8 - 0.8e1 * t633 * t528 * t634 + 0.8e1 * t754 * t528 * t739 + 0.128e3 * t800 * t119 * t80 * t19 * t106 + 0.8e1 * t674 * t739 * t13;
557   t812 = t649 * t80;
558   t817 = t83 * ZB;
559   t824 = t746 * ZB;
560   t828 = t800 * t14;
561   t855 = -0.64e2 * t812 * xc * t274 * t667 + 0.4e1 * t817 * t786 + 0.4e1 * t727 * ZA * t652 * ZB - 0.32e2 * t824 * t657 * t667 - 0.32e2 * t828 * t34 * t767 * xc - 0.8e1 * t633 * t15 * t34 * t634 - 0.8e1 * t674 * t739 * t15 + 0.32e2 * t768 * t64 * t275 + 0.4e1 * t708 * t14 * t307 + 0.2e1 * t708 * t206 * t652 + 0.8e1 * t632 * t35 * t13 * t34 * t634 * xc;
562   t858 = t35 * t19;
563   t873 = t2 * t8;
564   t878 = t61 * t1;
565   t901 = -0.16e2 * t632 * t2 * xc * t528 * t858 + 0.8e1 * t824 * t658 + 0.4e1 * t301 * t14 * t777 - 0.8e1 * t665 * t634 * t509 - 0.8e1 * t674 * t739 * t623 - 0.16e2 * t781 * t873 * t19 * xc + 0.8e1 * t878 * t14 * t127 + 0.8e1 * t878 * ZA * t51 * ZB + 0.8e1 * t686 * t764 + 0.8e1 * t665 * xc * t634 * t15 * t8 + 0.8e1 * t633 * t15 * t5 * t634 + 0.4e1 * t387 * t14 * t107 * ZA;
566   t903 = t739 * t790;
567   t923 = t737 * t80;
568   t924 = t667 * t790;
569   t927 = t780 * t2;
570   t937 = t15 * t19 * xc;
571   t943 = 0.32e2 * t738 * t903 + 0.16e2 * t781 * t873 * t19 + 0.8e1 * t754 * t15 * t652 * t5 + 0.16e2 * t666 * t858 * t623 + 0.64e2 * t828 * t25 * t767 * xc - 0.16e2 * t762 * t456 * t19 + 0.64e2 * t923 * t924 + 0.16e2 * t927 * t668 - 0.64e2 * t768 * ZA * t790 * t66 - 0.64e2 * t773 * t903 + 0.16e2 * t927 * t937 + 0.16e2 * t666 * t667 * t562;
572   t977 = 0.64e2 * t812 * t5 * t924 + 0.8e1 * t639 * t504 + 0.8e1 * t238 * t35 * t118 * t19 + 0.4e1 * t642 * t658 - 0.16e2 * t817 * t437 * t8 - 0.128e3 * t772 * ZB * t80 * t924 + 0.16e2 * t666 * t667 * t13 - 0.4e1 * t301 * t643 - 0.16e2 * t824 * t653 * t8 - 0.4e1 * t642 * t777 - 0.64e2 * t923 * t667 * t268 - 0.16e2 * t666 * t937;
573 
574   _PC3A = (t72 + t132 + t169 + t210 + t246 + t298 + t337 + t378 + t409 + t441 + t483 + t512 + t542 + t573 + t600 + t627) / (t682 + t717 + t758 + t809 + t855 + t901 + t943 + t977);
575 
576   _PC4A = 0;
577   t1    = nx * 0.3141592654e1;
578   t2    = t1 * xc;
579   t3    = PetscCosReal(t2);
580   t4    = nx * nx;
581   t6    = n * 0.3141592654e1;
582   t7    = t3 * t4 * t6;
583   t8    = ZA * ZB;
584   t9    = PetscExpReal(t6);
585   t10   = t9 * t9;
586   t11   = xc * n;
587   t13   = PetscExpReal(t11 * 0.3141592654e1);
588   t14   = t13 * t13;
589   t15   = t14 * t13;
590   t16   = t14 * t14;
591   t17   = t16 * t15;
592   t18   = t10 * t17;
593   t19   = t8 * t18;
594   t22   = PetscSinReal(t2);
595   t23   = nx * t22;
596   t24   = t23 * n;
597   t25   = ZB * ZB;
598   t30   = n * n;
599   t31   = t30 * n;
600   t32   = t31 * nx;
601   t33   = 0.3141592654e1 * 0.3141592654e1;
602   t35   = t32 * t22 * t33;
603   t36   = ZA * ZA;
604   t37   = t36 * xc;
605   t38   = t16 * t13;
606   t39   = t10 * t38;
607   t40   = t37 * t39;
608   t43   = PetscSinReal(t1);
609   t44   = nx * t43;
610   t45   = t30 * 0.3141592654e1;
611   t46   = t44 * t45;
612   t47   = ZA * xc;
613   t49   = ZB * t16 * t9;
614   t54   = t4 * nx * t43;
615   t55   = xc * xc;
616   t57   = t54 * t30 * t55;
617   t58   = t33 * 0.3141592654e1;
618   t59   = t58 * t25;
619   t60   = t16 * t9;
620   t61   = t59 * t60;
621   t64   = xc * t25;
622   t65   = t14 * t9;
623   t66   = t64 * t65;
624   t70   = t44 * t31 * t33;
625   t71   = t37 * t65;
626   t74   = t10 * t15;
627   t75   = t64 * t74;
628   t78   = t25 * t10;
629   t83   = t54 * n * t33;
630   t84   = t55 * t25;
631   t85   = t10 * t9;
632   t86   = t14 * t85;
633   t87   = t84 * t86;
634   t90   = t30 * t30;
635   t92   = t44 * t90 * t58;
636   t93   = t55 * xc;
637   t94   = t93 * t25;
638   t95   = t85 * t16;
639   t96   = t94 * t95;
640   t102  = t23 * t45;
641   t103  = t10 * t10;
642   t104  = ZB * t103;
643   t106  = t47 * t104 * t15;
644   t111  = t54 * 0.3141592654e1;
645   t112  = t25 * t85;
646   t113  = t112 * t16;
647   t115  = t8 * t39;
648   t118  = t16 * t14;
649   t119  = t85 * t118;
650   t120  = t37 * t119;
651   t123  = t16 * t16;
652   t124  = t36 * t123;
653   t125  = t124 * t9;
654   t127 = -0.8e1 * t7 * t19 + 0.2e1 * t24 * t25 * t13 * t10 - 0.16e2 * t35 * t40 - 0.16e2 * t46 * t47 * t49 - 0.8e1 * t57 * t61 + 0.4e1 * t46 * t66 + 0.2e1 * t70 * t71 - 0.16e2 * t35 * t75 + 0.6e1 * t24 * t78 * t38 - 0.2e1 * t83 * t87 - 0.8e1 * t92 * t96 - 0.8e1 * t46 * t37 * t95 - 0.12e2 * t102 * t106 + 0.2e1 * t83 * t71 + t111 * t113 + 0.8e1 * t7 * t115 + 0.2e1 * t83 * t120 + t111 * t125;
655   t128 = t37 * t74;
656   t131 = t44 * n;
657   t133 = t25 * t9 * t118;
658   t136 = t36 * t14;
659   t137 = t136 * t9;
660   t140 = t30 * t4;
661   t142 = t140 * t3 * t33;
662   t143 = t64 * t39;
663   t147 = t30 * nx * t43;
664   t148 = 0.3141592654e1 * t36;
665   t149 = t9 * t118;
666   t153 = t44 * t31 * ZA;
667   t154 = t33 * xc;
668   t155 = t154 * t49;
669   t160 = ZA * t17 * xc * ZB;
670   t163 = t103 * t13;
671   t164 = t64 * t163;
672   t170 = t44 * t90 * t55;
673   t171 = t58 * ZB;
674   t172 = ZA * t16;
675   t174 = t171 * t172 * t9;
676   t177 = t36 * t55;
677   t178 = t177 * t149;
678   t181 = t54 * t11;
679   t182 = t33 * t25;
680   t186 = t25 * t14;
681   t187 = t186 * t9;
682   t193 = t186 * t85;
683   t198 = ZB * t55;
684   t199 = ZA * t103;
685   t201 = t198 * t199 * t15;
686   t204 = 0.2e1 * t7 * t128 - 0.2e1 * t131 * t133 - 0.2e1 * t131 * t137 + 0.16e2 * t142 * t143 - t147 * t148 * t149 + 0.8e1 * t153 * t155 - 0.4e1 * t7 * t160 + 0.2e1 * t7 * t164 + 0.10e2 * t102 * t40 + 0.16e2 * t170 * t174 + 0.2e1 * t83 * t178 - 0.2e1 * t181 * t182 * t65 - t111 * t187 - 0.2e1 * t70 * t87 + 0.4e1 * t102 * t160 - 0.2e1 * t131 * t193 - 0.16e2 * t142 * t75 + 0.16e2 * t35 * t201;
687   t210 = t32 * t22;
688   t211 = t33 * t55;
689   t212 = t25 * t38;
690   t213 = t211 * t212;
691   t216 = n * nx;
692   t217 = t22 * t25;
693   t222 = ZB * t85 * t16;
694   t226 = t23 * t30;
695   t227 = t13 * t10;
696   t228 = t148 * t227;
697   t233 = t37 * t163;
698   t237 = n * t4 * t3;
699   t238 = t148 * t74;
700   t241 = t64 * t86;
701   t245 = t148 * xc * t15;
702   t248 = t112 * t118;
703   t250 = t22 * t36;
704   t256 = 0.3141592654e1 * t25;
705   t257 = t256 * t39;
706   t262 = t38 * t103;
707   t263 = t37 * t262;
708   t267 = t148 * t17 * xc;
709   t270 = -0.6e1 * t7 * t143 - 0.4e1 * t24 * t19 - 0.8e1 * t210 * t213 - 0.2e1 * t216 * t217 * t15 - 0.32e2 * t153 * t211 * t222 + 0.4e1 * t226 * t228 + 0.16e2 * t142 * t201 + 0.2e1 * t7 * t233 - 0.4e1 * t237 * t238 - 0.2e1 * t83 * t241 - 0.2e1 * t237 * t245 + t111 * t248 + 0.2e1 * t216 * t250 * t15 - 0.2e1 * t131 * t125 - 0.4e1 * t226 * t257 + t147 * t148 * t95 - 0.2e1 * t102 * t263 + 0.2e1 * t237 * t267;
710   t273 = t37 * t149;
711   t277 = t47 * t104 * t13;
712   t285 = t31 * t36;
713   t286 = t44 * t285;
714   t291 = t25 * t123 * t9;
715   t304 = 0.3141592654e1 * xc;
716   t305 = t304 * t212;
717   t312 = t256 * t18;
718   t315 = t8 * t60;
719   t319 = t54 * t30 * t58;
720   t323 = t90 * t36;
721   t324 = t44 * t323;
722   t325 = t55 * t58;
723   t326 = t325 * t60;
724   t329 = 0.2e1 * t102 * t164 + 0.2e1 * t83 * t273 - 0.4e1 * t102 * t277 - 0.2e1 * t7 * t263 + 0.4e1 * t24 * t8 * t17 - 0.4e1 * t286 * t154 * t60 - 0.2e1 * t131 * t291 - t147 * t148 * t119 + 0.2e1 * t24 * t78 * t17 + 0.2e1 * t54 * t85 * 0.3141592654e1 * ZA * ZB - 0.4e1 * t226 * t305 - 0.2e1 * t70 * t66 + t147 * t256 * t95 + 0.4e1 * t237 * t312 + 0.2e1 * t111 * t315 - 0.8e1 * t319 * t96 - t111 * t193 - 0.8e1 * t324 * t326;
725   t332 = t8 * t95;
726   t335 = t136 * t85;
727   t337 = t256 * t227;
728   t340 = t177 * t119;
729   t346 = t37 * t86;
730   t351 = t103 * t15;
731   t352 = t177 * t351;
732   t355 = t64 * t119;
733   t358 = t8 * t227;
734   t361 = t85 * 0.3141592654e1;
735   t365 = t84 * t39;
736   t372 = ZB * t10;
737   t373 = t372 * t38;
738   t374 = t47 * t373;
739   t379 = t177 * t39;
740   t384 = -0.2e1 * t46 * t332 + t111 * t335 + 0.4e1 * t237 * t337 - 0.2e1 * t83 * t340 + 0.16e2 * t286 * t211 * t95 + 0.2e1 * t70 * t346 - 0.8e1 * t170 * t61 - 0.8e1 * t142 * t352 - 0.2e1 * t83 * t355 - 0.4e1 * t24 * t358 + 0.2e1 * t147 * t361 * t8 + 0.8e1 * t35 * t365 - 0.2e1 * t226 * t267 + 0.8e1 * t102 * t115 - 0.12e2 * t102 * t374 + 0.16e2 * t142 * t40 - 0.8e1 * t142 * t379 + 0.4e1 * t237 * t228;
741   t386 = t54 * t30 * t93;
742   t387 = ZA * t85;
743   t389 = t171 * t387 * t16;
744   t394 = t64 * t60;
745   t398 = t304 * t25 * t15;
746   t401 = t361 * t25;
747   t405 = t84 * t65;
748   t410 = t148 * t18;
749   t414 = t25 * t16 * t9;
750   t417 = t84 * t74;
751   t422 = t177 * t86;
752   t428 = ZB * t38;
753   t429 = t47 * t428;
754   t432 = t148 * t39;
755   t439 = 0.16e2 * t386 * t389 - 0.16e2 * t386 * t174 + 0.8e1 * t46 * t394 + 0.2e1 * t237 * t398 - t147 * t401 + 0.4e1 * t7 * t374 + 0.2e1 * t83 * t405 - 0.4e1 * t46 * t241 - 0.4e1 * t226 * t410 + 0.2e1 * t131 * t414 + 0.8e1 * t35 * t417 - 0.8e1 * t142 * t365 + 0.2e1 * t70 * t422 - 0.4e1 * t181 * t182 * t60 + 0.12e2 * t102 * t429 - 0.4e1 * t226 * t432 + 0.32e2 * t35 * t374 - 0.4e1 * t7 * t106;
756   t442 = t36 * t9 * t118;
757   t444 = t123 * t9;
758   t445 = t8 * t444;
759   t448 = t361 * t36;
760   t451 = t47 * t372 * t17;
761   t454 = t94 * t60;
762   t457 = t25 * t103;
763   t465 = t47 * t372 * t15;
764   t468 = t36 * t85;
765   t469 = t468 * t16;
766   t474 = t43 * t85;
767   t478 = t8 * t74;
768   t484 = t256 * t74;
769   t489 = t198 * ZA * t10 * t15;
770   t501 = -t111 * t442 + 0.4e1 * t131 * t445 - t147 * t448 + 0.4e1 * t7 * t451 + 0.8e1 * t92 * t454 - 0.2e1 * t24 * t457 * t13 - 0.2e1 * t286 * t211 * t65 + 0.4e1 * t7 * t465 + t111 * t469 - 0.2e1 * t216 * t250 * t17 - 0.2e1 * t216 * t474 * t25 - 0.4e1 * t24 * t478 + 0.4e1 * t24 * t8 * t38 + 0.4e1 * t226 * t484 - 0.16e2 * t142 * t489 - 0.2e1 * t24 * t212 * t103 - 0.2e1 * t216 * t22 * t17 * t25 + 0.2e1 * t70 * t120;
771   t504 = t33 * t36 * t55 * t38;
772   t507 = t37 * t18;
773   t512 = t47 * ZB * t13 * t10;
774   t518 = t59 * t95;
775   t530 = t84 * t351;
776   t534 = t37 * t227;
777   t549 = -0.8e1 * t210 * t504 + 0.2e1 * t102 * t507 + 0.4e1 * t7 * t512 + t111 * t133 - 0.16e2 * t35 * t489 + 0.8e1 * t170 * t518 + 0.2e1 * t24 * t36 * t13 * t10 + 0.4e1 * t131 * t387 * ZB + 0.12e2 * t102 * t465 - 0.8e1 * t142 * t530 + t111 * t291 - 0.2e1 * t102 * t534 - 0.4e1 * t70 * t394 - 0.10e2 * t102 * t128 + 0.4e1 * t237 * t305 + 0.8e1 * t102 * t19 + 0.2e1 * t83 * t346 - 0.16e2 * t35 * t128;
778   t557 = t468 * t118;
779   t562 = t93 * t58;
780   t563 = t562 * t60;
781   t567 = t44 * t90 * t93;
782   t575 = ZA * t55;
783   t576 = t575 * t428;
784   t583 = t37 * t60;
785   t590 = t140 * t3;
786   t601 = -0.2e1 * t226 * t398 - 0.2e1 * t70 * t340 - 0.2e1 * t131 * t557 - 0.4e1 * t24 * t115 + 0.8e1 * t324 * t563 + 0.16e2 * t567 * t389 + 0.16e2 * t70 * t84 * t95 + 0.2e1 * t70 * t178 - 0.16e2 * t142 * t576 - 0.4e1 * t237 * t257 - 0.4e1 * t226 * t312 + 0.8e1 * t46 * t583 + 0.2e1 * t24 * t36 * t38 * t103 + 0.8e1 * t590 * t213 + 0.2e1 * t102 * t143 - 0.16e2 * t35 * t143 + 0.2e1 * t131 * t248 + 0.4e1 * t46 * t346;
787   t604 = n * t36;
788   t606 = t154 * t95;
789   t625 = t36 * t103;
790   t640 = t30 * t36;
791   t641 = t54 * t640;
792   t642 = t325 * t95;
793   t647 = -0.4e1 * t131 * t315 - 0.4e1 * t54 * t604 * t606 - t147 * t148 * t60 + 0.16e2 * t35 * t576 - 0.8e1 * t102 * t478 + 0.32e2 * t142 * t465 - 0.4e1 * t237 * t484 - 0.2e1 * t70 * t355 + 0.2e1 * t70 * t273 + 0.2e1 * t102 * t233 - 0.2e1 * t24 * t625 * t13 - 0.8e1 * t7 * t358 - 0.2e1 * t111 * t445 - 0.4e1 * t7 * t429 + 0.16e2 * t46 * t47 * t222 + 0.2e1 * t131 * t113 + 0.8e1 * t641 * t642 - 0.2e1 * t7 * t534;
794   t652 = t36 * t16;
795   t653 = t652 * t9;
796   t655 = t64 * t227;
797   t658 = t182 * t95;
798   t663 = t562 * t95;
799   t684 = t64 * t351;
800   t689 = t36 * t10;
801   t695 = t154 * t222;
802   t698 = -0.4e1 * t216 * t217 * t38 - t111 * t653 - 0.2e1 * t7 * t655 - 0.4e1 * t181 * t658 + 0.2e1 * t131 * t469 - 0.8e1 * t641 * t663 - 0.4e1 * t83 * t583 - 0.2e1 * t83 * t177 * t65 - 0.4e1 * t24 * t457 * t15 + 0.16e2 * t70 * t84 * t60 + 0.8e1 * t57 * t518 - 0.32e2 * t142 * t374 + 0.4e1 * t24 * t8 * t351 + 0.4e1 * t102 * t684 - t147 * t256 * t86 - 0.2e1 * t24 * t689 * t15 - 0.2e1 * t70 * t241 + 0.8e1 * t153 * t695;
803   t711 = t575 * t373;
804   t717 = t304 * t17 * t25;
805   t736 = t177 * t74;
806   t739 = 0.2e1 * t226 * t245 - 0.8e1 * t102 * t358 - 0.16e2 * t57 * t389 - 0.2e1 * t102 * t655 + 0.8e1 * t590 * t504 - 0.8e1 * t641 * t326 - 0.16e2 * t35 * t711 - t111 * t557 + t111 * t137 - 0.2e1 * t226 * t717 + 0.8e1 * t102 * t37 * t351 + 0.2e1 * t131 * t335 - 0.4e1 * t131 * t332 - 0.2e1 * t216 * t474 * t36 - 0.2e1 * t111 * t332 + 0.16e2 * t142 * t711 - t147 * t256 * t60 + 0.8e1 * t142 * t736;
807   t750 = t64 * t262;
808   t763 = t44 * t640;
809   t770 = t84 * t119;
810   t782 = 0.4e1 * t102 * t512 + 0.8e1 * t142 * t417 + 0.8e1 * t641 * t563 - 0.2e1 * t7 * t507 + 0.2e1 * t7 * t750 - 0.8e1 * t35 * t352 + 0.4e1 * t237 * t410 + 0.4e1 * t7 * t684 - 0.2e1 * t46 * t445 + t147 * t148 * t65 + 0.4e1 * t763 * t304 * t119 + 0.16e2 * t70 * t177 * t60 + 0.2e1 * t70 * t770 - t111 * t414 - 0.16e2 * t567 * t174 - 0.4e1 * t46 * t71 - 0.4e1 * t46 * t355 - 0.4e1 * t7 * t277;
811   t797 = t64 * t149;
812   t821 = -t54 * t448 + 0.2e1 * t131 * t442 + 0.8e1 * t7 * t478 + 0.8e1 * t35 * t379 - 0.2e1 * t181 * t182 * t149 + 0.2e1 * t70 * t405 + 0.2e1 * t83 * t770 - 0.2e1 * t70 * t797 - 0.6e1 * t7 * t75 - 0.4e1 * t286 * t606 - 0.4e1 * t237 * t432 + t147 * t256 * t149 - 0.4e1 * t763 * t304 * t149 - 0.2e1 * t102 * t75 + 0.2e1 * t237 * t717 + 0.8e1 * t324 * t642 - 0.16e2 * t170 * t389 + 0.2e1 * t83 * t422;
813   t827 = t84 * t149;
814   t846 = t54 * n * ZA;
815   t854 = t64 * t18;
816   t867 = -0.16e2 * t142 * t128 + 0.32e2 * t35 * t465 - 0.2e1 * t83 * t827 + 0.2e1 * t46 * t315 + t147 * t148 * t86 - 0.4e1 * t102 * t451 - 0.8e1 * t226 * t148 * xc * t38 - 0.2e1 * t24 * t689 * t38 + 0.2e1 * t131 * t187 + 0.8e1 * t846 * t155 + 0.8e1 * t35 * t736 + 0.2e1 * t24 * t689 * t17 - 0.2e1 * t7 * t854 + t147 * t256 * t119 + 0.2e1 * t102 * t854 - 0.8e1 * t35 * t530 + 0.4e1 * t46 * t797 + 0.2e1 * t102 * t750;
817   t909 = -0.8e1 * t324 * t663 + t147 * t256 * t444 - t147 * t256 * t65 + 0.4e1 * t226 * t238 + 0.2e1 * t7 * t40 - t54 * t401 + 0.16e2 * t57 * t174 + 0.4e1 * t226 * t337 + 0.4e1 * t24 * t8 * t163 + 0.8e1 * t846 * t695 + 0.8e1 * t319 * t454 + 0.2e1 * t131 * t653 - 0.8e1 * t46 * t64 * t95 + 0.6e1 * t24 * t78 * t15 - 0.4e1 * t44 * t31 * xc * t658 - 0.32e2 * t153 * t211 * t49 - 0.2e1 * t70 * t827 + t147 * t148 * t444;
818   t914 = t25 * ZB;
819   t915 = t33 * t914;
820   t919 = t4 * t4;
821   t920 = t16 * t919;
822   t929 = t123 * t90;
823   t932 = t919 * t103;
824   t935 = t33 * ZB;
825   t939 = t652 * t919;
826   t942 = t16 * t30;
827   t943 = t942 * t4;
828   t949 = t103 * t16;
829   t950 = t949 * t90;
830   t953 = -0.2e1 * t915 * t103 * t90 + 0.2e1 * t915 * t920 - 0.2e1 * t915 * t123 * t919 + 0.2e1 * t915 * t16 * t90 - 0.2e1 * t915 * t929 - 0.2e1 * t915 * t932 - 0.2e1 * t935 * t323 * t123 + 0.2e1 * t935 * t939 + 0.4e1 * t915 * t943 + 0.4e1 * t182 * t172 * t90 + 0.2e1 * t915 * t950;
831   t954  = t171 * t36;
832   t955  = t90 * n;
833   t956  = xc * t955;
834   t957  = t118 * t10;
835   t964  = t33 * t33;
836   t965  = t964 * ZB;
837   t966  = t965 * t640;
838   t967  = t10 * t919;
839   t968  = t55 * t16;
840   t969  = t967 * t968;
841   t972  = t935 * t36;
842   t974  = t103 * t30 * t4;
843   t977  = xc * t16;
844   t978  = t967 * t977;
845   t981  = t90 * t30;
846   t983  = t16 * t10;
847   t987  = t182 * ZA;
848   t988  = t4 * t10;
849   t992  = t171 * t604;
850   t993  = xc * t14;
851   t994  = t932 * t993;
852   t997  = t182 * t30;
853   t1005 = t171 * t285;
854   t1006 = t988 * t993;
855   t1009 = t58 * t914;
856   t1010 = t1009 * t31;
857   t1013 = 0.8e1 * t954 * t956 * t957 + 0.2e1 * t915 * t932 * t16 + 0.32e2 * t966 * t969 - 0.4e1 * t972 * t974 - 0.32e2 * t966 * t978 + 0.32e2 * t965 * t981 * t177 * t983 - 0.32e2 * t987 * t942 * t988 + 0.8e1 * t992 * t994 + 0.8e1 * t997 * t949 * ZA * t4 - 0.2e1 * t935 * t124 * t919 - 0.16e2 * t1005 * t1006 + 0.16e2 * t1010 * t1006;
858   t1015 = t964 * t25;
859   t1016 = ZA * t30;
860   t1017 = t1015 * t1016;
861   t1020 = t967 * t993;
862   t1031 = t1009 * t118;
863   t1032 = t31 * t10;
864   t1040 = t964 * t914;
865   t1041 = t1040 * t90;
866   t1044 = t55 * t10 * t4 * t16;
867   t1047 = t1040 * t30;
868   t1050 = t123 * ZA;
869   t1054 = t977 * t988;
870   t1057 = 0.64e2 * t1017 * t978 - 0.8e1 * t992 * t1020 + 0.2e1 * t972 * t950 + 0.4e1 * t182 * t929 * ZA + 0.4e1 * t182 * t199 * t90 - 0.16e2 * t1031 * t1032 * t4 * xc + 0.4e1 * t182 * t172 * t919 + 0.64e2 * t1041 * t1044 + 0.32e2 * t1047 * t969 + 0.4e1 * t182 * t1050 * t919 - 0.64e2 * t1041 * t1054;
871   t1058 = t1009 * n;
872   t1063 = t932 * ZA;
873   t1069 = t123 * t30 * t4;
874   t1080 = t993 * t103 * t4;
875   t1088 = t935 * t103;
876   t1094 = -0.8e1 * t1058 * t994 - 0.32e2 * t1047 * t978 + 0.4e1 * t182 * t1063 - 0.4e1 * t915 * t974 - 0.4e1 * t915 * t1069 - 0.2e1 * t935 * t625 * t90 - 0.8e1 * t1009 * t10 * t14 * t955 - 0.16e2 * t1010 * t1080 - 0.2e1 * t935 * t625 * t919 - 0.64e2 * t1017 * t969 + 0.2e1 * t1088 * t939 + 0.8e1 * t1009 * t957 * t955;
877   t1113 = t955 * t118 * xc;
878   t1120 = t4 * t118;
879   t1125 = t981 * xc;
880   t1133 = n * t10;
881   t1140 = -0.8e1 * t954 * t955 * t10 * t993 + 0.2e1 * t935 * t652 * t90 - 0.64e2 * t1015 * t981 * t575 * t983 + 0.8e1 * t182 * t103 * t1016 * t4 + 0.8e1 * t1009 * t1113 + 0.16e2 * t954 * t1032 * t4 * t14 - 0.16e2 * t954 * t1032 * t1120 + 0.64e2 * t1015 * t10 * t172 * t1125 + 0.8e1 * t171 * t103 * t136 * t956 - 0.8e1 * t1031 * t1133 * t919 * xc + 0.8e1 * t1058 * t1020;
882   t1153 = xc * t118;
883   t1165 = t182 * t16;
884   t1170 = t171 * t10;
885   t1178 = ZA * t90;
886   t1182 = 0.4e1 * t1088 * t652 * t140 + 0.8e1 * t954 * t1133 * t919 * t14 + 0.4e1 * t972 * t943 - 0.4e1 * t972 * t1069 - 0.16e2 * t954 * t31 * t4 * t1153 - 0.8e1 * t954 * n * t919 * t1153 - 0.8e1 * t954 * t1133 * t919 * t118 + 0.4e1 * t1165 * t1063 + 0.16e2 * t1005 * t1080 - 0.8e1 * t1170 * t118 * t36 * t955 - 0.16e2 * t987 * t920 * t10 - 0.16e2 * t1165 * t1178 * t10;
887   t1195 = t1040 * t981;
888   t1199 = t1009 * t955;
889   t1203 = t1009 * t10;
890   t1211 = t965 * t323;
891   t1225 = -0.32e2 * t965 * t10 * t652 * t1125 + 0.4e1 * t915 * t16 * t974 + 0.4e1 * t182 * t90 * t949 * ZA + 0.32e2 * t1195 * t968 * t10 - 0.8e1 * t1199 * t993 * t103 + 0.8e1 * t1203 * t118 * n * t919 + 0.8e1 * t1170 * t136 * t955 + 0.64e2 * t1211 * t1044 + 0.16e2 * t1031 * t1032 * t4 + 0.8e1 * t987 * t943 + 0.8e1 * t1199 * t993 * t10 + 0.8e1 * t997 * t1050 * t4;
892   t1263 = -0.128e3 * t1015 * t1178 * t1044 + 0.16e2 * t1005 * t988 * t1153 + 0.8e1 * t1058 * t1153 * t919 + 0.16e2 * t1010 * t1120 * xc - 0.8e1 * t954 * t1113 - 0.8e1 * t1203 * t14 * n * t919 - 0.16e2 * t1203 * t14 * t31 * t4 - 0.8e1 * t1203 * t1113 - 0.32e2 * t1195 * t977 * t10 - 0.64e2 * t1211 * t1054 + 0.8e1 * t992 * t967 * t1153 + 0.128e3 * t1015 * t983 * t90 * t4 * t47;
893 
894   _PC1B = (t127 + t204 + t270 + t329 + t384 + t439 + t501 + t549 + t601 + t647 + t698 + t739 + t782 + t821 + t867 + t909) / (t953 + t1013 + t1057 + t1094 + t1140 + t1182 + t1225 + t1263);
895 
896   t1   = n * n;
897   t2   = t1 * n;
898   t3   = nx * t2;
899   t4   = 0.3141592654e1 * ZA;
900   t5   = t3 * t4;
901   t6   = nx * 0.3141592654e1;
902   t7   = t6 * xc;
903   t8   = PetscSinReal(t7);
904   t9   = t8 * ZB;
905   t10  = n * 0.3141592654e1;
906   t11  = PetscExpReal(t10);
907   t12  = t11 * t11;
908   t15  = PetscExpReal(xc * n * 0.3141592654e1);
909   t16  = t15 * t15;
910   t17  = t16 * t16;
911   t18  = t17 * t15;
912   t19  = t12 * t18;
913   t23  = t1 * t1;
914   t24  = nx * t23;
915   t25  = ZB * ZB;
916   t27  = t18 * t8;
917   t28  = 0.3141592654e1 * 0.3141592654e1;
918   t29  = xc * xc;
919   t30  = t28 * t29;
920   t34  = t1 * xc;
921   t35  = 0.3141592654e1 * ZB;
922   t36  = t34 * t35;
923   t37  = PetscCosReal(t7);
924   t38  = ZA * t37;
925   t39  = nx * nx;
926   t40  = t39 * t12;
927   t41  = t16 * t15;
928   t43  = t38 * t40 * t41;
929   t46  = t25 * n;
930   t47  = t46 * 0.3141592654e1;
931   t48  = t39 * nx;
932   t49  = PetscSinReal(t6);
933   t50  = t48 * t49;
934   t51  = t12 * t11;
935   t52  = t51 * t17;
936   t53  = t50 * t52;
937   t56  = t34 * 0.3141592654e1 * t25;
938   t57  = t37 * t39;
939   t58  = t17 * t41;
940   t59  = t12 * t58;
941   t60  = t57 * t59;
942   t63  = t25 * t18;
943   t64  = t57 * n;
944   t67  = ZA * ZA;
945   t68  = t67 * n;
946   t69  = 0.3141592654e1 * t48;
947   t70  = t68 * t69;
948   t71  = t49 * xc;
949   t72  = t17 * t16;
950   t73  = t11 * t72;
951   t74  = t71 * t73;
952   t77  = t1 * t67;
953   t78  = t77 * 0.3141592654e1;
954   t81  = nx * t25;
955   t82  = t81 * t49;
956   t83  = t17 * t17;
957   t85  = t1 * t83 * t11;
958   t87  = nx * ZB;
959   t88  = t8 * t2;
960   t89  = t87 * t88;
961   t90  = 0.3141592654e1 * xc;
962   t91  = t12 * t12;
963   t92  = ZA * t91;
964   t97  = ZB * ZA;
965   t98  = t97 * t37;
966   t99  = t39 * n;
967   t100 = t12 * t41;
968   t104 = 0.8e1 * t5 * t9 * t19 + 0.8e1 * t24 * t25 * t27 * t30 + 0.12e2 * t36 * t43 - t47 * t53 - 0.2e1 * t56 * t60 - 0.4e1 * t63 * t64 + 0.6e1 * t70 * t74 + 0.4e1 * t78 * t60 - t82 * t85 + 0.4e1 * t89 * t90 * t92 * t41 + 0.4e1 * t98 * t99 * t100;
969   t105 = t67 * t48;
970   t106 = t49 * t51;
971   t107 = t106 * t72;
972   t109 = t1 * 0.3141592654e1;
973   t110 = t109 * xc;
974   t115 = nx * t67;
975   t116 = t115 * t49;
976   t117 = t1 * t16;
977   t118 = t117 * t11;
978   t120 = t2 * t25;
979   t121 = t28 * 0.3141592654e1;
980   t122 = t121 * t29;
981   t123 = t120 * t122;
982   t129 = t1 * ZB;
983   t130 = t129 * t4;
984   t131 = t57 * t100;
985   t134 = t12 * t16;
986   t136 = t109 * t39;
987   t139 = ZB * t18;
988   t141 = t39 * t1;
989   t142 = t141 * t90;
990   t145 = t77 * t90;
991   t146 = t91 * t41;
992   t147 = t57 * t146;
993   t151 = t25 * t39 * t1;
994   t152 = t72 * t12;
995   t156 = t49 * t2;
996   t158 = t83 * t11;
997   t162 = -t105 * t107 + 0.8e1 * t110 * t72 * t25 * t39 - t116 * t118 + 0.8e1 * t123 * t53 + 0.8e1 * t5 * t9 * t59 - 0.8e1 * t130 * t131 - 0.8e1 * t134 * t25 * t136 - 0.12e2 * t139 * t38 * t142 - 0.8e1 * t145 * t147 - 0.8e1 * t151 * t90 * t152 - 0.2e1 * t87 * t156 * t4 * t158;
998   t164 = t115 * t88;
999   t165 = t90 * t19;
1000   t168 = t25 * t48;
1001   t169 = t49 * t16;
1002   t170 = t169 * t11;
1003   t174 = ZA * n * t69;
1004   t175 = ZB * t51;
1005   t176 = t175 * t17;
1006   t177 = t71 * t176;
1007   t180 = t1 * t29;
1008   t181 = t28 * t25;
1009   t182 = t180 * t181;
1010   t183 = t50 * t73;
1011   t186 = ZA * t1;
1012   t187 = t28 * t48;
1013   t188 = t186 * t187;
1014   t189 = ZB * t17;
1015   t190 = t189 * t11;
1016   t191 = t71 * t190;
1017   t194 = t50 * t158;
1018   t196 = t115 * t156;
1019   t197 = t90 * t73;
1020   t201 = t49 * t17 * t11;
1021   t204 = t88 * t90;
1022   t207 = t68 * 0.3141592654e1;
1023   t208 = t17 * t11;
1024   t209 = t50 * t208;
1025   t211 = -0.2e1 * t164 * t165 - t168 * t170 + t168 * t107 + 0.8e1 * t174 * t177 + 0.2e1 * t182 * t183 + 0.8e1 * t188 * t191 + t47 * t194 - 0.6e1 * t196 * t197 - t168 * t201 - 0.4e1 * t81 * t18 * t204 - t207 * t209;
1026   t212 = t2 * 0.3141592654e1;
1027   t213 = t212 * t52;
1028   t215 = t81 * t8;
1029   t216 = t212 * t59;
1030   t219 = t3 * t90;
1031   t220 = t25 * t8;
1032   t221 = t18 * t91;
1033   t225 = t71 * t52;
1034   t231 = t16 * t51;
1035   t232 = t50 * t231;
1036   t237 = ZA * t12;
1037   t243 = t67 * t28;
1038   t244 = t24 * t243;
1039   t245 = t71 * t231;
1040   t249 = -t116 * t213 - 0.4e1 * t215 * t216 + 0.2e1 * t219 * t220 * t221 - 0.4e1 * t70 * t225 + 0.4e1 * t98 * t99 * t146 + t47 * t232 - 0.2e1 * t145 * t57 * t221 + 0.4e1 * t89 * t90 * t237 * t41 - t105 * t201 - 0.6e1 * t244 * t245 + t105 * t170;
1041   t252 = t25 * t37;
1042   t253 = t252 * t39;
1043   t255 = n * t15 * t12;
1044   t258 = t2 * t29;
1045   t259 = ZB * t28;
1046   t260 = t258 * t259;
1047   t263 = t106 * t17;
1048   t265 = xc * t25;
1049   t269 = t25 * t49;
1050   t270 = t269 * t52;
1051   t273 = t1 * t25;
1052   t274 = t273 * 0.3141592654e1;
1053   t275 = t57 * t19;
1054   t278 = t24 * t30;
1055   t288 = t1 * t11 * t72;
1056   t290 = t212 * t208;
1057   t292 = t2 * xc;
1058   t296 = 0.2e1 * t253 * t255 + 0.16e2 * t260 * t43 + t105 * t263 - 0.4e1 * t10 * t265 * t53 + 0.4e1 * t219 * t270 - 0.12e2 * t274 * t275 + 0.8e1 * t278 * t270 - 0.2e1 * ZB * n * t69 * t49 * ZA * t158 - t82 * t288 - t116 * t290 + 0.16e2 * t292 * t243 * t275;
1059   t301 = t50 * t176;
1060   t304 = t51 * t72;
1061   t305 = t71 * t304;
1062   t308 = t25 * t41;
1063   t311 = ZA * t48;
1064   t312 = t311 * t49;
1065   t317 = t91 * t15;
1066   t318 = t57 * t317;
1067   t321 = t81 * t88;
1068   t322 = t90 * t59;
1069   t325 = t212 * t231;
1070   t327 = t15 * t12;
1071   t328 = t57 * t327;
1072   t331 = t77 * t187;
1073   t334 = t2 * ZA;
1074   t335 = t334 * t122;
1075   t336 = t50 * t190;
1076   t339 = 0.8e1 * t151 * t90 * t134 + 0.16e2 * t186 * t30 * t301 - 0.2e1 * t70 * t305 + 0.2e1 * t308 * t64 - 0.2e1 * t312 * ZB * t83 * t11 + 0.2e1 * t56 * t318 + 0.2e1 * t321 * t322 - t116 * t325 - 0.4e1 * t274 * t328 + 0.2e1 * t331 * t305 - 0.16e2 * t335 * t336;
1077   t341 = t169 * t51;
1078   t344 = t49 * t11 * t72;
1079   t346 = t77 * t30;
1080   t347 = t50 * t304;
1081   t350 = t25 * t51;
1082   t352 = nx * ZA;
1083   t353 = t49 * t23;
1084   t354 = t352 * t353;
1085   t355 = t28 * xc;
1086   t362 = t25 * t91;
1087   t365 = t23 * n;
1088   t366 = nx * t365;
1089   t367 = t366 * t122;
1090   t368 = ZB * t49;
1091   t369 = ZA * t51;
1092   t370 = t369 * t17;
1093   t371 = t368 * t370;
1094   t374 = t115 * t353;
1095   t375 = t355 * t73;
1096   t381 = t105 * t341 - t105 * t344 - 0.2e1 * t346 * t347 - t350 * t50 - 0.8e1 * t354 * t355 * t176 - 0.4e1 * t98 * t99 * t317 - 0.2e1 * t362 * t99 - 0.16e2 * t367 * t371 + 0.6e1 * t374 * t375 - 0.8e1 * t182 * t53 - t82 * t290;
1097   t382 = t71 * t208;
1098   t394 = t2 * t67;
1099   t395 = t394 * t122;
1100   t398 = t352 * t156;
1101   t402 = t17 * t12;
1102   t403 = t39 * ZA;
1103   t404 = t402 * t403;
1104   t407 = t269 * t208;
1105   t411 = t49 * t83 * t11;
1106   t413 = t46 * t69;
1107   t419 = -0.4e1 * t331 * t382 + 0.2e1 * t115 * t58 * t204 - 0.2e1 * t145 * t60 + 0.12e2 * t274 * t131 + 0.2e1 * t346 * t232 + 0.8e1 * t395 * t53 - 0.8e1 * t398 * t90 * t176 - 0.64e2 * t260 * t404 + 0.4e1 * t219 * t407 + t168 * t411 - 0.6e1 * t413 * t74 - 0.2e1 * t110 * t308 * t57;
1108   t424 = t16 * t11;
1109   t425 = t212 * t424;
1110   t427 = t258 * t181;
1111   t430 = t67 * t29;
1112   t431 = t366 * t430;
1113   t432 = t121 * t49;
1114   t433 = t432 * t52;
1115   t436 = n * t12;
1116   t437 = t436 * t18;
1117   t440 = t29 * xc;
1118   t441 = t440 * t121;
1119   t442 = t394 * t441;
1120   t445 = t67 * t37;
1121   t446 = t445 * t39;
1122   t448 = n * t18 * t91;
1123   t453 = t352 * t49;
1124   t458 = t8 * t23;
1125   t462 = t81 * t458;
1126   t463 = t30 * t19;
1127   t466 = -t47 * t209 + t116 * t425 - 0.8e1 * t427 * t275 + 0.8e1 * t431 * t433 - 0.2e1 * t253 * t437 - 0.8e1 * t442 * t53 - 0.2e1 * t446 * t448 + 0.2e1 * t175 * t312 + 0.6e1 * t453 * t129 * t208 + 0.8e1 * t115 * t18 * t458 * t30 + 0.8e1 * t462 * t463;
1128   t470 = t436 * t58;
1129   t475 = t2 * t121 * t440 * t25;
1130   t485 = t212 * t73;
1131   t488 = t67 * t72 * t1;
1132   t490 = t39 * xc;
1133   t501 = 0.4e1 * t374 * t355 * t52 + 0.2e1 * t446 * t470 - 0.8e1 * t475 * t53 - 0.2e1 * t446 * t437 - 0.4e1 * t36 * t38 * t39 * t15 * t12 - t116 * t485 + 0.8e1 * t488 * 0.3141592654e1 * t12 * t490 - t207 * t183 - 0.2e1 * t182 * t232 - 0.6e1 * t413 * t245 - 0.4e1 * t413 * t382;
1134   t503 = t115 * t8;
1135   t510 = t355 * t19;
1136   t513 = t432 * t208;
1137   t525 = t38 * t40 * t18;
1138   t533 = -0.4e1 * t503 * t216 - 0.4e1 * t89 * t90 * t92 * t15 - 0.16e2 * t462 * t510 + 0.8e1 * t431 * t513 - 0.4e1 * t78 * t131 + t47 * t183 - 0.2e1 * t67 * t83 * t99 + 0.4e1 * t331 * t225 + 0.16e2 * t260 * t525 - 0.4e1 * t89 * t90 * t237 * t58 - t207 * t53;
1139   t536 = t28 * t37;
1140   t538 = t490 * t100;
1141   t541 = t334 * t441;
1142   t547 = t394 * t30;
1143   t550 = t212 * t19;
1144   t553 = t366 * t441;
1145   t556 = n * t17;
1146   t571 = -0.8e1 * t427 * t131 + 0.16e2 * t394 * t536 * t538 + 0.16e2 * t541 * t336 + 0.2e1 * t453 * t129 * t158 - 0.8e1 * t547 * t147 + 0.4e1 * t503 * t550 - 0.8e1 * t553 * t270 + 0.4e1 * t556 * ZB * t92 * t39 - 0.2e1 * t67 * t91 * t99 - t82 * t425 + 0.4e1 * t78 * t275 + 0.2e1 * t78 * xc * t41 * t57;
1147   t583 = t90 * t317;
1148   t594 = t212 * t158;
1149   t596 = t152 * t67;
1150   t602 = t67 * t17;
1151   t607 = 0.8e1 * t367 * t407 - 0.4e1 * t98 * t99 * t59 + 0.16e2 * t260 * t18 * ZA * t57 + 0.2e1 * t321 * t583 - 0.6e1 * t174 * t368 * t52 - 0.4e1 * t89 * t90 * ZA * t15 * t12 + t116 * t594 - 0.8e1 * t596 * t136 - 0.4e1 * t98 * t99 * t327 + 0.2e1 * t602 * t99 + 0.2e1 * t164 * t583;
1152   t613 = t83 * t25;
1153   t616 = t81 * t156;
1154   t627 = t90 * t231;
1155   t630 = t91 * t16;
1156   t638 = 0.4e1 * t196 * t90 * t208 - 0.8e1 * t130 * t60 - 0.2e1 * t613 * t99 + 0.6e1 * t616 * t197 - 0.8e1 * t547 * t131 + 0.8e1 * t67 * t18 * t37 * t142 + 0.2e1 * t145 * t328 - 0.6e1 * t196 * t627 + 0.8e1 * t630 * t67 * t142 - 0.8e1 * t547 * t275 + 0.8e1 * t395 * t209;
1157   t643 = t77 * t355;
1158   t648 = t115 * t458;
1159   t651 = t134 * t67;
1160   t657 = t30 * t304;
1161   t660 = t30 * t146;
1162   t665 = t25 * t17;
1163   t668 = t50 * t424;
1164   t671 = -0.4e1 * t321 * t90 * t146 - 0.6e1 * t643 * t232 + 0.8e1 * t182 * t209 - 0.16e2 * t648 * t510 + 0.8e1 * t651 * t136 + 0.8e1 * t89 * t4 * t100 - 0.2e1 * t374 * t657 - 0.8e1 * t648 * t660 + 0.8e1 * t130 * t328 + 0.2e1 * t665 * t99 + 0.2e1 * t346 * t668;
1165   t672 = t90 * t424;
1166   t676 = t120 * t536;
1167   t680 = t436 * t41;
1168   t688 = t366 * t67 * t440;
1169   t696 = xc * t12;
1170   t697 = t696 * t18;
1171   t701 = t252 * t141;
1172   t702 = t90 * t221;
1173   t705 = 0.2e1 * t196 * t672 - t47 * t347 + 0.16e2 * t676 * t538 - t116 * t85 - 0.2e1 * t253 * t680 + t207 * t194 + 0.4e1 * t98 * t99 * t19 - 0.8e1 * t688 * t433 + 0.16e2 * t541 * t301 - 0.6e1 * t312 * t190 + 0.4e1 * t352 * t88 * t35 * t697 + 0.2e1 * t701 * t702;
1174   t712 = t24 * t430;
1175   t713 = t28 * t49;
1176   t721 = t1 * t17 * t11;
1177   t726 = ZB * xc;
1178   t737 = n * t91;
1179   t741 = 0.8e1 * t346 * t209 + 0.2e1 * t712 * t713 * t424 + 0.8e1 * t130 * t275 - t47 * t668 + t116 * t721 - 0.8e1 * t688 * t513 + 0.4e1 * t352 * t27 * t212 * t726 + 0.8e1 * t648 * t463 + 0.4e1 * t274 * t60 - 0.4e1 * t374 * t355 * t208 - 0.4e1 * t253 * t737 * t41;
1180   t745 = t269 * t231;
1181   t749 = t1 * t28 * t265;
1182   t757 = t16 * t39;
1183   t758 = t696 * t757;
1184   t762 = t69 * t49;
1185   t772 = t355 * t100;
1186   t775 = t81 * t353;
1187   t778 = -0.8e1 * t398 * t90 * t190 - 0.2e1 * t278 * t745 + 0.4e1 * t749 * t53 + 0.32e2 * t394 * t29 * t28 * t17 * t40 - 0.8e1 * t78 * t758 + t350 * n * t762 - 0.6e1 * t87 * t49 * t186 * t52 - 0.8e1 * t553 * t407 - 0.4e1 * t749 * t209 + 0.16e2 * t648 * t772 - 0.6e1 * t775 * t375;
1188   t790 = t212 * t304;
1189   t793 = t156 * 0.3141592654e1;
1190   t795 = t355 * t304;
1191   t800 = t91 * t39;
1192   t801 = t800 * n;
1193   t807 = t2 * t28;
1194   t808 = t807 * t726;
1195   t811 = -0.2e1 * t616 * t672 - 0.2e1 * t446 * t680 - 0.2e1 * t78 * xc * t58 * t57 + 0.8e1 * t367 * t270 - t82 * t790 + t115 * t51 * t793 - 0.2e1 * t775 * t795 + 0.8e1 * t123 * t209 + 0.2e1 * t665 * t801 - 0.2e1 * t67 * t41 * t64 - 0.32e2 * t808 * t43;
1196   t812 = t117 * t51;
1197   t821 = t24 * t355;
1198   t827 = t90 * t304;
1199   t840 = t800 * t41;
1200   t844 = -t116 * t812 - 0.2e1 * t110 * t25 * t58 * t57 - 0.4e1 * t78 * t328 + t82 * t485 - 0.4e1 * t821 * t407 + 0.4e1 * t196 * t90 * t52 + 0.2e1 * t196 * t827 + t82 * t325 + 0.2e1 * t253 * t448 - 0.32e2 * t402 * t67 * t807 * t490 - t207 * t232 + 0.12e2 * t186 * t90 * ZB * t37 * t840;
1201   t849 = t1 * t51;
1202   t850 = t849 * t17;
1203   t860 = t269 * t424;
1204   t863 = t273 * t187;
1205   t874 = 0.16e2 * t462 * t772 - t116 * t850 + 0.16e2 * t553 * t371 + t116 * t288 - 0.12e2 * t97 * t57 * t109 * t697 + t82 * t594 - 0.2e1 * t278 * t860 - 0.2e1 * t863 * t305 - 0.16e2 * t180 * t259 * t311 * t201 - 0.6e1 * t863 * t74 + 0.8e1 * t174 * t191;
1206   t879 = xc * ZA;
1207   t888 = t67 * t51;
1208   t901 = ZA * t17;
1209   t903 = t368 * t901 * t11;
1210   t908 = -0.2e1 * t352 * t51 * t156 * t35 + 0.64e2 * t879 * t189 * t807 * t40 + 0.2e1 * t46 * t58 * t37 * t39 - t888 * t50 + t105 * t411 - 0.16e2 * t335 * t301 + 0.8e1 * t152 * t25 * t136 - 0.8e1 * t278 * t407 + 0.2e1 * t712 * t713 * t231 - 0.16e2 * t367 * t903 + 0.2e1 * t145 * t318;
1211   t923 = t71 * t424;
1212   t926 = t87 * t458;
1213   t927 = t28 * ZA;
1214   t944 = 0.8e1 * t354 * t355 * t190 - 0.8e1 * t110 * t16 * t25 * t800 - 0.2e1 * t374 * t30 * t73 - 0.16e2 * t354 * t30 * t176 - 0.2e1 * t244 * t923 - 0.32e2 * t926 * t927 * t696 * t41 - 0.32e2 * t808 * t525 + 0.6e1 * t749 * t232 - 0.8e1 * t188 * t177 + 0.4e1 * t36 * t58 * ZA * t57 + 0.4e1 * t821 * t270;
1215   t948 = t90 * t327;
1216   t961 = t30 * t100;
1217   t964 = t29 * t49;
1218   t981 = t106 * t1;
1219   t983 = -0.2e1 * t219 * t220 * t100 + 0.2e1 * t321 * t948 - 0.16e2 * t189 * ZA * t99 * t12 - 0.2e1 * t369 * n * t69 * t368 + 0.2e1 * t374 * t795 - 0.8e1 * t462 * t961 - 0.8e1 * t244 * t964 * t208 + 0.2e1 * t413 * t923 + 0.4e1 * t36 * t38 * t40 * t58 - 0.2e1 * t87 * t51 * t49 * t1 * ZA + t888 * n * t762 + t115 * t981;
1220   t1012 = 0.6e1 * t616 * t627 - t82 * t213 + 0.2e1 * t775 * t657 - 0.12e2 * t215 * t550 - 0.6e1 * t145 * t131 + 0.2e1 * t81 * t41 * t204 + 0.6e1 * ZB * t48 * t49 * t370 - 0.4e1 * t70 * t382 + 0.2e1 * t446 * t255 + 0.8e1 * t89 * t4 * t327 - 0.4e1 * t56 * t147;
1221   t1018 = t212 * t100;
1222   t1029 = t212 * t327;
1223   t1040 = 0.6e1 * t70 * t245 + 0.2e1 * t56 * t328 + t207 * t668 + 0.4e1 * t503 * t1018 + 0.2e1 * t253 * t470 - 0.6e1 * t398 * t35 * t208 - 0.8e1 * t331 * t964 * t52 - 0.4e1 * t503 * t1029 + 0.6e1 * t821 * t745 + 0.4e1 * t63 * t37 * t142 + 0.16e2 * t260 * t38 * t840;
1224   t1068 = t207 * t347 - 0.2e1 * t164 * t702 - 0.2e1 * t331 * t964 * t73 + 0.8e1 * t374 * t30 * t52 + 0.16e2 * t278 * t903 + 0.2e1 * t863 * t923 + 0.6e1 * t445 * t141 * t165 - 0.2e1 * t164 * t90 * t100 + 0.6e1 * t331 * t74 - 0.2e1 * t182 * t668 - 0.2e1 * t115 * t41 * t204;
1225   t1079 = t58 * t8;
1226   t1091 = t807 * t29;
1227   t1092 = t665 * t40;
1228   t1101 = ZB * t91;
1229   t1102 = t403 * n;
1230   t1105 = -0.4e1 * t58 * ZB * ZA * t64 - t82 * t850 + 0.2e1 * t821 * t860 + t81 * t51 * t793 + 0.2e1 * t3 * t25 * t1079 * t90 + t82 * t721 - 0.2e1 * t643 * t668 + 0.16e2 * t926 * t927 * t29 * t91 * t41 + 0.32e2 * t1091 * t1092 - 0.2e1 * t219 * t220 * t19 + 0.4e1 * t139 * ZA * t64 + 0.4e1 * t1101 * t1102;
1231   t1108 = t849 * t72;
1232   t1121 = t737 * t15;
1233   t1124 = t29 * t12;
1234   t1133 = t116 * t1108 - 0.8e1 * t475 * t209 - 0.32e2 * t807 * xc * t1092 + 0.2e1 * t278 * t269 * t73 + t82 * t812 - 0.6e1 * t56 * t131 + 0.2e1 * t253 * t1121 + 0.16e2 * t926 * t927 * t1124 * t41 + t168 * t263 - 0.2e1 * t616 * t827 + t81 * t981;
1235   t1134 = t394 * t28;
1236   t1159 = -0.8e1 * t1134 * t29 * t18 * t57 + t82 * t118 - 0.12e2 * t215 * t1018 + 0.2e1 * t602 * t801 - t168 * t341 + 0.2e1 * t67 * t58 * t64 + t168 * t344 - 0.6e1 * t174 * t368 * t208 + 0.16e2 * t553 * t903 + t116 * t790 - 0.4e1 * t36 * t38 * t800 * t15;
1237   t1161 = n * t83;
1238   t1173 = ZB * t12;
1239   t1196 = 0.4e1 * t1161 * ZB * t39 * ZA - 0.4e1 * t215 * t1029 - 0.8e1 * t488 * 0.3141592654e1 * t39 * xc + 0.32e2 * t821 * ZA * t8 * t1173 * t18 - 0.8e1 * t427 * t147 + 0.6e1 * t701 * t165 - 0.16e2 * t926 * t927 * t1124 * t18 - 0.8e1 * t1091 * t63 * t57 - 0.8e1 * t442 * t209 - 0.8e1 * t462 * t660 - 0.6e1 * t398 * t35 * t52;
1240   t1228 = 0.2e1 * t413 * t305 - 0.8e1 * t648 * t961 - 0.16e2 * t87 * t27 * t23 * t28 * ZA * t29 + 0.4e1 * t189 * t1102 - 0.4e1 * t87 * t1079 * t212 * t879 + 0.2e1 * t164 * t948 - 0.2e1 * t70 * t923 + 0.2e1 * t164 * t322 + 0.2e1 * t446 * t1121 + 0.2e1 * t863 * t964 * t304 - t82 * t1108 + 0.16e2 * t676 * t490 * t19;
1241   t1234 = t25 * ZB;
1242   t1235 = t1234 * t28;
1243   t1236 = t365 * t91;
1244   t1240 = ZB * t121;
1245   t1241 = t1240 * t77;
1246   t1242 = t39 * t39;
1247   t1243 = t12 * t1242;
1248   t1244 = xc * t72;
1249   t1245 = t1243 * t1244;
1250   t1248 = t365 * t25;
1251   t1252 = t243 * n;
1252   t1257 = t23 * t1;
1253   t1258 = t1240 * t1257;
1254   t1259 = t67 * t12;
1255   t1260 = xc * t16;
1256   t1268 = t1234 * t121;
1257   t1269 = t1268 * t23;
1258   t1272 = t1242 * t91;
1259   t1280 = t67 * xc;
1260   t1284 = t28 * t28;
1261   t1285 = t67 * t1284;
1262   t1287 = t1285 * t2 * ZB;
1263   t1288 = t17 * xc;
1264   t1289 = t1243 * t1288;
1265   t1292 = 0.2e1 * t1235 * t1236 * t17 + 0.8e1 * t1241 * t1245 + 0.4e1 * t927 * t1248 * t91 - 0.2e1 * t1252 * ZB * t1242 * t91 - 0.8e1 * t1258 * t1259 * t1260 - 0.4e1 * t1235 * t2 * t83 * t39 + 0.16e2 * t1269 * t758 + 0.2e1 * t1252 * t189 * t1272 - 0.2e1 * t1252 * t83 * t1242 * ZB + 0.8e1 * t1258 * t630 * t1280 - 0.32e2 * t1287 * t1289;
1266   t1293 = t365 * t83;
1267   t1300 = ZA * t1284;
1268   t1304 = t17 * t1242 * t25 * t12;
1269   t1307 = t927 * t2;
1270   t1311 = t23 * t2;
1271   t1312 = t1300 * t1311;
1272   t1316 = t1234 * t1284;
1273   t1317 = t1316 * t1311;
1274   t1321 = t1240 * t23;
1275   t1331 = t1240 * t23 * t67;
1276   t1332 = t40 * t1244;
1277   t1338 = t1243 * t1260;
1278   t1344 = -0.2e1 * t1235 * t1293 - 0.16e2 * t181 * t365 * t901 * t12 - 0.64e2 * t1300 * t258 * t1304 + 0.8e1 * t1307 * t613 * t39 + 0.64e2 * t1312 * t265 * t402 - 0.32e2 * t1317 * t1288 * t12 - 0.16e2 * t1321 * t67 * t39 * t1244 + 0.2e1 * t1235 * n * t1272 * t17 + 0.16e2 * t1331 * t1332 + 0.64e2 * t1300 * t292 * t1304 - 0.8e1 * t1241 * t1338 - 0.2e1 * t243 * t1293 * ZB;
1279   t1346 = t1316 * t2;
1280   t1349 = t927 * n;
1281   t1350 = t25 * t1242;
1282   t1354 = t1268 * t1257;
1283   t1366 = t1268 * t1;
1284   t1370 = t29 * t17;
1285   t1371 = t1243 * t1370;
1286   t1386 = -0.32e2 * t1346 * t1289 + 0.4e1 * t1349 * t1350 * t91 + 0.8e1 * t1354 * t1260 * t12 - 0.16e2 * t181 * n * t901 * t1243 - 0.4e1 * t1235 * t2 * t91 * t39 + 0.8e1 * t1366 * t152 * t1242 + 0.32e2 * t1287 * t1371 + 0.8e1 * t1258 * t1280 * t152 - 0.8e1 * t1354 * t1260 * t91 + 0.128e3 * t1300 * t365 * xc * t1092 + 0.8e1 * t1366 * t1338;
1287   t1387 = t1257 * t12;
1288   t1391 = t1240 * t1;
1289   t1399 = t1272 * t1260;
1290   t1412 = t1285 * t1311;
1291   t1427 = -0.8e1 * t1268 * t1387 * t16 - 0.8e1 * t1391 * t67 * t1242 * t1244 - 0.4e1 * t1134 * t1101 * t39 + 0.8e1 * t1241 * t1399 - 0.8e1 * t1258 * t596 + 0.4e1 * t927 * t1293 * t25 - 0.16e2 * t1331 * t758 + 0.8e1 * t1307 * t665 * t39 + 0.32e2 * t1412 * t1370 * t1173 + 0.8e1 * t1307 * t665 * t800 + 0.8e1 * t1391 * t1259 * t1242 * t16 - 0.8e1 * t1391 * t1259 * t1242 * t72;
1292   t1456 = t365 * ZB;
1293   t1468 = 0.4e1 * t927 * t1248 * t17 - 0.2e1 * t1235 * n * t1242 * t91 + 0.8e1 * t1366 * t1244 * t1242 - 0.16e2 * t1269 * t134 * t39 + 0.8e1 * t1268 * t1257 * t72 * xc + 0.16e2 * t1321 * t1259 * t757 + 0.32e2 * t1317 * t1370 * t12 + 0.4e1 * t1349 * t613 * t1242 + 0.2e1 * t243 * t1456 * t17 - 0.64e2 * t1285 * t365 * t12 * t189 * t490 - 0.8e1 * t1354 * t152 * xc;
1294   t1472 = t1316 * t365;
1295   t1474 = t1124 * t39 * t17;
1296   t1478 = t17 * t91;
1297   t1504 = t72 * t39;
1298   t1511 = 0.4e1 * t1134 * t189 * t800 + 0.64e2 * t1472 * t1474 + 0.4e1 * t1235 * t2 * t1478 * t39 + 0.4e1 * t1349 * t665 * t1242 - 0.8e1 * t1258 * t1280 * t72 + 0.2e1 * t1252 * t189 * t1242 + 0.2e1 * t243 * t365 * t189 * t91 + 0.4e1 * t927 * t365 * t1478 * t25 - 0.128e3 * t1300 * t1248 * t1474 - 0.2e1 * t1235 * t1236 + 0.16e2 * t1269 * t1504 * xc + 0.2e1 * t1235 * t365 * t17;
1299   t1545 = -0.2e1 * t1235 * t1161 * t1242 + 0.4e1 * t1349 * t1350 * t1478 - 0.8e1 * t1366 * t1245 + 0.2e1 * t1235 * t556 * t1242 - 0.32e2 * t1412 * t402 * t726 - 0.8e1 * t1366 * t1399 + 0.8e1 * t1258 * t651 - 0.2e1 * t243 * t1456 * t91 + 0.8e1 * t1268 * t1387 * t72 - 0.16e2 * t1269 * t1332 + 0.4e1 * t1134 * t189 * t39 + 0.16e2 * t1269 * t152 * t39;
1300   t1564 = t1260 * t800;
1301   t1583 = 0.64e2 * t1285 * t1456 * t1474 - 0.64e2 * t1472 * t1288 * t40 - 0.8e1 * t1366 * t134 * t1242 + 0.8e1 * t1307 * t362 * t39 + 0.4e1 * t1235 * t2 * t17 * t39 + 0.32e2 * t1346 * t1371 - 0.16e2 * t1269 * t1564 - 0.16e2 * t1321 * t1259 * t1504 + 0.16e2 * t1331 * t1564 - 0.64e2 * t1312 * t29 * t25 * t402 - 0.4e1 * t1134 * t83 * t39 * ZB - 0.32e2 * t181 * t2 * t404;
1302 
1303   _PC2B = (t1133 + t1196 + t1068 + t811 + t466 + t1012 + t381 + t162 + t249 + t533 + t844 + t104 + t1159 + t571 + t211 + t874 + t607 + t339 + t296 + t638 + t908 + t671 + t419 + t983 + t705 + t1105 + t501 + t778 + t1040 + t1228 + t741 + t944) / (t1292 + t1344 + t1386 + t1427 + t1468 + t1511 + t1545 + t1583);
1304 
1305   t1   = n * n;
1306   t2   = t1 * n;
1307   t3   = t2 * nx;
1308   t4   = nx * 0.3141592654e1;
1309   t5   = t4 * xc;
1310   t6   = PetscSinReal(t5);
1311   t7   = 0.3141592654e1 * 0.3141592654e1;
1312   t9   = t3 * t6 * t7;
1313   t10  = xc * xc;
1314   t11  = ZA * ZA;
1315   t12  = t10 * t11;
1316   t13  = n * 0.3141592654e1;
1317   t14  = PetscExpReal(t13);
1318   t15  = t14 * t14;
1319   t16  = xc * n;
1320   t18  = PetscExpReal(t16 * 0.3141592654e1);
1321   t19  = t18 * t18;
1322   t20  = t19 * t18;
1323   t21  = t15 * t20;
1324   t22  = t12 * t21;
1325   t25  = nx * t6;
1326   t26  = t1 * 0.3141592654e1;
1327   t27  = t25 * t26;
1328   t28  = ZA * ZB;
1329   t29  = t18 * t15;
1330   t30  = t28 * t29;
1331   t33  = t25 * n;
1332   t34  = t11 * t15;
1333   t35  = t19 * t19;
1334   t36  = t35 * t18;
1335   t40  = t25 * t1;
1336   t41  = 0.3141592654e1 * t11;
1337   t42  = t15 * t36;
1338   t43  = t41 * t42;
1339   t46  = nx * nx;
1340   t47  = t1 * t46;
1341   t48  = t47 * t11;
1342   t49  = t7 * xc;
1343   t50  = t35 * t15;
1344   t51  = t49 * t50;
1345   t55  = PetscSinReal(t4);
1346   t56  = t46 * nx * t55;
1347   t58  = t56 * n * t7;
1348   t59  = ZB * ZB;
1349   t60  = t10 * t59;
1350   t61  = t15 * t14;
1351   t62  = t19 * t61;
1352   t63  = t60 * t62;
1353   t66  = t19 * t14;
1354   t67  = t60 * t66;
1355   t70  = t28 * t42;
1356   t73  = PetscCosReal(t5);
1357   t74  = t47 * t73;
1358   t75  = t7 * t11;
1359   t77  = t75 * t10 * t36;
1360   t80  = t73 * t46;
1361   t81  = t80 * n;
1362   t82  = 0.3141592654e1 * t59;
1363   t83  = t82 * t42;
1364   t87  = xc * t11;
1365   t88  = t87 * t62;
1366   t91  = n * nx;
1367   t92  = t55 * t61;
1368   t96  = nx * t55;
1369   t98  = t96 * t2 * t7;
1370   t101 = xc * t59;
1371   t102 = t101 * t62;
1372   t108 = t1 * t1;
1373   t109 = t108 * t7;
1374   t111 = t59 * t35;
1375   t112 = t111 * t15;
1376   t115 = t35 * t20;
1377   t123 = t1 * nx * t55;
1378   t124 = t61 * t35;
1379   t127 = t35 * t19;
1380   t128 = t61 * t127;
1381   t129 = t60 * t128;
1382   t132 = t56 * t16;
1383   t133 = t7 * t59;
1384   t134 = t133 * t124;
1385   t137 = 0.6e1 * t58 * t88 - 0.2e1 * t91 * t92 * t11 + 0.2e1 * t98 * t63 - 0.6e1 * t58 * t102 - 0.2e1 * t91 * t92 * t59 - 0.16e2 * t109 * xc * t112 - 0.2e1 * t91 * t6 * t115 * t59 + 0.12e2 * t40 * t83 + t123 * t41 * t124 - 0.2e1 * t58 * t129 + 0.4e1 * t132 * t134;
1386   t139 = t56 * 0.3141592654e1;
1387   t140 = t111 * t14;
1388   t144 = t49 * t124;
1389   t147 = t91 * t55;
1390   t148 = t61 * ZA;
1391   t154 = ZA * t115 * xc * ZB;
1392   t157 = t7 * 0.3141592654e1;
1393   t159 = t96 * t108 * t157;
1394   t160 = t10 * xc;
1395   t161 = t160 * t59;
1396   t162 = t35 * t14;
1397   t163 = t161 * t162;
1398   t166 = t28 * t162;
1399   t169 = t80 * t13;
1400   t170 = t101 * t42;
1401   t173 = t2 * t11;
1402   t174 = t96 * t173;
1403   t175 = t7 * t10;
1404   t179 = t59 * t15;
1405   t184 = t15 * t15;
1406   t193 = t139 * t140 + 0.4e1 * t56 * n * t11 * t144 + 0.4e1 * t147 * t148 * ZB + 0.4e1 * t27 * t154 + 0.8e1 * t159 * t163 - 0.12e2 * t147 * t166 + 0.2e1 * t169 * t170 - 0.16e2 * t174 * t175 * t124 + 0.2e1 * t33 * t179 * t20 - 0.2e1 * t33 * t11 * t36 * t184 + 0.2e1 * t56 * t61 * 0.3141592654e1 * ZA * ZB;
1407   t194 = t173 * 0.3141592654e1;
1408   t195 = xc * t15;
1409   t196 = t195 * t19;
1410   t202 = t15 * t115;
1411   t203 = t28 * t202;
1412   t206 = t96 * t26;
1413   t207 = t14 * t127;
1414   t208 = t101 * t207;
1415   t211 = t12 * t128;
1416   t218 = t11 * t61;
1417   t219 = t218 * t35;
1418   t221 = t108 * ZA;
1419   t223 = t7 * ZB;
1420   t224 = t223 * t50;
1421   t227 = ZA * xc;
1422   t228 = ZB * t15;
1423   t229 = t228 * t36;
1424   t230 = t227 * t229;
1425   t233 = t87 * t207;
1426   t236 = t6 * t11;
1427   t240 = -0.4e1 * t194 * t196 + 0.4e1 * t194 * t195 * t127 + 0.4e1 * t33 * t203 - 0.12e2 * t206 * t208 + 0.2e1 * t58 * t211 - 0.16e2 * t47 * t10 * t133 * t50 + t139 * t219 - 0.32e2 * t221 * t10 * t224 - 0.4e1 * t169 * t230 - 0.6e1 * t98 * t233 + 0.2e1 * t91 * t236 * t20;
1428   t244 = t227 * t228 * t20;
1429   t252 = t184 * t18;
1430   t253 = t101 * t252;
1431   t256 = t35 * t35;
1432   t257 = t256 * t14;
1433   t258 = t28 * t257;
1434   t261 = t108 * t11;
1435   t263 = t7 * t35;
1436   t268 = ZB * t61 * t35;
1437   t273 = t96 * t108 * t160;
1438   t274 = t157 * ZB;
1439   t276 = t274 * t148 * t35;
1440   t279 = t101 * t21;
1441   t282 = 0.3141592654e1 * xc;
1442   t283 = t59 * t36;
1443   t284 = t282 * t283;
1444   t289 = 0.4e1 * t169 * t244 - 0.4e1 * t132 * t133 * t162 - 0.2e1 * t147 * t140 - 0.2e1 * t27 * t253 + 0.2e1 * t139 * t258 + 0.16e2 * t261 * t10 * t263 * t15 - 0.16e2 * t206 * t227 * t268 - 0.16e2 * t273 * t276 - 0.6e1 * t27 * t279 - 0.4e1 * t40 * t284 - 0.32e2 * t9 * t230;
1445   t290 = t1 * t11;
1446   t291 = t96 * t290;
1447   t297 = t59 * t61;
1448   t298 = t297 * t127;
1449   t300 = ZB * t36;
1450   t301 = t227 * t300;
1451   t304 = t1 * t59;
1452   t305 = t184 * t35;
1453   t310 = t46 * ZB;
1454   t311 = t184 * ZA;
1455   t312 = t310 * t311;
1456   t314 = t60 * t21;
1457   t317 = t1 * ZA;
1458   t318 = ZB * t35;
1459   t321 = t1 * t256;
1460   t324 = t96 * t261;
1461   t325 = t10 * t157;
1462   t326 = t325 * t124;
1463   t329 = -0.4e1 * t291 * t282 * t128 + t123 * t82 * t62 - t139 * t298 + 0.12e2 * t27 * t301 + t304 * t305 - 0.2e1 * t58 * t12 * t66 - 0.2e1 * t312 + 0.8e1 * t9 * t314 + 0.2e1 * t317 * t318 + 0.2e1 * t321 * t28 - 0.8e1 * t324 * t326;
1464   t331 = t28 * t124;
1465   t334 = 0.3141592654e1 * t15;
1466   t335 = t334 * t127;
1467   t338 = t35 * ZA;
1468   t341 = t46 * t256;
1469   t344 = t46 * t11;
1470   t346 = t46 * t59;
1471   t348 = t297 * t35;
1472   t351 = ZA * t10;
1473   t352 = t351 * t300;
1474   t355 = t1 * ZB;
1475   t362 = 0.12e2 * t147 * t331 - 0.4e1 * t173 * t335 - 0.2e1 * t310 * t338 - 0.2e1 * t341 * t28 - t344 * t305 - t346 * t305 + 0.2e1 * t147 * t348 + 0.16e2 * t9 * t352 + 0.2e1 * t355 * t311 + t290 * t305 + 0.2e1 * t33 * t34 * t20;
1476   t363 = t36 * t184;
1477   t364 = t87 * t363;
1478   t368 = t47 * t73 * t7;
1479   t373 = t160 * t157;
1480   t374 = t373 * t124;
1481   t377 = t311 * t35;
1482   t380 = t12 * t62;
1483   t386 = ZB * t10 * ZA * t15 * t20;
1484   t389 = t87 * t66;
1485   t393 = t56 * t1 * t10;
1486   t401 = 0.2e1 * t27 * t364 - 0.16e2 * t368 * t279 - t123 * t41 * t257 + 0.8e1 * t324 * t374 + 0.2e1 * t355 * t377 - 0.2e1 * t98 * t380 - 0.16e2 * t9 * t386 + 0.2e1 * t58 * t389 + 0.16e2 * t393 * t276 + t123 * t82 * t162 - 0.2e1 * t33 * t179 * t36;
1487   t412 = t11 * t14 * t127;
1488   t416 = t11 * t19;
1489   t417 = t416 * t61;
1490   t421 = t96 * t2 * ZA;
1491   t426 = t56 * n * ZA;
1492   t427 = t318 * t14;
1493   t428 = t49 * t427;
1494   t431 = t82 * t29;
1495   t434 = t87 * t21;
1496   t442 = 0.2e1 * t33 * t11 * t184 * t18 + 0.4e1 * t81 * t284 - t139 * t412 + 0.2e1 * t147 * t219 - 0.2e1 * t147 * t417 + 0.32e2 * t421 * t175 * t268 + 0.8e1 * t426 * t428 + 0.4e1 * t81 * t431 - 0.2e1 * t169 * t434 - 0.2e1 * t98 * t129 - 0.32e2 * t47 * t28 * t51;
1497   t443 = t184 * t20;
1498   t447 = t61 * 0.3141592654e1;
1499   t448 = t447 * t11;
1500   t450 = t49 * t268;
1501   t453 = t60 * t42;
1502   t456 = t41 * t202;
1503   t463 = t101 * t443;
1504   t469 = t41 * xc * t20;
1505   t474 = -0.8e1 * t27 * t87 * t443 - t56 * t448 - 0.8e1 * t426 * t450 + 0.8e1 * t368 * t453 + 0.4e1 * t40 * t456 + 0.4e1 * t40 * t431 - 0.4e1 * t81 * t456 - 0.4e1 * t27 * t463 + 0.6e1 * t139 * t331 + 0.2e1 * t40 * t469 - 0.16e2 * t9 * t434;
1506   t482 = t108 * t10;
1507   t492 = n * t46;
1508   t493 = t492 * t11;
1509   t495 = t282 * t19 * t184;
1510   t498 = t56 * t290;
1511   t499 = t325 * t162;
1512   t502 = t416 * t14;
1513   t504 = t60 * t207;
1514   t507 = -t123 * t82 * t257 - 0.4e1 * t169 * t301 + t123 * t41 * t162 + 0.16e2 * t482 * t7 * t112 - 0.12e2 * t206 * t102 - t123 * t82 * t66 - 0.4e1 * t147 * t258 - 0.4e1 * t493 * t495 - 0.8e1 * t498 * t499 + t139 * t502 - 0.2e1 * t98 * t504;
1515   t508 = t101 * t162;
1516   t512 = t41 * t115 * xc;
1517   t515 = t87 * t42;
1518   t520 = ZB * t184;
1519   t522 = t227 * t520 * t18;
1520   t525 = t492 * t59;
1521   t528 = t6 * t59;
1522   t532 = t520 * t20;
1523   t533 = t351 * t532;
1524   t539 = t447 * t59;
1525   t544 = 0.8e1 * t206 * t508 - 0.2e1 * t40 * t512 - 0.16e2 * t368 * t515 + 0.12e2 * t206 * t88 + 0.4e1 * t27 * t522 + 0.4e1 * t525 * t495 - 0.4e1 * t91 * t528 * t36 - 0.16e2 * t368 * t533 - 0.16e2 * t206 * t227 * t427 - t56 * t539 - 0.2e1 * t132 * t133 * t66;
1526   t551 = t87 * t162;
1527   t554 = t351 * t229;
1528   t560 = t59 * t19;
1529   t561 = t560 * t14;
1530   t564 = t101 * t202;
1531   t567 = t87 * t252;
1532   t573 = t227 * t228 * t115;
1533   t578 = 0.4e1 * t33 * t70 + 0.4e1 * t493 * t335 - 0.4e1 * t58 * t551 + 0.16e2 * t9 * t554 - 0.4e1 * t33 * t28 * t252 + 0.2e1 * t147 * t561 + 0.2e1 * t169 * t564 - 0.2e1 * t27 * t567 - 0.8e1 * t324 * t499 - 0.4e1 * t169 * t573 + 0.12e2 * t27 * t244;
1534   t579 = t82 * t202;
1535   t591 = t282 * t115 * t59;
1536   t598 = t101 * t66;
1537   t606 = -0.4e1 * t81 * t579 - 0.2e1 * t169 * t567 - 0.6e1 * t27 * t170 + 0.8e1 * t169 * t203 + 0.2e1 * t98 * t67 + 0.2e1 * t81 * t591 + 0.32e2 * t368 * t244 - 0.2e1 * t27 * t564 + 0.4e1 * t206 * t598 + 0.16e2 * t9 * t170 + 0.2e1 * t33 * t283 * t184;
1538   t608 = t373 * t162;
1539   t611 = t59 * t184;
1540   t617 = t101 * t29;
1541   t624 = t227 * ZB * t18 * t15;
1542   t629 = t157 * t59;
1543   t630 = t629 * t124;
1544   t633 = t3 * t6;
1545   t634 = t175 * t283;
1546   t644 = 0.8e1 * t498 * t608 + 0.2e1 * t33 * t611 * t18 - 0.4e1 * t206 * t389 - 0.2e1 * t27 * t617 - 0.4e1 * t169 * t154 + 0.4e1 * t27 * t624 + 0.12e2 * t27 * t230 - 0.8e1 * t393 * t630 - 0.8e1 * t633 * t634 + 0.16e2 * t47 * t7 * t101 * t50 + 0.2e1 * t123 * t447 * t28;
1547   t645 = t41 * t29;
1548   t648 = t2 * 0.3141592654e1;
1549   t649 = t648 * xc;
1550   t650 = t560 * t184;
1551   t656 = t56 * t1 * t157;
1552   t659 = t87 * t128;
1553   t662 = t96 * t482;
1554   t663 = t629 * t162;
1555   t671 = t161 * t124;
1556   t674 = t218 * t127;
1557   t679 = 0.4e1 * t81 * t645 - 0.4e1 * t649 * t650 - 0.8e1 * t169 * t70 + 0.8e1 * t656 * t163 - 0.2e1 * t98 * t659 - 0.8e1 * t662 * t663 - 0.32e2 * t421 * t175 * t427 - 0.2e1 * t147 * t502 + 0.8e1 * t656 * t671 + 0.2e1 * t147 * t674 - 0.16e2 * t368 * t386;
1558   t714 = t334 * t19;
1559   t719 = t12 * t42;
1560   t722 = t304 * t35 - t346 * t35 + t341 * t59 - t344 * t35 + t344 * t256 + t346 * t184 - 0.16e2 * t368 * t554 - 0.16e2 * t48 * t175 * t50 + 0.4e1 * t525 * t714 - 0.2e1 * t58 * t659 + 0.8e1 * t368 * t719;
1561   t730 = xc * t19;
1562   t735 = t59 * t256 * t14;
1563   t752 = 0.4e1 * t173 * t714 - 0.6e1 * t27 * t515 - 0.16e2 * t9 * t279 + 0.4e1 * t194 * t730 * t184 - t139 * t735 - 0.4e1 * t492 * t127 * t82 * xc - 0.4e1 * t98 * t508 - t123 * t41 * t207 - 0.2e1 * t147 * t298 + 0.8e1 * t368 * t314 + 0.6e1 * t132 * t133 * t207;
1564   t755 = t28 * t21;
1565   t759 = t274 * t338 * t14;
1566   t767 = t11 * t35;
1567   t768 = t767 * t14;
1568   t778 = t560 * t61;
1569   t781 = -0.2e1 * t58 * t504 - 0.8e1 * t27 * t755 + 0.16e2 * t662 * t759 + 0.12e2 * t291 * t282 * t207 - 0.6e1 * t27 * t434 + t139 * t768 - 0.8e1 * t498 * t326 + 0.4e1 * t33 * t611 * t20 + 0.2e1 * t81 * t512 - t139 * t561 + 0.2e1 * t147 * t778;
1570   t786 = t12 * t443;
1571   t790 = t282 * t59 * t20;
1572   t796 = t59 * t14 * t127;
1573   t806 = t41 * t21;
1574   t811 = -0.8e1 * t393 * t663 + 0.8e1 * t368 * t786 + 0.2e1 * t81 * t790 + 0.4e1 * t169 * t624 + t139 * t796 + 0.2e1 * t206 * t258 - 0.2e1 * t40 * t591 - 0.8e1 * t662 * t630 - 0.4e1 * t33 * t30 - 0.4e1 * t40 * t806 + 0.8e1 * t9 * t786;
1575   t819 = t282 * t15 * t127;
1576   t822 = t101 * t363;
1577   t830 = t11 * t256 * t14;
1578   t835 = t227 * t532;
1579   t842 = 0.2e1 * t33 * t11 * t18 * t15 + t123 * t41 * t66 - 0.4e1 * t493 * t819 - 0.2e1 * t27 * t822 - 0.16e2 * t368 * t170 - 0.4e1 * t169 * t463 - t139 * t830 - 0.4e1 * t649 * t179 * t127 + 0.12e2 * t27 * t835 - 0.16e2 * t368 * t434 - 0.2e1 * t40 * t790;
1580   t845 = t87 * t202;
1581   t854 = t338 * t15;
1582   t859 = t12 * t207;
1583   t868 = t139 * t348 - 0.2e1 * t27 * t845 + 0.8e1 * t169 * t755 - 0.2e1 * t58 * t380 + 0.6e1 * t206 * t331 + 0.8e1 * t310 * t854 - 0.2e1 * t169 * t822 + 0.2e1 * t98 * t859 + 0.8e1 * t159 * t671 + 0.8e1 * t74 * t634 - 0.2e1 * t169 * t253;
1584   t880 = t60 * t443;
1585   t891 = t101 * t128;
1586   t894 = -t123 * t539 - 0.2e1 * t147 * t796 + 0.32e2 * t368 * t230 + t139 * t674 - 0.16e2 * t98 * t60 * t124 + 0.32e2 * t9 * t244 + 0.8e1 * t368 * t880 - 0.8e1 * t40 * t41 * xc * t36 - t123 * t82 * t128 - 0.6e1 * t58 * t233 + 0.2e1 * t58 * t891;
1587   t903 = t179 * t19;
1588   t920 = t56 * t1 * t160;
1589   t925 = -0.2e1 * t174 * t175 * t66 - 0.4e1 * t493 * t714 + 0.4e1 * t649 * t903 - 0.4e1 * t81 * t43 + t123 * t82 * t207 + 0.4e1 * t206 * t891 - 0.16e2 * t273 * t759 - 0.8e1 * t27 * t203 + 0.32e2 * t221 * ZB * t51 - 0.16e2 * t920 * t759 - 0.8e1 * t9 * t453;
1590   t932 = t87 * t29;
1591   t945 = t82 * t21;
1592   t953 = -0.16e2 * t920 * t276 - 0.8e1 * t169 * t30 - 0.8e1 * t633 * t77 - 0.2e1 * t27 * t932 - 0.4e1 * t174 * t49 * t162 + 0.8e1 * t206 * t87 * t124 - 0.2e1 * t147 * t768 + 0.4e1 * t169 * t522 - 0.12e2 * t81 * t945 + 0.4e1 * t33 * t28 * t115 + 0.4e1 * t525 * t819;
1593   t971  = t282 * t127;
1594   t978  = -0.6e1 * t98 * t102 + 0.2e1 * t169 * t515 - 0.2e1 * t310 * t377 + 0.2e1 * t147 * t830 + 0.8e1 * t368 * t22 - 0.2e1 * t169 * t617 + 0.16e2 * t662 * t276 - 0.8e1 * t355 * t854 + 0.4e1 * t493 * t971 - 0.16e2 * t9 * t533 - 0.2e1 * t169 * t279;
1595   t997  = xc * t127;
1596   t998  = t997 * t59;
1597   t1003 = 0.4e1 * t40 * t579 + 0.2e1 * t169 * t845 + 0.16e2 * t9 * t515 + 0.8e1 * t206 * t551 + t123 * t41 * t128 + 0.16e2 * t98 * t60 * t162 + 0.2e1 * t169 * t364 - 0.2e1 * t169 * t932 + t139 * t778 + 0.4e1 * t648 * t998 + 0.2e1 * t147 * t412;
1598   t1006 = t2 * t59;
1599   t1017 = xc * t35;
1600   t1033 = 0.4e1 * t1006 * t335 + 0.4e1 * t81 * t806 - 0.2e1 * t33 * t34 * t115 + 0.8e1 * t498 * t374 - 0.16e2 * t261 * t7 * t1017 * t15 + 0.8e1 * t206 * t101 * t124 - t123 * t448 + 0.2e1 * t147 * t735 + 0.6e1 * t98 * t208 + 0.6e1 * t98 * t88 - 0.4e1 * t33 * t755;
1601   t1055 = -0.4e1 * t173 * t971 + 0.2e1 * t98 * t891 + 0.8e1 * t9 * t880 + 0.4e1 * t169 * t835 - t304 * t184 + t344 * t184 - t123 * t41 * t62 - 0.2e1 * t98 * t598 + 0.2e1 * t58 * t859 + 0.32e2 * t47 * t351 * t224 + 0.2e1 * t98 * t389;
1602   t1070 = t15 * t19;
1603   t1089 = -0.16e2 * t368 * t352 - 0.8e1 * t9 * t719 + 0.4e1 * t96 * t2 * xc * t134 - 0.2e1 * t91 * t236 * t115 + 0.4e1 * t27 * t573 + 0.4e1 * t493 * t282 * t1070 + 0.2e1 * t33 * t59 * t18 * t15 + 0.12e2 * t40 * t945 - 0.4e1 * t492 * xc * t82 * t1070 - 0.2e1 * t91 * t528 * t20 + 0.8e1 * t324 * t608;
1604   t1113 = t123 * t82 * t124 + 0.8e1 * t421 * t428 - t139 * t417 + 0.4e1 * t40 * t645 + 0.16e2 * t393 * t759 - 0.2e1 * t33 * t179 * t115 - 0.4e1 * t525 * t335 + 0.4e1 * t33 * t28 * t36 - 0.4e1 * t1006 * t714 + 0.6e1 * t206 * t166 - 0.8e1 * t421 * t450;
1605   t1119 = t321 * t46;
1606   t1122 = t157 * t11;
1607   t1123 = t1122 * t2;
1608   t1124 = t184 * t46;
1609   t1128 = t108 * n;
1610   t1132 = t7 * t7;
1611   t1133 = t1132 * t11;
1612   t1134 = t1133 * t108;
1613   t1135 = t15 * t46;
1614   t1139 = t7 * ZA;
1615   t1140 = t1139 * ZB;
1616   t1141 = t1 * t35;
1617   t1145 = t629 * t2;
1618   t1146 = t1135 * t730;
1619   t1149 = t157 * t1128;
1620   t1150 = t1149 * xc;
1621   t1153 = t46 * xc;
1622   t1154 = t1153 * t127;
1623   t1158 = t184 * t1 * t46;
1624   t1161 = t46 * t46;
1625   t1162 = t35 * t1161;
1626   t1166 = t7 * t1;
1627   t1170 = -0.4e1 * t133 * t1119 + 0.16e2 * t1123 * t1124 * t730 - 0.8e1 * t1122 * t1128 * t196 - 0.64e2 * t1134 * t1135 * t1017 - 0.32e2 * t1140 * t1141 * t1135 + 0.16e2 * t1145 * t1146 - 0.8e1 * t1150 * t650 - 0.16e2 * t1123 * t1154 - 0.4e1 * t133 * t1158 - 0.16e2 * t1140 * t1162 * t15 + 0.8e1 * t1166 * t35 * t312;
1628   t1171 = t1161 * t184;
1629   t1175 = t1122 * n;
1630   t1176 = t15 * t1161;
1631   t1180 = t1132 * ZA;
1632   t1181 = t1180 * t355;
1633   t1182 = t1176 * t1017;
1634   t1185 = t1161 * xc;
1635   t1189 = t1133 * t1;
1636   t1192 = t108 * t1;
1637   t1193 = t1132 * t1192;
1638   t1195 = t10 * t35;
1639   t1199 = t157 * t15;
1640   t1203 = t1141 * t46;
1641   t1211 = t184 * t108;
1642   t1218 = 0.2e1 * t133 * t1171 * t35 + 0.8e1 * t1175 * t1176 * t997 + 0.64e2 * t1181 * t1182 - 0.8e1 * t1175 * t1185 * t127 - 0.32e2 * t1189 * t1182 - 0.64e2 * t1193 * ZA * t1195 * t228 + 0.8e1 * t1199 * t416 * t1128 + 0.8e1 * t1140 * t1203 - 0.4e1 * t75 * t1158 - 0.8e1 * t1199 * t560 * t1128 - 0.2e1 * t133 * t1211 - 0.8e1 * t1199 * t127 * t11 * t1128;
1643   t1221 = t256 * t1161;
1644   t1224 = t35 * t108;
1645   t1233 = t7 * t256;
1646   t1236 = -t75 * t1211 - t75 * t1221 - t133 * t1221 + t75 * t1224 - t75 * t1171 - t133 * t1171 + t133 * t1224 + t75 * t1162 - t75 * t108 * t256 + t133 * t1162 - t1233 * t59 * t108;
1647   t1240 = t1135 * t1195;
1648   t1252 = t629 * t127;
1649   t1263 = t1171 * ZA;
1650   t1280 = -0.128e3 * t1180 * ZB * t108 * t1240 + 0.32e2 * t1193 * t10 * t112 + 0.4e1 * t133 * t1203 + 0.4e1 * t109 * t256 * ZA * ZB - 0.8e1 * t1252 * n * t15 * t1185 + 0.8e1 * t1175 * t1171 * t730 - 0.8e1 * t1175 * t1176 * t127 + 0.4e1 * t223 * t1263 - 0.8e1 * t1175 * t1176 * t730 + 0.8e1 * t1166 * ZA * t341 * ZB + 0.64e2 * t1134 * t1240 + 0.8e1 * t1122 * xc * t1128 * t127 * t15;
1651   t1283 = t1199 * t19;
1652   t1287 = t1199 * t127;
1653   t1289 = t59 * n * t1161;
1654   t1293 = t157 * n * xc;
1655   t1304 = t1132 * t108;
1656   t1310 = t263 * ZB;
1657   t1316 = t2 * t15;
1658   t1323 = -0.16e2 * t1283 * t1006 * t46 + 0.8e1 * t1287 * t1289 + 0.8e1 * t1293 * t127 * t1161 * t59 + 0.16e2 * t1123 * t1135 * t19 + 0.8e1 * t1293 * t560 * t1176 + 0.64e2 * t1304 * t59 * t1240 + 0.4e1 * t75 * t1203 + 0.4e1 * t1310 * t1263 + 0.4e1 * t223 * t338 * t108 - 0.16e2 * t1252 * t1316 * t1153 - 0.16e2 * t1310 * t221 * t15;
1659   t1330 = t1132 * t15;
1660   t1336 = t1132 * t1;
1661   t1338 = t1162 * t179;
1662   t1370 = 0.8e1 * t1175 * t1176 * t19 + 0.4e1 * t1139 * t318 * t1161 + 0.128e3 * t1330 * t318 * t108 * t46 * t227 - 0.32e2 * t1336 * xc * t1338 + 0.4e1 * t1233 * ZA * t1161 * ZB - 0.8e1 * t1287 * t59 * t1128 * xc + 0.2e1 * t75 * t305 * t108 + 0.8e1 * t1199 * t127 * t59 * t1128 - 0.8e1 * t1283 * t1289 - 0.8e1 * t1293 * t560 * t1171 + 0.4e1 * t133 * t35 * t1158 + 0.8e1 * t157 * t184 * t19 * t11 * t1128 * xc;
1663   t1376 = t7 * t184;
1664   t1380 = t1176 * t1195;
1665   t1393 = t1330 * t35;
1666   t1411 = 0.16e2 * t1145 * t1154 + 0.8e1 * t1149 * t998 + 0.4e1 * t1376 * t35 * t48 + 0.32e2 * t1189 * t1380 + 0.32e2 * t1193 * t11 * t1195 * t15 - 0.64e2 * t1304 * xc * t111 * t1135 - 0.16e2 * t1123 * t1146 + 0.64e2 * t1393 * t28 * t1192 * xc - 0.16e2 * t1123 * t1135 * t127 - 0.8e1 * t1122 * xc * t1128 * t127 - 0.32e2 * t1193 * xc * t112 + 0.16e2 * t1252 * t1316 * t46;
1667   t1450 = 0.2e1 * t1376 * t767 * t1161 + 0.2e1 * t1376 * t111 * t108 + 0.4e1 * t223 * t311 * t108 + 0.4e1 * t109 * t35 * t520 * ZA + 0.16e2 * t1123 * t1135 * t997 - 0.64e2 * t1181 * t1380 + 0.8e1 * t1150 * t903 - 0.32e2 * t1393 * t11 * t1192 * xc - 0.16e2 * t157 * t2 * xc * t560 * t1124 + 0.8e1 * t223 * t184 * t317 * t46 + 0.32e2 * t1336 * t10 * t1338 - 0.4e1 * t75 * t1119;
1668   _PC3B = (t606 + t722 + t1089 + t781 + 0.16e2 * t48 * t51 + t978 + t868 + t507 - t304 * t256 + 0.8e1 * t9 * t22 + t752 + 0.4e1 * t174 * t144 - 0.2e1 * t81 * t469 + 0.6e1 * t139 * t166 + t362 + 0.2e1 * t98 * t211 + t925 + t137 - t290 * t184 + 0.12e2 * t81 * t83 + t842 + 0.8e1 * t74 * t77 + 0.16e2 * t98 * t12 * t162 - 0.4e1 * t33 * t28 * t443 - 0.8e1 * t27 * t70 - 0.2e1 * t33 * t34 * t36 - 0.8e1 * t27 * t30 + 0.2e1 * t58 * t67 - 0.4e1 * t40 * t43 + 0.2e1 * t58 * t63 + t1033 - t290 * t256 + t290 * t35 + t193 + t1113 + t578 + t442 + t474 + t544 + t329 + t679 + t401 + t953 + t811 + t644 + t894 + t289 + t240 + t1055 + t1003) / (t1170 + t1218 + 0.2e1 * t1236 + t1280 + t1323 + t1370 + t1411 + t1450);
1669 
1670   t1   = n * n;
1671   t2   = t1 * xc;
1672   t3   = ZB * ZB;
1673   t5   = t2 * 0.3141592654e1 * t3;
1674   t6   = nx * 0.3141592654e1;
1675   t7   = t6 * xc;
1676   t8   = PetscCosReal(t7);
1677   t9   = nx * nx;
1678   t10  = t8 * t9;
1679   t11  = n * 0.3141592654e1;
1680   t12  = PetscExpReal(t11);
1681   t13  = t12 * t12;
1682   t16  = PetscExpReal(xc * n * 0.3141592654e1);
1683   t17  = t16 * t16;
1684   t18  = t17 * t16;
1685   t19  = t17 * t17;
1686   t20  = t19 * t18;
1687   t21  = t13 * t20;
1688   t22  = t10 * t21;
1689   t25  = ZA * ZA;
1690   t26  = t1 * t25;
1691   t27  = xc * 0.3141592654e1;
1692   t28  = t26 * t27;
1693   t29  = t19 * t16;
1694   t30  = t13 * t13;
1695   t31  = t29 * t30;
1696   t35  = t9 * nx;
1697   t36  = t3 * t35;
1698   t37  = PetscSinReal(t6);
1699   t38  = t13 * t12;
1700   t39  = t37 * t38;
1701   t40  = t39 * t19;
1702   t42  = t1 * t1;
1703   t43  = nx * t42;
1704   t44  = xc * xc;
1705   t45  = t25 * t44;
1706   t46  = t43 * t45;
1707   t47  = 0.3141592654e1 * 0.3141592654e1;
1708   t48  = t47 * t37;
1709   t49  = t17 * t38;
1710   t54  = 0.3141592654e1 * t35;
1711   t55  = ZA * n * t54;
1712   t56  = t37 * ZB;
1713   t57  = t19 * t12;
1714   t61  = t25 * t8;
1715   t62  = t61 * t9;
1716   t63  = n * t30;
1717   t64  = t63 * t16;
1718   t67  = t1 * n;
1719   t69  = t47 * ZB;
1720   t70  = t67 * t44 * t69;
1721   t75  = nx * t3;
1722   t76  = t75 * t37;
1723   t77  = t67 * 0.3141592654e1;
1724   t78  = t19 * t19;
1725   t79  = t78 * t12;
1726   t80  = t77 * t79;
1727   t82  = t3 * t38;
1728   t84  = t54 * t37;
1729   t87  = PetscSinReal(t7);
1730   t88  = t29 * t87;
1731   t89  = t47 * t44;
1732   t93  = nx * t25;
1733   t94  = t87 * t42;
1734   t95  = t93 * t94;
1735   t96  = t47 * xc;
1736   t97  = t13 * t29;
1737   t98  = t96 * t97;
1738   t101 = t87 * t67;
1739   t102 = t93 * t101;
1740   t103 = t13 * t18;
1741   t107 = t47 * t35;
1742   t108 = t26 * t107;
1743   t109 = t37 * t44;
1744   t110 = t19 * t17;
1745   t111 = t12 * t110;
1746   t116 = t37 * t19 * t12;
1747   t118 = t37 * xc;
1748   t119 = ZB * t19;
1749   t120 = t119 * t12;
1750   t121 = t118 * t120;
1751   t125 = xc * t3;
1752   t126 = t1 * t47 * t125;
1753   t127 = t35 * t37;
1754   t128 = t38 * t19;
1755   t129 = t127 * t128;
1756   t132 = t26 * 0.3141592654e1;
1757   t133 = t16 * t13;
1758   t134 = t10 * t133;
1759   t137 = 0.3141592654e1 * ZB;
1760   t138 = t2 * t137;
1761   t139 = ZA * t8;
1762   t140 = t9 * t13;
1763   t145 = t30 * t18;
1764   t146 = t10 * t145;
1765   t149 = t3 * t8;
1766   t150 = t149 * t9;
1767   t153 = 0.2e1 * t5 * t22 + 0.2e1 * t28 * t10 * t31 + t36 * t40 - 0.2e1 * t46 * t48 * t49 - 0.2e1 * t55 * t56 * t57 - 0.2e1 * t62 * t64 + 0.16e2 * t70 * t29 * ZA * t10 - t76 * t80 + t82 * n * t84 + 0.8e1 * t43 * t3 * t88 * t89 + 0.16e2 * t95 * t98 + 0.2e1 * t102 * t27 * t103 - 0.2e1 * t108 * t109 * t111 + t36 * t116 + 0.8e1 * t55 * t121 - 0.4e1 * t126 * t129 - 0.4e1 * t132 * t134 - 0.4e1 * t138 * t139 * t140 * t20 + 0.8e1 * t28 * t146 - 0.2e1 * t150 * t64;
1768   t154 = t42 * n;
1769   t155 = nx * t154;
1770   t156 = t44 * xc;
1771   t157 = t47 * 0.3141592654e1;
1772   t158 = t156 * t157;
1773   t159 = t155 * t158;
1774   t162 = t56 * ZA * t19 * t12;
1775   t165 = t77 * t49;
1776   t167 = t1 * t3;
1777   t168 = t167 * t89;
1778   t169 = t127 * t49;
1779   t172 = t37 * t67;
1780   t173 = t75 * t172;
1781   t174 = t38 * t110;
1782   t175 = t27 * t174;
1783   t179 = t47 * t25;
1784   t181 = t10 * t97;
1785   t184 = t27 * t31;
1786   t187 = t67 * t47;
1787   t188 = t44 * t3;
1788   t189 = t187 * t188;
1789   t192 = t25 * t35;
1790   t193 = t37 * t17;
1791   t194 = t193 * t12;
1792   t196 = nx * ZA;
1793   t197 = t196 * t172;
1794   t198 = ZB * t38;
1795   t199 = t198 * t19;
1796   t204 = t1 * t12 * t110;
1797   t207 = nx * ZB;
1798   t209 = t1 * ZA;
1799   t215 = t67 * t3;
1800   t216 = t47 * t8;
1801   t217 = t215 * t216;
1802   t218 = t9 * xc;
1803   t222 = nx * t67;
1804   t223 = t222 * t27;
1805   t224 = t3 * t87;
1806   t228 = t167 * t107;
1807   t232 = t26 * t96;
1808   t235 = t207 * t94;
1809   t236 = t47 * ZA;
1810   t243 = xc * t13;
1811   t244 = t243 * t29;
1812   t248 = t25 * n;
1813   t249 = t248 * 0.3141592654e1;
1814   t253 = ZB * ZA;
1815   t254 = t253 * t8;
1816   t255 = t9 * n;
1817   t256 = t30 * t16;
1818   t260 = 0.2e1 * t207 * t37 * t209 * t128 + 0.2e1 * t5 * t134 - 0.16e2 * t217 * t218 * t97 - 0.2e1 * t223 * t224 * t31 - 0.2e1 * t228 * t109 * t174 - 0.2e1 * t232 * t169 - 0.16e2 * t235 * t236 * t44 * t30 * t18 - 0.4e1 * t196 * t101 * t137 * t244 + t249 * t169 + 0.8e1 * t168 * t129 + 0.4e1 * t254 * t255 * t256;
1819   t263 = t43 * t179;
1820   t267 = t3 * n;
1821   t268 = t267 * t54;
1822   t269 = t118 * t57;
1823   t272 = t39 * t1;
1824   t274 = t67 * t25;
1825   t275 = t274 * t158;
1826   t278 = t75 * t87;
1827   t279 = t77 * t103;
1828   t282 = t25 * t38;
1829   t285 = ZA * t38;
1830   t290 = t267 * 0.3141592654e1;
1831   t296 = t77 * t111;
1832   t298 = t196 * t37;
1833   t299 = t1 * ZB;
1834   t303 = t37 * t42;
1835   t304 = t196 * t303;
1836   t308 = t77 * t57;
1837   t310 = t26 * t89;
1838   t313 = t77 * t128;
1839   t316 = t101 * t27;
1840   t319 = t93 * t87;
1841   t320 = t77 * t97;
1842   t323 = t127 * t57;
1843   t326 = t10 * n;
1844   t329 = t118 * t174;
1845   t332 = -0.8e1 * t263 * t109 * t57 - 0.4e1 * t268 * t269 + t93 * t272 + 0.8e1 * t275 * t129 - 0.4e1 * t278 * t279 + t282 * n * t84 - 0.2e1 * t285 * n * t54 * t56 - t290 * t169 - 0.2e1 * t196 * t38 * t172 * t137 + t76 * t296 - 0.2e1 * t298 * t299 * t79 + 0.8e1 * t304 * t96 * t120 + t76 * t308 - 0.2e1 * t310 * t169 - t76 * t313 + 0.2e1 * t75 * t18 * t316 + 0.4e1 * t319 * t320 + t249 * t323 - 0.2e1 * t25 * t18 * t326 + 0.2e1 * t228 * t329;
1846   t335 = t75 * t101;
1847   t336 = t27 * t21;
1848   t342 = t77 * t133;
1849   t347 = t209 * t137;
1850   t350 = t9 * t1;
1851   t351 = t149 * t350;
1852   t355 = t37 * t78 * t12;
1853   t359 = t93 * t303;
1854   t367 = t172 * 0.3141592654e1;
1855   t369 = t96 * t103;
1856   t376 = t209 * t107;
1857   t379 = t10 * t103;
1858   t383 = t207 * t101;
1859   t389 = 0.3141592654e1 * ZA;
1860   t390 = t222 * t389;
1861   t391 = t87 * ZB;
1862   t398 = -0.2e1 * t102 * t336 + t93 * t38 * t367 + 0.16e2 * t95 * t369 - t82 * t127 - 0.8e1 * t197 * t27 * t120 + 0.8e1 * t376 * t121 - 0.8e1 * t189 * t379 - t249 * t129 - 0.4e1 * t383 * t27 * ZA * t16 * t13 - 0.8e1 * t390 * t391 * t21 - 0.2e1 * t197 * t137 * t57;
1863   t402 = t39 * t110;
1864   t404 = t193 * t38;
1865   t406 = t127 * t174;
1866   t408 = t167 * 0.3141592654e1;
1867   t411 = t44 * t157;
1868   t412 = t155 * t411;
1869   t413 = t285 * t19;
1870   t414 = t56 * t413;
1871   t417 = ZA * t30;
1872   t424 = t93 * t37;
1873   t426 = t248 * t54;
1874   t427 = t17 * t12;
1875   t428 = t118 * t427;
1876   t431 = t77 * t21;
1877   t438 = ZA * t13;
1878   t443 = t93 * t172;
1879   t444 = t27 * t427;
1880   t448 = t1 * t78 * t12;
1881   t455 = t274 * t89;
1882   t461 = t118 * t111;
1883   t464 = -t36 * t402 + t36 * t404 - t249 * t406 - 0.4e1 * t408 * t134 + 0.16e2 * t412 * t414 - 0.4e1 * t383 * t27 * t417 * t18 + 0.2e1 * t28 * t22 - t424 * t80 - 0.2e1 * t426 * t428 + 0.4e1 * t278 * t431 + 0.4e1 * t254 * t255 * t103 + t290 * t323 + 0.4e1 * t383 * t27 * t438 * t20 + 0.2e1 * t443 * t444 + t424 * t448 - t36 * t194 - 0.32e2 * t235 * t236 * t243 * t18 + 0.8e1 * t455 * t181 - 0.4e1 * t359 * t96 * t128 - 0.2e1 * t426 * t461;
1884   t469 = n * t16 * t13;
1885   t474 = t1 * t38;
1886   t475 = t474 * t19;
1887   t480 = t89 * t103;
1888   t483 = t67 * ZA;
1889   t484 = t483 * t411;
1890   t485 = t127 * t120;
1891   t488 = t127 * t111;
1892   t497 = t77 * t427;
1893   t502 = t27 * t97;
1894   t508 = t1 * t19 * t12;
1895   t511 = t155 * t25 * t156;
1896   t512 = t157 * t37;
1897   t513 = t512 * t128;
1898   t527 = t1 * t17;
1899   t528 = t527 * t38;
1900   t530 = -t76 * t497 - 0.4e1 * t254 * t255 * t97 - 0.2e1 * t102 * t502 - 0.4e1 * t108 * t269 - t76 * t508 + 0.8e1 * t511 * t513 + 0.4e1 * t150 * t63 * t18 + 0.4e1 * t383 * t27 * t438 * t18 + 0.4e1 * t132 * t379 + 0.2e1 * t168 * t488 - t76 * t528;
1901   t535 = t44 * t13;
1902   t542 = t527 * t12;
1903   t544 = n * t13;
1904   t545 = t544 * t20;
1905   t548 = t75 * t303;
1906   t549 = t96 * t111;
1907   t552 = ZA * t35;
1908   t553 = t552 * t37;
1909   t562 = t43 * t96;
1910   t563 = t3 * t37;
1911   t564 = t563 * t128;
1912   t579 = t474 * t110;
1913   t590 = t9 * t30;
1914   t591 = t590 * t18;
1915   t595 = t127 * t427;
1916   t598 = t77 * t174;
1917   t600 = 0.4e1 * t5 * t146 + 0.16e2 * t235 * t236 * t535 * t18 + 0.8e1 * t455 * t146 + t76 * t542 - 0.2e1 * t150 * t545 + 0.2e1 * t548 * t549 - 0.2e1 * t553 * t120 + t290 * t488 - 0.8e1 * t274 * t47 * t44 * t29 * t10 - 0.4e1 * t562 * t564 - 0.2e1 * t132 * xc * t20 * t10 - 0.32e2 * t562 * ZA * t87 * ZB * t13 * t29 - 0.8e1 * t347 * t379 + t76 * t579 - 0.4e1 * t359 * t96 * t57 + 0.4e1 * t408 * t181 - 0.4e1 * t223 * t564 - 0.12e2 * t209 * t27 * ZB * t8 * t591 + 0.2e1 * t310 * t595 + t76 * t598;
1918   t601 = t27 * t49;
1919   t604 = t127 * t79;
1920   t606 = ZB * t29;
1921   t616 = t139 * t140 * t18;
1922   t638 = t10 * t256;
1923   t643 = t118 * t199;
1924   t653 = t544 * t29;
1925   t658 = t3 * t29;
1926   t660 = t350 * t27;
1927   t663 = -0.4e1 * t254 * t255 * t145 + 0.2e1 * t267 * t20 * t8 * t9 - 0.4e1 * t138 * t139 * t9 * t16 * t13 - 0.2e1 * t5 * t638 + 0.2e1 * t126 * t169 + 0.8e1 * t376 * t643 + 0.4e1 * t335 * t27 * t145 + 0.16e2 * t235 * t236 * t535 * t29 + 0.6e1 * t150 * t653 - 0.4e1 * t426 * t269 + 0.4e1 * t658 * t8 * t660;
1928   t670 = t274 * t411;
1929   t673 = t118 * t49;
1930   t694 = t155 * t45;
1931   t713 = n * t29 * t30;
1932   t717 = t20 * t87;
1933   t723 = t512 * t57;
1934   t728 = -0.2e1 * t443 * t175 - 0.8e1 * t670 * t129 + 0.2e1 * t426 * t673 - 0.16e2 * t207 * t88 * t42 * t47 * ZA * t44 + 0.4e1 * t254 * t255 * t21 + t249 * t595 + 0.8e1 * t25 * t29 * t8 * t660 + 0.2e1 * t268 * t461 + 0.8e1 * t189 * t181 - 0.8e1 * t694 * t513 + 0.2e1 * t198 * t553 - 0.12e2 * t606 * t139 * t660 - 0.2e1 * t359 * t549 + 0.4e1 * t138 * t139 * t590 * t16 + 0.8e1 * t93 * t29 * t94 * t89 - 0.2e1 * t150 * t713 + 0.2e1 * t222 * t3 * t717 * t27 + 0.8e1 * t670 * t323 + 0.8e1 * t694 * t723 - 0.2e1 * t62 * t653;
1935   t734 = t43 * t89;
1936   t735 = t563 * t427;
1937   t740 = t75 * t94;
1938   t744 = ZB * xc;
1939   t750 = t563 * t57;
1940   t754 = t218 * t103;
1941   t771 = t127 * t199;
1942   t776 = t89 * t174;
1943   t791 = -0.4e1 * t207 * t717 * t77 * xc * ZA + 0.4e1 * t443 * t27 * t57 + t192 * t40 - 0.8e1 * t55 * t643 - 0.16e2 * t209 * t89 * t771 - 0.8e1 * t275 * t323 + 0.2e1 * t359 * t776 + 0.16e2 * t304 * t89 * t199 + 0.4e1 * t278 * t320 + 0.2e1 * t207 * t172 * t389 * t79 - 0.8e1 * t390 * t391 * t97;
1944   t794 = t483 * t158;
1945   t801 = t2 * 0.3141592654e1;
1946   t818 = t215 * t411;
1947   t827 = t96 * t174;
1948   t837 = t37 * t12 * t110;
1949   t845 = 0.16e2 * t794 * t485 + 0.8e1 * t159 * t564 - 0.8e1 * t455 * t379 - 0.2e1 * t801 * t3 * t20 * t10 - 0.4e1 * t132 * t22 - 0.8e1 * t734 * t564 - 0.8e1 * t187 * t44 * t658 * t10 - 0.8e1 * t412 * t564 + 0.4e1 * t132 * t181 - 0.8e1 * t818 * t129 + 0.2e1 * t46 * t48 * t427 - 0.4e1 * t75 * t29 * t316 - 0.2e1 * t359 * t827 - t290 * t595 + 0.16e2 * t217 * t754 - t424 * t542 - 0.8e1 * t734 * t750 - t192 * t837 - 0.4e1 * t254 * t255 * t133 + 0.8e1 * t304 * t96 * t199;
1950   t864 = t544 * t18;
1951   t867 = t3 * t18;
1952   t884 = t27 * t256;
1953   t891 = t187 * t744;
1954   t894 = t563 * t49;
1955   t900 = -0.2e1 * t263 * t428 + 0.2e1 * t228 * t428 - 0.6e1 * t223 * t224 * t103 - t192 * t404 + 0.2e1 * t268 * t428 - 0.2e1 * t335 * t884 - t424 * t296 + 0.2e1 * t93 * t20 * t316 - 0.32e2 * t891 * t616 + 0.2e1 * t562 * t894 - 0.2e1 * t801 * t867 * t10;
1956   t904 = t27 * t111;
1957   t907 = t118 * t128;
1958   t915 = t89 * t145;
1959   t947 = t139 * t140 * t29;
1960   t952 = -0.2e1 * t173 * t904 + 0.4e1 * t426 * t907 + 0.12e2 * t253 * t10 * t1 * 0.3141592654e1 * t244 + 0.8e1 * t95 * t915 - t36 * t355 - 0.16e2 * t794 * t771 - 0.8e1 * t511 * t723 + 0.16e2 * t734 * t162 + t36 * t837 + 0.2e1 * t298 * t299 * t57 - 0.2e1 * t28 * t638 - 0.2e1 * t62 * t545 + 0.2e1 * t310 * t406 + 0.12e2 * t138 * t616 + 0.4e1 * t223 * t750 + t424 * t497 + 0.2e1 * t734 * t894 + 0.2e1 * t132 * xc * t18 * t10 - 0.16e2 * t70 * t947 + 0.32e2 * t891 * t947;
1961   t969  = t67 * t157 * t156 * t3;
1962   t974  = t27 * t133;
1963   t1001 = -0.8e1 * t159 * t750 - 0.16e2 * t412 * t162 - t290 * t129 + 0.8e1 * t310 * t323 - 0.4e1 * t319 * t342 + t75 * t272 + t192 * t402 - 0.8e1 * t359 * t89 * t128 - 0.10e2 * t61 * t350 * t502 + 0.8e1 * t818 * t323 - 0.4e1 * t108 * t907;
1964   t1042 = t89 * t97;
1965   t1055 = -0.2e1 * t168 * t595 + 0.16e2 * t484 * t771 + 0.4e1 * t11 * t125 * t129 - 0.2e1 * t173 * t444 + 0.2e1 * ZB * n * t54 * t37 * ZA * t79 - t424 * t475 + 0.2e1 * t562 * t735 - 0.2e1 * t548 * t776 + t424 * t204 + 0.2e1 * t25 * t20 * t326 + 0.8e1 * t383 * t389 * t133 + t75 * t38 * t367 + 0.2e1 * t62 * t469 + 0.2e1 * t197 * t137 * t128 - 0.2e1 * t102 * t884 - 0.2e1 * t5 * t379 - 0.8e1 * t740 * t1042 - 0.16e2 * t159 * t414 - 0.2e1 * ZB * t35 * t37 * t413 + 0.2e1 * t553 * ZB * t78 * t12;
1966   t1096 = 0.2e1 * t443 * t904 - 0.2e1 * t268 * t329 - 0.2e1 * t443 * t601 + 0.2e1 * t102 * t974 - 0.2e1 * t263 * t673 + t424 * t165 + 0.2e1 * t62 * t713 + t424 * t308 - t424 * t313 + 0.8e1 * t347 * t22 - t424 * t598;
1967   t1103 = t42 * t1 * t157;
1968   t1104 = t1103 * t25;
1969   t1108 = t3 * t19;
1970   t1112 = n * t47;
1971   t1113 = t9 * t9;
1972   t1118 = t42 * t157;
1973   t1119 = t1118 * t9;
1974   t1120 = t25 * xc;
1975   t1121 = t13 * t110;
1976   t1122 = t1120 * t1121;
1977   t1125 = t47 * t47;
1978   t1126 = t67 * t1125;
1979   t1127 = t1113 * ZA;
1980   t1128 = t1126 * t1127;
1981   t1129 = t19 * t13;
1982   t1130 = t744 * t1129;
1983   t1133 = t154 * t1125;
1984   t1134 = t1133 * t9;
1985   t1135 = t45 * t1129;
1986   t1138 = t154 * t47;
1987   t1139 = t25 * t30;
1988   t1142 = t1126 * t1113;
1989   t1145 = t125 * t1129;
1990   t1148 = t1103 * xc;
1991   t1149 = t3 * t13;
1992   t1150 = t1149 * t17;
1993   t1153 = t25 * t78;
1994   t1156 = -0.8e1 * t1104 * t243 * t17 + 0.4e1 * t187 * t1108 * t9 - 0.2e1 * t1112 * t3 * t1113 * t30 + 0.16e2 * t1119 * t1122 + 0.64e2 * t1128 * t1130 + 0.64e2 * t1134 * t1135 - 0.2e1 * t1138 * t1139 + 0.32e2 * t1142 * t1135 - 0.32e2 * t1142 * t1145 + 0.8e1 * t1148 * t1150 - 0.2e1 * t1138 * t1153;
1995   t1157 = t25 * t13;
1996   t1158 = t1157 * t17;
1997   t1161 = t13 * t17;
1998   t1162 = t1120 * t1161;
1999   t1165 = t3 * t78;
2000   t1170 = t42 * t67 * t1125;
2001   t1172 = t1108 * t13;
2002   t1175 = t1 * t157;
2003   t1176 = t1175 * t1113;
2004   t1182 = t1120 * t1129;
2005   t1189 = t110 * t9 * xc;
2006   t1192 = t1149 * t110;
2007   t1201 = 0.8e1 * t1103 * t1158 - 0.16e2 * t1119 * t1162 - 0.2e1 * t1112 * t1165 * t1113 + 0.32e2 * t1170 * t44 * t1172 - 0.8e1 * t1176 * t1162 + 0.8e1 * t1104 * t243 * t110 - 0.64e2 * t1134 * t1182 - 0.64e2 * t1134 * t1145 + 0.16e2 * t1118 * t3 * t1189 + 0.16e2 * t1119 * t1192 - 0.4e1 * t187 * t1165 * t9 - 0.4e1 * t187 * t1139 * t9;
2008   t1209 = t17 * t30;
2009   t1210 = t125 * t1209;
2010   t1213 = t1138 * ZA;
2011   t1214 = ZB * t30;
2012   t1218 = t1157 * t110;
2013   t1226 = t3 * t30;
2014   t1237 = t1170 * t25;
2015   t1242 = 0.4e1 * t1112 * ZA * t119 * t1113 - 0.16e2 * t1119 * t1150 - 0.8e1 * t1176 * t1210 + 0.4e1 * t1213 * t1214 * t19 - 0.16e2 * t1119 * t1218 - 0.32e2 * t1142 * t1182 - 0.8e1 * t1103 * t1120 * t110 - 0.4e1 * t187 * t1226 * t9 + 0.8e1 * t1103 * t1192 + 0.4e1 * t1112 * ZB * t1113 * t30 * ZA - 0.32e2 * t1237 * xc * t19 * t13;
2016   t1251 = t125 * t1121;
2017   t1260 = t1120 * t1209;
2018   t1263 = t1139 * t19;
2019   t1282 = 0.8e1 * t1103 * t110 * t3 * xc + 0.8e1 * t1104 * xc * t17 * t30 - 0.8e1 * t1176 * t1251 + 0.16e2 * t1119 * t1158 + 0.4e1 * t1112 * t78 * t1127 * ZB + 0.16e2 * t1119 * t1260 + 0.2e1 * t1138 * t1263 - 0.32e2 * t1170 * xc * t1172 - 0.16e2 * t1213 * t119 * t13 + 0.4e1 * t1138 * t1214 * ZA + 0.32e2 * t1237 * t44 * t19 * t13 - 0.16e2 * t1118 * t25 * t1189;
2020   t1287 = t188 * t1129;
2021   t1292 = t25 * t19;
2022   t1296 = t187 * t9;
2023   t1297 = t1226 * t19;
2024   t1311 = t1112 * t1113;
2025   t1317 = -0.8e1 * t1176 * t1150 + 0.32e2 * t1142 * t1287 - 0.8e1 * t1103 * t1150 + 0.2e1 * t1112 * t1292 * t1113 + 0.4e1 * t1296 * t1297 + 0.8e1 * t1176 * t1192 + 0.4e1 * t1296 * t1263 + 0.8e1 * t1176 * t1158 - 0.8e1 * t1175 * t25 * t1113 * xc * t110 + 0.2e1 * t1311 * t1297 + 0.2e1 * t1112 * t1108 * t1113;
2026   t1320 = t253 * t1129;
2027   t1328 = t253 * t30 * t19;
2028   t1333 = t125 * t1161;
2029   t1343 = ZB * t44 * t1129;
2030   t1350 = -0.8e1 * t1176 * t1218 - 0.16e2 * t1311 * t1320 + 0.8e1 * t1176 * t1260 - 0.16e2 * t1119 * t1210 + 0.4e1 * t1311 * t1328 + 0.2e1 * t1311 * t1263 + 0.8e1 * t1176 * t1333 + 0.8e1 * t187 * ZB * t417 * t9 - 0.2e1 * t1138 * t1165 - 0.64e2 * t1128 * t1343 + 0.64e2 * t1134 * t1287 + 0.2e1 * t1138 * t1108;
2031   t1369 = t1133 * t9 * ZA;
2032   t1378 = t187 * ZA;
2033   t1383 = t1170 * ZA;
2034   t1388 = 0.2e1 * t1138 * t1297 - 0.8e1 * t1148 * t1192 + 0.2e1 * t1138 * t1292 - 0.16e2 * t1119 * t1251 + 0.8e1 * t1175 * xc * t110 * t1113 * t3 - 0.2e1 * t1112 * t1153 * t1113 + 0.128e3 * t1369 * t1130 + 0.16e2 * t1119 * t1333 + 0.4e1 * t1138 * t78 * ZA * ZB + 0.8e1 * t1378 * t78 * t9 * ZB - 0.64e2 * t1383 * t1343 + 0.64e2 * t1383 * t1130;
2035   t1420 = 0.4e1 * t1138 * t119 * ZA - 0.128e3 * t1369 * t1343 - 0.4e1 * t187 * t1153 * t9 - 0.2e1 * t1138 * t1226 + 0.8e1 * t1296 * t1328 - 0.2e1 * t1112 * t1139 * t1113 - 0.8e1 * t1148 * t3 * t17 * t30 - 0.32e2 * t1296 * t1320 + 0.8e1 * t1176 * t1122 + 0.4e1 * t187 * t1292 * t9 + 0.8e1 * t1378 * t119 * t9 - 0.8e1 * t1103 * t1218;
2036 
2037   _PC4B = (-t424 * t508 + 0.8e1 * t412 * t750 - 0.2e1 * t232 * t595 - 0.4e1 * t126 * t323 + t1096 - t76 * t204 + t728 + 0.2e1 * t548 * t827 + 0.2e1 * t150 * t469 + t398 + 0.8e1 * t189 * t146 + t260 - 0.2e1 * t351 * t184 - 0.2e1 * t268 * t673 - 0.4e1 * t319 * t279 + t464 - 0.2e1 * t108 * t461 + 0.16e2 * t740 * t369 + 0.16e2 * t274 * t216 * t754 - 0.16e2 * t70 * t139 * t591 + 0.2e1 * t55 * t56 * t128 - 0.2e1 * t359 * t89 * t111 + 0.2e1 * t734 * t563 * t111 + 0.6e1 * t223 * t224 * t97 + 0.8e1 * t383 * t389 * t103 + 0.4e1 * t606 * ZA * t326 - 0.2e1 * t93 * t18 * t316 - 0.4e1 * t443 * t27 * t128 + 0.8e1 * t197 * t27 * t199 + 0.8e1 * t108 * t109 * t128 - t249 * t604 + 0.16e2 * t70 * t616 - 0.8e1 * t969 * t323 + t845 - t424 * t579 + 0.16e2 * t159 * t162 + t290 * t406 - 0.6e1 * t150 * t864 + t192 * t116 + 0.2e1 * t867 * t326 - 0.4e1 * t658 * t326 - 0.2e1 * t351 * t502 - t76 * t165 + t900 + 0.8e1 * t168 * t323 + t791 + 0.8e1 * t740 * t915 - 0.4e1 * t562 * t750 - 0.4e1 * t278 * t342 + 0.4e1 * t319 * t431 + 0.2e1 * t173 * t175 + t424 * t528 + 0.8e1 * t969 * t129 - 0.8e1 * t347 * t181 + t332 + t530 - 0.2e1 * t108 * t329 - 0.2e1 * t207 * t38 * t37 * t1 * ZA + t1001 + 0.4e1 * t408 * t379 + t76 * t448 + 0.2e1 * t102 * t184 + 0.2e1 * t426 * t329 + 0.16e2 * t740 * t98 - t282 * t127 - 0.16e2 * t1 * t44 * t69 * t552 * t116 + 0.2e1 * t168 * t169 + 0.2e1 * t28 * t134 - t290 * t604 - 0.16e2 * t484 * t485 - 0.8e1 * t740 * t480 + 0.2e1 * t173 * t601 - 0.2e1 * t335 * t336 + t600 + 0.2e1 * t62 * t864 + t952 + 0.8e1 * t347 * t134 - t192 * t355 + t192 * t194 + 0.2e1 * t228 * t461 + t663 + 0.4e1 * t383 * t27 * t417 * t16 + 0.4e1 * t138 * t20 * ZA * t10 - 0.4e1 * t20 * ZB * ZA * t326 + 0.4e1 * t196 * t88 * t77 * t744 - 0.16e2 * t67 * xc * t179 * t181 - 0.8e1 * t95 * t480 - t249 * t488 - t76 * t475 + t1055 - 0.4e1 * t408 * t22 - 0.10e2 * t28 * t379 + 0.2e1 * t335 * t974 + t153 - 0.8e1 * t95 * t1042 - 0.2e1 * t734 * t735) / (t1156 + t1201 + t1242 + t1282 + t1317 + t1350 + t1388 + t1420);
2038 
2039   if (x > xc) {
2040     _PC1 = _PC1B;
2041     _PC2 = _PC2B;
2042     _PC3 = _PC3B;
2043     _PC4 = _PC4B;
2044     Z    = ZB;
2045   } else {
2046     _PC1 = _PC1A;
2047     _PC2 = _PC2A;
2048     _PC3 = _PC3A;
2049     _PC4 = _PC4A;
2050     Z    = ZA;
2051   }
2052 
2053   t1  = n * n;
2054   t2  = t1 * t1;
2055   t3  = t2 * n;
2056   t4  = x * t3;
2057   t5  = 0.3141592654e1 * 0.3141592654e1;
2058   t6  = t5 * 0.3141592654e1;
2059   t11 = _PC3 * t6;
2060   t12 = x * n;
2061   t13 = nx * nx;
2062   t14 = t13 * t13;
2063   t15 = t12 * t14;
2064   t19 = PetscExpReal(t12 * 0.3141592654e1);
2065   t20 = t19 * t19;
2066   t21 = t4 * t20;
2067   t24 = _PC1 * t5;
2068   t25 = Z * t20;
2069   t29 = _PC1 * t6;
2070   t30 = t29 * Z;
2071   t31 = t1 * n;
2072   t32 = x * t31;
2073   t33 = t32 * t13;
2074   t36 = t11 * x;
2075   t41 = n * t20;
2076   t45 = t6 * _PC4;
2077   t49 = t20 * t1;
2078   t51 = _PC2 * Z;
2079   t55 = -0.2e1 * t4 * t6 * _PC2 * Z - 0.2e1 * t11 * t15 - 0.2e1 * t11 * t21 + 0.2e1 * t24 * t25 * t14 - t13 + 0.4e1 * t30 * t33 - 0.4e1 * t36 * t31 * t20 * t13 - 0.2e1 * t36 * t41 * t14 - 0.2e1 * t4 * t45 * t20 - t49 - 0.2e1 * t4 * t6 * t51 * t20;
2080   t58 = t32 * t6;
2081   t59 = _PC4 * t20;
2082   t63 = t20 * t13;
2083   t67 = t12 * t6;
2084   t68 = t20 * t14;
2085   t87 = t49 * t13;
2086   t90 = -0.4e1 * t11 * t33 - 0.4e1 * t58 * t59 * t13 - 0.4e1 * t58 * t51 * t63 - 0.2e1 * t67 * t51 * t68 + 0.4e1 * t32 * t45 * t13 - 0.2e1 * t67 * t59 * t14 - 0.2e1 * t30 * t21 + t1 + 0.2e1 * t24 * t25 * t2 + 0.2e1 * t12 * t45 * t14 + 0.4e1 * t24 * Z * t87;
2087   t106 = _PC3 * t5;
2088   t120 = -0.4e1 * t30 * t32 * t63 + t63 + 0.4e1 * t24 * Z * t1 * t13 + 0.2e1 * t29 * Z * x * t3 - 0.4e1 * t58 * t51 * t13 - 0.2e1 * t106 * t2 + t32 * 0.3141592654e1 - 0.2e1 * t106 * t14 - 0.2e1 * t30 * t12 * t68 - 0.2e1 * t67 * t51 * t14 + 0.4e1 * t106 * t87;
2089   t129 = PetscSinReal(nx * 0.3141592654e1 * x);
2090   t155 = 0.2e1 * t30 * t15 + x * 0.3141592654e1 * t41 * t13 - 0.4e1 * t19 * nx * t129 * n + t32 * 0.3141592654e1 * t20 + 0.2e1 * t106 * t68 + 0.2e1 * t106 * t20 * t2 - 0.4e1 * t106 * t1 * t13 - 0.2e1 * t11 * t4 + 0.2e1 * t4 * t45 + 0.2e1 * t24 * Z * t2 + 0.2e1 * t24 * Z * t14 + t12 * 0.3141592654e1 * t13;
2091   t158 = t5 * Z;
2092 
2093   u1 = (t55 + t90 + t120 + t155) / (0.4e1 * t158 * t19 * t2 + 0.8e1 * t158 * t19 * t1 * t13 + 0.4e1 * t158 * t19 * t14);
2094 
2095   t1  = n * n;
2096   t2  = t1 * n;
2097   t3  = x * t2;
2098   t4  = 0.3141592654e1 * 0.3141592654e1;
2099   t5  = t4 * 0.3141592654e1;
2100   t6  = t3 * t5;
2101   t7  = _PC2 * Z;
2102   t8  = nx * nx;
2103   t12 = t1 * t1;
2104   t13 = t12 * n;
2105   t14 = x * t13;
2106   t15 = t5 * _PC4;
2107   t16 = x * n;
2108   t18 = PetscExpReal(t16 * 0.3141592654e1);
2109   t19 = t18 * t18;
2110   t23 = t16 * t5;
2111   t24 = t8 * t8;
2112   t28 = _PC3 * t5;
2113   t29 = t14 * t19;
2114   t32 = _PC1 * t5;
2115   t33 = t32 * Z;
2116   t34 = t16 * t24;
2117   t37 = _PC4 * t19;
2118   t45 = _PC2 * t4;
2119   t53 = t19 * t8;
2120   t58 = _PC4 * t4;
2121   t60 = t1 * t19 * t8;
2122   t63 = t19 * t24;
2123   t67 = t3 * t8;
2124   t73 = n * t19;
2125   t86 = t28 * x;
2126   t91 = 0.4e1 * t58 * t60 + 0.2e1 * t33 * t16 * t63 + 0.4e1 * t33 * t67 + 0.2e1 * t33 * t29 - x * 0.3141592654e1 * t73 * t8 - 0.2e1 * t53 + 0.2e1 * t32 * Z * x * t13 - 0.2e1 * t58 * t12 - 0.2e1 * t58 * t24 + t3 * 0.3141592654e1 + 0.4e1 * t86 * t2 * t19 * t8;
2127   t94 = Z * t12;
2128   t121 = -0.2e1 * t8 + 0.2e1 * t45 * t94 * t19 + 0.2e1 * t14 * t5 * t7 * t19 + 0.4e1 * t6 * t7 * t53 + 0.2e1 * t23 * t7 * t63 - 0.4e1 * t28 * t67 + 0.2e1 * t45 * t94 + 0.2e1 * t58 * t12 * t19 + t16 * 0.3141592654e1 * t8 + 0.2e1 * t14 * t15 - 0.2e1 * t28 * t14;
2129   t146 = PetscCosReal(nx * 0.3141592654e1 * x);
2130   t156 = -t3 * 0.3141592654e1 * t19 + 0.2e1 * t58 * t63 - 0.4e1 * t58 * t1 * t8 + 0.4e1 * t45 * Z * t1 * t8 - 0.2e1 * t28 * t34 + 0.2e1 * t86 * t73 * t24 + 0.4e1 * t3 * t15 * t8 + 0.4e1 * t45 * Z * t60 + 0.4e1 * t18 * t146 * t8 + 0.2e1 * t45 * Z * t24 + 0.2e1 * t16 * t15 * t24;
2131   t159 = t4 * Z;
2132 
2133   u2 = (-0.4e1 * t6 * t7 * t8 + 0.2e1 * t14 * t15 * t19 - 0.2e1 * t23 * t7 * t24 + 0.2e1 * t28 * t29 + 0.2e1 * t33 * t34 + 0.4e1 * t6 * t37 * t8 - 0.2e1 * t14 * t5 * _PC2 * Z + 0.2e1 * t45 * Z * t19 * t24 + 0.2e1 * t23 * t37 * t24 + 0.4e1 * t33 * t3 * t53 + t91 + t121 + t156) / (0.4e1 * t159 * t18 * t12 + 0.8e1 * t159 * t18 * t1 * t8 + 0.4e1 * t159 * t18 * t24);
2134 
2135   t1  = 0.3141592654e1 * 0.3141592654e1;
2136   t2  = t1 * 0.3141592654e1;
2137   t3  = _PC1 * t2;
2138   t4  = t3 * Z;
2139   t5  = n * n;
2140   t6  = t5 * t5;
2141   t7  = t6 * n;
2142   t8  = x * t7;
2143   t9  = x * n;
2144   t11 = PetscExpReal(t9 * 0.3141592654e1);
2145   t12 = t11 * t11;
2146   t13 = t8 * t12;
2147   t16 = t5 * n;
2148   t17 = x * t16;
2149   t18 = t17 * t2;
2150   t19 = _PC4 * t12;
2151   t20 = nx * nx;
2152   t24 = t2 * _PC4;
2153   t28 = _PC3 * t2;
2154   t29 = t28 * x;
2155   t30 = t12 * n;
2156   t31 = t20 * t20;
2157   t40 = _PC2 * Z;
2158   t44 = t9 * t2;
2159   t48 = t12 * t20;
2160   t52 = t17 * t20;
2161   t57 = -0.2e1 * t4 * t13 - 0.4e1 * t18 * t19 * t20 - 0.2e1 * t8 * t24 * t12 - 0.2e1 * t29 * t30 * t31 + 0.2e1 * t8 * t2 * _PC2 * Z - 0.2e1 * t8 * t2 * t40 * t12 - 0.2e1 * t44 * t19 * t31 - 0.4e1 * t18 * t40 * t48 + t20 + 0.4e1 * t28 * t52 + t17 * 0.3141592654e1 * t12;
2162   t58 = t9 * t31;
2163   t61 = _PC3 * t1;
2164   t62 = t12 * t31;
2165   t73 = t5 * t20;
2166   t78 = _PC1 * t1;
2167   t90 = Z * t12;
2168   t94 = 0.2e1 * t28 * t58 + 0.2e1 * t61 * t62 + 0.2e1 * t61 * t12 * t6 - 0.4e1 * t4 * t17 * t48 + 0.2e1 * t28 * t8 + 0.4e1 * t61 * t73 - 0.2e1 * t8 * t24 - 0.2e1 * t78 * Z * t6 - 0.2e1 * t44 * t40 * t62 - 0.2e1 * t78 * Z * t31 - t9 * 0.3141592654e1 * t20 + 0.2e1 * t78 * t90 * t6;
2169   t101 = PetscCosReal(nx * 0.3141592654e1 * x);
2170   t102 = t11 * t101;
2171   t109 = t12 * t5;
2172   t110 = t109 * t20;
2173   t128 = 0.2e1 * t61 * t6 - t17 * 0.3141592654e1 + 0.2e1 * t102 * t5 - 0.4e1 * t17 * t24 * t20 + 0.4e1 * t78 * Z * t110 - 0.2e1 * t9 * t24 * t31 - 0.4e1 * t4 * t52 - 0.2e1 * t4 * t9 * t62 + x * 0.3141592654e1 * t30 * t20 - t5 - 0.4e1 * t78 * Z * t5 * t20;
2174   t156 = 0.2e1 * t78 * t90 * t31 - 0.2e1 * t3 * Z * x * t7 + t48 + 0.4e1 * t61 * t110 + 0.4e1 * t18 * t40 * t20 - 0.2e1 * t102 * t20 + 0.2e1 * t61 * t31 + 0.2e1 * t44 * t40 * t31 - t109 - 0.2e1 * t4 * t58 - 0.2e1 * t28 * t13 - 0.4e1 * t29 * t16 * t12 * t20;
2175   t159 = t1 * t11;
2176 
2177   u3  = (t57 + t94 + t128 + t156) / (0.4e1 * t159 * t6 + 0.8e1 * t159 * t73 + 0.4e1 * t159 * t31);
2178   t1  = _PC2 * Z;
2179   t2  = 0.3141592654e1 * 0.3141592654e1;
2180   t3  = t2 * 0.3141592654e1;
2181   t4  = n * n;
2182   t5  = t4 * t4;
2183   t6  = t5 * t4;
2184   t8  = t3 * t6 * x;
2185   t11 = x * t4;
2186   t12 = t11 * t3;
2187   t15 = PetscExpReal(x * n * 0.3141592654e1);
2188   t16 = t15 * t15;
2189   t17 = _PC3 * t16;
2190   t18 = nx * nx;
2191   t19 = t18 * t18;
2192   t23 = t5 * n;
2193   t24 = t2 * t23;
2194   t28 = t1 * t3;
2195   t29 = t6 * x;
2196   t30 = t29 * t16;
2197   t33 = _PC4 * t3;
2198   t34 = t5 * x;
2199   t35 = t34 * t18;
2200   t41 = PetscSinReal(nx * 0.3141592654e1 * x);
2201   t47 = t11 * t19;
2202   t54 = t3 * _PC3;
2203   t57 = 0.2e1 * t1 * t8 + 0.2e1 * t12 * t17 * t19 + 0.2e1 * t1 * t24 * t16 + 0.2e1 * t28 * t30 - 0.4e1 * t33 * t35 + 0.2e1 * t15 * nx * t41 * t4 + 0.4e1 * t28 * t35 - 0.2e1 * t33 * t47 - 0.2e1 * t1 * t24 - 0.2e1 * t33 * t29 + 0.2e1 * t29 * t54;
2204   t58 = 0.3141592654e1 * t16;
2205   t60 = t2 * _PC4;
2206   t69 = t4 * n;
2207   t73 = t1 * t2;
2208   t75 = t69 * t16 * t18;
2209   t79 = x * t16;
2210   t83 = n * t16;
2211   t84 = t83 * t19;
2212   t95 = -t34 * t58 + 0.2e1 * t60 * t23 * t16 + 0.2e1 * t60 * n * t19 - t11 * 0.3141592654e1 * t18 + 0.4e1 * t60 * t69 * t18 + 0.4e1 * t73 * t75 + 0.4e1 * t33 * t5 * t79 * t18 + 0.2e1 * t73 * t84 + 0.2e1 * t60 * t84 + 0.2e1 * t33 * t4 * t79 * t19 + 0.4e1 * t60 * t75;
2213   t97  = t34 * t3;
2214   t101 = Z * _PC1;
2215   t102 = t16 * t19;
2216   t106 = t16 * t18;
2217   t127 = t2 * t69;
2218   t131 = t2 * n;
2219   t135 = 0.4e1 * t97 * t17 * t18 + 0.2e1 * t12 * t101 * t102 + 0.4e1 * t28 * t34 * t106 + 0.2e1 * t28 * t11 * t102 - 0.2e1 * t29 * t3 * Z * _PC1 - 0.4e1 * t97 * t101 * t18 - 0.2e1 * t12 * t101 * t19 + 0.2e1 * t60 * t23 - 0.2e1 * t83 * t18 - 0.4e1 * t1 * t127 * t18 - 0.2e1 * t1 * t131 * t19;
2220   t164 = 0.2e1 * t28 * t47 + 0.2e1 * t11 * t54 * t19 + 0.2e1 * t8 * t101 * t16 + 0.2e1 * t33 * t30 - t11 * t58 * t18 + 0.2e1 * t29 * t54 * t16 + 0.4e1 * t34 * t54 * t18 + 0.4e1 * t97 * t101 * t106 - 0.2e1 * t15 * t18 * nx * t41 - t34 * 0.3141592654e1 + 0.2e1 * n * t18;
2221 
2222   u4 = (t57 + t95 + t135 + t164) / (0.4e1 * t24 * t15 + 0.8e1 * t127 * t15 * t18 + 0.4e1 * t131 * t15 * t19);
2223 
2224   /****************************************************************************************/
2225   /****************************************************************************************/
2226 
2227   u5 = (PetscReal)(-2 * Z * n * PETSC_PI * u2 - u3 * 2 * n * PETSC_PI) * PetscCosReal(n * PETSC_PI * z); /* pressure */
2228 
2229   u6 = (PetscReal)(u3 * 2 * n * PETSC_PI + 4 * Z * n * PETSC_PI * u2) * PetscCosReal(n * PETSC_PI * z); /* zz stress */
2230   sum5 += u5;
2231   sum6 += u6;
2232 
2233   u1 *= PetscCosReal(n * PETSC_PI * z); /* x velocity */
2234   sum1 += u1;
2235   u2 *= PetscSinReal(n * PETSC_PI * z); /* z velocity */
2236   sum2 += u2;
2237   u3 *= 2 * n * PETSC_PI * PetscCosReal(n * PETSC_PI * z); /* xx stress */
2238   sum3 += u3;
2239   u4 *= 2 * n * PETSC_PI * PetscSinReal(n * PETSC_PI * z); /* zx stress */
2240   sum4 += u4;
2241 
2242   /* Output */
2243   if (vel) {
2244     vel[0] = sum1;
2245     vel[1] = sum2;
2246   }
2247   if (pressure) (*pressure) = sum5;
2248 
2249   if (total_stress) {
2250     total_stress[0] = sum3;
2251     total_stress[1] = sum6;
2252     total_stress[2] = sum4;
2253   }
2254   if (strain_rate) {
2255     if (x > xc) Z = ZB;
2256     else Z = ZA;
2257 
2258     strain_rate[0] = (sum3 + sum5) / (2.0 * Z);
2259     strain_rate[1] = (sum6 + sum5) / (2.0 * Z);
2260     strain_rate[2] = (sum4 + sum5) / (2.0 * Z);
2261   }
2262 }
2263