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