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