xref: /petsc/src/snes/tutorials/ex69.c (revision 53134ebed6238e4143779c3e5c7db66beceae3ba)
1 static char help[] = "The variable-viscosity Stokes Problem in 2d with finite elements.\n\
2 We solve the Stokes problem in a square domain\n\
3 and compare against exact solutions from Mirko Velic.\n\n\n";
4 
5 /* We discretize the variable-viscosity Stokes problem using the finite element method on an unstructured mesh. The weak form equations are
6 \begin{align*}
7   (\nabla v, \mu (\nabla u + {\nabla u}^T)) - (\nabla\cdot v, p) + (v, f) &= 0 \\
8   (q, \nabla\cdot u)                                                      &= 0
9 \end{align*}
10 Free slip conditions for velocity are enforced on every wall. The pressure is constrained to have zero integral over the domain.
11 
12 To produce nice output, use
13 
14   -dm_refine 3 -show_error -dm_view hdf5:sol1.h5 -error_vec_view hdf5:sol1.h5::append -sol_vec_view hdf5:sol1.h5::append -exact_vec_view hdf5:sol1.h5::append
15 */
16 
17 #include <petscdmplex.h>
18 #include <petscsnes.h>
19 #include <petscds.h>
20 #include <petscbag.h>
21 
22 typedef enum {SOLKX, SOLCX, NUM_SOL_TYPES} SolutionType;
23 const char *solTypes[NUM_SOL_TYPES+1] = {"solkx", "solcx", "unknown"};
24 
25 typedef struct {
26   PetscInt  n, m;       /* x- and y-wavelengths for variation across the domain */
27   /* SolKx */
28   PetscReal B;          /* Exponential scale for viscosity variation */
29   /* SolCx */
30   PetscReal etaA, etaB; /* Two viscosities for discontinuous change */
31   PetscReal xc;         /* The location of viscosity jump */
32 } Parameter;
33 
34 typedef struct {
35   SolutionType solType; /* The type of exact solution */
36   PetscBag     bag;     /* Holds problem parameters */
37 } AppCtx;
38 
39 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
40 {
41   PetscInt c;
42   for (c = 0; c < Nc; ++c) u[c] = 0.0;
43   return 0;
44 }
45 static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
46 {
47   PetscInt c;
48   for (c = 0; c < Nc; ++c) u[c] = 1.0;
49   return 0;
50 }
51 
52 static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
53                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
54                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
55                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
56 {
57   f0[0] = 0.0;
58   f0[1] = -PetscSinScalar(constants[1]*PETSC_PI*x[1])*PetscCosScalar(constants[0]*PETSC_PI*x[0]);
59 }
60 
61 static void stokes_momentum_kx(PetscInt dim, PetscInt Nf, PetscInt NfAux,
62                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
63                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
64                                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
65 {
66   const PetscReal mu = PetscExpReal(2.0*PetscRealPart(constants[2])*x[0]);
67   PetscInt        c, d;
68 
69   for (c = 0; c < dim; ++c) {
70     for (d = 0; d < dim; ++d) {
71       f1[c*dim+d] = mu * (u_x[c*dim+d] + u_x[d*dim+c]);
72     }
73     f1[c*dim+c] -= u[dim];
74   }
75 }
76 
77 static void stokes_momentum_cx(PetscInt dim, PetscInt Nf, PetscInt NfAux,
78                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
79                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
80                                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
81 {
82   const PetscReal mu = x[0] < PetscRealPart(constants[4]) ? PetscRealPart(constants[2]) : PetscRealPart(constants[3]);
83   PetscInt        c, d;
84 
85   for (c = 0; c < dim; ++c) {
86     for (d = 0; d < dim; ++d) {
87       f1[c*dim+d] = mu * (u_x[c*dim+d] + u_x[d*dim+c]);
88     }
89     f1[c*dim+c] -= u[dim];
90   }
91 }
92 
93 static void stokes_mass(PetscInt dim, PetscInt Nf, PetscInt NfAux,
94                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
95                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
97 {
98   PetscInt d;
99   f0[0] = 0.0;
100   for (d = 0; d < dim; ++d) f0[0] -= u_x[d*dim+d];
101 }
102 
103 static void f1_zero(PetscInt dim, PetscInt Nf, PetscInt NfAux,
104                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
105                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
106                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
107 {
108   PetscInt d;
109   for (d = 0; d < dim*dim; ++d) f1[d] = 0.0;
110 }
111 
112 /* < q, \nabla\cdot u >, J_{pu} */
113 static void stokes_mass_J(PetscInt dim, PetscInt Nf, PetscInt NfAux,
114                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
115                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
116                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
117 {
118   PetscInt d;
119   for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
120 }
121 
122 /* -< \nabla\cdot v, p >, J_{up} */
123 static void stokes_momentum_pres_J(PetscInt dim, PetscInt Nf, PetscInt NfAux,
124                                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
125                                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
126                                    PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
127 {
128   PetscInt d;
129   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
130 }
131 
132 /* < \nabla v, \nabla u + {\nabla u}^T >, J_{uu}
133    This just gives \nabla u, give the perdiagonal for the transpose */
134 static void stokes_momentum_vel_J_kx(PetscInt dim, PetscInt Nf, PetscInt NfAux,
135                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
136                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
137                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
138 {
139   const PetscReal mu  = PetscExpReal(2.0*PetscRealPart(constants[2])*x[0]);
140   PetscInt        cI, d;
141 
142   for (cI = 0; cI < dim; ++cI) {
143     for (d = 0; d < dim; ++d) {
144       g3[((cI*dim+cI)*dim+d)*dim+d] += mu; /*g3[cI, cI, d, d]*/
145       g3[((cI*dim+d)*dim+d)*dim+cI] += mu; /*g3[cI, d, d, cI]*/
146     }
147   }
148 }
149 static void stokes_momentum_vel_J_cx(PetscInt dim, PetscInt Nf, PetscInt NfAux,
150                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
151                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
152                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
153 {
154   const PetscReal mu = x[0] < PetscRealPart(constants[4]) ? PetscRealPart(constants[2]) : PetscRealPart(constants[3]);
155   PetscInt        cI, d;
156 
157   for (cI = 0; cI < dim; ++cI) {
158     for (d = 0; d < dim; ++d) {
159       g3[((cI*dim+cI)*dim+d)*dim+d] += mu; /*g3[cI, cI, d, d]*/
160       g3[((cI*dim+d)*dim+d)*dim+cI] += mu; /*g3[cI, d, d, cI]*/
161     }
162   }
163 }
164 
165 /* 1/mu < q, I q >, Jp_{pp} */
166 static void stokes_identity_J_kx(PetscInt dim, PetscInt Nf, PetscInt NfAux,
167                                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
168                                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
169                                  PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
170 {
171   const PetscReal mu = PetscExpReal(2.0*PetscRealPart(constants[2])*x[0]);
172   g0[0] = 1.0/mu;
173 }
174 
175 static void stokes_identity_J_cx(PetscInt dim, PetscInt Nf, PetscInt NfAux,
176                                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
177                                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
178                                  PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
179 {
180   const PetscReal mu = x[0] < PetscRealPart(constants[4]) ? PetscRealPart(constants[2]) : PetscRealPart(constants[3]);
181   g0[0] = 1.0/mu;
182 }
183 
184 /*
185   SolKxSolution - Exact Stokes solutions for exponentially varying viscosity
186 
187  Input Parameters:
188 + pos   - The (x,z) coordinate at which to evaluate the solution
189 . n     - The constant defining the x-dependence of the forcing function
190 . m     - The constant defining the z-dependence of the forcing function
191 - B     - The viscosity coefficient
192 
193   Output Parameters:
194 + vel   - The (x,z)-velocity at (x,z), or NULL
195 . p     - The pressure at (x,z), or NULL
196 . s     - The total stress (sigma_xx, sigma_xz, sigma_zz) at (x,z), or NULL
197 . gamma - The strain rate, or NULL
198 - mu    - The viscosity at (x,z), or NULL
199 
200   Note:
201 $  The domain is the square 0 <= x,z <= 1. We solve the Stokes equation for incompressible flow with free-slip boundary
202 $  conditions everywhere. The forcing term f is given by
203 $
204 $    fx = 0
205 $    fz = sigma*sin(km*z)*cos(kn*x)
206 $
207 $  where
208 $
209 $    km = m*Pi (m may be non-integral)
210 $    kn = n*Pi
211 $
212 $  meaning that the density rho is -sigma*sin(km*z)*cos(kn*x). Here we set sigma = 1.
213 $  The viscosity eta is exp(2*B*x).
214 */
215 static PetscErrorCode SolKxSolution(const PetscReal pos[], PetscReal m, PetscInt n, PetscReal B,
216                                     PetscScalar vel[], PetscScalar *p, PetscScalar s[], PetscScalar gamma[], PetscScalar *mu)
217 {
218   PetscReal sigma = 1.0;
219   PetscReal Z;
220   PetscReal u1,u2,u3,u4,u5,u6;
221   PetscReal sum1,sum2,sum3,sum4,sum5,sum6;
222   PetscReal kn,km,x,z;
223   PetscReal _PC1,_PC2,_PC3,_PC4;
224   PetscReal Rp, UU, VV;
225   PetscReal a,b,r,_aa,_bb,AA,BB,Rm;
226   PetscReal num1,num2,num3,num4,den1;
227 
228   PetscReal t1,t2,t3,t4,t5,t6,t7,t8,t9,t10;
229   PetscReal t11,t12,t13,t14,t15,t16,t17,t18,t19,t20;
230   PetscReal t21,t22,t23,t24,t25,t26,t27,t28,t29,t30;
231   PetscReal t31,t32,t33,t34,t35,t36,t37,t38,t39,t40;
232   PetscReal t41,t42,t43,t44,t45,t46,t47,t49,t51,t53;
233   PetscReal t56,t58,t61,t62,t63,t64,t65,t66,t67,t68;
234   PetscReal t69,t70,t71,t72,t73,t74,t75,t76,t77,t78;
235   PetscReal t79,t80,t81,t82,t83,t84,t85,t86,t87,t88;
236   PetscReal t89,t90,t91,t92,t93,t94,t95,t96,t97,t99;
237   PetscReal t100,t101,t103,t104,t105,t106,t107,t108,t109,t110;
238   PetscReal t111,t112,t113,t114,t115,t116,t117,t118,t119,t120;
239   PetscReal t121,t124,t125,t126,t127,t129,t130,t132,t133,t135;
240   PetscReal t136,t138,t140,t141,t142,t143,t152,t160,t162;
241 
242   PetscFunctionBegin;
243   /*************************************************************************/
244   /*************************************************************************/
245   /* rho = -sin(km*z)*cos(kn*x) */
246   x = pos[0];
247   z = pos[1];
248   Z = PetscExpReal( 2.0 * B * x);
249   km = m*PETSC_PI; /* solution valid for km not zero -- should get trivial solution if km=0 */
250   kn = (PetscReal) n*PETSC_PI;
251   /*************************************************************************/
252   /*************************************************************************/
253   a = B*B + km*km;
254   b = 2.0*km*B;
255   r = PetscSqrtReal(a*a + b*b);
256   Rp = PetscSqrtReal( (r+a)/2.0);
257   Rm  = PetscSqrtReal( (r-a)/2.0);
258   UU  = Rp - B;
259   VV = Rp + B;
260 
261   sum1=0.0;
262   sum2=0.0;
263   sum3=0.0;
264   sum4=0.0;
265   sum5=0.0;
266   sum6=0.0;
267   /*sum7=0.0;*/
268 
269   /*******************************************/
270   /*         calculate the constants         */
271   /*******************************************/
272 
273   t1 = kn * kn;
274   t4 = km * km;
275   t5 = t4 + t1;
276   t6 = t5 * t5;
277   t8 = pow(km + kn, 0.2e1);
278   t9 = B * B;
279   t16 = pow(km - kn, 0.2e1);
280   _aa = -0.4e1 * B * t1 * sigma * t5 / (t6 + 0.4e1 * t8 * t9) / (t6 + 0.4e1 * t16 * t9);
281 
282   t2 = km * km;
283   t3 = kn * kn;
284   t5 = pow(t2 + t3, 0.2e1);
285   t6 = km - kn;
286   t7 = km + kn;
287   t9 = B * B;
288   t13 = t7 * t7;
289   t19 = t6 * t6;
290   _bb = sigma * kn * (t5 + 0.4e1 * t6 * t7 * t9) / (t5 + 0.4e1 * t13 * t9) / (t5 + 0.4e1 * t19 * t9);
291 
292   AA = _aa;
293   BB = _bb;
294 
295   /*******************************************/
296   /*       calculate the velocities etc      */
297   /*******************************************/
298   t1 = Rm * Rm;
299   t2 = B - Rp;
300   t4 = Rp + B;
301   t6 = UU * x;
302   t9 = PetscExpReal(t6 - 0.4e1 * Rp);
303   t13 = B * B;
304   t16 = Rp * t1;
305   t18 = Rp * Rp;
306   t19 = B * t18;
307   t20 = t13 * Rp;
308   t22 = kn * kn;
309   t24 = B * t1;
310   t32 = 0.8e1 * t13 * BB * kn * Rp;
311   t34 = 0.2e1 * Rm;
312   t35 = PetscCosReal(t34);
313   t37 = Rp * Rm;
314   t49 = PetscSinReal(t34);
315   t63 = PetscExpReal(t6 - 0.2e1 * Rp);
316   t65 = Rm * t2;
317   t67 = 0.2e1 * B * kn;
318   t68 = B * Rm;
319   t69 = t67 + t68 + t37;
320   t73 = 0.3e1 * t13;
321   t75 = 0.2e1 * B * Rp;
322   t76 = t73 - t75 + t1 - t22 - t18;
323   t78 = t65 * t76 * BB;
324   t80 = Rm - kn;
325   t81 = PetscCosReal(t80);
326   t83 = -t67 + t68 + t37;
327   t88 = Rm + kn;
328   t89 = PetscCosReal(t88);
329   t92 = t65 * t76 * AA;
330   t97 = PetscSinReal(t80);
331   t103 = PetscSinReal(t88);
332   t108 = PetscExpReal(t6 - 0.3e1 * Rp - B);
333   t110 = Rm * t4;
334   t111 = t67 + t68 - t37;
335   t115 = t73 + t75 + t1 - t22 - t18;
336   t117 = t110 * t115 * BB;
337   t120 = -t67 + t68 - t37;
338   t127 = t110 * t115 * AA;
339   t140 = PetscExpReal(t6 - Rp - B);
340   num1 = -0.4e1 * t1 * t2 * t4 * AA * t9 + ((0.2e1 * Rp * (0.3e1 * t13 * B - 0.2e1 * t16 - t19 - 0.2e1 * t20 - B * t22 - t24) * AA - t32) * t35 + (0.2e1 * t37 * (t1 + 0.5e1 * t13 - t22 - t18) * AA - 0.8e1 * B * BB * kn * Rm * Rp) * t49 - 0.2e1 * B * (0.3e1 * t20 - Rp * t22 - t18 * Rp - 0.2e1 * t19 - t16 - 0.2e1 * t24) * AA + t32) * t63 + ((0.2e1 * t65 * t69 * AA + t78) * t81 + (0.2e1 * t65 * t83 * AA - t78) * t89 + (t92 - 0.2e1 * t65 * t69 * BB) * t97 + (t92 + 0.2e1 * t65 * t83 * BB) * t103) * t108 + ((-0.2e1 * t110 * t111 * AA - t117) * t81 + (-0.2e1 * t110 * t120 * AA + t117) * t89 + (-t127 + 0.2e1 * t110 * t111 * BB) * t97 + (-t127 - 0.2e1 * t110 * t120 * BB) * t103) * t140;
341 
342   t1 = Rp + B;
343   t2 = Rm * t1;
344   t3 = B * B;
345   t4 = 0.3e1 * t3;
346   t5 = B * Rp;
347   t7 = Rm * Rm;
348   t8 = kn * kn;
349   t9 = Rp * Rp;
350   t10 = t4 + 0.2e1 * t5 + t7 - t8 - t9;
351   t12 = t2 * t10 * AA;
352   t14 = B * Rm;
353   t20 = UU * x;
354   t23 = PetscExpReal(t20 - 0.4e1 * Rp);
355   t25 = Rp * Rm;
356   t32 = Rm * kn;
357   t37 = 0.2e1 * Rm;
358   t38 = PetscCosReal(t37);
359   t40 = t3 * B;
360   t44 = B * t9;
361   t45 = t3 * Rp;
362   t53 = t3 * BB;
363   t58 = PetscSinReal(t37);
364   t69 = PetscExpReal(t20 - 0.2e1 * Rp);
365   t72 = 0.3e1 * t40 * Rm;
366   t73 = t9 * Rp;
367   t74 = t73 * Rm;
368   t75 = t7 * Rm;
369   t76 = B * t75;
370   t77 = t14 * t8;
371   t78 = Rp * t75;
372   t80 = 0.8e1 * t45 * kn;
373   t81 = t25 * t8;
374   t83 = 0.5e1 * t45 * Rm;
375   t84 = t44 * Rm;
376   t85 = t72 - t74 + t76 - t77 + t78 + t80 - t81 + t83 + t84;
377   t88 = 0.2e1 * t9 * t3;
378   t90 = 0.3e1 * t40 * Rp;
379   t91 = t7 * t3;
380   t93 = 0.2e1 * t5 * t32;
381   t94 = t5 * t7;
382   t95 = t5 * t8;
383   t96 = B * t73;
384   t97 = t7 * t9;
385   t100 = 0.2e1 * t3 * Rm * kn;
386   t101 = -t88 + t90 - t91 - t93 - t94 - t95 - t96 - t97 - t100;
387   t105 = Rm - kn;
388   t106 = PetscCosReal(t105);
389   t108 = t72 - t80 + t83 + t76 + t84 - t81 - t74 + t78 - t77;
390   t110 = -t97 - t96 - t88 + t100 + t90 - t95 + t93 - t91 - t94;
391   t114 = Rm + kn;
392   t115 = PetscCosReal(t114);
393   t121 = PetscSinReal(t105);
394   t127 = PetscSinReal(t114);
395   t132 = PetscExpReal(t20 - 0.3e1 * Rp - B);
396   t135 = 0.2e1 * B * kn;
397   t136 = t135 + t14 - t25;
398   t142 = -t135 + t14 - t25;
399   t152 = t2 * t10 * BB;
400   t162 = PetscExpReal(t20 - Rp - B);
401   num2 = (0.2e1 * t12 - 0.8e1 * t14 * kn * t1 * BB) * t23 + ((-0.2e1 * t25 * (t7 + 0.5e1 * t3 - t8 - t9) * AA + 0.8e1 * B * BB * t32 * Rp) * t38 + (0.2e1 * Rp * (0.3e1 * t40 - 0.2e1 * Rp * t7 - t44 - 0.2e1 * t45 - B * t8 - B * t7) * AA - 0.8e1 * t53 * kn * Rp) * t58 - 0.2e1 * t14 * (-t8 + t9 + t4 + t7) * AA + 0.8e1 * t53 * t32) * t69 + ((-t85 * AA - 0.2e1 * t101 * BB) * t106 + (-t108 * AA + 0.2e1 * t110 * BB) * t115 + (-0.2e1 * t101 * AA + t85 * BB) * t121 + (-0.2e1 * t110 * AA - t108 * BB) * t127) * t132 + ((t12 - 0.2e1 * t2 * t136 * BB) * t106 + (t12 + 0.2e1 * t2 * t142 * BB) * t115 + (-0.2e1 * t2 * t136 * AA - t152) * t121 + (-0.2e1 * t2 * t142 * AA + t152) * t127) * t162;
402 
403   t1 = Rm * Rm;
404   t2 = B - Rp;
405   t4 = Rp + B;
406   t6 = VV * x;
407   t7 = PetscExpReal(-t6);
408   t11 = kn * kn;
409   t13 = B * t1;
410   t14 = Rp * Rp;
411   t15 = B * t14;
412   t16 = B * B;
413   t17 = t16 * Rp;
414   t21 = Rp * t1;
415   t30 = 0.8e1 * t16 * BB * kn * Rp;
416   t32 = 0.2e1 * Rm;
417   t33 = PetscCosReal(t32);
418   t35 = Rp * Rm;
419   t47 = PetscSinReal(t32);
420   t61 = PetscExpReal(-t6 - 0.2e1 * Rp);
421   t63 = Rm * t2;
422   t65 = 0.2e1 * B * kn;
423   t66 = B * Rm;
424   t67 = t65 + t66 + t35;
425   t71 = 0.3e1 * t16;
426   t73 = 0.2e1 * B * Rp;
427   t74 = t71 - t73 + t1 - t11 - t14;
428   t76 = t63 * t74 * BB;
429   t78 = Rm - kn;
430   t79 = PetscCosReal(t78);
431   t81 = -t65 + t66 + t35;
432   t86 = Rm + kn;
433   t87 = PetscCosReal(t86);
434   t90 = t63 * t74 * AA;
435   t95 = PetscSinReal(t78);
436   t101 = PetscSinReal(t86);
437   t106 = PetscExpReal(-t6 - 0.3e1 * Rp - B);
438   t108 = Rm * t4;
439   t109 = t65 + t66 - t35;
440   t113 = t71 + t73 + t1 - t11 - t14;
441   t115 = t108 * t113 * BB;
442   t118 = -t65 + t66 - t35;
443   t125 = t108 * t113 * AA;
444   t138 = PetscExpReal(-t6 - Rp - B);
445   num3 = -0.4e1 * t1 * t2 * t4 * AA * t7 + ((-0.2e1 * Rp * (-B * t11 - t13 - t15 + 0.2e1 * t17 + 0.3e1 * t16 * B + 0.2e1 * t21) * AA + t30) * t33 + (-0.2e1 * t35 * (t1 + 0.5e1 * t16 - t11 - t14) * AA + 0.8e1 * B * BB * kn * Rm * Rp) * t47 + 0.2e1 * B * (0.3e1 * t17 - t21 + 0.2e1 * t15 + 0.2e1 * t13 - Rp * t11 - t14 * Rp) * AA - t30) * t61 + ((-0.2e1 * t63 * t67 * AA - t76) * t79 + (-0.2e1 * t63 * t81 * AA + t76) * t87 + (-t90 + 0.2e1 * t63 * t67 * BB) * t95 + (-t90 - 0.2e1 * t63 * t81 * BB) * t101) * t106 + ((0.2e1 * t108 * t109 * AA + t115) * t79 + (0.2e1 * t108 * t118 * AA - t115) * t87 + (t125 - 0.2e1 * t108 * t109 * BB) * t95 + (t125 + 0.2e1 * t108 * t118 * BB) * t101) * t138;
446 
447   t1 = B - Rp;
448   t2 = Rm * t1;
449   t3 = B * B;
450   t4 = 0.3e1 * t3;
451   t5 = B * Rp;
452   t7 = Rm * Rm;
453   t8 = kn * kn;
454   t9 = Rp * Rp;
455   t10 = t4 - 0.2e1 * t5 + t7 - t8 - t9;
456   t12 = t2 * t10 * AA;
457   t14 = B * Rm;
458   t20 = VV * x;
459   t21 = PetscExpReal(-t20);
460   t23 = Rp * Rm;
461   t30 = Rm * kn;
462   t35 = 0.2e1 * Rm;
463   t36 = PetscCosReal(t35);
464   t40 = B * t9;
465   t41 = t3 * Rp;
466   t43 = t3 * B;
467   t51 = t3 * BB;
468   t56 = PetscSinReal(t35);
469   t67 = PetscExpReal(-t20 - 0.2e1 * Rp);
470   t70 = 0.2e1 * B * kn;
471   t71 = t70 + t14 + t23;
472   t76 = Rm - kn;
473   t77 = PetscCosReal(t76);
474   t79 = -t70 + t14 + t23;
475   t84 = Rm + kn;
476   t85 = PetscCosReal(t84);
477   t91 = t2 * t10 * BB;
478   t93 = PetscSinReal(t76);
479   t99 = PetscSinReal(t84);
480   t104 = PetscExpReal(-t20 - 0.3e1 * Rp - B);
481   t107 = 0.3e1 * t43 * Rm;
482   t108 = t9 * Rp;
483   t109 = t108 * Rm;
484   t110 = t7 * Rm;
485   t111 = B * t110;
486   t112 = t14 * t8;
487   t113 = Rp * t110;
488   t115 = 0.8e1 * t41 * kn;
489   t116 = t23 * t8;
490   t118 = 0.5e1 * t41 * Rm;
491   t119 = t40 * Rm;
492   t120 = t107 + t109 + t111 - t112 - t113 - t115 + t116 - t118 + t119;
493   t124 = 0.2e1 * t3 * Rm * kn;
494   t125 = t5 * t8;
495   t126 = B * t108;
496   t127 = t7 * t9;
497   t129 = 0.2e1 * t9 * t3;
498   t130 = t5 * t7;
499   t132 = 0.3e1 * t43 * Rp;
500   t133 = t7 * t3;
501   t135 = 0.2e1 * t5 * t30;
502   t136 = t124 - t125 - t126 + t127 + t129 - t130 + t132 + t133 - t135;
503   t141 = t107 + t115 - t118 + t111 + t119 + t116 + t109 - t113 - t112;
504   t143 = t132 + t129 - t125 + t133 + t127 - t124 - t130 - t126 + t135;
505   t160 = PetscExpReal(-t20 - Rp - B);
506   num4 = (0.2e1 * t12 - 0.8e1 * t14 * kn * t1 * BB) * t21 + ((0.2e1 * t23 * (t7 + 0.5e1 * t3 - t8 - t9) * AA - 0.8e1 * B * BB * t30 * Rp) * t36 + (-0.2e1 * Rp * (-B * t8 - B * t7 - t40 + 0.2e1 * t41 + 0.3e1 * t43 + 0.2e1 * Rp * t7) * AA + 0.8e1 * t51 * kn * Rp) * t56 - 0.2e1 * t14 * (-t8 + t9 + t4 + t7) * AA + 0.8e1 * t51 * t30) * t67 + ((t12 - 0.2e1 * t2 * t71 * BB) * t77 + (t12 + 0.2e1 * t2 * t79 * BB) * t85 + (-0.2e1 * t2 * t71 * AA - t91) * t93 + (-0.2e1 * t2 * t79 * AA + t91) * t99) * t104 + ((-t120 * AA + 0.2e1 * t136 * BB) * t77 + (-t141 * AA - 0.2e1 * t143 * BB) * t85 + (0.2e1 * t136 * AA + t120 * BB) * t93 + (0.2e1 * t143 * AA - t141 * BB) * t99) * t160;
507 
508   t1 = Rm * Rm;
509   t2 = Rp * Rp;
510   t3 = t1 * t2;
511   t4 = B * B;
512   t5 = t1 * t4;
513   t9 = PetscExpReal(-0.4e1 * Rp);
514   t15 = PetscCosReal(0.2e1 * Rm);
515   t22 = PetscExpReal(-0.2e1 * Rp);
516   den1 = (-0.4e1 * t3 + 0.4e1 * t5) * t9 + ((0.8e1 * t1 + 0.8e1 * t4) * t2 * t15 - 0.8e1 * t5 - 0.8e1 * t2 * t4) * t22 - 0.4e1 * t3 + 0.4e1 * t5;
517 
518   _PC1=num1/den1; _PC2=num2/den1; _PC3=num3/den1; _PC4=num4/den1;
519 
520   t1 = Rm * x;
521   t2 = PetscCosReal(t1);
522   t4 = PetscSinReal(t1);
523   t10 = PetscExpReal(-0.2e1 * x * B);
524   t12 = kn * x;
525   t13 = PetscCosReal(t12);
526   t16 = PetscSinReal(t12);
527   u1 = -km * (_PC1 * t2 + _PC2 * t4 + _PC3 * t2 + _PC4 * t4 + t10 * AA * t13 + t10 * BB * t16);
528 
529   t2 = Rm * x;
530   t3 = PetscCosReal(t2);
531   t6 = PetscSinReal(t2);
532   t22 = PetscExpReal(-0.2e1 * x * B);
533   t23 = B * t22;
534   t24 = kn * x;
535   t25 = PetscCosReal(t24);
536   t29 = PetscSinReal(t24);
537   u2 = UU * _PC1 * t3 + UU * _PC2 * t6 - _PC1 * t6 * Rm + _PC2 * t3 * Rm - VV * _PC3 * t3 - VV * _PC4 * t6 - _PC3 * t6 * Rm + _PC4 * t3 * Rm - 0.2e1 * t23 * AA * t25 - 0.2e1 * t23 * BB * t29 - t22 * AA * t29 * kn + t22 * BB * t25 * kn;
538 
539   t3 = PetscExpReal(0.2e1 * x * B);
540   t4 = t3 * B;
541   t8 = km * km;
542   t9 = t3 * t8;
543   t11 = 0.3e1 * t9 * Rm;
544   t12 = Rm * Rm;
545   t14 = t3 * t12 * Rm;
546   t15 = UU * UU;
547   t19 = 0.4e1 * t4 * UU * Rm - t11 - t14 + 0.3e1 * t3 * t15 * Rm;
548   t20 = Rm * x;
549   t21 = PetscSinReal(t20);
550   t27 = 0.2e1 * B * t9;
551   t33 = 0.2e1 * t4 * t12;
552   t36 = 0.3e1 * t3 * UU * t12 - t27 - 0.2e1 * t4 * t15 + 0.3e1 * t9 * UU + t33 - t3 * t15 * UU;
553   t37 = PetscCosReal(t20);
554   t49 = VV * VV;
555   t53 = -0.4e1 * t4 * VV * Rm - t11 + 0.3e1 * t3 * t49 * Rm - t14;
556   t64 = t3 * t49 * VV + t33 - 0.3e1 * t9 * VV - 0.2e1 * t4 * t49 - t27 - 0.3e1 * t3 * VV * t12;
557   t76 = B * t8;
558   t80 = kn * kn;
559   t83 = B * B;
560   t87 = t80 * kn;
561   t90 = kn * x;
562   t91 = PetscSinReal(t90);
563   t106 = PetscCosReal(t90);
564   u3 = -((t19 * t21 + t36 * t37) * _PC1 + (t36 * t21 - t19 * t37) * _PC2 + (t53 * t21 + t64 * t37) * _PC3 + (t64 * t21 - t53 * t37) * _PC4 + (-0.3e1 * t8 * AA * kn - 0.8e1 * t76 * BB - 0.4e1 * BB * B * t80 + 0.4e1 * AA * t83 * kn - AA * t87) * t91 + (-0.4e1 * AA * t80 * B - 0.4e1 * t83 * BB * kn + 0.3e1 * t8 * BB * kn - sigma + BB * t87 - 0.8e1 * t76 * AA) * t106) / km;
565 
566   t3 = PetscExpReal(0.2e1 * x * B);
567   t4 = km * km;
568   t5 = t3 * t4;
569   t6 = Rm * x;
570   t7 = PetscCosReal(t6);
571   t8 = _PC1 * t7;
572   t10 = PetscSinReal(t6);
573   t11 = _PC2 * t10;
574   t13 = _PC3 * t7;
575   t15 = _PC4 * t10;
576   t18 = kn * x;
577   t19 = PetscCosReal(t18);
578   t22 = PetscSinReal(t18);
579   t24 = UU * UU;
580   t25 = t3 * t24;
581   t28 = t3 * UU;
582   t38 = Rm * Rm;
583   t39 = t7 * t38;
584   t42 = t10 * t38;
585   t44 = t5 * t8 + t5 * t11 + t5 * t13 + t5 * t15 + t4 * AA * t19 + t4 * BB * t22 + t25 * t8 + t25 * t11 - 0.2e1 * t28 * _PC1 * t10 * Rm + 0.2e1 * t28 * _PC2 * t7 * Rm - t3 * _PC1 * t39 - t3 * _PC2 * t42;
586   t45 = VV * VV;
587   t46 = t3 * t45;
588   t49 = t3 * VV;
589   t62 = B * B;
590   t78 = kn * kn;
591   t82 = t46 * t13 + t46 * t15 + 0.2e1 * t49 * _PC3 * t10 * Rm - 0.2e1 * t49 * _PC4 * t7 * Rm - t3 * _PC3 * t39 - t3 * _PC4 * t42 + 0.4e1 * t62 * AA * t19 + 0.4e1 * t62 * BB * t22 + 0.4e1 * B * AA * t22 * kn - 0.4e1 * B * BB * t19 * kn - AA * t19 * t78 - BB * t22 * t78;
592   u4 = t44 + t82;
593 
594   t3 = PetscExpReal(0.2e1 * x * B);
595   t4 = t3 * B;
596   t8 = km * km;
597   t9 = t3 * t8;
598   t10 = t9 * Rm;
599   t11 = Rm * Rm;
600   t13 = t3 * t11 * Rm;
601   t14 = UU * UU;
602   t18 = 0.4e1 * t4 * UU * Rm - t10 - t13 + 0.3e1 * t3 * t14 * Rm;
603   t19 = Rm * x;
604   t20 = PetscSinReal(t19);
605   t26 = 0.2e1 * B * t9;
606   t31 = 0.2e1 * t4 * t11;
607   t34 = 0.3e1 * t3 * UU * t11 - t26 - 0.2e1 * t4 * t14 + t9 * UU + t31 - t3 * t14 * UU;
608   t35 = PetscCosReal(t19);
609   t47 = VV * VV;
610   t51 = -0.4e1 * t4 * VV * Rm - t10 + 0.3e1 * t3 * t47 * Rm - t13;
611   t61 = t3 * t47 * VV + t31 - t9 * VV - 0.2e1 * t4 * t47 - t26 - 0.3e1 * t3 * VV * t11;
612   t72 = B * t8;
613   t76 = kn * kn;
614   t79 = B * B;
615   t83 = t76 * kn;
616   t86 = kn * x;
617   t87 = PetscSinReal(t86);
618   t101 = PetscCosReal(t86);
619   u5 = ((t18 * t20 + t34 * t35) * _PC1 + (t34 * t20 - t18 * t35) * _PC2 + (t51 * t20 + t61 * t35) * _PC3 + (t61 * t20 - t51 * t35) * _PC4 + (-t8 * AA * kn - 0.4e1 * t72 * BB - 0.4e1 * BB * B * t76 + 0.4e1 * AA * t79 * kn - AA * t83) * t87 + (-0.4e1 * AA * t76 * B - 0.4e1 * t79 * BB * kn + t8 * BB * kn - sigma + BB * t83 - 0.4e1 * t72 * AA) * t101) / km;
620 
621   t3 = PetscExpReal(0.2e1 * x * B);
622   t4 = UU * UU;
623   t8 = km * km;
624   t9 = t3 * t8;
625   t10 = t9 * Rm;
626   t11 = Rm * Rm;
627   t13 = t3 * t11 * Rm;
628   t14 = t3 * B;
629   t18 = 0.3e1 * t3 * t4 * Rm + t10 - t13 + 0.4e1 * t14 * UU * Rm;
630   t19 = Rm * x;
631   t20 = PetscSinReal(t19);
632   t28 = 0.2e1 * B * t9;
633   t33 = 0.2e1 * t14 * t11;
634   t34 = -0.2e1 * t4 * t14 + 0.3e1 * t3 * UU * t11 - t28 - t3 * t4 * UU - t9 * UU + t33;
635   t35 = PetscCosReal(t19);
636   t47 = VV * VV;
637   t51 = -0.4e1 * t14 * VV * Rm - t13 + t10 + 0.3e1 * t3 * t47 * Rm;
638   t61 = -0.3e1 * t3 * VV * t11 + t33 + t3 * t47 * VV + t9 * VV - 0.2e1 * t14 * t47 - t28;
639   t71 = kn * kn;
640   t74 = B * B;
641   t80 = t71 * kn;
642   t83 = kn * x;
643   t84 = PetscSinReal(t83);
644   t96 = PetscCosReal(t83);
645   u6 = -((t18 * t20 + t34 * t35) * _PC1 + (t34 * t20 - t18 * t35) * _PC2 + (t51 * t20 + t61 * t35) * _PC3 + (t61 * t20 - t51 * t35) * _PC4 + (-0.4e1 * BB * B * t71 + 0.4e1 * AA * t74 * kn + t8 * AA * kn - AA * t80) * t84 + (-0.4e1 * AA * t71 * B - t8 * BB * kn - 0.4e1 * t74 * BB * kn - sigma + BB * t80) * t96) / km;
646 
647   /*SS = sin(km*z)*(exp(UU*x)*(_PC1*cos(Rm*x)+_PC2*sin(Rm*x)) + exp(-VV*x)*(_PC3*cos(Rm*x)+_PC4*sin(Rm*x)) + exp(-2*x*B)*(AA*cos(kn*x)+BB*sin(kn*x)));*/
648 
649   /* u1 = Vx, u2 = Vz, u3 = txx, u4 = tzx, u5 = pressure, u6 = tzz */
650 
651   sum5 += u5*PetscCosReal(km*z);  /* pressure */
652   sum6 += u6*PetscCosReal(km*z);  /* zz total stress */
653 
654   u1 *= PetscCosReal(km*z); /* x velocity */
655   sum1 += u1;
656   u2 *= PetscSinReal(km*z); /* z velocity */
657   sum2 += u2;
658 
659   u3 *= PetscCosReal(km*z); /* xx total stress */
660   sum3 += u3;
661   u4 *= PetscSinReal(km*z); /* zx stress */
662   sum4 += u4;
663 
664   /* rho = -sigma*sin(km*z)*cos(kn*x); */ /* density */
665   /* sum7 += rho; */
666 
667   /* Output */
668   if (mu) {
669     *mu = Z;
670   }
671   if (vel) {
672     vel[0] = sum1;
673     vel[1] = sum2;
674   }
675   if (p) {
676     (*p) = sum5;
677   }
678   if (s) {
679     s[0] = sum3;
680     s[1] = sum4;
681     s[2] = sum6;
682   }
683   if (gamma) {
684     /* sigma = tau - p, tau = sigma + p, tau[] = 2*eta*gamma[] */
685     gamma[0] = (sum3+sum5)/(2.0*Z);
686     gamma[1] = (sum4)/(2.0*Z);
687     gamma[2] = (sum6+sum5)/(2.0*Z);
688   }
689   PetscFunctionReturn(0);
690 }
691 
692 static PetscErrorCode SolKxSolutionVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar v[], void *ctx)
693 {
694   Parameter     *s = (Parameter *) ctx;
695 
696   PetscFunctionBegin;
697   PetscCall(SolKxSolution(x, s->m, s->n, s->B, v, NULL, NULL, NULL, NULL));
698   PetscFunctionReturn(0);
699 }
700 
701 static PetscErrorCode SolKxSolutionPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar p[], void *ctx)
702 {
703   Parameter     *s = (Parameter *) ctx;
704 
705   PetscFunctionBegin;
706   PetscCall(SolKxSolution(x, s->m, s->n, s->B, NULL, p, NULL, NULL, NULL));
707   PetscFunctionReturn(0);
708 }
709 
710 /*
711   SolCxSolution - Exact Stokes solutions for discontinuous viscosity
712 
713  Input Parameters:
714 + pos   - The (x,z) coordinate at which to evaluate the solution
715 . n     - The constant defining the x-dependence of the forcing function
716 . m     - The constant defining the z-dependence of the forcing function
717 . xc    - The x coordinate at which the viscosity is discontinuous
718 . etaA  - The viscosity coefficient for x < xc
719 - etaB  - The viscosity coefficient for x > xc
720 
721   Output Parameters:
722 + vel   - The (x,z)-velocity at (x,z), or NULL
723 . p     - The pressure at (x,z), or NULL
724 . s     - The total stress (sigma_xx, sigma_xz, sigma_zz) at (x,z), or NULL
725 . gamma - The strain rate, or NULL
726 - mu    - The viscosity at (x,z), or NULL
727 
728   Note:
729 $  The domain is the square 0 <= x,z <= 1. We solve the Stokes equation for incompressible flow with free-slip boundary
730 $  conditions everywhere. The forcing term f is given by
731 $
732 $    fx = 0
733 $    fz = sigma*sin(km*z)*cos(kn*x)
734 $
735 $  where
736 $
737 $    km = m*Pi (m may be non-integral)
738 $    kn = n*Pi
739 $
740 $  meaning that the density rho is -sigma*sin(km*z)*cos(kn*x). Here we set sigma = 1.
741 $  The viscosity eta jumps from etaA to etaB at x = xc.
742 */
743 static PetscErrorCode SolCxSolution(const PetscReal pos[], PetscReal m, PetscInt n, PetscReal xc, PetscReal etaA, PetscReal etaB,
744                                     PetscScalar vel[], PetscScalar *p, PetscScalar s[], PetscScalar gamma[], PetscScalar *mu)
745 {
746   PetscReal _PC1A,_PC2A,_PC3A,_PC4A,_PC1B,_PC2B,_PC3B,_PC4B,_PC1,_PC2,_PC3,_PC4;
747   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;
748   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;
749   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;
750   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;
751   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;
752   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;
753   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;
754   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;
755   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;
756   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;
757   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;
758   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;
759   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;
760   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;
761   PetscReal t561,t562,t563,t564,t567,t568,t571,t573,t575,t576,t578,t579,t583,t590,t591,t594,t595,t596,t597,t598,t600;
762   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;
763   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;
764   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;
765   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;
766   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;
767   PetscReal t801,t806,t807,t808,t809,t811,t812,t817,t818,t819,t821,t822,t824,t827,t828,t830,t834,t835,t837,t840;
768   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;
769   PetscReal t884,t888,t891,t894,t900,t901,t903,t904,t907,t908,t909,t911,t914,t915,t916,t919,t920;
770   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;
771   PetscReal t961,t964,t965,t966,t967,t968,t969,t971,t972,t974,t977,t978,t981,t983,t987,t988,t992,t993,t994,t997,t998;
772   PetscReal t1001,t1003,t1005,t1006,t1009,t1010,t1012,t1013,t1015,t1016,t1017,t1018,t1020,t1021,t1029,t1031,t1032,t1033,t1040;
773   PetscReal t1041,t1042,t1044,t1047,t1050,t1054,t1055,t1057,t1058,t1063,t1068,t1069,t1070,t1079,t1080;
774   PetscReal t1088,t1089,t1091,t1092,t1094,t1096,t1101,t1102,t1103,t1104,t1105,t1108,t1112,t1113,t1118,t1119,t1120;
775   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;
776   PetscReal t1161,t1162,t1165,t1166,t1170,t1171,t1172,t1173,t1175,t1176,t1178,t1180,t1181,t1182,t1185,t1189,t1192,t1193,t1195,t1196,t1199;
777   PetscReal t1201,t1203,t1209,t1210,t1211,t1213,t1214,t1218,t1221,t1224,t1225,t1226,t1228,t1233,t1234,t1235,t1236,t1237,t1240;
778   PetscReal t1241,t1242,t1243,t1244,t1245,t1248,t1251,t1252,t1257,t1258,t1259,t1260,t1263,t1268,t1269,t1272,t1280;
779   PetscReal t1282,t1283,t1284,t1285,t1287,t1288,t1289,t1292,t1293,t1296,t1297,t1300,t1304,t1307,t1310,t1311,t1312,t1316,t1317,t1320;
780   PetscReal t1321,t1323,t1328,t1330,t1331,t1332,t1333,t1336,t1338,t1343,t1344,t1346,t1349,t1350,t1354;
781   PetscReal t1366,t1369,t1370,t1371,t1376,t1378,t1380,t1383,t1386,t1387,t1388,t1391,t1393,t1399;
782   PetscReal t1411,t1412,t1420,t1427;
783   PetscReal t1450,t1456,t1468,t1472,t1474,t1478;
784   PetscReal t1504,t1511;
785   PetscReal t1545;
786   PetscReal t1564,t1583;
787 
788   PetscReal      sum1 = 0.0, sum2 = 0.0, sum3 = 0.0, sum4 = 0.0, sum5 = 0.0, sum6 = 0.0;
789   PetscReal      ZA = etaA, ZB = etaB;
790   PetscInt       nz  = m,    nx = n;
791   PetscReal      u1, u2, u3, u4, u5, u6, Z, x = pos[0], z = pos[1];
792 
793   PetscFunctionBegin;
794   /* Note that there is no Fourier sum here. */
795   /****************************************************************************************/
796   _PC1A = 0;
797   /****************************************************************************************/
798   t1 = nx * 0.3141592654e1;
799   t2 = PetscSinReal(t1);
800   t3 = nx * t2;
801   t4 = nz * nz;
802   t5 = t4 * t4;
803   t6 = 0.3141592654e1 * 0.3141592654e1;
804   t8 = t3 * t5 * t6;
805   t9 = ZA * xc;
806   t12 = PetscExpReal(xc * nz * 0.3141592654e1);
807   t13 = t12 * t12;
808   t15 = nz * 0.3141592654e1;
809   t16 = PetscExpReal(t15);
810   t17 = t16 * t16;
811   t18 = t17 * t16;
812   t19 = ZB * t13 * t18;
813   t20 = t9 * t19;
814   t23 = ZA * ZA;
815   t24 = nx * nx;
816   t25 = t24 * nx;
817   t26 = t23 * t25;
818   t28 = t13 * t13;
819   t29 = t28 * t13;
820   t33 = nx * ZB;
821   t34 = t1 * xc;
822   t35 = PetscSinReal(t34);
823   t36 = t4 * nz;
824   t37 = t35 * t36;
825   t38 = t33 * t37;
826   t39 = 0.3141592654e1 * ZA;
827   t40 = t13 * t12;
828   t41 = t17 * t40;
829   t45 = ZB * ZB;
830   t46 = t45 * t24;
831   t47 = t46 * t4;
832   t48 = 0.3141592654e1 * xc;
833   t49 = t13 * t17;
834   t53 = xc * xc;
835   t54 = t36 * t53;
836   t56 = t54 * t6 * t45;
837   t57 = PetscCosReal(t34);
838   t58 = t57 * t24;
839   t59 = t28 * t12;
840   t60 = t17 * t59;
841   t61 = t58 * t60;
842   t64 = t25 * t2;
843   t65 = t64 * t15;
844   t72 = nx * t23;
845   t74 = t72 * t2 * t5;
846   t75 = t6 * t53;
847   t76 = t16 * t29;
848   t80 = t23 * nz;
849   t81 = t80 * 0.3141592654e1;
850   t82 = t18 * t28;
851   t86 = nx * t5;
852   t87 = t23 * t6;
853   t89 = xc * t2;
854   t90 = t13 * t18;
855   t91 = t89 * t90;
856   t94 = t28 * t28;
857   t96 = t24 * nz;
858   t98 = t4 * t45;
859   t99 = t98 * 0.3141592654e1;
860   t100 = t58 * t41;
861   t104 = 0.3141592654e1 * t25;
862   t105 = ZA * nz * t104;
863   t106 = t2 * ZB;
864   t110 = t17 * t17;
865   t111 = ZA * t110;
866   t116 = nz * t28;
867   t122 = t64 * t4 * t6;
868   t126 = t23 * t29 * t4;
869   t128 = t24 * xc;
870   t132 = t36 * t23;
871   t133 = t6 * t57;
872   t135 = t128 * t41;
873   t138 = t6 * xc;
874   t142 = t72 * t2;
875   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;
876   t149 = ZA * t25 * t2;
877   t150 = ZB * t28;
878   t154 = t35 * t5;
879   t155 = t72 * t154;
880   t156 = t75 * t41;
881   t159 = nx * ZA;
882   t160 = t2 * t36;
883   t161 = t159 * t160;
884   t162 = 0.3141592654e1 * ZB;
885   t163 = t28 * t16;
886   t167 = t23 * t57;
887   t168 = t167 * t24;
888   t169 = nz * t110;
889   t170 = t169 * t40;
890   t173 = ZA * ZB;
891   t174 = t173 * t90;
892   t177 = t36 * 0.3141592654e1;
893   t181 = t80 * t104;
894   t184 = nz * t17;
895   t188 = t17 * t29;
896   t190 = t4 * 0.3141592654e1;
897   t191 = t190 * t24;
898   t206 = t138 * t60;
899   t209 = t23 * t4;
900   t211 = t209 * t6 * t25;
901   t212 = t89 * t76;
902   t216 = ZB * t16 * t29;
903   t217 = t9 * t216;
904   t220 = ZB * t110;
905   t221 = ZA * t24;
906   t222 = t221 * nz;
907   t225 = t132 * t75;
908   t232 = t45 * t28;
909   t233 = t110 * t24;
910   t234 = t233 * nz;
911   t236 = t209 * 0.3141592654e1;
912   t237 = t17 * xc;
913   t239 = t237 * t13 * t24;
914   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;
915   t244 = nx * t45;
916   t245 = t244 * t37;
917   t246 = t110 * t40;
918   t251 = t237 * t59;
919   t256 = t64 * t90;
920   t260 = t36 * t45 * t133;
921   t263 = t45 * t57;
922   t264 = t263 * t24;
923   t265 = t169 * t12;
924   t269 = t6 * t36;
925   t270 = t17 * t24;
926   t274 = t110 * t13;
927   t276 = t190 * t128;
928   t279 = nx * t36;
929   t281 = t28 * t40;
930   t282 = t281 * t35;
931   t286 = t138 * t41;
932   t289 = t75 * t60;
933   t296 = t190 * t173;
934   t305 = t86 * t45 * t35;
935   t312 = t33 * t154;
936   t313 = t6 * ZA;
937   t324 = t232 * t270;
938   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;
939   t328 = t3 * t4;
940   t331 = t72 * t37;
941   t332 = t48 * t60;
942   t335 = nz * t94;
943   t345 = t72 * t35;
944   t349 = t173 * t57;
945   t355 = t53 * t17;
946   t364 = t54 * t6 * ZB;
947   t365 = t28 * t17;
948   t369 = xc * ZB;
949   t370 = t269 * t369;
950   t371 = ZA * t57;
951   t373 = t371 * t270 * t40;
952   t385 = nx * t35;
953   t396 = t4 * xc;
954   t397 = t396 * t162;
955   t415 = t37 * t48;
956   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;
957   t427 = t371 * t270 * t59;
958   t439 = t209 * t48;
959   t440 = t110 * t12;
960   t441 = t58 * t440;
961   t447 = t36 * xc;
962   t455 = t48 * t440;
963   t471 = ZB * t17;
964   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;
965   t499 = ZA * t53 * t19;
966   t505 = t396 * 0.3141592654e1;
967   t518 = t173 * t53 * t16 * t29;
968   t533 = t23 * t28;
969   t535 = t188 * t45;
970   t538 = t24 * t4;
971   t545 = t3 * t177;
972   t546 = t173 * t76;
973   t555 = t45 * t110;
974   t557 = t72 * t160;
975   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;
976   t567 = t396 * 0.3141592654e1 * t45;
977   t568 = t58 * t246;
978   t597 = t58 * nz;
979   t615 = t13 * t45;
980   t616 = t615 * t233;
981   t619 = t94 * t45;
982   t621 = t45 * t59;
983   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 * nz * 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;
984   t626 = t49 * t45;
985   t660 = t29 * t45;
986   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;
987   t690 = t6 * 0.3141592654e1;
988   t691 = ZA * t690;
989   t693 = t24 * t24;
990   t694 = t693 * xc;
991   t695 = t188 * t694;
992   t698 = t23 * ZA;
993   t699 = t698 * t690;
994   t700 = t699 * t5;
995   t704 = t5 * t4;
996   t705 = t691 * t704;
997   t709 = t691 * t5;
998   t713 = t5 * nz;
999   t714 = t713 * ZB;
1000   t718 = t698 * t6;
1001   t719 = t713 * t28;
1002   t722 = t699 * t704;
1003   t726 = t713 * t94;
1004   t733 = t713 * t45;
1005   t736 = t87 * t36;
1006   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;
1007   t746 = t313 * t36;
1008   t752 = t6 * t6;
1009   t753 = t23 * t752;
1010   t759 = t698 * t752;
1011   t760 = t759 * t36;
1012   t761 = t17 * t693;
1013   t762 = xc * t28;
1014   t763 = t761 * t762;
1015   t766 = t87 * t713;
1016   t773 = t699 * t4;
1017   t774 = t110 * t693;
1018   t775 = xc * t13;
1019   t785 = t704 * t17;
1020   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;
1021   t791 = t691 * t4;
1022   t792 = t45 * t693;
1023   t793 = t49 * t792;
1024   t796 = t759 * t713;
1025   t797 = t53 * t28;
1026   t798 = t270 * t797;
1027   t801 = t87 * nz;
1028   t818 = t5 * t36;
1029   t819 = t753 * t818;
1030   t827 = t753 * t36 * ZB;
1031   t830 = xc * t45;
1032   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;
1033   t842 = xc * t29;
1034   t843 = t270 * t842;
1035   t849 = t759 * t818;
1036   t853 = t691 * t396;
1037   t857 = t691 * t5 * t45;
1038   t869 = t313 * nz;
1039   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;
1040   t877 = ZA * t752;
1041   t878 = t877 * t818;
1042   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;
1043   t916 = t110 * t28;
1044   t945 = t28 * t693 * t45 * t17;
1045   t948 = 0.32e2 * t877 * t733 * t798 + 0.2e1 * t718 * t36 * t916 * t24 - 0.4e1 * t705 * t626 + t718 * nz * 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;
1046   t974 = t761 * t797;
1047   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;
1048   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;
1049 
1050   _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);
1051   /****************************************************************************************/
1052   t1 = nz * nz;
1053   t2 = t1 * nz;
1054   t3 = t2 * 0.3141592654e1;
1055   t4 = t3 * xc;
1056   t5 = ZB * ZB;
1057   t7 = PetscExpReal(nz * 0.3141592654e1);
1058   t8 = t7 * t7;
1059   t9 = t5 * t8;
1060   t12 = PetscExpReal(xc * nz * 0.3141592654e1);
1061   t13 = t12 * t12;
1062   t14 = t13 * t13;
1063   t15 = t14 * t13;
1064   t19 = nx * nx;
1065   t21 = nx * 0.3141592654e1;
1066   t22 = PetscSinReal(t21);
1067   t23 = t19 * nx * t22;
1068   t24 = t23 * 0.3141592654e1;
1069   t25 = ZA * ZB;
1070   t26 = t7 * t15;
1071   t27 = t25 * t26;
1072   t30 = t21 * xc;
1073   t31 = PetscSinReal(t30);
1074   t32 = t31 * nx;
1075   t33 = t32 * nz;
1076   t34 = ZA * ZA;
1077   t35 = t8 * t8;
1078   t36 = t34 * t35;
1079   t40 = t2 * t34;
1080   t41 = 0.3141592654e1 * t8;
1081   t42 = t41 * t15;
1082   t45 = t1 * t5;
1083   t46 = t14 * t14;
1084   t49 = t19 * t5;
1085   t51 = t19 * t46;
1086   t53 = t19 * t34;
1087   t55 = t8 * t7;
1088   t56 = t13 * t55;
1089   t57 = t25 * t56;
1090   t60 = t2 * nx;
1091   t61 = 0.3141592654e1 * 0.3141592654e1;
1092   t63 = t60 * t31 * t61;
1093   t64 = xc * xc;
1094   t65 = ZA * t64;
1095   t66 = ZB * t8;
1096   t67 = t14 * t12;
1097   t68 = t66 * t67;
1098   t69 = t65 * t68;
1099   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;
1100   t73 = t1 * nx;
1101   t75 = t73 * t31 * 0.3141592654e1;
1102   t76 = t8 * t67;
1103   t77 = t25 * t76;
1104   t80 = t1 * t1;
1105   t81 = t80 * t34;
1106   t83 = t61 * t14;
1107   t87 = t1 * t19;
1108   t88 = PetscCosReal(t30);
1109   t90 = t87 * t88 * t61;
1110   t91 = t5 * t64;
1111   t92 = t13 * t12;
1112   t93 = t8 * t92;
1113   t94 = t91 * t93;
1114   t100 = ZB * t64 * ZA * t8 * t92;
1115   t103 = nz * t19;
1116   t105 = t103 * t88 * 0.3141592654e1;
1117   t106 = ZA * xc;
1118   t107 = ZB * t35;
1119   t109 = t106 * t107 * t12;
1120   t112 = t34 * xc;
1121   t113 = t112 * t93;
1122   t116 = t35 * t14;
1123   t118 = t1 * ZA;
1124   t119 = ZB * t14;
1125   t122 = t1 * t46;
1126   t125 = t19 * ZB;
1127   t126 = t35 * ZA;
1128   t127 = t125 * t126;
1129   t129 = t1 * ZB;
1130   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;
1131   t134 = t1 * t34;
1132   t136 = t34 * t64;
1133   t137 = t136 * t76;
1134   t141 = t91 * t76;
1135   t145 = t103 * t34;
1136   t146 = 0.3141592654e1 * xc;
1137   t147 = t8 * t13;
1138   t153 = t14 * ZA;
1139   t156 = xc * t5;
1140   t157 = t156 * t93;
1141   t160 = t103 * t5;
1142   t162 = t146 * t8 * t15;
1143   t166 = t34 * t7 * t15;
1144   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;
1145   t170 = t106 * t68;
1146   t177 = t35 * t92;
1147   t178 = t112 * t177;
1148   t181 = t156 * t76;
1149   t186 = t35 * t12;
1150   t187 = t112 * t186;
1151   t193 = t5 * 0.3141592654e1;
1152   t206 = t34 * t14;
1153   t207 = t206 * t7;
1154   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 * nz * t34 * t67 + 0.4e1 * t24 * t207;
1155   t217 = t106 * t66 * t92;
1156   t226 = t88 * t19 * nz;
1157   t227 = 0.3141592654e1 * t34;
1158   t229 = t227 * xc * t67;
1159   t232 = t73 * t31;
1160   t234 = t146 * t5 * t67;
1161   t238 = t61 * ZB;
1162   t239 = t14 * t8;
1163   t240 = t238 * t239;
1164   t243 = t136 * t93;
1165   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;
1166   t247 = t14 * t92;
1167   t249 = t227 * t247 * xc;
1168   t254 = t73 * t22;
1169   t259 = t60 * t22 * t61;
1170   t260 = t112 * t26;
1171   t264 = t146 * t247 * t5;
1172   t268 = xc * t14;
1173   t274 = t5 * t14;
1174   t275 = t274 * t8;
1175   t280 = nz * nx;
1176   t281 = t280 * t22;
1177   t282 = t55 * t14;
1178   t283 = t25 * t282;
1179   t290 = ZA * t247 * xc * ZB;
1180   t295 = t22 * nx * t1 * 0.3141592654e1;
1181   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;
1182   t301 = t61 * t5;
1183   t307 = t87 * t34;
1184   t312 = t61 * xc;
1185   t313 = t312 * t239;
1186   t317 = t34 * t55 * t14;
1187   t329 = ZB * t13 * t55;
1188   t330 = t65 * t329;
1189   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;
1190   t340 = t5 * t35;
1191   t344 = t25 * t93;
1192   t356 = t41 * t13;
1193   t360 = t23 * nz * t61;
1194   t363 = t25 * t64 * t7 * t15;
1195   t366 = t156 * t177;
1196   t369 = t14 * t7;
1197   t370 = t25 * t369;
1198   t373 = t156 * t186;
1199   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;
1200   t382 = t112 * t76;
1201   t387 = t80 * t61;
1202   t391 = t136 * t26;
1203   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;
1204   t423 = t153 * t8;
1205   t426 = t34 * t13;
1206   t427 = t426 * t55;
1207   t430 = t34 * t8;
1208   t437 = t80 * ZA;
1209   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;
1210   t453 = t106 * ZB * t7 * t15;
1211   t456 = t2 * t5;
1212   t459 = t112 * t56;
1213   t462 = t126 * t14;
1214   t474 = t40 * 0.3141592654e1;
1215   t475 = xc * t8;
1216   t480 = t146 * t13 * t35;
1217   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;
1218   t497 = t136 * t56;
1219   t504 = t9 * t13;
1220   t509 = t475 * t13;
1221   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;
1222   t528 = t5 * t13;
1223   t529 = t528 * t35;
1224   t532 = t106 * t329;
1225   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;
1226   t557 = t146 * t15;
1227   t562 = xc * t15;
1228   t563 = t562 * t5;
1229   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;
1230   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;
1231   t623 = xc * t13;
1232   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;
1233   t632 = t61 * 0.3141592654e1;
1234   t633 = t632 * t8;
1235   t634 = t80 * nz;
1236   t638 = t632 * t634;
1237   t639 = t638 * xc;
1238   t642 = t61 * t34;
1239   t643 = t122 * t19;
1240   t649 = t61 * t61;
1241   t650 = t649 * t1;
1242   t652 = t19 * t19;
1243   t653 = t14 * t652;
1244   t654 = t653 * t9;
1245   t657 = t14 * t1;
1246   t658 = t657 * t19;
1247   t665 = t632 * t34;
1248   t666 = t665 * t2;
1249   t667 = t8 * t19;
1250   t668 = t667 * t623;
1251   t674 = t665 * nz;
1252   t675 = t652 * xc;
1253   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;
1254   t683 = t46 * t652;
1255   t686 = t633 * t15;
1256   t691 = t35 * t80;
1257   t698 = t35 * t652;
1258   t705 = t14 * t80;
1259   t708 = t61 * t35;
1260   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;
1261   t727 = t61 * t46;
1262   t737 = t649 * t34;
1263   t738 = t737 * t1;
1264   t739 = t8 * t652;
1265   t740 = t739 * t268;
1266   t746 = t61 * ZA;
1267   t754 = t632 * nz * xc;
1268   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;
1269   t762 = t633 * t13;
1270   t764 = t5 * nz * t652;
1271   t767 = t80 * t1;
1272   t768 = t649 * t767;
1273   t772 = t649 * ZA;
1274   t773 = t772 * t129;
1275   t777 = t35 * t1 * t19;
1276   t780 = t632 * t5;
1277   t781 = t780 * t15;
1278   t786 = t698 * ZA;
1279   t790 = t64 * t14;
1280   t800 = t649 * t8;
1281   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 * nz * 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;
1282   t812 = t649 * t80;
1283   t817 = t83 * ZB;
1284   t824 = t746 * ZB;
1285   t828 = t800 * t14;
1286   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;
1287   t858 = t35 * t19;
1288   t873 = t2 * t8;
1289   t878 = t61 * t1;
1290   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;
1291   t903 = t739 * t790;
1292   t923 = t737 * t80;
1293   t924 = t667 * t790;
1294   t927 = t780 * t2;
1295   t937 = t15 * t19 * xc;
1296   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;
1297   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;
1298 
1299   _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);
1300   /****************************************************************************************/
1301   _PC4A = 0;
1302   /****************************************************************************************/
1303   t1 = nx * 0.3141592654e1;
1304   t2 = t1 * xc;
1305   t3 = PetscCosReal(t2);
1306   t4 = nx * nx;
1307   t6 = nz * 0.3141592654e1;
1308   t7 = t3 * t4 * t6;
1309   t8 = ZA * ZB;
1310   t9 = PetscExpReal(t6);
1311   t10 = t9 * t9;
1312   t11 = xc * nz;
1313   t13 = PetscExpReal(t11 * 0.3141592654e1);
1314   t14 = t13 * t13;
1315   t15 = t14 * t13;
1316   t16 = t14 * t14;
1317   t17 = t16 * t15;
1318   t18 = t10 * t17;
1319   t19 = t8 * t18;
1320   t22 = PetscSinReal(t2);
1321   t23 = nx * t22;
1322   t24 = t23 * nz;
1323   t25 = ZB * ZB;
1324   t30 = nz * nz;
1325   t31 = t30 * nz;
1326   t32 = t31 * nx;
1327   t33 = 0.3141592654e1 * 0.3141592654e1;
1328   t35 = t32 * t22 * t33;
1329   t36 = ZA * ZA;
1330   t37 = t36 * xc;
1331   t38 = t16 * t13;
1332   t39 = t10 * t38;
1333   t40 = t37 * t39;
1334   t43 = PetscSinReal(t1);
1335   t44 = nx * t43;
1336   t45 = t30 * 0.3141592654e1;
1337   t46 = t44 * t45;
1338   t47 = ZA * xc;
1339   t49 = ZB * t16 * t9;
1340   t54 = t4 * nx * t43;
1341   t55 = xc * xc;
1342   t57 = t54 * t30 * t55;
1343   t58 = t33 * 0.3141592654e1;
1344   t59 = t58 * t25;
1345   t60 = t16 * t9;
1346   t61 = t59 * t60;
1347   t64 = xc * t25;
1348   t65 = t14 * t9;
1349   t66 = t64 * t65;
1350   t70 = t44 * t31 * t33;
1351   t71 = t37 * t65;
1352   t74 = t10 * t15;
1353   t75 = t64 * t74;
1354   t78 = t25 * t10;
1355   t83 = t54 * nz * t33;
1356   t84 = t55 * t25;
1357   t85 = t10 * t9;
1358   t86 = t14 * t85;
1359   t87 = t84 * t86;
1360   t90 = t30 * t30;
1361   t92 = t44 * t90 * t58;
1362   t93 = t55 * xc;
1363   t94 = t93 * t25;
1364   t95 = t85 * t16;
1365   t96 = t94 * t95;
1366   t102 = t23 * t45;
1367   t103 = t10 * t10;
1368   t104 = ZB * t103;
1369   t106 = t47 * t104 * t15;
1370   t111 = t54 * 0.3141592654e1;
1371   t112 = t25 * t85;
1372   t113 = t112 * t16;
1373   t115 = t8 * t39;
1374   t118 = t16 * t14;
1375   t119 = t85 * t118;
1376   t120 = t37 * t119;
1377   t123 = t16 * t16;
1378   t124 = t36 * t123;
1379   t125 = t124 * t9;
1380   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;
1381   t128 = t37 * t74;
1382   t131 = t44 * nz;
1383   t133 = t25 * t9 * t118;
1384   t136 = t36 * t14;
1385   t137 = t136 * t9;
1386   t140 = t30 * t4;
1387   t142 = t140 * t3 * t33;
1388   t143 = t64 * t39;
1389   t147 = t30 * nx * t43;
1390   t148 = 0.3141592654e1 * t36;
1391   t149 = t9 * t118;
1392   t153 = t44 * t31 * ZA;
1393   t154 = t33 * xc;
1394   t155 = t154 * t49;
1395   t160 = ZA * t17 * xc * ZB;
1396   t163 = t103 * t13;
1397   t164 = t64 * t163;
1398   t170 = t44 * t90 * t55;
1399   t171 = t58 * ZB;
1400   t172 = ZA * t16;
1401   t174 = t171 * t172 * t9;
1402   t177 = t36 * t55;
1403   t178 = t177 * t149;
1404   t181 = t54 * t11;
1405   t182 = t33 * t25;
1406   t186 = t25 * t14;
1407   t187 = t186 * t9;
1408   t193 = t186 * t85;
1409   t198 = ZB * t55;
1410   t199 = ZA * t103;
1411   t201 = t198 * t199 * t15;
1412   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;
1413   t210 = t32 * t22;
1414   t211 = t33 * t55;
1415   t212 = t25 * t38;
1416   t213 = t211 * t212;
1417   t216 = nz * nx;
1418   t217 = t22 * t25;
1419   t222 = ZB * t85 * t16;
1420   t226 = t23 * t30;
1421   t227 = t13 * t10;
1422   t228 = t148 * t227;
1423   t233 = t37 * t163;
1424   t237 = nz * t4 * t3;
1425   t238 = t148 * t74;
1426   t241 = t64 * t86;
1427   t245 = t148 * xc * t15;
1428   t248 = t112 * t118;
1429   t250 = t22 * t36;
1430   t256 = 0.3141592654e1 * t25;
1431   t257 = t256 * t39;
1432   t262 = t38 * t103;
1433   t263 = t37 * t262;
1434   t267 = t148 * t17 * xc;
1435   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;
1436   t273 = t37 * t149;
1437   t277 = t47 * t104 * t13;
1438   t285 = t31 * t36;
1439   t286 = t44 * t285;
1440   t291 = t25 * t123 * t9;
1441   t304 = 0.3141592654e1 * xc;
1442   t305 = t304 * t212;
1443   t312 = t256 * t18;
1444   t315 = t8 * t60;
1445   t319 = t54 * t30 * t58;
1446   t323 = t90 * t36;
1447   t324 = t44 * t323;
1448   t325 = t55 * t58;
1449   t326 = t325 * t60;
1450   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;
1451   t332 = t8 * t95;
1452   t335 = t136 * t85;
1453   t337 = t256 * t227;
1454   t340 = t177 * t119;
1455   t346 = t37 * t86;
1456   t351 = t103 * t15;
1457   t352 = t177 * t351;
1458   t355 = t64 * t119;
1459   t358 = t8 * t227;
1460   t361 = t85 * 0.3141592654e1;
1461   t365 = t84 * t39;
1462   t372 = ZB * t10;
1463   t373 = t372 * t38;
1464   t374 = t47 * t373;
1465   t379 = t177 * t39;
1466   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;
1467   t386 = t54 * t30 * t93;
1468   t387 = ZA * t85;
1469   t389 = t171 * t387 * t16;
1470   t394 = t64 * t60;
1471   t398 = t304 * t25 * t15;
1472   t401 = t361 * t25;
1473   t405 = t84 * t65;
1474   t410 = t148 * t18;
1475   t414 = t25 * t16 * t9;
1476   t417 = t84 * t74;
1477   t422 = t177 * t86;
1478   t428 = ZB * t38;
1479   t429 = t47 * t428;
1480   t432 = t148 * t39;
1481   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;
1482   t442 = t36 * t9 * t118;
1483   t444 = t123 * t9;
1484   t445 = t8 * t444;
1485   t448 = t361 * t36;
1486   t451 = t47 * t372 * t17;
1487   t454 = t94 * t60;
1488   t457 = t25 * t103;
1489   t465 = t47 * t372 * t15;
1490   t468 = t36 * t85;
1491   t469 = t468 * t16;
1492   t474 = t43 * t85;
1493   t478 = t8 * t74;
1494   t484 = t256 * t74;
1495   t489 = t198 * ZA * t10 * t15;
1496   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;
1497   t504 = t33 * t36 * t55 * t38;
1498   t507 = t37 * t18;
1499   t512 = t47 * ZB * t13 * t10;
1500   t518 = t59 * t95;
1501   t530 = t84 * t351;
1502   t534 = t37 * t227;
1503   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;
1504   t557 = t468 * t118;
1505   t562 = t93 * t58;
1506   t563 = t562 * t60;
1507   t567 = t44 * t90 * t93;
1508   t575 = ZA * t55;
1509   t576 = t575 * t428;
1510   t583 = t37 * t60;
1511   t590 = t140 * t3;
1512   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;
1513   t604 = nz * t36;
1514   t606 = t154 * t95;
1515   t625 = t36 * t103;
1516   t640 = t30 * t36;
1517   t641 = t54 * t640;
1518   t642 = t325 * t95;
1519   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;
1520   t652 = t36 * t16;
1521   t653 = t652 * t9;
1522   t655 = t64 * t227;
1523   t658 = t182 * t95;
1524   t663 = t562 * t95;
1525   t684 = t64 * t351;
1526   t689 = t36 * t10;
1527   t695 = t154 * t222;
1528   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;
1529   t711 = t575 * t373;
1530   t717 = t304 * t17 * t25;
1531   t736 = t177 * t74;
1532   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;
1533   t750 = t64 * t262;
1534   t763 = t44 * t640;
1535   t770 = t84 * t119;
1536   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;
1537   t797 = t64 * t149;
1538   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;
1539   t827 = t84 * t149;
1540   t846 = t54 * nz * ZA;
1541   t854 = t64 * t18;
1542   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;
1543   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;
1544   t914 = t25 * ZB;
1545   t915 = t33 * t914;
1546   t919 = t4 * t4;
1547   t920 = t16 * t919;
1548   t929 = t123 * t90;
1549   t932 = t919 * t103;
1550   t935 = t33 * ZB;
1551   t939 = t652 * t919;
1552   t942 = t16 * t30;
1553   t943 = t942 * t4;
1554   t949 = t103 * t16;
1555   t950 = t949 * t90;
1556   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;
1557   t954 = t171 * t36;
1558   t955 = t90 * nz;
1559   t956 = xc * t955;
1560   t957 = t118 * t10;
1561   t964 = t33 * t33;
1562   t965 = t964 * ZB;
1563   t966 = t965 * t640;
1564   t967 = t10 * t919;
1565   t968 = t55 * t16;
1566   t969 = t967 * t968;
1567   t972 = t935 * t36;
1568   t974 = t103 * t30 * t4;
1569   t977 = xc * t16;
1570   t978 = t967 * t977;
1571   t981 = t90 * t30;
1572   t983 = t16 * t10;
1573   t987 = t182 * ZA;
1574   t988 = t4 * t10;
1575   t992 = t171 * t604;
1576   t993 = xc * t14;
1577   t994 = t932 * t993;
1578   t997 = t182 * t30;
1579   t1005 = t171 * t285;
1580   t1006 = t988 * t993;
1581   t1009 = t58 * t914;
1582   t1010 = t1009 * t31;
1583   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;
1584   t1015 = t964 * t25;
1585   t1016 = ZA * t30;
1586   t1017 = t1015 * t1016;
1587   t1020 = t967 * t993;
1588   t1031 = t1009 * t118;
1589   t1032 = t31 * t10;
1590   t1040 = t964 * t914;
1591   t1041 = t1040 * t90;
1592   t1044 = t55 * t10 * t4 * t16;
1593   t1047 = t1040 * t30;
1594   t1050 = t123 * ZA;
1595   t1054 = t977 * t988;
1596   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;
1597   t1058 = t1009 * nz;
1598   t1063 = t932 * ZA;
1599   t1069 = t123 * t30 * t4;
1600   t1080 = t993 * t103 * t4;
1601   t1088 = t935 * t103;
1602   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;
1603   t1113 = t955 * t118 * xc;
1604   t1120 = t4 * t118;
1605   t1125 = t981 * xc;
1606   t1133 = nz * t10;
1607   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;
1608   t1153 = xc * t118;
1609   t1165 = t182 * t16;
1610   t1170 = t171 * t10;
1611   t1178 = ZA * t90;
1612   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 * nz * 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;
1613   t1195 = t1040 * t981;
1614   t1199 = t1009 * t955;
1615   t1203 = t1009 * t10;
1616   t1211 = t965 * t323;
1617   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 * nz * 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;
1618   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 * nz * 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;
1619 
1620   _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);
1621   /****************************************************************************************/
1622   t1 = nz * nz;
1623   t2 = t1 * nz;
1624   t3 = nx * t2;
1625   t4 = 0.3141592654e1 * ZA;
1626   t5 = t3 * t4;
1627   t6 = nx * 0.3141592654e1;
1628   t7 = t6 * xc;
1629   t8 = PetscSinReal(t7);
1630   t9 = t8 * ZB;
1631   t10 = nz * 0.3141592654e1;
1632   t11 = PetscExpReal(t10);
1633   t12 = t11 * t11;
1634   t15 = PetscExpReal(xc * nz * 0.3141592654e1);
1635   t16 = t15 * t15;
1636   t17 = t16 * t16;
1637   t18 = t17 * t15;
1638   t19 = t12 * t18;
1639   t23 = t1 * t1;
1640   t24 = nx * t23;
1641   t25 = ZB * ZB;
1642   t27 = t18 * t8;
1643   t28 = 0.3141592654e1 * 0.3141592654e1;
1644   t29 = xc * xc;
1645   t30 = t28 * t29;
1646   t34 = t1 * xc;
1647   t35 = 0.3141592654e1 * ZB;
1648   t36 = t34 * t35;
1649   t37 = PetscCosReal(t7);
1650   t38 = ZA * t37;
1651   t39 = nx * nx;
1652   t40 = t39 * t12;
1653   t41 = t16 * t15;
1654   t43 = t38 * t40 * t41;
1655   t46 = t25 * nz;
1656   t47 = t46 * 0.3141592654e1;
1657   t48 = t39 * nx;
1658   t49 = PetscSinReal(t6);
1659   t50 = t48 * t49;
1660   t51 = t12 * t11;
1661   t52 = t51 * t17;
1662   t53 = t50 * t52;
1663   t56 = t34 * 0.3141592654e1 * t25;
1664   t57 = t37 * t39;
1665   t58 = t17 * t41;
1666   t59 = t12 * t58;
1667   t60 = t57 * t59;
1668   t63 = t25 * t18;
1669   t64 = t57 * nz;
1670   t67 = ZA * ZA;
1671   t68 = t67 * nz;
1672   t69 = 0.3141592654e1 * t48;
1673   t70 = t68 * t69;
1674   t71 = t49 * xc;
1675   t72 = t17 * t16;
1676   t73 = t11 * t72;
1677   t74 = t71 * t73;
1678   t77 = t1 * t67;
1679   t78 = t77 * 0.3141592654e1;
1680   t81 = nx * t25;
1681   t82 = t81 * t49;
1682   t83 = t17 * t17;
1683   t85 = t1 * t83 * t11;
1684   t87 = nx * ZB;
1685   t88 = t8 * t2;
1686   t89 = t87 * t88;
1687   t90 = 0.3141592654e1 * xc;
1688   t91 = t12 * t12;
1689   t92 = ZA * t91;
1690   t97 = ZB * ZA;
1691   t98 = t97 * t37;
1692   t99 = t39 * nz;
1693   t100 = t12 * t41;
1694   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;
1695   t105 = t67 * t48;
1696   t106 = t49 * t51;
1697   t107 = t106 * t72;
1698   t109 = t1 * 0.3141592654e1;
1699   t110 = t109 * xc;
1700   t115 = nx * t67;
1701   t116 = t115 * t49;
1702   t117 = t1 * t16;
1703   t118 = t117 * t11;
1704   t120 = t2 * t25;
1705   t121 = t28 * 0.3141592654e1;
1706   t122 = t121 * t29;
1707   t123 = t120 * t122;
1708   t129 = t1 * ZB;
1709   t130 = t129 * t4;
1710   t131 = t57 * t100;
1711   t134 = t12 * t16;
1712   t136 = t109 * t39;
1713   t139 = ZB * t18;
1714   t141 = t39 * t1;
1715   t142 = t141 * t90;
1716   t145 = t77 * t90;
1717   t146 = t91 * t41;
1718   t147 = t57 * t146;
1719   t151 = t25 * t39 * t1;
1720   t152 = t72 * t12;
1721   t156 = t49 * t2;
1722   t158 = t83 * t11;
1723   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;
1724   t164 = t115 * t88;
1725   t165 = t90 * t19;
1726   t168 = t25 * t48;
1727   t169 = t49 * t16;
1728   t170 = t169 * t11;
1729   t174 = ZA * nz * t69;
1730   t175 = ZB * t51;
1731   t176 = t175 * t17;
1732   t177 = t71 * t176;
1733   t180 = t1 * t29;
1734   t181 = t28 * t25;
1735   t182 = t180 * t181;
1736   t183 = t50 * t73;
1737   t186 = ZA * t1;
1738   t187 = t28 * t48;
1739   t188 = t186 * t187;
1740   t189 = ZB * t17;
1741   t190 = t189 * t11;
1742   t191 = t71 * t190;
1743   t194 = t50 * t158;
1744   t196 = t115 * t156;
1745   t197 = t90 * t73;
1746   t201 = t49 * t17 * t11;
1747   t204 = t88 * t90;
1748   t207 = t68 * 0.3141592654e1;
1749   t208 = t17 * t11;
1750   t209 = t50 * t208;
1751   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;
1752   t212 = t2 * 0.3141592654e1;
1753   t213 = t212 * t52;
1754   t215 = t81 * t8;
1755   t216 = t212 * t59;
1756   t219 = t3 * t90;
1757   t220 = t25 * t8;
1758   t221 = t18 * t91;
1759   t225 = t71 * t52;
1760   t231 = t16 * t51;
1761   t232 = t50 * t231;
1762   t237 = ZA * t12;
1763   t243 = t67 * t28;
1764   t244 = t24 * t243;
1765   t245 = t71 * t231;
1766   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;
1767   t252 = t25 * t37;
1768   t253 = t252 * t39;
1769   t255 = nz * t15 * t12;
1770   t258 = t2 * t29;
1771   t259 = ZB * t28;
1772   t260 = t258 * t259;
1773   t263 = t106 * t17;
1774   t265 = xc * t25;
1775   t269 = t25 * t49;
1776   t270 = t269 * t52;
1777   t273 = t1 * t25;
1778   t274 = t273 * 0.3141592654e1;
1779   t275 = t57 * t19;
1780   t278 = t24 * t30;
1781   t288 = t1 * t11 * t72;
1782   t290 = t212 * t208;
1783   t292 = t2 * xc;
1784   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 * nz * t69 * t49 * ZA * t158 - t82 * t288 - t116 * t290 + 0.16e2 * t292 * t243 * t275;
1785   t301 = t50 * t176;
1786   t304 = t51 * t72;
1787   t305 = t71 * t304;
1788   t308 = t25 * t41;
1789   t311 = ZA * t48;
1790   t312 = t311 * t49;
1791   t317 = t91 * t15;
1792   t318 = t57 * t317;
1793   t321 = t81 * t88;
1794   t322 = t90 * t59;
1795   t325 = t212 * t231;
1796   t327 = t15 * t12;
1797   t328 = t57 * t327;
1798   t331 = t77 * t187;
1799   t334 = t2 * ZA;
1800   t335 = t334 * t122;
1801   t336 = t50 * t190;
1802   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;
1803   t341 = t169 * t51;
1804   t344 = t49 * t11 * t72;
1805   t346 = t77 * t30;
1806   t347 = t50 * t304;
1807   t350 = t25 * t51;
1808   t352 = nx * ZA;
1809   t353 = t49 * t23;
1810   t354 = t352 * t353;
1811   t355 = t28 * xc;
1812   t362 = t25 * t91;
1813   t365 = t23 * nz;
1814   t366 = nx * t365;
1815   t367 = t366 * t122;
1816   t368 = ZB * t49;
1817   t369 = ZA * t51;
1818   t370 = t369 * t17;
1819   t371 = t368 * t370;
1820   t374 = t115 * t353;
1821   t375 = t355 * t73;
1822   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;
1823   t382 = t71 * t208;
1824   t394 = t2 * t67;
1825   t395 = t394 * t122;
1826   t398 = t352 * t156;
1827   t402 = t17 * t12;
1828   t403 = t39 * ZA;
1829   t404 = t402 * t403;
1830   t407 = t269 * t208;
1831   t411 = t49 * t83 * t11;
1832   t413 = t46 * t69;
1833   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;
1834   t424 = t16 * t11;
1835   t425 = t212 * t424;
1836   t427 = t258 * t181;
1837   t430 = t67 * t29;
1838   t431 = t366 * t430;
1839   t432 = t121 * t49;
1840   t433 = t432 * t52;
1841   t436 = nz * t12;
1842   t437 = t436 * t18;
1843   t440 = t29 * xc;
1844   t441 = t440 * t121;
1845   t442 = t394 * t441;
1846   t445 = t67 * t37;
1847   t446 = t445 * t39;
1848   t448 = nz * t18 * t91;
1849   t453 = t352 * t49;
1850   t458 = t8 * t23;
1851   t462 = t81 * t458;
1852   t463 = t30 * t19;
1853   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;
1854   t470 = t436 * t58;
1855   t475 = t2 * t121 * t440 * t25;
1856   t485 = t212 * t73;
1857   t488 = t67 * t72 * t1;
1858   t490 = t39 * xc;
1859   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;
1860   t503 = t115 * t8;
1861   t510 = t355 * t19;
1862   t513 = t432 * t208;
1863   t525 = t38 * t40 * t18;
1864   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;
1865   t536 = t28 * t37;
1866   t538 = t490 * t100;
1867   t541 = t334 * t441;
1868   t547 = t394 * t30;
1869   t550 = t212 * t19;
1870   t553 = t366 * t441;
1871   t556 = nz * t17;
1872   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;
1873   t583 = t90 * t317;
1874   t594 = t212 * t158;
1875   t596 = t152 * t67;
1876   t602 = t67 * t17;
1877   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;
1878   t613 = t83 * t25;
1879   t616 = t81 * t156;
1880   t627 = t90 * t231;
1881   t630 = t91 * t16;
1882   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;
1883   t643 = t77 * t355;
1884   t648 = t115 * t458;
1885   t651 = t134 * t67;
1886   t657 = t30 * t304;
1887   t660 = t30 * t146;
1888   t665 = t25 * t17;
1889   t668 = t50 * t424;
1890   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;
1891   t672 = t90 * t424;
1892   t676 = t120 * t536;
1893   t680 = t436 * t41;
1894   t688 = t366 * t67 * t440;
1895   t696 = xc * t12;
1896   t697 = t696 * t18;
1897   t701 = t252 * t141;
1898   t702 = t90 * t221;
1899   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;
1900   t712 = t24 * t430;
1901   t713 = t28 * t49;
1902   t721 = t1 * t17 * t11;
1903   t726 = ZB * xc;
1904   t737 = nz * t91;
1905   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;
1906   t745 = t269 * t231;
1907   t749 = t1 * t28 * t265;
1908   t757 = t16 * t39;
1909   t758 = t696 * t757;
1910   t762 = t69 * t49;
1911   t772 = t355 * t100;
1912   t775 = t81 * t353;
1913   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 * nz * t762 - 0.6e1 * t87 * t49 * t186 * t52 - 0.8e1 * t553 * t407 - 0.4e1 * t749 * t209 + 0.16e2 * t648 * t772 - 0.6e1 * t775 * t375;
1914   t790 = t212 * t304;
1915   t793 = t156 * 0.3141592654e1;
1916   t795 = t355 * t304;
1917   t800 = t91 * t39;
1918   t801 = t800 * nz;
1919   t807 = t2 * t28;
1920   t808 = t807 * t726;
1921   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;
1922   t812 = t117 * t51;
1923   t821 = t24 * t355;
1924   t827 = t90 * t304;
1925   t840 = t800 * t41;
1926   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;
1927   t849 = t1 * t51;
1928   t850 = t849 * t17;
1929   t860 = t269 * t424;
1930   t863 = t273 * t187;
1931   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;
1932   t879 = xc * ZA;
1933   t888 = t67 * t51;
1934   t901 = ZA * t17;
1935   t903 = t368 * t901 * t11;
1936   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;
1937   t923 = t71 * t424;
1938   t926 = t87 * t458;
1939   t927 = t28 * ZA;
1940   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;
1941   t948 = t90 * t327;
1942   t961 = t30 * t100;
1943   t964 = t29 * t49;
1944   t981 = t106 * t1;
1945   t983 = -0.2e1 * t219 * t220 * t100 + 0.2e1 * t321 * t948 - 0.16e2 * t189 * ZA * t99 * t12 - 0.2e1 * t369 * nz * 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 * nz * t762 + t115 * t981;
1946   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;
1947   t1018 = t212 * t100;
1948   t1029 = t212 * t327;
1949   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;
1950   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;
1951   t1079 = t58 * t8;
1952   t1091 = t807 * t29;
1953   t1092 = t665 * t40;
1954   t1101 = ZB * t91;
1955   t1102 = t403 * nz;
1956   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;
1957   t1108 = t849 * t72;
1958   t1121 = t737 * t15;
1959   t1124 = t29 * t12;
1960   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;
1961   t1134 = t394 * t28;
1962   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;
1963   t1161 = nz * t83;
1964   t1173 = ZB * t12;
1965   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;
1966   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;
1967   t1234 = t25 * ZB;
1968   t1235 = t1234 * t28;
1969   t1236 = t365 * t91;
1970   t1240 = ZB * t121;
1971   t1241 = t1240 * t77;
1972   t1242 = t39 * t39;
1973   t1243 = t12 * t1242;
1974   t1244 = xc * t72;
1975   t1245 = t1243 * t1244;
1976   t1248 = t365 * t25;
1977   t1252 = t243 * nz;
1978   t1257 = t23 * t1;
1979   t1258 = t1240 * t1257;
1980   t1259 = t67 * t12;
1981   t1260 = xc * t16;
1982   t1268 = t1234 * t121;
1983   t1269 = t1268 * t23;
1984   t1272 = t1242 * t91;
1985   t1280 = t67 * xc;
1986   t1284 = t28 * t28;
1987   t1285 = t67 * t1284;
1988   t1287 = t1285 * t2 * ZB;
1989   t1288 = t17 * xc;
1990   t1289 = t1243 * t1288;
1991   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;
1992   t1293 = t365 * t83;
1993   t1300 = ZA * t1284;
1994   t1304 = t17 * t1242 * t25 * t12;
1995   t1307 = t927 * t2;
1996   t1311 = t23 * t2;
1997   t1312 = t1300 * t1311;
1998   t1316 = t1234 * t1284;
1999   t1317 = t1316 * t1311;
2000   t1321 = t1240 * t23;
2001   t1331 = t1240 * t23 * t67;
2002   t1332 = t40 * t1244;
2003   t1338 = t1243 * t1260;
2004   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 * nz * t1272 * t17 + 0.16e2 * t1331 * t1332 + 0.64e2 * t1300 * t292 * t1304 - 0.8e1 * t1241 * t1338 - 0.2e1 * t243 * t1293 * ZB;
2005   t1346 = t1316 * t2;
2006   t1349 = t927 * nz;
2007   t1350 = t25 * t1242;
2008   t1354 = t1268 * t1257;
2009   t1366 = t1268 * t1;
2010   t1370 = t29 * t17;
2011   t1371 = t1243 * t1370;
2012   t1386 = -0.32e2 * t1346 * t1289 + 0.4e1 * t1349 * t1350 * t91 + 0.8e1 * t1354 * t1260 * t12 - 0.16e2 * t181 * nz * 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;
2013   t1387 = t1257 * t12;
2014   t1391 = t1240 * t1;
2015   t1399 = t1272 * t1260;
2016   t1412 = t1285 * t1311;
2017   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;
2018   t1456 = t365 * ZB;
2019   t1468 = 0.4e1 * t927 * t1248 * t17 - 0.2e1 * t1235 * nz * 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;
2020   t1472 = t1316 * t365;
2021   t1474 = t1124 * t39 * t17;
2022   t1478 = t17 * t91;
2023   t1504 = t72 * t39;
2024   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;
2025   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;
2026   t1564 = t1260 * t800;
2027   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;
2028 
2029   _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);
2030   /****************************************************************************************/
2031   t1 = nz * nz;
2032   t2 = t1 * nz;
2033   t3 = t2 * nx;
2034   t4 = nx * 0.3141592654e1;
2035   t5 = t4 * xc;
2036   t6 = PetscSinReal(t5);
2037   t7 = 0.3141592654e1 * 0.3141592654e1;
2038   t9 = t3 * t6 * t7;
2039   t10 = xc * xc;
2040   t11 = ZA * ZA;
2041   t12 = t10 * t11;
2042   t13 = nz * 0.3141592654e1;
2043   t14 = PetscExpReal(t13);
2044   t15 = t14 * t14;
2045   t16 = xc * nz;
2046   t18 = PetscExpReal(t16 * 0.3141592654e1);
2047   t19 = t18 * t18;
2048   t20 = t19 * t18;
2049   t21 = t15 * t20;
2050   t22 = t12 * t21;
2051   t25 = nx * t6;
2052   t26 = t1 * 0.3141592654e1;
2053   t27 = t25 * t26;
2054   t28 = ZA * ZB;
2055   t29 = t18 * t15;
2056   t30 = t28 * t29;
2057   t33 = t25 * nz;
2058   t34 = t11 * t15;
2059   t35 = t19 * t19;
2060   t36 = t35 * t18;
2061   t40 = t25 * t1;
2062   t41 = 0.3141592654e1 * t11;
2063   t42 = t15 * t36;
2064   t43 = t41 * t42;
2065   t46 = nx * nx;
2066   t47 = t1 * t46;
2067   t48 = t47 * t11;
2068   t49 = t7 * xc;
2069   t50 = t35 * t15;
2070   t51 = t49 * t50;
2071   t55 = PetscSinReal(t4);
2072   t56 = t46 * nx * t55;
2073   t58 = t56 * nz * t7;
2074   t59 = ZB * ZB;
2075   t60 = t10 * t59;
2076   t61 = t15 * t14;
2077   t62 = t19 * t61;
2078   t63 = t60 * t62;
2079   t66 = t19 * t14;
2080   t67 = t60 * t66;
2081   t70 = t28 * t42;
2082   t73 = PetscCosReal(t5);
2083   t74 = t47 * t73;
2084   t75 = t7 * t11;
2085   t77 = t75 * t10 * t36;
2086   t80 = t73 * t46;
2087   t81 = t80 * nz;
2088   t82 = 0.3141592654e1 * t59;
2089   t83 = t82 * t42;
2090   t87 = xc * t11;
2091   t88 = t87 * t62;
2092   t91 = nz * nx;
2093   t92 = t55 * t61;
2094   t96 = nx * t55;
2095   t98 = t96 * t2 * t7;
2096   t101 = xc * t59;
2097   t102 = t101 * t62;
2098   t108 = t1 * t1;
2099   t109 = t108 * t7;
2100   t111 = t59 * t35;
2101   t112 = t111 * t15;
2102   t115 = t35 * t20;
2103   t123 = t1 * nx * t55;
2104   t124 = t61 * t35;
2105   t127 = t35 * t19;
2106   t128 = t61 * t127;
2107   t129 = t60 * t128;
2108   t132 = t56 * t16;
2109   t133 = t7 * t59;
2110   t134 = t133 * t124;
2111   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;
2112   t139 = t56 * 0.3141592654e1;
2113   t140 = t111 * t14;
2114   t144 = t49 * t124;
2115   t147 = t91 * t55;
2116   t148 = t61 * ZA;
2117   t154 = ZA * t115 * xc * ZB;
2118   t157 = t7 * 0.3141592654e1;
2119   t159 = t96 * t108 * t157;
2120   t160 = t10 * xc;
2121   t161 = t160 * t59;
2122   t162 = t35 * t14;
2123   t163 = t161 * t162;
2124   t166 = t28 * t162;
2125   t169 = t80 * t13;
2126   t170 = t101 * t42;
2127   t173 = t2 * t11;
2128   t174 = t96 * t173;
2129   t175 = t7 * t10;
2130   t179 = t59 * t15;
2131   t184 = t15 * t15;
2132   t193 = t139 * t140 + 0.4e1 * t56 * nz * 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;
2133   t194 = t173 * 0.3141592654e1;
2134   t195 = xc * t15;
2135   t196 = t195 * t19;
2136   t202 = t15 * t115;
2137   t203 = t28 * t202;
2138   t206 = t96 * t26;
2139   t207 = t14 * t127;
2140   t208 = t101 * t207;
2141   t211 = t12 * t128;
2142   t218 = t11 * t61;
2143   t219 = t218 * t35;
2144   t221 = t108 * ZA;
2145   t223 = t7 * ZB;
2146   t224 = t223 * t50;
2147   t227 = ZA * xc;
2148   t228 = ZB * t15;
2149   t229 = t228 * t36;
2150   t230 = t227 * t229;
2151   t233 = t87 * t207;
2152   t236 = t6 * t11;
2153   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;
2154   t244 = t227 * t228 * t20;
2155   t252 = t184 * t18;
2156   t253 = t101 * t252;
2157   t256 = t35 * t35;
2158   t257 = t256 * t14;
2159   t258 = t28 * t257;
2160   t261 = t108 * t11;
2161   t263 = t7 * t35;
2162   t268 = ZB * t61 * t35;
2163   t273 = t96 * t108 * t160;
2164   t274 = t157 * ZB;
2165   t276 = t274 * t148 * t35;
2166   t279 = t101 * t21;
2167   t282 = 0.3141592654e1 * xc;
2168   t283 = t59 * t36;
2169   t284 = t282 * t283;
2170   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;
2171   t290 = t1 * t11;
2172   t291 = t96 * t290;
2173   t297 = t59 * t61;
2174   t298 = t297 * t127;
2175   t300 = ZB * t36;
2176   t301 = t227 * t300;
2177   t304 = t1 * t59;
2178   t305 = t184 * t35;
2179   t310 = t46 * ZB;
2180   t311 = t184 * ZA;
2181   t312 = t310 * t311;
2182   t314 = t60 * t21;
2183   t317 = t1 * ZA;
2184   t318 = ZB * t35;
2185   t321 = t1 * t256;
2186   t324 = t96 * t261;
2187   t325 = t10 * t157;
2188   t326 = t325 * t124;
2189   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;
2190   t331 = t28 * t124;
2191   t334 = 0.3141592654e1 * t15;
2192   t335 = t334 * t127;
2193   t338 = t35 * ZA;
2194   t341 = t46 * t256;
2195   t344 = t46 * t11;
2196   t346 = t46 * t59;
2197   t348 = t297 * t35;
2198   t351 = ZA * t10;
2199   t352 = t351 * t300;
2200   t355 = t1 * ZB;
2201   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;
2202   t363 = t36 * t184;
2203   t364 = t87 * t363;
2204   t368 = t47 * t73 * t7;
2205   t373 = t160 * t157;
2206   t374 = t373 * t124;
2207   t377 = t311 * t35;
2208   t380 = t12 * t62;
2209   t386 = ZB * t10 * ZA * t15 * t20;
2210   t389 = t87 * t66;
2211   t393 = t56 * t1 * t10;
2212   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;
2213   t412 = t11 * t14 * t127;
2214   t416 = t11 * t19;
2215   t417 = t416 * t61;
2216   t421 = t96 * t2 * ZA;
2217   t426 = t56 * nz * ZA;
2218   t427 = t318 * t14;
2219   t428 = t49 * t427;
2220   t431 = t82 * t29;
2221   t434 = t87 * t21;
2222   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;
2223   t443 = t184 * t20;
2224   t447 = t61 * 0.3141592654e1;
2225   t448 = t447 * t11;
2226   t450 = t49 * t268;
2227   t453 = t60 * t42;
2228   t456 = t41 * t202;
2229   t463 = t101 * t443;
2230   t469 = t41 * xc * t20;
2231   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;
2232   t482 = t108 * t10;
2233   t492 = nz * t46;
2234   t493 = t492 * t11;
2235   t495 = t282 * t19 * t184;
2236   t498 = t56 * t290;
2237   t499 = t325 * t162;
2238   t502 = t416 * t14;
2239   t504 = t60 * t207;
2240   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;
2241   t508 = t101 * t162;
2242   t512 = t41 * t115 * xc;
2243   t515 = t87 * t42;
2244   t520 = ZB * t184;
2245   t522 = t227 * t520 * t18;
2246   t525 = t492 * t59;
2247   t528 = t6 * t59;
2248   t532 = t520 * t20;
2249   t533 = t351 * t532;
2250   t539 = t447 * t59;
2251   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;
2252   t551 = t87 * t162;
2253   t554 = t351 * t229;
2254   t560 = t59 * t19;
2255   t561 = t560 * t14;
2256   t564 = t101 * t202;
2257   t567 = t87 * t252;
2258   t573 = t227 * t228 * t115;
2259   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;
2260   t579 = t82 * t202;
2261   t591 = t282 * t115 * t59;
2262   t598 = t101 * t66;
2263   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;
2264   t608 = t373 * t162;
2265   t611 = t59 * t184;
2266   t617 = t101 * t29;
2267   t624 = t227 * ZB * t18 * t15;
2268   t629 = t157 * t59;
2269   t630 = t629 * t124;
2270   t633 = t3 * t6;
2271   t634 = t175 * t283;
2272   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;
2273   t645 = t41 * t29;
2274   t648 = t2 * 0.3141592654e1;
2275   t649 = t648 * xc;
2276   t650 = t560 * t184;
2277   t656 = t56 * t1 * t157;
2278   t659 = t87 * t128;
2279   t662 = t96 * t482;
2280   t663 = t629 * t162;
2281   t671 = t161 * t124;
2282   t674 = t218 * t127;
2283   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;
2284   t714 = t334 * t19;
2285   t719 = t12 * t42;
2286   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;
2287   t730 = xc * t19;
2288   t735 = t59 * t256 * t14;
2289   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;
2290   t755 = t28 * t21;
2291   t759 = t274 * t338 * t14;
2292   t767 = t11 * t35;
2293   t768 = t767 * t14;
2294   t778 = t560 * t61;
2295   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;
2296   t786 = t12 * t443;
2297   t790 = t282 * t59 * t20;
2298   t796 = t59 * t14 * t127;
2299   t806 = t41 * t21;
2300   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;
2301   t819 = t282 * t15 * t127;
2302   t822 = t101 * t363;
2303   t830 = t11 * t256 * t14;
2304   t835 = t227 * t532;
2305   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;
2306   t845 = t87 * t202;
2307   t854 = t338 * t15;
2308   t859 = t12 * t207;
2309   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;
2310   t880 = t60 * t443;
2311   t891 = t101 * t128;
2312   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;
2313   t903 = t179 * t19;
2314   t920 = t56 * t1 * t160;
2315   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;
2316   t932 = t87 * t29;
2317   t945 = t82 * t21;
2318   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;
2319   t971 = t282 * t127;
2320   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;
2321   t997 = xc * t127;
2322   t998 = t997 * t59;
2323   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;
2324   t1006 = t2 * t59;
2325   t1017 = xc * t35;
2326   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;
2327   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;
2328   t1070 = t15 * t19;
2329   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;
2330   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;
2331   t1119 = t321 * t46;
2332   t1122 = t157 * t11;
2333   t1123 = t1122 * t2;
2334   t1124 = t184 * t46;
2335   t1128 = t108 * nz;
2336   t1132 = t7 * t7;
2337   t1133 = t1132 * t11;
2338   t1134 = t1133 * t108;
2339   t1135 = t15 * t46;
2340   t1139 = t7 * ZA;
2341   t1140 = t1139 * ZB;
2342   t1141 = t1 * t35;
2343   t1145 = t629 * t2;
2344   t1146 = t1135 * t730;
2345   t1149 = t157 * t1128;
2346   t1150 = t1149 * xc;
2347   t1153 = t46 * xc;
2348   t1154 = t1153 * t127;
2349   t1158 = t184 * t1 * t46;
2350   t1161 = t46 * t46;
2351   t1162 = t35 * t1161;
2352   t1166 = t7 * t1;
2353   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;
2354   t1171 = t1161 * t184;
2355   t1175 = t1122 * nz;
2356   t1176 = t15 * t1161;
2357   t1180 = t1132 * ZA;
2358   t1181 = t1180 * t355;
2359   t1182 = t1176 * t1017;
2360   t1185 = t1161 * xc;
2361   t1189 = t1133 * t1;
2362   t1192 = t108 * t1;
2363   t1193 = t1132 * t1192;
2364   t1195 = t10 * t35;
2365   t1199 = t157 * t15;
2366   t1203 = t1141 * t46;
2367   t1211 = t184 * t108;
2368   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;
2369   t1221 = t256 * t1161;
2370   t1224 = t35 * t108;
2371   t1233 = t7 * t256;
2372   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;
2373   t1240 = t1135 * t1195;
2374   t1252 = t629 * t127;
2375   t1263 = t1171 * ZA;
2376   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 * nz * 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;
2377   t1283 = t1199 * t19;
2378   t1287 = t1199 * t127;
2379   t1289 = t59 * nz * t1161;
2380   t1293 = t157 * nz * xc;
2381   t1304 = t1132 * t108;
2382   t1310 = t263 * ZB;
2383   t1316 = t2 * t15;
2384   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;
2385   t1330 = t1132 * t15;
2386   t1336 = t1132 * t1;
2387   t1338 = t1162 * t179;
2388   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;
2389   t1376 = t7 * t184;
2390   t1380 = t1176 * t1195;
2391   t1393 = t1330 * t35;
2392   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;
2393   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;
2394   _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);
2395   /****************************************************************************************/
2396   t1 = nz * nz;
2397   t2 = t1 * xc;
2398   t3 = ZB * ZB;
2399   t5 = t2 * 0.3141592654e1 * t3;
2400   t6 = nx * 0.3141592654e1;
2401   t7 = t6 * xc;
2402   t8 = PetscCosReal(t7);
2403   t9 = nx * nx;
2404   t10 = t8 * t9;
2405   t11 = nz * 0.3141592654e1;
2406   t12 = PetscExpReal(t11);
2407   t13 = t12 * t12;
2408   t16 = PetscExpReal(xc * nz * 0.3141592654e1);
2409   t17 = t16 * t16;
2410   t18 = t17 * t16;
2411   t19 = t17 * t17;
2412   t20 = t19 * t18;
2413   t21 = t13 * t20;
2414   t22 = t10 * t21;
2415   t25 = ZA * ZA;
2416   t26 = t1 * t25;
2417   t27 = xc * 0.3141592654e1;
2418   t28 = t26 * t27;
2419   t29 = t19 * t16;
2420   t30 = t13 * t13;
2421   t31 = t29 * t30;
2422   t35 = t9 * nx;
2423   t36 = t3 * t35;
2424   t37 = PetscSinReal(t6);
2425   t38 = t13 * t12;
2426   t39 = t37 * t38;
2427   t40 = t39 * t19;
2428   t42 = t1 * t1;
2429   t43 = nx * t42;
2430   t44 = xc * xc;
2431   t45 = t25 * t44;
2432   t46 = t43 * t45;
2433   t47 = 0.3141592654e1 * 0.3141592654e1;
2434   t48 = t47 * t37;
2435   t49 = t17 * t38;
2436   t54 = 0.3141592654e1 * t35;
2437   t55 = ZA * nz * t54;
2438   t56 = t37 * ZB;
2439   t57 = t19 * t12;
2440   t61 = t25 * t8;
2441   t62 = t61 * t9;
2442   t63 = nz * t30;
2443   t64 = t63 * t16;
2444   t67 = t1 * nz;
2445   t69 = t47 * ZB;
2446   t70 = t67 * t44 * t69;
2447   t75 = nx * t3;
2448   t76 = t75 * t37;
2449   t77 = t67 * 0.3141592654e1;
2450   t78 = t19 * t19;
2451   t79 = t78 * t12;
2452   t80 = t77 * t79;
2453   t82 = t3 * t38;
2454   t84 = t54 * t37;
2455   t87 = PetscSinReal(t7);
2456   t88 = t29 * t87;
2457   t89 = t47 * t44;
2458   t93 = nx * t25;
2459   t94 = t87 * t42;
2460   t95 = t93 * t94;
2461   t96 = t47 * xc;
2462   t97 = t13 * t29;
2463   t98 = t96 * t97;
2464   t101 = t87 * t67;
2465   t102 = t93 * t101;
2466   t103 = t13 * t18;
2467   t107 = t47 * t35;
2468   t108 = t26 * t107;
2469   t109 = t37 * t44;
2470   t110 = t19 * t17;
2471   t111 = t12 * t110;
2472   t116 = t37 * t19 * t12;
2473   t118 = t37 * xc;
2474   t119 = ZB * t19;
2475   t120 = t119 * t12;
2476   t121 = t118 * t120;
2477   t125 = xc * t3;
2478   t126 = t1 * t47 * t125;
2479   t127 = t35 * t37;
2480   t128 = t38 * t19;
2481   t129 = t127 * t128;
2482   t132 = t26 * 0.3141592654e1;
2483   t133 = t16 * t13;
2484   t134 = t10 * t133;
2485   t137 = 0.3141592654e1 * ZB;
2486   t138 = t2 * t137;
2487   t139 = ZA * t8;
2488   t140 = t9 * t13;
2489   t145 = t30 * t18;
2490   t146 = t10 * t145;
2491   t149 = t3 * t8;
2492   t150 = t149 * t9;
2493   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 * nz * 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;
2494   t154 = t42 * nz;
2495   t155 = nx * t154;
2496   t156 = t44 * xc;
2497   t157 = t47 * 0.3141592654e1;
2498   t158 = t156 * t157;
2499   t159 = t155 * t158;
2500   t162 = t56 * ZA * t19 * t12;
2501   t165 = t77 * t49;
2502   t167 = t1 * t3;
2503   t168 = t167 * t89;
2504   t169 = t127 * t49;
2505   t172 = t37 * t67;
2506   t173 = t75 * t172;
2507   t174 = t38 * t110;
2508   t175 = t27 * t174;
2509   t179 = t47 * t25;
2510   t181 = t10 * t97;
2511   t184 = t27 * t31;
2512   t187 = t67 * t47;
2513   t188 = t44 * t3;
2514   t189 = t187 * t188;
2515   t192 = t25 * t35;
2516   t193 = t37 * t17;
2517   t194 = t193 * t12;
2518   t196 = nx * ZA;
2519   t197 = t196 * t172;
2520   t198 = ZB * t38;
2521   t199 = t198 * t19;
2522   t204 = t1 * t12 * t110;
2523   t207 = nx * ZB;
2524   t209 = t1 * ZA;
2525   t215 = t67 * t3;
2526   t216 = t47 * t8;
2527   t217 = t215 * t216;
2528   t218 = t9 * xc;
2529   t222 = nx * t67;
2530   t223 = t222 * t27;
2531   t224 = t3 * t87;
2532   t228 = t167 * t107;
2533   t232 = t26 * t96;
2534   t235 = t207 * t94;
2535   t236 = t47 * ZA;
2536   t243 = xc * t13;
2537   t244 = t243 * t29;
2538   t248 = t25 * nz;
2539   t249 = t248 * 0.3141592654e1;
2540   t253 = ZB * ZA;
2541   t254 = t253 * t8;
2542   t255 = t9 * nz;
2543   t256 = t30 * t16;
2544   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;
2545   t263 = t43 * t179;
2546   t267 = t3 * nz;
2547   t268 = t267 * t54;
2548   t269 = t118 * t57;
2549   t272 = t39 * t1;
2550   t274 = t67 * t25;
2551   t275 = t274 * t158;
2552   t278 = t75 * t87;
2553   t279 = t77 * t103;
2554   t282 = t25 * t38;
2555   t285 = ZA * t38;
2556   t290 = t267 * 0.3141592654e1;
2557   t296 = t77 * t111;
2558   t298 = t196 * t37;
2559   t299 = t1 * ZB;
2560   t303 = t37 * t42;
2561   t304 = t196 * t303;
2562   t308 = t77 * t57;
2563   t310 = t26 * t89;
2564   t313 = t77 * t128;
2565   t316 = t101 * t27;
2566   t319 = t93 * t87;
2567   t320 = t77 * t97;
2568   t323 = t127 * t57;
2569   t326 = t10 * nz;
2570   t329 = t118 * t174;
2571   t332 = -0.8e1 * t263 * t109 * t57 - 0.4e1 * t268 * t269 + t93 * t272 + 0.8e1 * t275 * t129 - 0.4e1 * t278 * t279 + t282 * nz * t84 - 0.2e1 * t285 * nz * 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;
2572   t335 = t75 * t101;
2573   t336 = t27 * t21;
2574   t342 = t77 * t133;
2575   t347 = t209 * t137;
2576   t350 = t9 * t1;
2577   t351 = t149 * t350;
2578   t355 = t37 * t78 * t12;
2579   t359 = t93 * t303;
2580   t367 = t172 * 0.3141592654e1;
2581   t369 = t96 * t103;
2582   t376 = t209 * t107;
2583   t379 = t10 * t103;
2584   t383 = t207 * t101;
2585   t389 = 0.3141592654e1 * ZA;
2586   t390 = t222 * t389;
2587   t391 = t87 * ZB;
2588   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;
2589   t402 = t39 * t110;
2590   t404 = t193 * t38;
2591   t406 = t127 * t174;
2592   t408 = t167 * 0.3141592654e1;
2593   t411 = t44 * t157;
2594   t412 = t155 * t411;
2595   t413 = t285 * t19;
2596   t414 = t56 * t413;
2597   t417 = ZA * t30;
2598   t424 = t93 * t37;
2599   t426 = t248 * t54;
2600   t427 = t17 * t12;
2601   t428 = t118 * t427;
2602   t431 = t77 * t21;
2603   t438 = ZA * t13;
2604   t443 = t93 * t172;
2605   t444 = t27 * t427;
2606   t448 = t1 * t78 * t12;
2607   t455 = t274 * t89;
2608   t461 = t118 * t111;
2609   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;
2610   t469 = nz * t16 * t13;
2611   t474 = t1 * t38;
2612   t475 = t474 * t19;
2613   t480 = t89 * t103;
2614   t483 = t67 * ZA;
2615   t484 = t483 * t411;
2616   t485 = t127 * t120;
2617   t488 = t127 * t111;
2618   t497 = t77 * t427;
2619   t502 = t27 * t97;
2620   t508 = t1 * t19 * t12;
2621   t511 = t155 * t25 * t156;
2622   t512 = t157 * t37;
2623   t513 = t512 * t128;
2624   t527 = t1 * t17;
2625   t528 = t527 * t38;
2626   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;
2627   t535 = t44 * t13;
2628   t542 = t527 * t12;
2629   t544 = nz * t13;
2630   t545 = t544 * t20;
2631   t548 = t75 * t303;
2632   t549 = t96 * t111;
2633   t552 = ZA * t35;
2634   t553 = t552 * t37;
2635   t562 = t43 * t96;
2636   t563 = t3 * t37;
2637   t564 = t563 * t128;
2638   t579 = t474 * t110;
2639   t590 = t9 * t30;
2640   t591 = t590 * t18;
2641   t595 = t127 * t427;
2642   t598 = t77 * t174;
2643   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;
2644   t601 = t27 * t49;
2645   t604 = t127 * t79;
2646   t606 = ZB * t29;
2647   t616 = t139 * t140 * t18;
2648   t638 = t10 * t256;
2649   t643 = t118 * t199;
2650   t653 = t544 * t29;
2651   t658 = t3 * t29;
2652   t660 = t350 * t27;
2653   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;
2654   t670 = t274 * t411;
2655   t673 = t118 * t49;
2656   t694 = t155 * t45;
2657   t713 = nz * t29 * t30;
2658   t717 = t20 * t87;
2659   t723 = t512 * t57;
2660   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;
2661   t734 = t43 * t89;
2662   t735 = t563 * t427;
2663   t740 = t75 * t94;
2664   t744 = ZB * xc;
2665   t750 = t563 * t57;
2666   t754 = t218 * t103;
2667   t771 = t127 * t199;
2668   t776 = t89 * t174;
2669   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;
2670   t794 = t483 * t158;
2671   t801 = t2 * 0.3141592654e1;
2672   t818 = t215 * t411;
2673   t827 = t96 * t174;
2674   t837 = t37 * t12 * t110;
2675   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;
2676   t864 = t544 * t18;
2677   t867 = t3 * t18;
2678   t884 = t27 * t256;
2679   t891 = t187 * t744;
2680   t894 = t563 * t49;
2681   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;
2682   t904 = t27 * t111;
2683   t907 = t118 * t128;
2684   t915 = t89 * t145;
2685   t947 = t139 * t140 * t29;
2686   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;
2687   t969 = t67 * t157 * t156 * t3;
2688   t974 = t27 * t133;
2689   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;
2690   t1042 = t89 * t97;
2691   t1055 = -0.2e1 * t168 * t595 + 0.16e2 * t484 * t771 + 0.4e1 * t11 * t125 * t129 - 0.2e1 * t173 * t444 + 0.2e1 * ZB * nz * 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;
2692   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;
2693   t1103 = t42 * t1 * t157;
2694   t1104 = t1103 * t25;
2695   t1108 = t3 * t19;
2696   t1112 = nz * t47;
2697   t1113 = t9 * t9;
2698   t1118 = t42 * t157;
2699   t1119 = t1118 * t9;
2700   t1120 = t25 * xc;
2701   t1121 = t13 * t110;
2702   t1122 = t1120 * t1121;
2703   t1125 = t47 * t47;
2704   t1126 = t67 * t1125;
2705   t1127 = t1113 * ZA;
2706   t1128 = t1126 * t1127;
2707   t1129 = t19 * t13;
2708   t1130 = t744 * t1129;
2709   t1133 = t154 * t1125;
2710   t1134 = t1133 * t9;
2711   t1135 = t45 * t1129;
2712   t1138 = t154 * t47;
2713   t1139 = t25 * t30;
2714   t1142 = t1126 * t1113;
2715   t1145 = t125 * t1129;
2716   t1148 = t1103 * xc;
2717   t1149 = t3 * t13;
2718   t1150 = t1149 * t17;
2719   t1153 = t25 * t78;
2720   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;
2721   t1157 = t25 * t13;
2722   t1158 = t1157 * t17;
2723   t1161 = t13 * t17;
2724   t1162 = t1120 * t1161;
2725   t1165 = t3 * t78;
2726   t1170 = t42 * t67 * t1125;
2727   t1172 = t1108 * t13;
2728   t1175 = t1 * t157;
2729   t1176 = t1175 * t1113;
2730   t1182 = t1120 * t1129;
2731   t1189 = t110 * t9 * xc;
2732   t1192 = t1149 * t110;
2733   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;
2734   t1209 = t17 * t30;
2735   t1210 = t125 * t1209;
2736   t1213 = t1138 * ZA;
2737   t1214 = ZB * t30;
2738   t1218 = t1157 * t110;
2739   t1226 = t3 * t30;
2740   t1237 = t1170 * t25;
2741   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;
2742   t1251 = t125 * t1121;
2743   t1260 = t1120 * t1209;
2744   t1263 = t1139 * t19;
2745   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;
2746   t1287 = t188 * t1129;
2747   t1292 = t25 * t19;
2748   t1296 = t187 * t9;
2749   t1297 = t1226 * t19;
2750   t1311 = t1112 * t1113;
2751   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;
2752   t1320 = t253 * t1129;
2753   t1328 = t253 * t30 * t19;
2754   t1333 = t125 * t1161;
2755   t1343 = ZB * t44 * t1129;
2756   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;
2757   t1369 = t1133 * t9 * ZA;
2758   t1378 = t187 * ZA;
2759   t1383 = t1170 * ZA;
2760   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;
2761   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;
2762 
2763   _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);
2764   /****************************************************************************************/
2765   /****************************************************************************************/
2766 
2767   if (x>xc) {
2768     _PC1=_PC1B; _PC2=_PC2B; _PC3=_PC3B; _PC4=_PC4B; Z=ZB;
2769   }
2770   else {
2771     _PC1=_PC1A; _PC2=_PC2A; _PC3=_PC3A; _PC4=_PC4A; Z=ZA;
2772   }
2773   /****************************************************************************************/
2774   /****************************************************************************************/
2775   t1 = nz * nz;
2776   t2 = t1 * t1;
2777   t3 = t2 * nz;
2778   t4 = x * t3;
2779   t5 = 0.3141592654e1 * 0.3141592654e1;
2780   t6 = t5 * 0.3141592654e1;
2781   t11 = _PC3 * t6;
2782   t12 = x * nz;
2783   t13 = nx * nx;
2784   t14 = t13 * t13;
2785   t15 = t12 * t14;
2786   t19 = PetscExpReal(t12 * 0.3141592654e1);
2787   t20 = t19 * t19;
2788   t21 = t4 * t20;
2789   t24 = _PC1 * t5;
2790   t25 = Z * t20;
2791   t29 = _PC1 * t6;
2792   t30 = t29 * Z;
2793   t31 = t1 * nz;
2794   t32 = x * t31;
2795   t33 = t32 * t13;
2796   t36 = t11 * x;
2797   t41 = nz * t20;
2798   t45 = t6 * _PC4;
2799   t49 = t20 * t1;
2800   t51 = _PC2 * Z;
2801   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;
2802   t58 = t32 * t6;
2803   t59 = _PC4 * t20;
2804   t63 = t20 * t13;
2805   t67 = t12 * t6;
2806   t68 = t20 * t14;
2807   t87 = t49 * t13;
2808   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;
2809   t106 = _PC3 * t5;
2810   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;
2811   t129 = PetscSinReal(nx * 0.3141592654e1 * x);
2812   t155 = 0.2e1 * t30 * t15 + x * 0.3141592654e1 * t41 * t13 - 0.4e1 * t19 * nx * t129 * nz + 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;
2813   t158 = t5 * Z;
2814 
2815   u1 = (t55 + t90 + t120 + t155) / (0.4e1 * t158 * t19 * t2 + 0.8e1 * t158 * t19 * t1 * t13 + 0.4e1 * t158 * t19 * t14);
2816   /****************************************************************************************/
2817   /****************************************************************************************/
2818   t1 = nz * nz;
2819   t2 = t1 * nz;
2820   t3 = x * t2;
2821   t4 = 0.3141592654e1 * 0.3141592654e1;
2822   t5 = t4 * 0.3141592654e1;
2823   t6 = t3 * t5;
2824   t7 = _PC2 * Z;
2825   t8 = nx * nx;
2826   t12 = t1 * t1;
2827   t13 = t12 * nz;
2828   t14 = x * t13;
2829   t15 = t5 * _PC4;
2830   t16 = x * nz;
2831   t18 = PetscExpReal(t16 * 0.3141592654e1);
2832   t19 = t18 * t18;
2833   t23 = t16 * t5;
2834   t24 = t8 * t8;
2835   t28 = _PC3 * t5;
2836   t29 = t14 * t19;
2837   t32 = _PC1 * t5;
2838   t33 = t32 * Z;
2839   t34 = t16 * t24;
2840   t37 = _PC4 * t19;
2841   t45 = _PC2 * t4;
2842   t53 = t19 * t8;
2843   t58 = _PC4 * t4;
2844   t60 = t1 * t19 * t8;
2845   t63 = t19 * t24;
2846   t67 = t3 * t8;
2847   t73 = nz * t19;
2848   t86 = t28 * x;
2849   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;
2850   t94 = Z * t12;
2851   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;
2852   t146 = PetscCosReal(nx * 0.3141592654e1 * x);
2853   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;
2854   t159 = t4 * Z;
2855 
2856   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);
2857   /****************************************************************************************/
2858   /****************************************************************************************/
2859   t1 = 0.3141592654e1 * 0.3141592654e1;
2860   t2 = t1 * 0.3141592654e1;
2861   t3 = _PC1 * t2;
2862   t4 = t3 * Z;
2863   t5 = nz * nz;
2864   t6 = t5 * t5;
2865   t7 = t6 * nz;
2866   t8 = x * t7;
2867   t9 = x * nz;
2868   t11 = PetscExpReal(t9 * 0.3141592654e1);
2869   t12 = t11 * t11;
2870   t13 = t8 * t12;
2871   t16 = t5 * nz;
2872   t17 = x * t16;
2873   t18 = t17 * t2;
2874   t19 = _PC4 * t12;
2875   t20 = nx * nx;
2876   t24 = t2 * _PC4;
2877   t28 = _PC3 * t2;
2878   t29 = t28 * x;
2879   t30 = t12 * nz;
2880   t31 = t20 * t20;
2881   t40 = _PC2 * Z;
2882   t44 = t9 * t2;
2883   t48 = t12 * t20;
2884   t52 = t17 * t20;
2885   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;
2886   t58 = t9 * t31;
2887   t61 = _PC3 * t1;
2888   t62 = t12 * t31;
2889   t73 = t5 * t20;
2890   t78 = _PC1 * t1;
2891   t90 = Z * t12;
2892   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;
2893   t101 = PetscCosReal(nx * 0.3141592654e1 * x);
2894   t102 = t11 * t101;
2895   t109 = t12 * t5;
2896   t110 = t109 * t20;
2897   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;
2898   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;
2899   t159 = t1 * t11;
2900 
2901   u3 = (t57 + t94 + t128 + t156) / (0.4e1 * t159 * t6 + 0.8e1 * t159 * t73 + 0.4e1 * t159 * t31);
2902   /****************************************************************************************/
2903   /****************************************************************************************/
2904   t1 = _PC2 * Z;
2905   t2 = 0.3141592654e1 * 0.3141592654e1;
2906   t3 = t2 * 0.3141592654e1;
2907   t4 = nz * nz;
2908   t5 = t4 * t4;
2909   t6 = t5 * t4;
2910   t8 = t3 * t6 * x;
2911   t11 = x * t4;
2912   t12 = t11 * t3;
2913   t15 = PetscExpReal(x * nz * 0.3141592654e1);
2914   t16 = t15 * t15;
2915   t17 = _PC3 * t16;
2916   t18 = nx * nx;
2917   t19 = t18 * t18;
2918   t23 = t5 * nz;
2919   t24 = t2 * t23;
2920   t28 = t1 * t3;
2921   t29 = t6 * x;
2922   t30 = t29 * t16;
2923   t33 = _PC4 * t3;
2924   t34 = t5 * x;
2925   t35 = t34 * t18;
2926   t41 = PetscSinReal(nx * 0.3141592654e1 * x);
2927   t47 = t11 * t19;
2928   t54 = t3 * _PC3;
2929   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;
2930   t58 = 0.3141592654e1 * t16;
2931   t60 = t2 * _PC4;
2932   t69 = t4 * nz;
2933   t73 = t1 * t2;
2934   t75 = t69 * t16 * t18;
2935   t79 = x * t16;
2936   t83 = nz * t16;
2937   t84 = t83 * t19;
2938   t95 = -t34 * t58 + 0.2e1 * t60 * t23 * t16 + 0.2e1 * t60 * nz * 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;
2939   t97 = t34 * t3;
2940   t101 = Z * _PC1;
2941   t102 = t16 * t19;
2942   t106 = t16 * t18;
2943   t127 = t2 * t69;
2944   t131 = t2 * nz;
2945   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;
2946   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 * nz * t18;
2947 
2948   u4 = (t57 + t95 + t135 + t164) / (0.4e1 * t24 * t15 + 0.8e1 * t127 * t15 * t18 + 0.4e1 * t131 * t15 * t19);
2949 
2950   /****************************************************************************************/
2951   /****************************************************************************************/
2952 
2953   u5 = (PetscReal)(-2*Z*nz*PETSC_PI*u2-u3*2*nz*PETSC_PI)*PetscCosReal(nz*PETSC_PI*z); /* pressure */
2954 
2955   u6 = (PetscReal)(u3*2*nz*PETSC_PI + 4*Z*nz*PETSC_PI*u2)*PetscCosReal(nz*PETSC_PI*z); /* zz stress */
2956   sum5 +=u5;
2957   sum6 +=u6;
2958 
2959   u1 *= PetscCosReal(nz*PETSC_PI*z); /* x velocity */
2960   sum1 += u1;
2961   u2 *= PetscSinReal(nz*PETSC_PI*z); /* z velocity */
2962   sum2 += u2;
2963   u3 *= 2*nz*PETSC_PI*PetscCosReal(nz*PETSC_PI*z); /* xx stress */
2964   sum3 += u3;
2965   u4 *= 2*nz*PETSC_PI*PetscSinReal(nz*PETSC_PI*z); /* zx stress */
2966   sum4 += u4;
2967 
2968   /* Output */
2969   if (mu) {
2970     *mu = Z;
2971   }
2972   if (vel) {
2973     vel[0] = sum1;
2974     vel[1] = sum2;
2975   }
2976   if (p) {
2977     (*p) = sum5;
2978   }
2979   if (s) {
2980     s[0] = sum3;
2981     s[1] = sum4;
2982     s[2] = sum6;
2983   }
2984   if (gamma) {
2985     /* sigma = tau - p, tau = sigma + p, tau[] = 2*eta*gamma[] */
2986     gamma[0] = (sum3+sum5)/(2.0*Z);
2987     gamma[1] = (sum4)/(2.0*Z);
2988     gamma[2] = (sum6+sum5)/(2.0*Z);
2989   }
2990   PetscFunctionReturn(0);
2991 }
2992 
2993 static PetscErrorCode SolCxSolutionVelocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar v[], void *ctx)
2994 {
2995   Parameter     *s = (Parameter *) ctx;
2996 
2997   PetscFunctionBegin;
2998   PetscCall(SolCxSolution(x, s->m, s->n, s->xc, s->etaA, s->etaB, v, NULL, NULL, NULL, NULL));
2999   PetscFunctionReturn(0);
3000 }
3001 
3002 static PetscErrorCode SolCxSolutionPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar p[], void *ctx)
3003 {
3004   Parameter     *s = (Parameter *) ctx;
3005 
3006   PetscFunctionBegin;
3007   PetscCall(SolCxSolution(x, s->m, s->n, s->xc, s->etaA, s->etaB, NULL, p, NULL, NULL, NULL));
3008   PetscFunctionReturn(0);
3009 }
3010 
3011 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
3012 {
3013   PetscInt       sol;
3014 
3015   PetscFunctionBeginUser;
3016   options->solType = SOLKX;
3017   PetscOptionsBegin(comm, "", "Variable-Viscosity Stokes Problem Options", "DMPLEX");
3018   sol  = options->solType;
3019   PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex69.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL));
3020   options->solType = (SolutionType) sol;
3021   PetscOptionsEnd();
3022   PetscFunctionReturn(0);
3023 }
3024 
3025 static PetscErrorCode SetUpParameters(AppCtx *user)
3026 {
3027   PetscBag       bag;
3028   Parameter     *p;
3029 
3030   PetscFunctionBeginUser;
3031   /* setup PETSc parameter bag */
3032   PetscCall(PetscBagGetData(user->bag, (void **) &p));
3033   PetscCall(PetscBagSetName(user->bag, "par", "Problem parameters"));
3034   bag  = user->bag;
3035   switch (user->solType) {
3036   case SOLKX:
3037     PetscCall(PetscBagRegisterInt(bag,  &p->n, 1,   "n", "x-wavelength for forcing variation"));
3038     PetscCall(PetscBagRegisterInt(bag,  &p->m, 1,   "m", "z-wavelength for forcing variation"));
3039     PetscCall(PetscBagRegisterReal(bag, &p->B, 1.0, "B", "Exponential scale for viscosity variation"));
3040     break;
3041   case SOLCX:
3042     PetscCall(PetscBagRegisterInt(bag,  &p->n,    1,   "n",    "x-wavelength for forcing variation"));
3043     PetscCall(PetscBagRegisterInt(bag,  &p->m,    1,   "m",    "z-wavelength for forcing variation"));
3044     PetscCall(PetscBagRegisterReal(bag, &p->etaA, 1.0, "etaA", "Viscosity for x < xc"));
3045     PetscCall(PetscBagRegisterReal(bag, &p->etaB, 1.0, "etaB", "Viscosity for x > xc"));
3046     PetscCall(PetscBagRegisterReal(bag, &p->xc,   0.5, "xc",   "x-coordinate of the viscosity jump"));
3047     break;
3048   default:
3049     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid solution type %d (%s)", user->solType, solTypes[PetscMin(user->solType, NUM_SOL_TYPES)]);
3050   }
3051   PetscCall(PetscBagSetFromOptions(bag));
3052   PetscCall(PetscBagViewFromOptions(bag, NULL, "-param_view"));
3053   PetscFunctionReturn(0);
3054 }
3055 
3056 /* Make split labels so that we can have corners in multiple labels */
3057 static PetscErrorCode CreateSplitLabels(DM dm)
3058 {
3059   const char    *names[4] = {"markerBottom", "markerRight", "markerTop", "markerLeft"};
3060   PetscInt       ids[4]   = {1, 2, 3, 4};
3061   DMLabel        label;
3062   IS             is;
3063   PetscInt       f;
3064 
3065   PetscFunctionBeginUser;
3066   for (f = 0; f < 4; ++f) {
3067     PetscCall(DMCreateLabel(dm, names[f]));
3068     PetscCall(DMGetStratumIS(dm, "marker", ids[f],  &is));
3069     if (!is) continue;
3070     PetscCall(DMGetLabel(dm, names[f], &label));
3071     PetscCall(DMLabelInsertIS(label, is, 1));
3072     PetscCall(ISDestroy(&is));
3073   }
3074   PetscFunctionReturn(0);
3075 }
3076 
3077 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
3078 {
3079   DM             cdm;
3080 
3081   PetscFunctionBeginUser;
3082   PetscCall(DMCreate(comm, dm));
3083   PetscCall(DMSetType(*dm, DMPLEX));
3084   PetscCall(DMSetFromOptions(*dm));
3085   cdm  = *dm;
3086   while (cdm) {
3087     PetscCall(CreateSplitLabels(cdm));
3088     PetscCall(DMGetCoarseDM(cdm, &cdm));
3089   }
3090   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
3091   PetscFunctionReturn(0);
3092 }
3093 
3094 static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
3095 {
3096   PetscErrorCode (*exactFunc)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
3097   PetscDS        prob;
3098   DMLabel        label;
3099   const PetscInt id  = 1;
3100   PetscInt       dim, comp;
3101   Parameter     *ctx;
3102   void          *data;
3103 
3104   PetscFunctionBeginUser;
3105   PetscCall(DMGetDimension(dm, &dim));
3106   PetscCall(DMGetDS(dm, &prob));
3107   switch (user->solType) {
3108   case SOLKX:
3109     PetscCall(PetscDSSetResidual(prob, 0, f0_u, stokes_momentum_kx));
3110     PetscCall(PetscDSSetResidual(prob, 1, stokes_mass, f1_zero));
3111     PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL,  NULL,  stokes_momentum_vel_J_kx));
3112     PetscCall(PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  stokes_momentum_pres_J, NULL));
3113     PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, stokes_mass_J, NULL,  NULL));
3114     PetscCall(PetscDSSetJacobianPreconditioner(prob, 0, 0, NULL, NULL, NULL, stokes_momentum_vel_J_kx));
3115     PetscCall(PetscDSSetJacobianPreconditioner(prob, 1, 1, stokes_identity_J_kx, NULL, NULL, NULL));
3116     break;
3117   case SOLCX:
3118     PetscCall(PetscDSSetResidual(prob, 0, f0_u, stokes_momentum_cx));
3119     PetscCall(PetscDSSetResidual(prob, 1, stokes_mass, f1_zero));
3120     PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL,  NULL,  stokes_momentum_vel_J_cx));
3121     PetscCall(PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  stokes_momentum_pres_J, NULL));
3122     PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, stokes_mass_J, NULL,  NULL));
3123     PetscCall(PetscDSSetJacobianPreconditioner(prob, 0, 0, NULL, NULL, NULL, stokes_momentum_vel_J_kx));
3124     PetscCall(PetscDSSetJacobianPreconditioner(prob, 1, 1, stokes_identity_J_cx, NULL, NULL, NULL));
3125     break;
3126   default:
3127     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid solution type %d (%s)", user->solType, solTypes[PetscMin(user->solType, NUM_SOL_TYPES)]);
3128   }
3129   PetscCall(PetscBagGetData(user->bag, &data));
3130   switch (dim) {
3131   case 2:
3132     switch (user->solType) {
3133     case SOLKX:
3134       PetscCall(PetscDSSetExactSolution(prob, 0, SolKxSolutionVelocity, data));
3135       PetscCall(PetscDSSetExactSolution(prob, 1, SolKxSolutionPressure, data));
3136       break;
3137     case SOLCX:
3138       PetscCall(PetscDSSetExactSolution(prob, 0, SolCxSolutionVelocity, data));
3139       PetscCall(PetscDSSetExactSolution(prob, 1, SolCxSolutionPressure, data));
3140       break;
3141     default:
3142       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid solution type %d (%s)", user->solType, solTypes[PetscMin(user->solType, NUM_SOL_TYPES)]);
3143     }
3144     break;
3145   default:
3146     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
3147   }
3148   /* Setup constants */
3149   {
3150     Parameter *param;
3151 
3152     PetscCall(PetscBagGetData(user->bag, (void **) &param));
3153     switch (user->solType) {
3154     case SOLKX:
3155     {
3156       PetscScalar constants[3];
3157 
3158       constants[0] = param->m;
3159       constants[1] = param->n;
3160       constants[2] = param->B;
3161       PetscCall(PetscDSSetConstants(prob, 3, constants));
3162     }
3163     break;
3164     case SOLCX:
3165     {
3166       PetscScalar constants[5];
3167 
3168       constants[0] = param->m;
3169       constants[1] = param->n;
3170       constants[2] = param->etaA;
3171       constants[3] = param->etaB;
3172       constants[4] = param->xc;
3173       PetscCall(PetscDSSetConstants(prob, 5, constants));
3174     }
3175     break;
3176     default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "No parameter information for solution type %d", user->solType);
3177     }
3178   }
3179   /* Setup Boundary Conditions */
3180   PetscCall(PetscDSGetExactSolution(prob, 0, &exactFunc, (void **) &ctx));
3181   comp = 1;
3182   PetscCall(DMGetLabel(dm, "markerBottom", &label));
3183   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wallB", label, 1, &id, 0, 1, &comp, (void (*)(void)) exactFunc, NULL, ctx, NULL));
3184   comp = 0;
3185   PetscCall(DMGetLabel(dm, "markerRight", &label));
3186   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wallR", label, 1, &id, 0, 1, &comp, (void (*)(void)) exactFunc, NULL, ctx, NULL));
3187   comp = 1;
3188   PetscCall(DMGetLabel(dm, "markerTop", &label));
3189   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wallT", label, 1, &id, 0, 1, &comp, (void (*)(void)) exactFunc, NULL, ctx, NULL));
3190   comp = 0;
3191   PetscCall(DMGetLabel(dm, "markerLeft", &label));
3192   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wallL", label, 1, &id, 0, 1, &comp, (void (*)(void)) exactFunc, NULL, ctx, NULL));
3193   PetscFunctionReturn(0);
3194 }
3195 
3196 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace)
3197 {
3198   Vec              vec;
3199   PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero, one};
3200 
3201   PetscFunctionBeginUser;
3202   PetscCheck(origField == 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Field %" PetscInt_FMT " should be 1 for pressure", origField);
3203   funcs[field] = one;
3204   {
3205     PetscDS ds;
3206     PetscCall(DMGetDS(dm, &ds));
3207     PetscCall(PetscObjectViewFromOptions((PetscObject) ds, NULL, "-ds_view"));
3208   }
3209   PetscCall(DMCreateGlobalVector(dm, &vec));
3210   PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
3211   PetscCall(VecNormalize(vec, NULL));
3212   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace));
3213   PetscCall(VecDestroy(&vec));
3214   /* New style for field null spaces */
3215   {
3216     PetscObject  pressure;
3217     MatNullSpace nullspacePres;
3218 
3219     PetscCall(DMGetField(dm, field, NULL, &pressure));
3220     PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres));
3221     PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres));
3222     PetscCall(MatNullSpaceDestroy(&nullspacePres));
3223   }
3224   PetscFunctionReturn(0);
3225 }
3226 
3227 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
3228 {
3229   DM             cdm = dm;
3230   PetscFE        fe[2];
3231   DMPolytopeType ct;
3232   PetscInt       dim, cStart;
3233   PetscBool      simplex;
3234 
3235   PetscFunctionBeginUser;
3236   PetscCall(DMGetDimension(dm, &dim));
3237   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
3238   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
3239   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
3240   /* Create discretization of solution fields */
3241   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
3242   PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity"));
3243   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
3244   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
3245   PetscCall(PetscObjectSetName((PetscObject) fe[1], "pressure"));
3246   /* Set discretization and boundary conditions for each mesh */
3247   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0]));
3248   PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1]));
3249   PetscCall(DMCreateDS(dm));
3250   PetscCall(SetupProblem(dm, user));
3251   while (cdm) {
3252     PetscCall(DMCopyDisc(dm, cdm));
3253     PetscCall(DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace));
3254     PetscCall(DMGetCoarseDM(cdm, &cdm));
3255   }
3256   PetscCall(PetscFEDestroy(&fe[0]));
3257   PetscCall(PetscFEDestroy(&fe[1]));
3258   {
3259     PetscObject  pressure;
3260     MatNullSpace nullSpacePres;
3261 
3262     PetscCall(DMGetField(dm, 1, NULL, &pressure));
3263     PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullSpacePres));
3264     PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject) nullSpacePres));
3265     PetscCall(MatNullSpaceDestroy(&nullSpacePres));
3266   }
3267   PetscFunctionReturn(0);
3268 }
3269 
3270 /* Add a vector in the nullspace to make the continuum integral 0.
3271 
3272    If int(u) = a and int(n) = b, then int(u - a/b n) = a - a/b b = 0
3273 */
3274 static void pressure(PetscInt dim, PetscInt Nf, PetscInt NfAux,
3275                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
3276                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
3277                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[])
3278 {
3279   p[0] = u[uOff[1]];
3280 }
3281 static PetscErrorCode CorrectDiscretePressure(DM dm, MatNullSpace nullspace, Vec u, AppCtx *user)
3282 {
3283   PetscDS        ds;
3284   const Vec     *nullvecs;
3285   PetscScalar    pintd, intc[2], intn[2];
3286   MPI_Comm       comm;
3287 
3288   PetscFunctionBeginUser;
3289   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
3290   PetscCall(DMGetDS(dm, &ds));
3291   PetscCall(PetscDSSetObjective(ds, 1, pressure));
3292   PetscCall(MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs));
3293   PetscCall(VecDot(nullvecs[0], u, &pintd));
3294   PetscCheck(PetscAbsScalar(pintd) <= PETSC_SMALL,comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g", (double) PetscRealPart(pintd));
3295   PetscCall(DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, user));
3296   PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, user));
3297   PetscCall(VecAXPY(u, -intc[1]/intn[1], nullvecs[0]));
3298   PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, user));
3299   PetscCheck(PetscAbsScalar(intc[1]) <= PETSC_SMALL,comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g", (double) PetscRealPart(intc[1]));
3300   PetscFunctionReturn(0);
3301 }
3302 
3303 static PetscErrorCode SNESConvergenceCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *user)
3304 {
3305   PetscFunctionBeginUser;
3306   PetscCall(SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, user));
3307   if (*reason > 0) {
3308     DM           dm;
3309     Mat          J;
3310     Vec          u;
3311     MatNullSpace nullspace;
3312 
3313     PetscCall(SNESGetDM(snes, &dm));
3314     PetscCall(SNESGetSolution(snes, &u));
3315     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
3316     PetscCall(MatGetNullSpace(J, &nullspace));
3317     PetscCheck(nullspace,PetscObjectComm((PetscObject) snes), PETSC_ERR_ARG_WRONG, "SNES Jacobian has no attached null space");
3318     PetscCall(CorrectDiscretePressure(dm, nullspace, u, (AppCtx *) user));
3319   }
3320   PetscFunctionReturn(0);
3321 }
3322 
3323 int main(int argc, char **argv)
3324 {
3325   SNES            snes;      /* nonlinear solver */
3326   DM              dm;        /* problem definition */
3327   Vec             u,r;       /* solution, residual vectors */
3328   Mat             J, M;      /* Jacobian and preconditiong matrix */
3329   MatNullSpace    nullSpace; /* May be necessary for pressure */
3330   AppCtx          user;      /* user-defined work context */
3331   PetscErrorCode  (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero, zero};
3332 
3333   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
3334   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3335   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
3336   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
3337   PetscCall(SNESSetDM(snes, dm));
3338   PetscCall(DMSetApplicationContext(dm, &user));
3339   /* Setup problem parameters */
3340   PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag));
3341   PetscCall(SetUpParameters(&user));
3342   /* Setup problem */
3343   PetscCall(SetupDiscretization(dm, &user));
3344   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
3345 
3346   PetscCall(DMCreateGlobalVector(dm, &u));
3347   PetscCall(VecDuplicate(u, &r));
3348 
3349   PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user));
3350   PetscCall(CreatePressureNullSpace(dm, 1, 1, &nullSpace));
3351 
3352   { /* set tolerances */
3353     KSP ksp;
3354 
3355     PetscCall(SNESGetKSP(snes, &ksp));
3356     PetscCall(KSPSetTolerances(ksp,1.e-2*PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,PETSC_DEFAULT));
3357   }
3358 
3359   /* There should be a way to express this using the DM */
3360   PetscCall(SNESSetFromOptions(snes));
3361   PetscCall(SNESSetUp(snes));
3362   PetscCall(SNESGetJacobian(snes, &J, &M, NULL, NULL));
3363   PetscCall(MatSetNullSpace(J, nullSpace));
3364   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) M, "prec_"));
3365   PetscCall(MatSetFromOptions(M));
3366   PetscCall(SNESSetConvergenceTest(snes, SNESConvergenceCorrectPressure, &user, NULL));
3367 
3368   PetscCall(DMSNESCheckFromOptions(snes, u));
3369   PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
3370   PetscCall(PetscObjectSetName((PetscObject) u, "Solution"));
3371   PetscCall(SNESSolve(snes, NULL, u));
3372   {
3373     PetscErrorCode (*exacts[2])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
3374     void            *ectxs[2];
3375     PetscDS          ds;
3376     Vec              e;
3377 
3378     PetscCall(DMGetDS(dm, &ds));
3379     PetscCall(PetscDSGetExactSolution(ds, 0, &exacts[0], &ectxs[0]));
3380     PetscCall(PetscDSGetExactSolution(ds, 1, &exacts[1], &ectxs[1]));
3381 
3382     PetscCall(DMGetGlobalVector(dm, &e));
3383     PetscCall(PetscObjectCompose((PetscObject) e, "__Vec_bc_zero__", (PetscObject) dm));
3384     PetscCall(DMPlexComputeL2DiffVec(dm, 0.0, exacts, ectxs, u, e));
3385     PetscCall(PetscObjectSetName((PetscObject) e, "Solution Error"));
3386     PetscCall(VecViewFromOptions(e, NULL, "-error_vec_view"));
3387     PetscCall(PetscObjectCompose((PetscObject) e, "__Vec_bc_zero__", NULL));
3388     PetscCall(DMRestoreGlobalVector(dm, &e));
3389   }
3390   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
3391 
3392   PetscCall(MatNullSpaceDestroy(&nullSpace));
3393   PetscCall(VecDestroy(&u));
3394   PetscCall(VecDestroy(&r));
3395   PetscCall(SNESDestroy(&snes));
3396   PetscCall(DMDestroy(&dm));
3397   PetscCall(PetscBagDestroy(&user.bag));
3398   PetscCall(PetscFinalize());
3399   return 0;
3400 }
3401 
3402 /*TEST
3403 
3404   # 2D serial discretization tests
3405   test:
3406     suffix: p2p1
3407     requires: triangle
3408     args: -dm_plex_separate_marker -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3409       -snes_error_if_not_converged -dmsnes_check .001 \
3410       -ksp_rtol 1.e-9 -ksp_error_if_not_converged -pc_use_amat \
3411       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition a11 \
3412         -fieldsplit_velocity_pc_type lu \
3413         -fieldsplit_pressure_ksp_rtol 1.e-9 -fieldsplit_pressure_pc_type lu
3414   test:
3415     suffix: p2p1_gmg
3416     requires: triangle
3417     args: -dm_plex_separate_marker -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3418       -snes_error_if_not_converged -dmsnes_check .001 \
3419       -ksp_type fgmres -ksp_rtol 1.e-9 -ksp_error_if_not_converged -pc_use_amat \
3420       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition a11 \
3421         -fieldsplit_velocity_pc_type mg \
3422         -fieldsplit_pressure_ksp_rtol 1.e-9 -fieldsplit_pressure_pc_type lu
3423   test:
3424     suffix: p2p1_conv
3425     requires: triangle
3426     # -dm_refine 2 gives L_2 convergence rate: [3.0, 2.2]
3427     args: -dm_plex_separate_marker -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3428       -snes_error_if_not_converged -snes_convergence_estimate -convest_num_refine 2 \
3429       -ksp_rtol 1.e-9 -ksp_error_if_not_converged -pc_use_amat \
3430       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition a11 \
3431         -fieldsplit_velocity_pc_type lu \
3432         -fieldsplit_pressure_ksp_rtol 1.e-9 -fieldsplit_pressure_pc_type lu
3433   test:
3434     suffix: q2q1_conv
3435     # -dm_refine 2 gives L_2 convergence rate: [3.0, 2.1]
3436     args: -dm_plex_simplex 0 -dm_plex_separate_marker -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3437       -snes_error_if_not_converged -snes_convergence_estimate -convest_num_refine 2 \
3438       -ksp_rtol 1.e-9 -ksp_error_if_not_converged -pc_use_amat \
3439       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition a11 \
3440         -fieldsplit_velocity_pc_type lu \
3441         -fieldsplit_pressure_ksp_rtol 1.e-9 -fieldsplit_pressure_pc_type lu
3442   test:
3443     suffix: q1p0_conv
3444     # -dm_refine 2 gives L_2 convergence rate: [2.0, 1.0]
3445     args: -dm_plex_simplex 0 -dm_plex_separate_marker -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
3446       -snes_error_if_not_converged -snes_convergence_estimate -convest_num_refine 2 \
3447       -ksp_rtol 1.e-9 -ksp_error_if_not_converged -pc_use_amat \
3448       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition a11 \
3449         -fieldsplit_velocity_pc_type lu \
3450         -fieldsplit_pressure_ksp_rtol 1.e-9 -fieldsplit_pressure_pc_type lu
3451   test:
3452     suffix: q2p1_conv
3453     # -dm_refine 2 gives L_2 convergence rate: [3.0, 2.0]
3454     args: -dm_plex_simplex 0 -dm_plex_separate_marker \
3455       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pres_petscspace_poly_tensor 0 -pres_petscdualspace_lagrange_continuity 0 \
3456       -snes_error_if_not_converged -snes_convergence_estimate -convest_num_refine 2 \
3457       -ksp_rtol 1.e-9 -ksp_error_if_not_converged -pc_use_amat \
3458       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition a11 \
3459         -fieldsplit_velocity_pc_type lu \
3460         -fieldsplit_pressure_ksp_rtol 1.e-9 -fieldsplit_pressure_pc_type lu
3461 
3462   # FETI-DP tests
3463   testset:
3464     output_file: output/ex69_q2p1fetidp.out
3465     suffix: q2p1fetidp
3466     requires: !single
3467     nsize: 5
3468     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine_pre 1 -dm_mat_type is -dm_view -petscpartitioner_type simple \
3469       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pres_petscspace_poly_tensor 0 -pres_petscdualspace_lagrange_continuity 0 -pres_petscdualspace_lagrange_node_endpoints 0 \
3470       -snes_error_if_not_converged \
3471       -ksp_error_if_not_converged \
3472       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural \
3473       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type svd -fetidp_bddc_pc_bddc_dirichlet_pc_type svd -fetidp_bddc_pc_bddc_neumann_pc_type svd
3474     test:
3475       suffix: aij
3476       args: -matis_localmat_type aij
3477     test:
3478       requires: viennacl !CUDA_VERSION_11PLUS broken
3479       suffix: aijviennacl
3480       args: -matis_localmat_type aijviennacl
3481     test:
3482       requires: cuda
3483       suffix: aijcusparse
3484       args: -matis_localmat_type aijcusparse
3485 
3486   testset:
3487     suffix: q2p1fetidp_deluxe
3488     output_file: output/ex69_q2p1fetidp_deluxe.out
3489     requires: mumps double
3490     nsize: 5
3491     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine_pre 1 -dm_mat_type is -dm_view -petscpartitioner_type simple \
3492       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pres_petscspace_poly_tensor 0 -pres_petscdualspace_lagrange_continuity 0 -pres_petscdualspace_lagrange_node_endpoints 0 \
3493       -snes_error_if_not_converged \
3494       -ksp_error_if_not_converged \
3495       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural \
3496       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_deluxe_singlemat -fetidp_bddc_pc_bddc_deluxe_zerorows \
3497       -fetidp_bddc_sub_schurs_mat_solver_type mumps -fetidp_bddc_sub_schurs_mat_mumps_icntl_14 500 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type svd
3498     test:
3499       suffix: aij
3500       args: -matis_localmat_type aij
3501     test:
3502       suffix: aij_seqdense
3503       args: -matis_localmat_type aij -fetidp_bddc_sub_schurs_schur_mat_type seqdense
3504     test:
3505       requires: viennacl !CUDA_VERSION_11PLUS
3506       suffix: aijviennacl
3507       args: -matis_localmat_type aijviennacl
3508     test:
3509       requires: cuda
3510       suffix: aijcusparse
3511       args: -matis_localmat_type aijcusparse
3512 
3513   testset:
3514     suffix: q2p1fetidp_deluxe_adaptive
3515     output_file: output/ex69_q2p1fetidp_deluxe_adaptive.out
3516     requires: mumps double
3517     nsize: 5
3518     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine_pre 1 -dm_mat_type is -dm_view -petscpartitioner_type simple \
3519       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pres_petscspace_poly_tensor 0 -pres_petscdualspace_lagrange_continuity 0 -pres_petscdualspace_lagrange_node_endpoints 0 \
3520       -snes_error_if_not_converged \
3521       -ksp_error_if_not_converged \
3522       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural \
3523       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_deluxe_singlemat -fetidp_bddc_pc_bddc_adaptive_userdefined -fetidp_bddc_pc_bddc_adaptive_threshold 1.3 \
3524       -fetidp_bddc_sub_schurs_mat_solver_type mumps -fetidp_bddc_sub_schurs_mat_mumps_icntl_14 500 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type svd
3525     test:
3526       suffix: aij
3527       args: -matis_localmat_type aij
3528     test:
3529       suffix: aij_seqdense
3530       args: -matis_localmat_type aij -fetidp_bddc_sub_schurs_schur_mat_type seqdense
3531     test:
3532       requires: viennacl !CUDA_VERSION_11PLUS
3533       suffix: aijviennacl
3534       args: -matis_localmat_type aijviennacl
3535     test:
3536       requires: cuda
3537       suffix: aijcusparse
3538       args: -matis_localmat_type aijcusparse
3539 
3540   test:
3541     suffix: p2p1fetidp
3542     requires: triangle
3543     nsize: 4
3544     args: -dm_plex_separate_marker -dm_refine_pre 2 -dm_refine_uniform_pre -dm_mat_type is -dm_view -petscpartitioner_type simple \
3545       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3546       -snes_error_if_not_converged \
3547       -ksp_error_if_not_converged \
3548       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural -fetidp_pc_fieldsplit_schur_fact_type diag \
3549       -fetidp_fieldsplit_p_pc_type jacobi -fetidp_fieldsplit_p_ksp_type preonly \
3550       -fetidp_fieldsplit_lag_ksp_type preonly \
3551       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type cholesky -fetidp_bddc_pc_bddc_dirichlet_pc_type svd -fetidp_bddc_pc_bddc_neumann_pc_type svd
3552 
3553   test:
3554     suffix: p2p1fetidp_allp
3555     requires: triangle
3556     nsize: 4
3557     args: -dm_plex_separate_marker -dm_refine_pre 2 -dm_refine_uniform_pre -dm_mat_type is -dm_view -petscpartitioner_type simple \
3558       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3559       -snes_error_if_not_converged \
3560       -ksp_error_if_not_converged \
3561       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural -fetidp_pc_fieldsplit_schur_fact_type diag -ksp_fetidp_pressure_all \
3562       -fetidp_fieldsplit_p_pc_type jacobi -fetidp_fieldsplit_p_ksp_type preonly \
3563       -fetidp_fieldsplit_lag_ksp_type preonly \
3564       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type cholesky -fetidp_bddc_pc_bddc_dirichlet_pc_type svd -fetidp_bddc_pc_bddc_neumann_pc_type svd
3565 
3566   test:
3567     suffix: p2p1fetidp_discharm
3568     requires: triangle
3569     nsize: 4
3570     args: -dm_plex_separate_marker -dm_refine_pre 2 -dm_refine_uniform_pre -dm_mat_type is -dm_view -petscpartitioner_type simple \
3571       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3572       -snes_error_if_not_converged \
3573       -ksp_error_if_not_converged \
3574       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural -fetidp_pc_fieldsplit_schur_fact_type diag -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_type cholesky \
3575       -fetidp_fieldsplit_p_pc_type jacobi -fetidp_fieldsplit_p_ksp_type preonly \
3576       -fetidp_fieldsplit_lag_ksp_type preonly \
3577       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type cholesky -fetidp_bddc_pc_bddc_dirichlet_pc_type none -fetidp_bddc_pc_bddc_neumann_pc_type svd
3578 
3579   test:
3580     suffix: p2p1fetidp_lumped
3581     requires: triangle
3582     nsize: 4
3583     args: -dm_plex_separate_marker -dm_refine_pre 2 -dm_refine_uniform_pre -dm_mat_type is -dm_view -petscpartitioner_type simple \
3584       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3585       -snes_error_if_not_converged \
3586       -ksp_error_if_not_converged \
3587       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural -fetidp_pc_fieldsplit_schur_fact_type diag -fetidp_pc_lumped \
3588       -fetidp_fieldsplit_p_pc_type jacobi -fetidp_fieldsplit_p_ksp_type preonly \
3589       -fetidp_fieldsplit_lag_ksp_type preonly \
3590       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type cholesky -fetidp_bddc_pc_bddc_dirichlet_pc_type none -fetidp_bddc_pc_bddc_neumann_pc_type svd
3591 
3592   test:
3593     suffix: p2p1fetidp_deluxe
3594     requires: triangle mumps
3595     nsize: 4
3596     args: -dm_plex_separate_marker -dm_refine_pre 2 -dm_refine_uniform_pre -dm_mat_type is -dm_view -petscpartitioner_type simple \
3597       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3598       -snes_error_if_not_converged \
3599       -ksp_error_if_not_converged \
3600       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural -fetidp_pc_fieldsplit_schur_fact_type diag \
3601       -fetidp_fieldsplit_p_pc_type jacobi -fetidp_fieldsplit_p_ksp_type preonly \
3602       -fetidp_fieldsplit_lag_ksp_type preonly \
3603       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type cholesky -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_deluxe_singlemat \
3604       -fetidp_bddc_sub_schurs_mat_solver_type mumps -fetidp_bddc_sub_schurs_mat_mumps_icntl_14 500 -fetidp_bddc_sub_schurs_posdef 0
3605 
3606   test:
3607     suffix: p2p1fetidp_deluxe_discharm
3608     requires: triangle mumps
3609     nsize: 4
3610     args: -dm_plex_separate_marker -dm_refine_pre 2 -dm_refine_uniform_pre -dm_mat_type is -dm_view -petscpartitioner_type simple \
3611       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3612       -snes_error_if_not_converged \
3613       -ksp_error_if_not_converged \
3614       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural -fetidp_pc_fieldsplit_schur_fact_type diag \
3615       -fetidp_fieldsplit_p_pc_type jacobi -fetidp_fieldsplit_p_ksp_type preonly \
3616       -fetidp_fieldsplit_lag_ksp_type preonly \
3617       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type cholesky -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_deluxe_singlemat \
3618        -fetidp_bddc_sub_schurs_mat_solver_type mumps -fetidp_bddc_sub_schurs_mat_mumps_icntl_14 500 -fetidp_bddc_sub_schurs_posdef 0  -fetidp_bddc_sub_schurs_discrete_harmonic
3619 
3620   testset:
3621     nsize: 4
3622     requires: triangle
3623     output_file: output/ex69_p2p1fetidp_olof.out
3624     args: -dm_plex_separate_marker -dm_refine_pre 2 -dm_refine_uniform_pre -dm_mat_type is -dm_view -petscpartitioner_type simple \
3625       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
3626       -snes_error_if_not_converged \
3627       -ksp_error_if_not_converged \
3628       -ksp_type fetidp -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_ksp_norm_type natural -fetidp_pc_discrete_harmonic 1 -fetidp_harmonic_pc_type cholesky -ksp_fetidp_pressure_schur \
3629       -fetidp_fieldsplit_p_pc_type bddc -fetidp_fieldsplit_p_pc_bddc_dirichlet_pc_type none -fetidp_fieldsplit_p_ksp_type preonly \
3630       -fetidp_fieldsplit_lag_ksp_error_if_not_converged 0 -fetidp_fieldsplit_lag_ksp_type chebyshev -fetidp_fieldsplit_lag_ksp_max_it 2 \
3631       -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_graph_maxcount 2 -fetidp_bddc_pc_bddc_coarse_redundant_pc_type cholesky -fetidp_bddc_pc_bddc_dirichlet_pc_type none -fetidp_bddc_pc_bddc_neumann_pc_type svd
3632     test:
3633       suffix: p2p1fetidp_olof_full
3634       args: -fetidp_pc_fieldsplit_schur_fact_type full
3635     test:
3636       suffix: p2p1fetidp_olof_diag
3637       args: -fetidp_pc_fieldsplit_schur_fact_type diag
3638     test:
3639       suffix: p2p1fetidp_olof_additive
3640       args: -fetidp_pc_fieldsplit_type additive
3641 
3642   #BDDC with benign trick
3643   testset:
3644     suffix: q2p1bddc
3645     output_file: output/ex69_q2p1fetidp_deluxe.out
3646     nsize: 5
3647     args: -dm_plex_simplex 0 -dm_plex_separate_marker -dm_refine_pre 1 -dm_mat_type is -dm_view -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -pres_petscspace_poly_tensor 0 -pres_petscdualspace_lagrange_continuity 0 -pres_petscdualspace_lagrange_node_endpoints 0 -petscds_jac_pre 0 -snes_error_if_not_converged -ksp_error_if_not_converged -ksp_type cg -ksp_norm_type natural -pc_type bddc -pc_bddc_benign_trick -pc_bddc_nonetflux -pc_bddc_detect_disconnected -pc_bddc_vertex_size 2 -pc_bddc_coarse_redundant_pc_type svd -pc_bddc_use_qr_single
3648     test:
3649       requires: double
3650       suffix: benign_card
3651       # no native support for saddle point factorizations from PETSc
3652       args: -pc_bddc_dirichlet_pc_type svd -pc_bddc_neumann_pc_type svd
3653     test:
3654       requires: mumps double
3655       suffix: benign_deluxe_mumps
3656       args: -pc_bddc_use_deluxe_scaling -pc_bddc_deluxe_zerorows -sub_schurs_mat_solver_type mumps -sub_schurs_mat_mumps_icntl_14 1000
3657     test:
3658       requires: mumps double
3659       suffix: benign_deluxe_adaptive_mumps
3660       args: -pc_bddc_adaptive_threshold 1.7 -pc_bddc_use_deluxe_scaling -pc_bddc_deluxe_zerorows -sub_schurs_mat_solver_type mumps -sub_schurs_mat_mumps_icntl_14 1000
3661     test:
3662       TODO: broken (INDEFINITE PC)
3663       requires: mkl_pardiso double !complex
3664       suffix: benign_deluxe_mkl
3665       args: -pc_bddc_use_deluxe_scaling -pc_bddc_deluxe_zerorows -sub_schurs_mat_solver_type mkl_pardiso -snes_rtol 1.e-7
3666     test:
3667       TODO: broken (INDEFINITE PC)
3668       requires: mkl_pardiso double !complex
3669       suffix: benign_deluxe_adaptive_mkl
3670       args: -pc_bddc_adaptive_threshold 1.7 -pc_bddc_use_deluxe_scaling -pc_bddc_deluxe_zerorows -sub_schurs_mat_solver_type mkl_pardiso -snes_rtol 1.e-7
3671 
3672 TEST*/
3673