xref: /petsc/src/snes/tutorials/ex30.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
1 static const char help[] = "Steady-state 2D subduction flow, pressure and temperature solver.\n\
2        The flow is driven by the subducting slab.\n\
3 ---------------------------------ex30 help---------------------------------\n\
4   -OPTION <DEFAULT> = (UNITS) DESCRIPTION.\n\n\
5   -width <320> = (km) width of domain.\n\
6   -depth <300> = (km) depth of domain.\n\
7   -slab_dip <45> = (degrees) dip angle of the slab (determines the grid aspect ratio).\n\
8   -lid_depth <35> = (km) depth of the static conductive lid.\n\
9   -fault_depth <35> = (km) depth of slab-wedge mechanical coupling\n\
10      (fault dept >= lid depth).\n\
11 \n\
12   -ni <82> = grid cells in x-direction. (nj adjusts to accommodate\n\
13       the slab dip & depth). DO NOT USE -da_grid_x option!!!\n\
14   -ivisc <3> = rheology option.\n\
15       0 --- constant viscosity.\n\
16       1 --- olivine diffusion creep rheology (T&P-dependent, newtonian).\n\
17       2 --- olivine dislocation creep rheology (T&P-dependent, non-newtonian).\n\
18       3 --- Full mantle rheology, combination of 1 & 2.\n\
19 \n\
20   -slab_velocity <5> = (cm/year) convergence rate of slab into subduction zone.\n\
21   -slab_age <50> = (million yrs) age of slab for thermal profile boundary condition.\n\
22   -lid_age <50> = (million yrs) age of lid for thermal profile boundary condition.\n\
23 \n\
24   FOR OTHER PARAMETER OPTIONS AND THEIR DEFAULT VALUES, see SetParams() in ex30.c.\n\
25 ---------------------------------ex30 help---------------------------------\n";
26 
27 /*F-----------------------------------------------------------------------
28 
29     This PETSc 2.2.0 example by Richard F. Katz
30     http://www.ldeo.columbia.edu/~katz/
31 
32     The problem is modeled by the partial differential equation system
33 
34 \begin{eqnarray}
35          -\nabla P + \nabla \cdot [\eta (\nabla v + \nabla v^T)] & = & 0  \\
36                                            \nabla \cdot v & = & 0   \\
37                     dT/dt + \nabla \cdot (vT) - 1/Pe \triangle^2(T) & = & 0  \\
38 \end{eqnarray}
39 
40  \begin{eqnarray}
41         \eta(T,Eps\_dot) &  = & \hbox{constant                        }    \hbox{if ivisc} ==0  \\
42                       &  = & \hbox{diffusion creep (T,P-dependent)    }     \hbox{if ivisc} ==1  \\
43                       &  = & \hbox{dislocation creep (T,P,v-dependent)}  \hbox{if ivisc} ==2  \\
44                       &  = & \hbox{mantle viscosity (difn and disl)   }  \hbox{if ivisc} ==3
45 \end{eqnarray}
46 
47     which is uniformly discretized on a staggered mesh:
48                       -------$w_{ij}$------
49                   $u_{i-1j}$    $P,T_{ij}$   $u_{ij}$
50                       ------$w_{ij-1}$-----
51 
52   ------------------------------------------------------------------------F*/
53 
54 #include <petscsnes.h>
55 #include <petscdm.h>
56 #include <petscdmda.h>
57 
58 #define VISC_CONST   0
59 #define VISC_DIFN    1
60 #define VISC_DISL    2
61 #define VISC_FULL    3
62 #define CELL_CENTER  0
63 #define CELL_CORNER  1
64 #define BC_ANALYTIC  0
65 #define BC_NOSTRESS  1
66 #define BC_EXPERMNT  2
67 #define ADVECT_FV    0
68 #define ADVECT_FROMM 1
69 #define PLATE_SLAB   0
70 #define PLATE_LID    1
71 #define EPS_ZERO     0.00000001
72 
73 typedef struct { /* holds the variables to be solved for */
74   PetscScalar u, w, p, T;
75 } Field;
76 
77 typedef struct { /* parameters needed to compute viscosity */
78   PetscReal A, n, Estar, Vstar;
79 } ViscParam;
80 
81 typedef struct { /* physical and miscelaneous parameters */
82   PetscReal width, depth, scaled_width, scaled_depth, peclet, potentialT;
83   PetscReal slab_dip, slab_age, slab_velocity, kappa, z_scale;
84   PetscReal c, d, sb, cb, skt, visc_cutoff, lid_age, eta0, continuation;
85   PetscReal L, V, lid_depth, fault_depth;
86   ViscParam diffusion, dislocation;
87   PetscInt  ivisc, adv_scheme, ibound, output_ivisc;
88   PetscBool quiet, param_test, output_to_file, pv_analytic;
89   PetscBool interrupted, stop_solve, toggle_kspmon, kspmon;
90   char      filename[PETSC_MAX_PATH_LEN];
91 } Parameter;
92 
93 typedef struct { /* grid parameters */
94   DMBoundaryType  bx, by;
95   DMDAStencilType stencil;
96   PetscInt        corner, ni, nj, jlid, jfault, inose;
97   PetscInt        dof, stencil_width, mglevels;
98   PetscReal       dx, dz;
99 } GridInfo;
100 
101 typedef struct { /* application context */
102   Vec        x, Xguess;
103   Parameter *param;
104   GridInfo  *grid;
105 } AppCtx;
106 
107 /* Callback functions (static interface) */
108 extern PetscErrorCode FormFunctionLocal(DMDALocalInfo *, Field **, Field **, void *);
109 
110 /* Main routines */
111 extern PetscErrorCode SetParams(Parameter *, GridInfo *);
112 extern PetscErrorCode ReportParams(Parameter *, GridInfo *);
113 extern PetscErrorCode Initialize(DM);
114 extern PetscErrorCode UpdateSolution(SNES, AppCtx *, PetscInt *);
115 extern PetscErrorCode DoOutput(SNES, PetscInt);
116 
117 /* Post-processing & misc */
118 extern PetscErrorCode ViscosityField(DM, Vec, Vec);
119 extern PetscErrorCode StressField(DM);
120 extern PetscErrorCode SNESConverged_Interactive(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
121 extern PetscErrorCode InteractiveHandler(int, void *);
122 
123 int main(int argc, char **argv) {
124   SNES      snes;
125   AppCtx   *user; /* user-defined work context */
126   Parameter param;
127   GridInfo  grid;
128   PetscInt  nits;
129   MPI_Comm  comm;
130   DM        da;
131 
132   PetscFunctionBeginUser;
133   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
134   PetscOptionsSetValue(NULL, "-file", "ex30_output");
135   PetscOptionsSetValue(NULL, "-snes_monitor_short", NULL);
136   PetscOptionsSetValue(NULL, "-snes_max_it", "20");
137   PetscOptionsSetValue(NULL, "-ksp_max_it", "1500");
138   PetscOptionsSetValue(NULL, "-ksp_gmres_restart", "300");
139   PetscOptionsInsert(NULL, &argc, &argv, NULL);
140 
141   comm = PETSC_COMM_WORLD;
142 
143   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144      Set up the problem parameters.
145      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146   PetscCall(SetParams(&param, &grid));
147   PetscCall(ReportParams(&param, &grid));
148 
149   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151   PetscCall(SNESCreate(comm, &snes));
152   PetscCall(DMDACreate2d(comm, grid.bx, grid.by, grid.stencil, grid.ni, grid.nj, PETSC_DECIDE, PETSC_DECIDE, grid.dof, grid.stencil_width, 0, 0, &da));
153   PetscCall(DMSetFromOptions(da));
154   PetscCall(DMSetUp(da));
155   PetscCall(SNESSetDM(snes, da));
156   PetscCall(DMDASetFieldName(da, 0, "x-velocity"));
157   PetscCall(DMDASetFieldName(da, 1, "y-velocity"));
158   PetscCall(DMDASetFieldName(da, 2, "pressure"));
159   PetscCall(DMDASetFieldName(da, 3, "temperature"));
160 
161   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162      Create user context, set problem data, create vector data structures.
163      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164   PetscCall(PetscNew(&user));
165   user->param = &param;
166   user->grid  = &grid;
167   PetscCall(DMSetApplicationContext(da, user));
168   PetscCall(DMCreateGlobalVector(da, &(user->Xguess)));
169 
170   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171      Set up the SNES solver with callback functions.
172      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, (void *)user));
174   PetscCall(SNESSetFromOptions(snes));
175 
176   PetscCall(SNESSetConvergenceTest(snes, SNESConverged_Interactive, (void *)user, NULL));
177   PetscCall(PetscPushSignalHandler(InteractiveHandler, (void *)user));
178 
179   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180      Initialize and solve the nonlinear system
181      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182   PetscCall(Initialize(da));
183   PetscCall(UpdateSolution(snes, user, &nits));
184 
185   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186      Output variables.
187      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
188   PetscCall(DoOutput(snes, nits));
189 
190   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191      Free work space.
192      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193   PetscCall(VecDestroy(&user->Xguess));
194   PetscCall(VecDestroy(&user->x));
195   PetscCall(PetscFree(user));
196   PetscCall(SNESDestroy(&snes));
197   PetscCall(DMDestroy(&da));
198   PetscCall(PetscPopSignalHandler());
199   PetscCall(PetscFinalize());
200   return 0;
201 }
202 
203 /*=====================================================================
204   PETSc INTERACTION FUNCTIONS (initialize & call SNESSolve)
205   =====================================================================*/
206 
207 /*  manages solve: adaptive continuation method  */
208 PetscErrorCode UpdateSolution(SNES snes, AppCtx *user, PetscInt *nits) {
209   KSP                 ksp;
210   PC                  pc;
211   SNESConvergedReason reason    = SNES_CONVERGED_ITERATING;
212   Parameter          *param     = user->param;
213   PetscReal           cont_incr = 0.3;
214   PetscInt            its;
215   PetscBool           q = PETSC_FALSE;
216   DM                  dm;
217 
218   PetscFunctionBeginUser;
219   PetscCall(SNESGetDM(snes, &dm));
220   PetscCall(DMCreateGlobalVector(dm, &user->x));
221   PetscCall(SNESGetKSP(snes, &ksp));
222   PetscCall(KSPGetPC(ksp, &pc));
223   PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
224 
225   *nits = 0;
226 
227   /* Isoviscous solve */
228   if (param->ivisc == VISC_CONST && !param->stop_solve) {
229     param->ivisc = VISC_CONST;
230 
231     PetscCall(SNESSolve(snes, 0, user->x));
232     PetscCall(SNESGetConvergedReason(snes, &reason));
233     PetscCall(SNESGetIterationNumber(snes, &its));
234     *nits += its;
235     PetscCall(VecCopy(user->x, user->Xguess));
236     if (param->stop_solve) goto done;
237   }
238 
239   /* Olivine diffusion creep */
240   if (param->ivisc >= VISC_DIFN && !param->stop_solve) {
241     if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Computing Variable Viscosity Solution\n"));
242 
243     /* continuation method on viscosity cutoff */
244     for (param->continuation = 0.0;; param->continuation += cont_incr) {
245       if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Continuation parameter = %g\n", (double)param->continuation));
246 
247       /* solve the non-linear system */
248       PetscCall(VecCopy(user->Xguess, user->x));
249       PetscCall(SNESSolve(snes, 0, user->x));
250       PetscCall(SNESGetConvergedReason(snes, &reason));
251       PetscCall(SNESGetIterationNumber(snes, &its));
252       *nits += its;
253       if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SNES iterations: %" PetscInt_FMT ", Cumulative: %" PetscInt_FMT "\n", its, *nits));
254       if (param->stop_solve) goto done;
255 
256       if (reason < 0) {
257         /* NOT converged */
258         cont_incr = -PetscAbsReal(cont_incr) / 2.0;
259         if (PetscAbsReal(cont_incr) < 0.01) goto done;
260 
261       } else {
262         /* converged */
263         PetscCall(VecCopy(user->x, user->Xguess));
264         if (param->continuation >= 1.0) goto done;
265         if (its <= 3) cont_incr = 0.30001;
266         else if (its <= 8) cont_incr = 0.15001;
267         else cont_incr = 0.10001;
268 
269         if (param->continuation + cont_incr > 1.0) cont_incr = 1.0 - param->continuation;
270       } /* endif reason<0 */
271     }
272   }
273 done:
274   if (param->stop_solve && !q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: stopping solve.\n"));
275   if (reason < 0 && !q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "FAILED TO CONVERGE: stopping solve.\n"));
276   PetscFunctionReturn(0);
277 }
278 
279 /*=====================================================================
280   PHYSICS FUNCTIONS (compute the discrete residual)
281   =====================================================================*/
282 
283 static inline PetscScalar UInterp(Field **x, PetscInt i, PetscInt j) {
284   return 0.25 * (x[j][i].u + x[j + 1][i].u + x[j][i + 1].u + x[j + 1][i + 1].u);
285 }
286 
287 static inline PetscScalar WInterp(Field **x, PetscInt i, PetscInt j) {
288   return 0.25 * (x[j][i].w + x[j + 1][i].w + x[j][i + 1].w + x[j + 1][i + 1].w);
289 }
290 
291 static inline PetscScalar PInterp(Field **x, PetscInt i, PetscInt j) {
292   return 0.25 * (x[j][i].p + x[j + 1][i].p + x[j][i + 1].p + x[j + 1][i + 1].p);
293 }
294 
295 static inline PetscScalar TInterp(Field **x, PetscInt i, PetscInt j) {
296   return 0.25 * (x[j][i].T + x[j + 1][i].T + x[j][i + 1].T + x[j + 1][i + 1].T);
297 }
298 
299 /*  isoviscous analytic solution for IC */
300 static inline PetscScalar HorizVelocity(PetscInt i, PetscInt j, AppCtx *user) {
301   Parameter  *param = user->param;
302   GridInfo   *grid  = user->grid;
303   PetscScalar st, ct, th, c = param->c, d = param->d;
304   PetscReal   x, z, r;
305 
306   x  = (i - grid->jlid) * grid->dx;
307   z  = (j - grid->jlid - 0.5) * grid->dz;
308   r  = PetscSqrtReal(x * x + z * z);
309   st = z / r;
310   ct = x / r;
311   th = PetscAtanReal(z / x);
312   return ct * (c * th * st + d * (st + th * ct)) + st * (c * (st - th * ct) + d * th * st);
313 }
314 
315 /*  isoviscous analytic solution for IC */
316 static inline PetscScalar VertVelocity(PetscInt i, PetscInt j, AppCtx *user)
317 
318 {
319   Parameter  *param = user->param;
320   GridInfo   *grid  = user->grid;
321   PetscScalar st, ct, th, c = param->c, d = param->d;
322   PetscReal   x, z, r;
323 
324   x  = (i - grid->jlid - 0.5) * grid->dx;
325   z  = (j - grid->jlid) * grid->dz;
326   r  = PetscSqrtReal(x * x + z * z);
327   st = z / r;
328   ct = x / r;
329   th = PetscAtanReal(z / x);
330   return st * (c * th * st + d * (st + th * ct)) - ct * (c * (st - th * ct) + d * th * st);
331 }
332 
333 /*  isoviscous analytic solution for IC */
334 static inline PetscScalar Pressure(PetscInt i, PetscInt j, AppCtx *user) {
335   Parameter  *param = user->param;
336   GridInfo   *grid  = user->grid;
337   PetscScalar x, z, r, st, ct, c = param->c, d = param->d;
338 
339   x  = (i - grid->jlid - 0.5) * grid->dx;
340   z  = (j - grid->jlid - 0.5) * grid->dz;
341   r  = PetscSqrtReal(x * x + z * z);
342   st = z / r;
343   ct = x / r;
344   return (-2.0 * (c * ct - d * st) / r);
345 }
346 
347 /*  computes the second invariant of the strain rate tensor */
348 static inline PetscScalar CalcSecInv(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user) {
349   Parameter  *param = user->param;
350   GridInfo   *grid  = user->grid;
351   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1;
352   PetscScalar uN, uS, uE, uW, wN, wS, wE, wW;
353   PetscScalar eps11, eps12, eps22;
354 
355   if (i < j) return EPS_ZERO;
356   if (i == ilim) i--;
357   if (j == jlim) j--;
358 
359   if (ipos == CELL_CENTER) { /* on cell center */
360     if (j <= grid->jlid) return EPS_ZERO;
361 
362     uE = x[j][i].u;
363     uW = x[j][i - 1].u;
364     wN = x[j][i].w;
365     wS = x[j - 1][i].w;
366     wE = WInterp(x, i, j - 1);
367     if (i == j) {
368       uN = param->cb;
369       wW = param->sb;
370     } else {
371       uN = UInterp(x, i - 1, j);
372       wW = WInterp(x, i - 1, j - 1);
373     }
374 
375     if (j == grid->jlid + 1) uS = 0.0;
376     else uS = UInterp(x, i - 1, j - 1);
377 
378   } else { /* on CELL_CORNER */
379     if (j < grid->jlid) return EPS_ZERO;
380 
381     uN = x[j + 1][i].u;
382     uS = x[j][i].u;
383     wE = x[j][i + 1].w;
384     wW = x[j][i].w;
385     if (i == j) {
386       wN = param->sb;
387       uW = param->cb;
388     } else {
389       wN = WInterp(x, i, j);
390       uW = UInterp(x, i - 1, j);
391     }
392 
393     if (j == grid->jlid) {
394       uE = 0.0;
395       uW = 0.0;
396       uS = -uN;
397       wS = -wN;
398     } else {
399       uE = UInterp(x, i, j);
400       wS = WInterp(x, i, j - 1);
401     }
402   }
403 
404   eps11 = (uE - uW) / grid->dx;
405   eps22 = (wN - wS) / grid->dz;
406   eps12 = 0.5 * ((uN - uS) / grid->dz + (wE - wW) / grid->dx);
407 
408   return PetscSqrtReal(0.5 * (eps11 * eps11 + 2.0 * eps12 * eps12 + eps22 * eps22));
409 }
410 
411 /*  computes the shear viscosity */
412 static inline PetscScalar Viscosity(PetscScalar T, PetscScalar eps, PetscScalar z, Parameter *param) {
413   PetscReal   result = 0.0;
414   ViscParam   difn = param->diffusion, disl = param->dislocation;
415   PetscInt    iVisc     = param->ivisc;
416   PetscScalar eps_scale = param->V / (param->L * 1000.0);
417   PetscScalar strain_power, v1, v2, P;
418   PetscScalar rho_g = 32340.0, R = 8.3144;
419 
420   P = rho_g * (z * param->L * 1000.0); /* Pa */
421 
422   if (iVisc == VISC_CONST) {
423     /* constant viscosity */
424     return 1.0;
425   } else if (iVisc == VISC_DIFN) {
426     /* diffusion creep rheology */
427     result = PetscRealPart((difn.A * PetscExpScalar((difn.Estar + P * difn.Vstar) / R / (T + 273.0)) / param->eta0));
428   } else if (iVisc == VISC_DISL) {
429     /* dislocation creep rheology */
430     strain_power = PetscPowScalar(eps * eps_scale, (1.0 - disl.n) / disl.n);
431 
432     result = PetscRealPart(disl.A * PetscExpScalar((disl.Estar + P * disl.Vstar) / disl.n / R / (T + 273.0)) * strain_power / param->eta0);
433   } else if (iVisc == VISC_FULL) {
434     /* dislocation/diffusion creep rheology */
435     strain_power = PetscPowScalar(eps * eps_scale, (1.0 - disl.n) / disl.n);
436 
437     v1 = difn.A * PetscExpScalar((difn.Estar + P * difn.Vstar) / R / (T + 273.0)) / param->eta0;
438     v2 = disl.A * PetscExpScalar((disl.Estar + P * disl.Vstar) / disl.n / R / (T + 273.0)) * strain_power / param->eta0;
439 
440     result = PetscRealPart(1.0 / (1.0 / v1 + 1.0 / v2));
441   }
442 
443   /* max viscosity is param->eta0 */
444   result = PetscMin(result, 1.0);
445   /* min viscosity is param->visc_cutoff */
446   result = PetscMax(result, param->visc_cutoff);
447   /* continuation method */
448   result = PetscPowReal(result, param->continuation);
449   return result;
450 }
451 
452 /*  computes the residual of the x-component of eqn (1) above */
453 static inline PetscScalar XMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user) {
454   Parameter  *param = user->param;
455   GridInfo   *grid  = user->grid;
456   PetscScalar dx = grid->dx, dz = grid->dz;
457   PetscScalar etaN, etaS, etaE, etaW, epsN = 0.0, epsS = 0.0, epsE = 0.0, epsW = 0.0;
458   PetscScalar TE = 0.0, TN = 0.0, TS = 0.0, TW = 0.0, dPdx, residual, z_scale;
459   PetscScalar dudxW, dudxE, dudzN, dudzS, dwdxN, dwdxS;
460   PetscInt    jlim = grid->nj - 1;
461 
462   z_scale = param->z_scale;
463 
464   if (param->ivisc == VISC_DIFN || param->ivisc >= VISC_DISL) { /* viscosity is T-dependent */
465     TS = param->potentialT * TInterp(x, i, j - 1) * PetscExpScalar((j - 1.0) * dz * z_scale);
466     if (j == jlim) TN = TS;
467     else TN = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
468     TW = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
469     TE = param->potentialT * x[j][i + 1].T * PetscExpScalar((j - 0.5) * dz * z_scale);
470     if (param->ivisc >= VISC_DISL) { /* olivine dislocation creep */
471       epsN = CalcSecInv(x, i, j, CELL_CORNER, user);
472       epsS = CalcSecInv(x, i, j - 1, CELL_CORNER, user);
473       epsE = CalcSecInv(x, i + 1, j, CELL_CENTER, user);
474       epsW = CalcSecInv(x, i, j, CELL_CENTER, user);
475     }
476   }
477   etaN = Viscosity(TN, epsN, dz * (j + 0.5), param);
478   etaS = Viscosity(TS, epsS, dz * (j - 0.5), param);
479   etaW = Viscosity(TW, epsW, dz * j, param);
480   etaE = Viscosity(TE, epsE, dz * j, param);
481 
482   dPdx = (x[j][i + 1].p - x[j][i].p) / dx;
483   if (j == jlim) dudzN = etaN * (x[j][i].w - x[j][i + 1].w) / dx;
484   else dudzN = etaN * (x[j + 1][i].u - x[j][i].u) / dz;
485   dudzS = etaS * (x[j][i].u - x[j - 1][i].u) / dz;
486   dudxE = etaE * (x[j][i + 1].u - x[j][i].u) / dx;
487   dudxW = etaW * (x[j][i].u - x[j][i - 1].u) / dx;
488 
489   residual = -dPdx /* X-MOMENTUM EQUATION*/
490            + (dudxE - dudxW) / dx + (dudzN - dudzS) / dz;
491 
492   if (param->ivisc != VISC_CONST) {
493     dwdxN = etaN * (x[j][i + 1].w - x[j][i].w) / dx;
494     dwdxS = etaS * (x[j - 1][i + 1].w - x[j - 1][i].w) / dx;
495 
496     residual += (dudxE - dudxW) / dx + (dwdxN - dwdxS) / dz;
497   }
498 
499   return residual;
500 }
501 
502 /*  computes the residual of the z-component of eqn (1) above */
503 static inline PetscScalar ZMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
504 
505 {
506   Parameter  *param = user->param;
507   GridInfo   *grid  = user->grid;
508   PetscScalar dx = grid->dx, dz = grid->dz;
509   PetscScalar etaN = 0.0, etaS = 0.0, etaE = 0.0, etaW = 0.0, epsN = 0.0, epsS = 0.0, epsE = 0.0, epsW = 0.0;
510   PetscScalar TE = 0.0, TN = 0.0, TS = 0.0, TW = 0.0, dPdz, residual, z_scale;
511   PetscScalar dudzE, dudzW, dwdxW, dwdxE, dwdzN, dwdzS;
512   PetscInt    ilim = grid->ni - 1;
513 
514   /* geometric and other parameters */
515   z_scale = param->z_scale;
516 
517   /* viscosity */
518   if (param->ivisc == VISC_DIFN || param->ivisc >= VISC_DISL) { /* viscosity is T-dependent */
519     TN = param->potentialT * x[j + 1][i].T * PetscExpScalar((j + 0.5) * dz * z_scale);
520     TS = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
521     TW = param->potentialT * TInterp(x, i - 1, j) * PetscExpScalar(j * dz * z_scale);
522     if (i == ilim) TE = TW;
523     else TE = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
524     if (param->ivisc >= VISC_DISL) { /* olivine dislocation creep */
525       epsN = CalcSecInv(x, i, j + 1, CELL_CENTER, user);
526       epsS = CalcSecInv(x, i, j, CELL_CENTER, user);
527       epsE = CalcSecInv(x, i, j, CELL_CORNER, user);
528       epsW = CalcSecInv(x, i - 1, j, CELL_CORNER, user);
529     }
530   }
531   etaN = Viscosity(TN, epsN, dz * (j + 1.0), param);
532   etaS = Viscosity(TS, epsS, dz * (j + 0.0), param);
533   etaW = Viscosity(TW, epsW, dz * (j + 0.5), param);
534   etaE = Viscosity(TE, epsE, dz * (j + 0.5), param);
535 
536   dPdz  = (x[j + 1][i].p - x[j][i].p) / dz;
537   dwdzN = etaN * (x[j + 1][i].w - x[j][i].w) / dz;
538   dwdzS = etaS * (x[j][i].w - x[j - 1][i].w) / dz;
539   if (i == ilim) dwdxE = etaE * (x[j][i].u - x[j + 1][i].u) / dz;
540   else dwdxE = etaE * (x[j][i + 1].w - x[j][i].w) / dx;
541   dwdxW = 2.0 * etaW * (x[j][i].w - x[j][i - 1].w) / dx;
542 
543   /* Z-MOMENTUM */
544   residual = -dPdz /* constant viscosity terms */
545            + (dwdzN - dwdzS) / dz + (dwdxE - dwdxW) / dx;
546 
547   if (param->ivisc != VISC_CONST) {
548     dudzE = etaE * (x[j + 1][i].u - x[j][i].u) / dz;
549     dudzW = etaW * (x[j + 1][i - 1].u - x[j][i - 1].u) / dz;
550 
551     residual += (dwdzN - dwdzS) / dz + (dudzE - dudzW) / dx;
552   }
553 
554   return residual;
555 }
556 
557 /*  computes the residual of eqn (2) above */
558 static inline PetscScalar ContinuityResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user) {
559   GridInfo   *grid = user->grid;
560   PetscScalar uE, uW, wN, wS, dudx, dwdz;
561 
562   uW   = x[j][i - 1].u;
563   uE   = x[j][i].u;
564   dudx = (uE - uW) / grid->dx;
565   wS   = x[j - 1][i].w;
566   wN   = x[j][i].w;
567   dwdz = (wN - wS) / grid->dz;
568 
569   return dudx + dwdz;
570 }
571 
572 /*  computes the residual of eqn (3) above */
573 static inline PetscScalar EnergyResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user) {
574   Parameter  *param = user->param;
575   GridInfo   *grid  = user->grid;
576   PetscScalar dx = grid->dx, dz = grid->dz;
577   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1, jlid = grid->jlid;
578   PetscScalar TE, TN, TS, TW, residual;
579   PetscScalar uE, uW, wN, wS;
580   PetscScalar fN, fS, fE, fW, dTdxW, dTdxE, dTdzN, dTdzS;
581 
582   dTdzN = (x[j + 1][i].T - x[j][i].T) / dz;
583   dTdzS = (x[j][i].T - x[j - 1][i].T) / dz;
584   dTdxE = (x[j][i + 1].T - x[j][i].T) / dx;
585   dTdxW = (x[j][i].T - x[j][i - 1].T) / dx;
586 
587   residual = ((dTdzN - dTdzS) / dz + /* diffusion term */
588               (dTdxE - dTdxW) / dx) *
589              dx * dz / param->peclet;
590 
591   if (j <= jlid && i >= j) {
592     /* don't advect in the lid */
593     return residual;
594   } else if (i < j) {
595     /* beneath the slab sfc */
596     uW = uE = param->cb;
597     wS = wN = param->sb;
598   } else {
599     /* advect in the slab and wedge */
600     uW = x[j][i - 1].u;
601     uE = x[j][i].u;
602     wS = x[j - 1][i].w;
603     wN = x[j][i].w;
604   }
605 
606   if (param->adv_scheme == ADVECT_FV || i == ilim - 1 || j == jlim - 1 || i == 1 || j == 1) {
607     /* finite volume advection */
608     TS = (x[j][i].T + x[j - 1][i].T) / 2.0;
609     TN = (x[j][i].T + x[j + 1][i].T) / 2.0;
610     TE = (x[j][i].T + x[j][i + 1].T) / 2.0;
611     TW = (x[j][i].T + x[j][i - 1].T) / 2.0;
612     fN = wN * TN * dx;
613     fS = wS * TS * dx;
614     fE = uE * TE * dz;
615     fW = uW * TW * dz;
616 
617   } else {
618     /* Fromm advection scheme */
619     fE = (uE * (-x[j][i + 2].T + 5.0 * (x[j][i + 1].T + x[j][i].T) - x[j][i - 1].T) / 8.0 - PetscAbsScalar(uE) * (-x[j][i + 2].T + 3.0 * (x[j][i + 1].T - x[j][i].T) + x[j][i - 1].T) / 8.0) * dz;
620     fW = (uW * (-x[j][i + 1].T + 5.0 * (x[j][i].T + x[j][i - 1].T) - x[j][i - 2].T) / 8.0 - PetscAbsScalar(uW) * (-x[j][i + 1].T + 3.0 * (x[j][i].T - x[j][i - 1].T) + x[j][i - 2].T) / 8.0) * dz;
621     fN = (wN * (-x[j + 2][i].T + 5.0 * (x[j + 1][i].T + x[j][i].T) - x[j - 1][i].T) / 8.0 - PetscAbsScalar(wN) * (-x[j + 2][i].T + 3.0 * (x[j + 1][i].T - x[j][i].T) + x[j - 1][i].T) / 8.0) * dx;
622     fS = (wS * (-x[j + 1][i].T + 5.0 * (x[j][i].T + x[j - 1][i].T) - x[j - 2][i].T) / 8.0 - PetscAbsScalar(wS) * (-x[j + 1][i].T + 3.0 * (x[j][i].T - x[j - 1][i].T) + x[j - 2][i].T) / 8.0) * dx;
623   }
624 
625   residual -= (fE - fW + fN - fS);
626 
627   return residual;
628 }
629 
630 /*  computes the shear stress---used on the boundaries */
631 static inline PetscScalar ShearStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user) {
632   Parameter  *param = user->param;
633   GridInfo   *grid  = user->grid;
634   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1;
635   PetscScalar uN, uS, wE, wW;
636 
637   if (j <= grid->jlid || i < j || i == ilim || j == jlim) return EPS_ZERO;
638 
639   if (ipos == CELL_CENTER) { /* on cell center */
640 
641     wE = WInterp(x, i, j - 1);
642     if (i == j) {
643       wW = param->sb;
644       uN = param->cb;
645     } else {
646       wW = WInterp(x, i - 1, j - 1);
647       uN = UInterp(x, i - 1, j);
648     }
649     if (j == grid->jlid + 1) uS = 0.0;
650     else uS = UInterp(x, i - 1, j - 1);
651 
652   } else { /* on cell corner */
653 
654     uN = x[j + 1][i].u;
655     uS = x[j][i].u;
656     wW = x[j][i].w;
657     wE = x[j][i + 1].w;
658   }
659 
660   return (uN - uS) / grid->dz + (wE - wW) / grid->dx;
661 }
662 
663 /*  computes the normal stress---used on the boundaries */
664 static inline PetscScalar XNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user) {
665   Parameter  *param = user->param;
666   GridInfo   *grid  = user->grid;
667   PetscScalar dx = grid->dx, dz = grid->dz;
668   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1, ivisc;
669   PetscScalar epsC = 0.0, etaC, TC, uE, uW, pC, z_scale;
670   if (i < j || j <= grid->jlid) return EPS_ZERO;
671 
672   ivisc   = param->ivisc;
673   z_scale = param->z_scale;
674 
675   if (ipos == CELL_CENTER) { /* on cell center */
676 
677     TC = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
678     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CENTER, user);
679     etaC = Viscosity(TC, epsC, dz * j, param);
680 
681     uW = x[j][i - 1].u;
682     uE = x[j][i].u;
683     pC = x[j][i].p;
684 
685   } else { /* on cell corner */
686     if (i == ilim || j == jlim) return EPS_ZERO;
687 
688     TC = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
689     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CORNER, user);
690     etaC = Viscosity(TC, epsC, dz * (j + 0.5), param);
691 
692     if (i == j) uW = param->sb;
693     else uW = UInterp(x, i - 1, j);
694     uE = UInterp(x, i, j);
695     pC = PInterp(x, i, j);
696   }
697 
698   return 2.0 * etaC * (uE - uW) / dx - pC;
699 }
700 
701 /*  computes the normal stress---used on the boundaries */
702 static inline PetscScalar ZNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user) {
703   Parameter  *param = user->param;
704   GridInfo   *grid  = user->grid;
705   PetscScalar dz    = grid->dz;
706   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1, ivisc;
707   PetscScalar epsC = 0.0, etaC, TC;
708   PetscScalar pC, wN, wS, z_scale;
709   if (i < j || j <= grid->jlid) return EPS_ZERO;
710 
711   ivisc   = param->ivisc;
712   z_scale = param->z_scale;
713 
714   if (ipos == CELL_CENTER) { /* on cell center */
715 
716     TC = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
717     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CENTER, user);
718     etaC = Viscosity(TC, epsC, dz * j, param);
719     wN   = x[j][i].w;
720     wS   = x[j - 1][i].w;
721     pC   = x[j][i].p;
722 
723   } else { /* on cell corner */
724     if ((i == ilim) || (j == jlim)) return EPS_ZERO;
725 
726     TC = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
727     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CORNER, user);
728     etaC = Viscosity(TC, epsC, dz * (j + 0.5), param);
729     if (i == j) wN = param->sb;
730     else wN = WInterp(x, i, j);
731     wS = WInterp(x, i, j - 1);
732     pC = PInterp(x, i, j);
733   }
734 
735   return 2.0 * etaC * (wN - wS) / dz - pC;
736 }
737 
738 /*=====================================================================
739   INITIALIZATION, POST-PROCESSING AND OUTPUT FUNCTIONS
740   =====================================================================*/
741 
742 /* initializes the problem parameters and checks for
743    command line changes */
744 PetscErrorCode SetParams(Parameter *param, GridInfo *grid) {
745   PetscReal SEC_PER_YR                     = 3600.00 * 24.00 * 365.2500;
746   PetscReal alpha_g_on_cp_units_inverse_km = 4.0e-5 * 9.8;
747 
748   /* domain geometry */
749   param->slab_dip    = 45.0;
750   param->width       = 320.0; /* km */
751   param->depth       = 300.0; /* km */
752   param->lid_depth   = 35.0;  /* km */
753   param->fault_depth = 35.0;  /* km */
754 
755   PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_dip", &(param->slab_dip), NULL));
756   PetscCall(PetscOptionsGetReal(NULL, NULL, "-width", &(param->width), NULL));
757   PetscCall(PetscOptionsGetReal(NULL, NULL, "-depth", &(param->depth), NULL));
758   PetscCall(PetscOptionsGetReal(NULL, NULL, "-lid_depth", &(param->lid_depth), NULL));
759   PetscCall(PetscOptionsGetReal(NULL, NULL, "-fault_depth", &(param->fault_depth), NULL));
760 
761   param->slab_dip = param->slab_dip * PETSC_PI / 180.0; /* radians */
762 
763   /* grid information */
764   PetscCall(PetscOptionsGetInt(NULL, NULL, "-jfault", &(grid->jfault), NULL));
765   grid->ni = 82;
766   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ni", &(grid->ni), NULL));
767 
768   grid->dx     = param->width / ((PetscReal)(grid->ni - 2)); /* km */
769   grid->dz     = grid->dx * PetscTanReal(param->slab_dip);   /* km */
770   grid->nj     = (PetscInt)(param->depth / grid->dz + 3.0);  /* gridpoints*/
771   param->depth = grid->dz * (grid->nj - 2);                  /* km */
772   grid->inose  = 0;                                          /* gridpoints*/
773   PetscCall(PetscOptionsGetInt(NULL, NULL, "-inose", &(grid->inose), NULL));
774   grid->bx            = DM_BOUNDARY_NONE;
775   grid->by            = DM_BOUNDARY_NONE;
776   grid->stencil       = DMDA_STENCIL_BOX;
777   grid->dof           = 4;
778   grid->stencil_width = 2;
779   grid->mglevels      = 1;
780 
781   /* boundary conditions */
782   param->pv_analytic = PETSC_FALSE;
783   param->ibound      = BC_NOSTRESS;
784   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ibound", &(param->ibound), NULL));
785 
786   /* physical constants */
787   param->slab_velocity = 5.0;       /* cm/yr */
788   param->slab_age      = 50.0;      /* Ma */
789   param->lid_age       = 50.0;      /* Ma */
790   param->kappa         = 0.7272e-6; /* m^2/sec */
791   param->potentialT    = 1300.0;    /* degrees C */
792 
793   PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_velocity", &(param->slab_velocity), NULL));
794   PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_age", &(param->slab_age), NULL));
795   PetscCall(PetscOptionsGetReal(NULL, NULL, "-lid_age", &(param->lid_age), NULL));
796   PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", &(param->kappa), NULL));
797   PetscCall(PetscOptionsGetReal(NULL, NULL, "-potentialT", &(param->potentialT), NULL));
798 
799   /* viscosity */
800   param->ivisc        = 3;    /* 0=isovisc, 1=difn creep, 2=disl creep, 3=full */
801   param->eta0         = 1e24; /* Pa-s */
802   param->visc_cutoff  = 0.0;  /* factor of eta_0 */
803   param->continuation = 1.0;
804 
805   /* constants for diffusion creep */
806   param->diffusion.A     = 1.8e7; /* Pa-s */
807   param->diffusion.n     = 1.0;   /* dim'less */
808   param->diffusion.Estar = 375e3; /* J/mol */
809   param->diffusion.Vstar = 5e-6;  /* m^3/mol */
810 
811   /* constants for param->dislocationocation creep */
812   param->dislocation.A     = 2.8969e4; /* Pa-s */
813   param->dislocation.n     = 3.5;      /* dim'less */
814   param->dislocation.Estar = 530e3;    /* J/mol */
815   param->dislocation.Vstar = 14e-6;    /* m^3/mol */
816 
817   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ivisc", &(param->ivisc), NULL));
818   PetscCall(PetscOptionsGetReal(NULL, NULL, "-visc_cutoff", &(param->visc_cutoff), NULL));
819 
820   param->output_ivisc = param->ivisc;
821 
822   PetscCall(PetscOptionsGetInt(NULL, NULL, "-output_ivisc", &(param->output_ivisc), NULL));
823   PetscCall(PetscOptionsGetReal(NULL, NULL, "-vstar", &(param->dislocation.Vstar), NULL));
824 
825   /* output options */
826   param->quiet      = PETSC_FALSE;
827   param->param_test = PETSC_FALSE;
828 
829   PetscCall(PetscOptionsHasName(NULL, NULL, "-quiet", &(param->quiet)));
830   PetscCall(PetscOptionsHasName(NULL, NULL, "-test", &(param->param_test)));
831   PetscCall(PetscOptionsGetString(NULL, NULL, "-file", param->filename, sizeof(param->filename), &(param->output_to_file)));
832 
833   /* advection */
834   param->adv_scheme = ADVECT_FROMM; /* advection scheme: 0=finite vol, 1=Fromm */
835 
836   PetscCall(PetscOptionsGetInt(NULL, NULL, "-adv_scheme", &(param->adv_scheme), NULL));
837 
838   /* misc. flags */
839   param->stop_solve    = PETSC_FALSE;
840   param->interrupted   = PETSC_FALSE;
841   param->kspmon        = PETSC_FALSE;
842   param->toggle_kspmon = PETSC_FALSE;
843 
844   /* derived parameters for slab angle */
845   param->sb = PetscSinReal(param->slab_dip);
846   param->cb = PetscCosReal(param->slab_dip);
847   param->c  = param->slab_dip * param->sb / (param->slab_dip * param->slab_dip - param->sb * param->sb);
848   param->d  = (param->slab_dip * param->cb - param->sb) / (param->slab_dip * param->slab_dip - param->sb * param->sb);
849 
850   /* length, velocity and time scale for non-dimensionalization */
851   param->L = PetscMin(param->width, param->depth);      /* km */
852   param->V = param->slab_velocity / 100.0 / SEC_PER_YR; /* m/sec */
853 
854   /* other unit conversions and derived parameters */
855   param->scaled_width = param->width / param->L;                   /* dim'less */
856   param->scaled_depth = param->depth / param->L;                   /* dim'less */
857   param->lid_depth    = param->lid_depth / param->L;               /* dim'less */
858   param->fault_depth  = param->fault_depth / param->L;             /* dim'less */
859   grid->dx            = grid->dx / param->L;                       /* dim'less */
860   grid->dz            = grid->dz / param->L;                       /* dim'less */
861   grid->jlid          = (PetscInt)(param->lid_depth / grid->dz);   /* gridcells */
862   grid->jfault        = (PetscInt)(param->fault_depth / grid->dz); /* gridcells */
863   param->lid_depth    = grid->jlid * grid->dz;                     /* dim'less */
864   param->fault_depth  = grid->jfault * grid->dz;                   /* dim'less */
865   grid->corner        = grid->jlid + 1;                            /* gridcells */
866   param->peclet       = param->V                                   /* m/sec */
867                 * param->L * 1000.0                                /* m */
868                 / param->kappa;                                    /* m^2/sec */
869   param->z_scale = param->L * alpha_g_on_cp_units_inverse_km;
870   param->skt     = PetscSqrtReal(param->kappa * param->slab_age * SEC_PER_YR);
871   PetscCall(PetscOptionsGetReal(NULL, NULL, "-peclet", &(param->peclet), NULL));
872 
873   return 0;
874 }
875 
876 /*  prints a report of the problem parameters to stdout */
877 PetscErrorCode ReportParams(Parameter *param, GridInfo *grid) {
878   char date[30];
879 
880   PetscCall(PetscGetDate(date, 30));
881 
882   if (!(param->quiet)) {
883     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---------------------BEGIN ex30 PARAM REPORT-------------------\n"));
884     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Domain: \n"));
885     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Width = %g km,         Depth = %g km\n", (double)param->width, (double)param->depth));
886     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Slab dip = %g degrees,  Slab velocity = %g cm/yr\n", (double)(param->slab_dip * 180.0 / PETSC_PI), (double)param->slab_velocity));
887     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Lid depth = %5.2f km,   Fault depth = %5.2f km\n", (double)(param->lid_depth * param->L), (double)(param->fault_depth * param->L)));
888 
889     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nGrid: \n"));
890     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  [ni,nj] = %" PetscInt_FMT ", %" PetscInt_FMT "       [dx,dz] = %g, %g km\n", grid->ni, grid->nj, (double)(grid->dx * param->L), (double)(grid->dz * param->L)));
891     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  jlid = %3" PetscInt_FMT "              jfault = %3" PetscInt_FMT " \n", grid->jlid, grid->jfault));
892     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Pe = %g\n", (double)param->peclet));
893 
894     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nRheology:"));
895     if (param->ivisc == VISC_CONST) {
896       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Isoviscous \n"));
897       if (param->pv_analytic) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Pressure and Velocity prescribed! \n"));
898     } else if (param->ivisc == VISC_DIFN) {
899       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Diffusion Creep (T-Dependent Newtonian) \n"));
900       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
901     } else if (param->ivisc == VISC_DISL) {
902       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Dislocation Creep (T-Dependent Non-Newtonian) \n"));
903       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
904     } else if (param->ivisc == VISC_FULL) {
905       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Full Rheology \n"));
906       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
907     } else {
908       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Invalid! \n"));
909       return 1;
910     }
911 
912     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Boundary condition:"));
913     if (param->ibound == BC_ANALYTIC) {
914       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Isoviscous Analytic Dirichlet \n"));
915     } else if (param->ibound == BC_NOSTRESS) {
916       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Stress-Free (normal & shear stress)\n"));
917     } else if (param->ibound == BC_EXPERMNT) {
918       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Experimental boundary condition \n"));
919     } else {
920       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Invalid! \n"));
921       return 1;
922     }
923 
924     if (param->output_to_file) {
925 #if defined(PETSC_HAVE_MATLAB_ENGINE)
926       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output Destination:       Mat file \"%s\"\n", param->filename));
927 #else
928       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output Destination:       PETSc binary file \"%s\"\n", param->filename));
929 #endif
930     }
931     if (param->output_ivisc != param->ivisc) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Output viscosity: -ivisc %" PetscInt_FMT "\n", param->output_ivisc));
932 
933     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---------------------END ex30 PARAM REPORT---------------------\n"));
934   }
935   if (param->param_test) PetscEnd();
936   return 0;
937 }
938 
939 /* ------------------------------------------------------------------- */
940 /*  generates an initial guess using the analytic solution for isoviscous
941     corner flow */
942 PetscErrorCode Initialize(DM da)
943 /* ------------------------------------------------------------------- */
944 {
945   AppCtx    *user;
946   Parameter *param;
947   GridInfo  *grid;
948   PetscInt   i, j, is, js, im, jm;
949   Field    **x;
950   Vec        Xguess;
951 
952   /* Get the fine grid */
953   PetscCall(DMGetApplicationContext(da, &user));
954   Xguess = user->Xguess;
955   param  = user->param;
956   grid   = user->grid;
957   PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
958   PetscCall(DMDAVecGetArray(da, Xguess, (void **)&x));
959 
960   /* Compute initial guess */
961   for (j = js; j < js + jm; j++) {
962     for (i = is; i < is + im; i++) {
963       if (i < j) x[j][i].u = param->cb;
964       else if (j <= grid->jlid) x[j][i].u = 0.0;
965       else x[j][i].u = HorizVelocity(i, j, user);
966 
967       if (i <= j) x[j][i].w = param->sb;
968       else if (j <= grid->jlid) x[j][i].w = 0.0;
969       else x[j][i].w = VertVelocity(i, j, user);
970 
971       if (i < j || j <= grid->jlid) x[j][i].p = 0.0;
972       else x[j][i].p = Pressure(i, j, user);
973 
974       x[j][i].T = PetscMin(grid->dz * (j - 0.5), 1.0);
975     }
976   }
977 
978   /* Restore x to Xguess */
979   PetscCall(DMDAVecRestoreArray(da, Xguess, (void **)&x));
980 
981   return 0;
982 }
983 
984 /*  controls output to a file */
985 PetscErrorCode DoOutput(SNES snes, PetscInt its) {
986   AppCtx     *user;
987   Parameter  *param;
988   GridInfo   *grid;
989   PetscInt    ivt;
990   PetscMPIInt rank;
991   PetscViewer viewer;
992   Vec         res, pars;
993   MPI_Comm    comm;
994   DM          da;
995 
996   PetscCall(SNESGetDM(snes, &da));
997   PetscCall(DMGetApplicationContext(da, &user));
998   param = user->param;
999   grid  = user->grid;
1000   ivt   = param->ivisc;
1001 
1002   param->ivisc = param->output_ivisc;
1003 
1004   /* compute final residual and final viscosity/strain rate fields */
1005   PetscCall(SNESGetFunction(snes, &res, NULL, NULL));
1006   PetscCall(ViscosityField(da, user->x, user->Xguess));
1007 
1008   /* get the communicator and the rank of the processor */
1009   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
1010   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1011 
1012   if (param->output_to_file) { /* send output to binary file */
1013     PetscCall(VecCreate(comm, &pars));
1014     if (rank == 0) { /* on processor 0 */
1015       PetscCall(VecSetSizes(pars, 20, PETSC_DETERMINE));
1016       PetscCall(VecSetFromOptions(pars));
1017       PetscCall(VecSetValue(pars, 0, (PetscScalar)(grid->ni), INSERT_VALUES));
1018       PetscCall(VecSetValue(pars, 1, (PetscScalar)(grid->nj), INSERT_VALUES));
1019       PetscCall(VecSetValue(pars, 2, (PetscScalar)(grid->dx), INSERT_VALUES));
1020       PetscCall(VecSetValue(pars, 3, (PetscScalar)(grid->dz), INSERT_VALUES));
1021       PetscCall(VecSetValue(pars, 4, (PetscScalar)(param->L), INSERT_VALUES));
1022       PetscCall(VecSetValue(pars, 5, (PetscScalar)(param->V), INSERT_VALUES));
1023       /* skipped 6 intentionally */
1024       PetscCall(VecSetValue(pars, 7, (PetscScalar)(param->slab_dip), INSERT_VALUES));
1025       PetscCall(VecSetValue(pars, 8, (PetscScalar)(grid->jlid), INSERT_VALUES));
1026       PetscCall(VecSetValue(pars, 9, (PetscScalar)(param->lid_depth), INSERT_VALUES));
1027       PetscCall(VecSetValue(pars, 10, (PetscScalar)(grid->jfault), INSERT_VALUES));
1028       PetscCall(VecSetValue(pars, 11, (PetscScalar)(param->fault_depth), INSERT_VALUES));
1029       PetscCall(VecSetValue(pars, 12, (PetscScalar)(param->potentialT), INSERT_VALUES));
1030       PetscCall(VecSetValue(pars, 13, (PetscScalar)(param->ivisc), INSERT_VALUES));
1031       PetscCall(VecSetValue(pars, 14, (PetscScalar)(param->visc_cutoff), INSERT_VALUES));
1032       PetscCall(VecSetValue(pars, 15, (PetscScalar)(param->ibound), INSERT_VALUES));
1033       PetscCall(VecSetValue(pars, 16, (PetscScalar)(its), INSERT_VALUES));
1034     } else { /* on some other processor */
1035       PetscCall(VecSetSizes(pars, 0, PETSC_DETERMINE));
1036       PetscCall(VecSetFromOptions(pars));
1037     }
1038     PetscCall(VecAssemblyBegin(pars));
1039     PetscCall(VecAssemblyEnd(pars));
1040 
1041     /* create viewer */
1042 #if defined(PETSC_HAVE_MATLAB_ENGINE)
1043     PetscCall(PetscViewerMatlabOpen(PETSC_COMM_WORLD, param->filename, FILE_MODE_WRITE, &viewer));
1044 #else
1045     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, param->filename, FILE_MODE_WRITE, &viewer));
1046 #endif
1047 
1048     /* send vectors to viewer */
1049     PetscCall(PetscObjectSetName((PetscObject)res, "res"));
1050     PetscCall(VecView(res, viewer));
1051     PetscCall(PetscObjectSetName((PetscObject)user->x, "out"));
1052     PetscCall(VecView(user->x, viewer));
1053     PetscCall(PetscObjectSetName((PetscObject)(user->Xguess), "aux"));
1054     PetscCall(VecView(user->Xguess, viewer));
1055     PetscCall(StressField(da)); /* compute stress fields */
1056     PetscCall(PetscObjectSetName((PetscObject)(user->Xguess), "str"));
1057     PetscCall(VecView(user->Xguess, viewer));
1058     PetscCall(PetscObjectSetName((PetscObject)pars, "par"));
1059     PetscCall(VecView(pars, viewer));
1060 
1061     /* destroy viewer and vector */
1062     PetscCall(PetscViewerDestroy(&viewer));
1063     PetscCall(VecDestroy(&pars));
1064   }
1065 
1066   param->ivisc = ivt;
1067   return 0;
1068 }
1069 
1070 /* ------------------------------------------------------------------- */
1071 /* Compute both the second invariant of the strain rate tensor and the viscosity, at both cell centers and cell corners */
1072 PetscErrorCode ViscosityField(DM da, Vec X, Vec V)
1073 /* ------------------------------------------------------------------- */
1074 {
1075   AppCtx    *user;
1076   Parameter *param;
1077   GridInfo  *grid;
1078   Vec        localX;
1079   Field    **v, **x;
1080   PetscReal  eps, /* dx,*/ dz, T, epsC, TC;
1081   PetscInt   i, j, is, js, im, jm, ilim, jlim, ivt;
1082 
1083   PetscFunctionBeginUser;
1084   PetscCall(DMGetApplicationContext(da, &user));
1085   param        = user->param;
1086   grid         = user->grid;
1087   ivt          = param->ivisc;
1088   param->ivisc = param->output_ivisc;
1089 
1090   PetscCall(DMGetLocalVector(da, &localX));
1091   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
1092   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
1093   PetscCall(DMDAVecGetArray(da, localX, (void **)&x));
1094   PetscCall(DMDAVecGetArray(da, V, (void **)&v));
1095 
1096   /* Parameters */
1097   /* dx = grid->dx; */ dz = grid->dz;
1098 
1099   ilim = grid->ni - 1;
1100   jlim = grid->nj - 1;
1101 
1102   /* Compute real temperature, strain rate and viscosity */
1103   PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
1104   for (j = js; j < js + jm; j++) {
1105     for (i = is; i < is + im; i++) {
1106       T = PetscRealPart(param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * param->z_scale));
1107       if (i < ilim && j < jlim) {
1108         TC = PetscRealPart(param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * param->z_scale));
1109       } else {
1110         TC = T;
1111       }
1112       eps  = PetscRealPart((CalcSecInv(x, i, j, CELL_CENTER, user)));
1113       epsC = PetscRealPart(CalcSecInv(x, i, j, CELL_CORNER, user));
1114 
1115       v[j][i].u = eps;
1116       v[j][i].w = epsC;
1117       v[j][i].p = Viscosity(T, eps, dz * (j - 0.5), param);
1118       v[j][i].T = Viscosity(TC, epsC, dz * j, param);
1119     }
1120   }
1121   PetscCall(DMDAVecRestoreArray(da, V, (void **)&v));
1122   PetscCall(DMDAVecRestoreArray(da, localX, (void **)&x));
1123   PetscCall(DMRestoreLocalVector(da, &localX));
1124 
1125   param->ivisc = ivt;
1126   PetscFunctionReturn(0);
1127 }
1128 
1129 /* ------------------------------------------------------------------- */
1130 /* post-processing: compute stress everywhere */
1131 PetscErrorCode StressField(DM da)
1132 /* ------------------------------------------------------------------- */
1133 {
1134   AppCtx  *user;
1135   PetscInt i, j, is, js, im, jm;
1136   Vec      locVec;
1137   Field  **x, **y;
1138 
1139   PetscCall(DMGetApplicationContext(da, &user));
1140 
1141   /* Get the fine grid of Xguess and X */
1142   PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
1143   PetscCall(DMDAVecGetArray(da, user->Xguess, (void **)&x));
1144 
1145   PetscCall(DMGetLocalVector(da, &locVec));
1146   PetscCall(DMGlobalToLocalBegin(da, user->x, INSERT_VALUES, locVec));
1147   PetscCall(DMGlobalToLocalEnd(da, user->x, INSERT_VALUES, locVec));
1148   PetscCall(DMDAVecGetArray(da, locVec, (void **)&y));
1149 
1150   /* Compute stress on the corner points */
1151   for (j = js; j < js + jm; j++) {
1152     for (i = is; i < is + im; i++) {
1153       x[j][i].u = ShearStress(y, i, j, CELL_CENTER, user);
1154       x[j][i].w = ShearStress(y, i, j, CELL_CORNER, user);
1155       x[j][i].p = XNormalStress(y, i, j, CELL_CENTER, user);
1156       x[j][i].T = ZNormalStress(y, i, j, CELL_CENTER, user);
1157     }
1158   }
1159 
1160   /* Restore the fine grid of Xguess and X */
1161   PetscCall(DMDAVecRestoreArray(da, user->Xguess, (void **)&x));
1162   PetscCall(DMDAVecRestoreArray(da, locVec, (void **)&y));
1163   PetscCall(DMRestoreLocalVector(da, &locVec));
1164   return 0;
1165 }
1166 
1167 /*=====================================================================
1168   UTILITY FUNCTIONS
1169   =====================================================================*/
1170 
1171 /* returns the velocity of the subducting slab and handles fault nodes for BC */
1172 static inline PetscScalar SlabVel(char c, PetscInt i, PetscInt j, AppCtx *user) {
1173   Parameter *param = user->param;
1174   GridInfo  *grid  = user->grid;
1175 
1176   if (c == 'U' || c == 'u') {
1177     if (i < j - 1) return param->cb;
1178     else if (j <= grid->jfault) return 0.0;
1179     else return param->cb;
1180 
1181   } else {
1182     if (i < j) return param->sb;
1183     else if (j <= grid->jfault) return 0.0;
1184     else return param->sb;
1185   }
1186 }
1187 
1188 /*  solution to diffusive half-space cooling model for BC */
1189 static inline PetscScalar PlateModel(PetscInt j, PetscInt plate, AppCtx *user) {
1190   Parameter  *param = user->param;
1191   PetscScalar z;
1192   if (plate == PLATE_LID) z = (j - 0.5) * user->grid->dz;
1193   else z = (j - 0.5) * user->grid->dz * param->cb; /* PLATE_SLAB */
1194 #if defined(PETSC_HAVE_ERF)
1195   return (PetscReal)(erf((double)PetscRealPart(z * param->L / 2.0 / param->skt)));
1196 #else
1197   (*PetscErrorPrintf)("erf() not available on this machine\n");
1198   MPI_Abort(PETSC_COMM_SELF, 1);
1199 #endif
1200 }
1201 
1202 /*=====================================================================
1203   INTERACTIVE SIGNAL HANDLING
1204   =====================================================================*/
1205 
1206 /* ------------------------------------------------------------------- */
1207 PetscErrorCode SNESConverged_Interactive(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *ctx)
1208 /* ------------------------------------------------------------------- */
1209 {
1210   AppCtx    *user  = (AppCtx *)ctx;
1211   Parameter *param = user->param;
1212   KSP        ksp;
1213 
1214   PetscFunctionBeginUser;
1215   if (param->interrupted) {
1216     param->interrupted = PETSC_FALSE;
1217     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: exiting SNES solve. \n"));
1218     *reason = SNES_CONVERGED_FNORM_ABS;
1219     PetscFunctionReturn(0);
1220   } else if (param->toggle_kspmon) {
1221     param->toggle_kspmon = PETSC_FALSE;
1222 
1223     PetscCall(SNESGetKSP(snes, &ksp));
1224 
1225     if (param->kspmon) {
1226       PetscCall(KSPMonitorCancel(ksp));
1227 
1228       param->kspmon = PETSC_FALSE;
1229       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: deactivating ksp singular value monitor. \n"));
1230     } else {
1231       PetscViewerAndFormat *vf;
1232       PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf));
1233       PetscCall(KSPMonitorSet(ksp, (PetscErrorCode(*)(KSP, PetscInt, PetscReal, void *))KSPMonitorSingularValue, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
1234 
1235       param->kspmon = PETSC_TRUE;
1236       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: activating ksp singular value monitor. \n"));
1237     }
1238   }
1239   PetscCall(SNESConvergedDefault(snes, it, xnorm, snorm, fnorm, reason, ctx));
1240   PetscFunctionReturn(0);
1241 }
1242 
1243 /* ------------------------------------------------------------------- */
1244 #include <signal.h>
1245 PetscErrorCode InteractiveHandler(int signum, void *ctx)
1246 /* ------------------------------------------------------------------- */
1247 {
1248   AppCtx    *user  = (AppCtx *)ctx;
1249   Parameter *param = user->param;
1250 
1251   if (signum == SIGILL) {
1252     param->toggle_kspmon = PETSC_TRUE;
1253 #if !defined(PETSC_MISSING_SIGCONT)
1254   } else if (signum == SIGCONT) {
1255     param->interrupted = PETSC_TRUE;
1256 #endif
1257 #if !defined(PETSC_MISSING_SIGURG)
1258   } else if (signum == SIGURG) {
1259     param->stop_solve = PETSC_TRUE;
1260 #endif
1261   }
1262   return 0;
1263 }
1264 
1265 /*  main call-back function that computes the processor-local piece of the residual */
1266 PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field **x, Field **f, void *ptr) {
1267   AppCtx     *user  = (AppCtx *)ptr;
1268   Parameter  *param = user->param;
1269   GridInfo   *grid  = user->grid;
1270   PetscScalar mag_w, mag_u;
1271   PetscInt    i, j, mx, mz, ilim, jlim;
1272   PetscInt    is, ie, js, je, ibound; /* ,ivisc */
1273 
1274   PetscFunctionBeginUser;
1275   /* Define global and local grid parameters */
1276   mx   = info->mx;
1277   mz   = info->my;
1278   ilim = mx - 1;
1279   jlim = mz - 1;
1280   is   = info->xs;
1281   ie   = info->xs + info->xm;
1282   js   = info->ys;
1283   je   = info->ys + info->ym;
1284 
1285   /* Define geometric and numeric parameters */
1286   /* ivisc = param->ivisc; */ ibound = param->ibound;
1287 
1288   for (j = js; j < je; j++) {
1289     for (i = is; i < ie; i++) {
1290       /************* X-MOMENTUM/VELOCITY *************/
1291       if (i < j) f[j][i].u = x[j][i].u - SlabVel('U', i, j, user);
1292       else if (j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1293         /* in the lithospheric lid */
1294         f[j][i].u = x[j][i].u - 0.0;
1295       } else if (i == ilim) {
1296         /* on the right side boundary */
1297         if (ibound == BC_ANALYTIC) {
1298           f[j][i].u = x[j][i].u - HorizVelocity(i, j, user);
1299         } else {
1300           f[j][i].u = XNormalStress(x, i, j, CELL_CENTER, user) - EPS_ZERO;
1301         }
1302 
1303       } else if (j == jlim) {
1304         /* on the bottom boundary */
1305         if (ibound == BC_ANALYTIC) {
1306           f[j][i].u = x[j][i].u - HorizVelocity(i, j, user);
1307         } else if (ibound == BC_NOSTRESS) {
1308           f[j][i].u = XMomentumResidual(x, i, j, user);
1309         } else {
1310           /* experimental boundary condition */
1311         }
1312 
1313       } else {
1314         /* in the mantle wedge */
1315         f[j][i].u = XMomentumResidual(x, i, j, user);
1316       }
1317 
1318       /************* Z-MOMENTUM/VELOCITY *************/
1319       if (i <= j) {
1320         f[j][i].w = x[j][i].w - SlabVel('W', i, j, user);
1321 
1322       } else if (j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1323         /* in the lithospheric lid */
1324         f[j][i].w = x[j][i].w - 0.0;
1325 
1326       } else if (j == jlim) {
1327         /* on the bottom boundary */
1328         if (ibound == BC_ANALYTIC) {
1329           f[j][i].w = x[j][i].w - VertVelocity(i, j, user);
1330         } else {
1331           f[j][i].w = ZNormalStress(x, i, j, CELL_CENTER, user) - EPS_ZERO;
1332         }
1333 
1334       } else if (i == ilim) {
1335         /* on the right side boundary */
1336         if (ibound == BC_ANALYTIC) {
1337           f[j][i].w = x[j][i].w - VertVelocity(i, j, user);
1338         } else if (ibound == BC_NOSTRESS) {
1339           f[j][i].w = ZMomentumResidual(x, i, j, user);
1340         } else {
1341           /* experimental boundary condition */
1342         }
1343 
1344       } else {
1345         /* in the mantle wedge */
1346         f[j][i].w = ZMomentumResidual(x, i, j, user);
1347       }
1348 
1349       /************* CONTINUITY/PRESSURE *************/
1350       if (i < j || j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1351         /* in the lid or slab */
1352         f[j][i].p = x[j][i].p;
1353 
1354       } else if ((i == ilim || j == jlim) && ibound == BC_ANALYTIC) {
1355         /* on an analytic boundary */
1356         f[j][i].p = x[j][i].p - Pressure(i, j, user);
1357 
1358       } else {
1359         /* in the mantle wedge */
1360         f[j][i].p = ContinuityResidual(x, i, j, user);
1361       }
1362 
1363       /************* TEMPERATURE *************/
1364       if (j == 0) {
1365         /* on the surface */
1366         f[j][i].T = x[j][i].T + x[j + 1][i].T + PetscMax(PetscRealPart(x[j][i].T), 0.0);
1367 
1368       } else if (i == 0) {
1369         /* slab inflow boundary */
1370         f[j][i].T = x[j][i].T - PlateModel(j, PLATE_SLAB, user);
1371 
1372       } else if (i == ilim) {
1373         /* right side boundary */
1374         mag_u     = 1.0 - PetscPowRealInt((1.0 - PetscMax(PetscMin(PetscRealPart(x[j][i - 1].u) / param->cb, 1.0), 0.0)), 5);
1375         f[j][i].T = x[j][i].T - mag_u * x[j - 1][i - 1].T - (1.0 - mag_u) * PlateModel(j, PLATE_LID, user);
1376 
1377       } else if (j == jlim) {
1378         /* bottom boundary */
1379         mag_w     = 1.0 - PetscPowRealInt((1.0 - PetscMax(PetscMin(PetscRealPart(x[j - 1][i].w) / param->sb, 1.0), 0.0)), 5);
1380         f[j][i].T = x[j][i].T - mag_w * x[j - 1][i - 1].T - (1.0 - mag_w);
1381 
1382       } else {
1383         /* in the mantle wedge */
1384         f[j][i].T = EnergyResidual(x, i, j, user);
1385       }
1386     }
1387   }
1388   PetscFunctionReturn(0);
1389 }
1390 
1391 /*TEST
1392 
1393    build:
1394       requires: !complex erf
1395 
1396    test:
1397       args: -ni 18
1398       filter: grep -v Destination
1399       requires: !single
1400 
1401 TEST*/
1402