xref: /petsc/src/snes/tutorials/ex30.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static const char help[] = "Steady-state 2D subduction flow, pressure and temperature solver.\n\
2c4762a1bSJed Brown        The flow is driven by the subducting slab.\n\
3c4762a1bSJed Brown ---------------------------------ex30 help---------------------------------\n\
4c4762a1bSJed Brown   -OPTION <DEFAULT> = (UNITS) DESCRIPTION.\n\n\
5c4762a1bSJed Brown   -width <320> = (km) width of domain.\n\
6c4762a1bSJed Brown   -depth <300> = (km) depth of domain.\n\
7c4762a1bSJed Brown   -slab_dip <45> = (degrees) dip angle of the slab (determines the grid aspect ratio).\n\
8c4762a1bSJed Brown   -lid_depth <35> = (km) depth of the static conductive lid.\n\
9c4762a1bSJed Brown   -fault_depth <35> = (km) depth of slab-wedge mechanical coupling\n\
10c4762a1bSJed Brown      (fault dept >= lid depth).\n\
11c4762a1bSJed Brown \n\
12c4762a1bSJed Brown   -ni <82> = grid cells in x-direction. (nj adjusts to accommodate\n\
13c4762a1bSJed Brown       the slab dip & depth). DO NOT USE -da_grid_x option!!!\n\
14c4762a1bSJed Brown   -ivisc <3> = rheology option.\n\
15c4762a1bSJed Brown       0 --- constant viscosity.\n\
16c4762a1bSJed Brown       1 --- olivine diffusion creep rheology (T&P-dependent, newtonian).\n\
17c4762a1bSJed Brown       2 --- olivine dislocation creep rheology (T&P-dependent, non-newtonian).\n\
18c4762a1bSJed Brown       3 --- Full mantle rheology, combination of 1 & 2.\n\
19c4762a1bSJed Brown \n\
20c4762a1bSJed Brown   -slab_velocity <5> = (cm/year) convergence rate of slab into subduction zone.\n\
21c4762a1bSJed Brown   -slab_age <50> = (million yrs) age of slab for thermal profile boundary condition.\n\
22c4762a1bSJed Brown   -lid_age <50> = (million yrs) age of lid for thermal profile boundary condition.\n\
23c4762a1bSJed Brown \n\
24c4762a1bSJed Brown   FOR OTHER PARAMETER OPTIONS AND THEIR DEFAULT VALUES, see SetParams() in ex30.c.\n\
25c4762a1bSJed Brown ---------------------------------ex30 help---------------------------------\n";
26c4762a1bSJed Brown 
27c4762a1bSJed Brown /*F-----------------------------------------------------------------------
28c4762a1bSJed Brown 
29c4762a1bSJed Brown     This PETSc 2.2.0 example by Richard F. Katz
30c4762a1bSJed Brown     http://www.ldeo.columbia.edu/~katz/
31c4762a1bSJed Brown 
32c4762a1bSJed Brown     The problem is modeled by the partial differential equation system
33c4762a1bSJed Brown 
34c4762a1bSJed Brown \begin{eqnarray}
35c4762a1bSJed Brown          -\nabla P + \nabla \cdot [\eta (\nabla v + \nabla v^T)] & = & 0  \\
36c4762a1bSJed Brown                                            \nabla \cdot v & = & 0   \\
37c4762a1bSJed Brown                     dT/dt + \nabla \cdot (vT) - 1/Pe \triangle^2(T) & = & 0  \\
38c4762a1bSJed Brown \end{eqnarray}
39c4762a1bSJed Brown 
40c4762a1bSJed Brown  \begin{eqnarray}
41c4762a1bSJed Brown         \eta(T,Eps\_dot) &  = & \hbox{constant                        }    \hbox{if ivisc} ==0  \\
42c4762a1bSJed Brown                       &  = & \hbox{diffusion creep (T,P-dependent)    }     \hbox{if ivisc} ==1  \\
43c4762a1bSJed Brown                       &  = & \hbox{dislocation creep (T,P,v-dependent)}  \hbox{if ivisc} ==2  \\
44c4762a1bSJed Brown                       &  = & \hbox{mantle viscosity (difn and disl)   }  \hbox{if ivisc} ==3
45c4762a1bSJed Brown \end{eqnarray}
46c4762a1bSJed Brown 
47c4762a1bSJed Brown     which is uniformly discretized on a staggered mesh:
48c4762a1bSJed Brown                       -------$w_{ij}$------
49c4762a1bSJed Brown                   $u_{i-1j}$    $P,T_{ij}$   $u_{ij}$
50c4762a1bSJed Brown                       ------$w_{ij-1}$-----
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   ------------------------------------------------------------------------F*/
53c4762a1bSJed Brown 
54c4762a1bSJed Brown #include <petscsnes.h>
55c4762a1bSJed Brown #include <petscdm.h>
56c4762a1bSJed Brown #include <petscdmda.h>
57c4762a1bSJed Brown 
58c4762a1bSJed Brown #define VISC_CONST   0
59c4762a1bSJed Brown #define VISC_DIFN    1
60c4762a1bSJed Brown #define VISC_DISL    2
61c4762a1bSJed Brown #define VISC_FULL    3
62c4762a1bSJed Brown #define CELL_CENTER  0
63c4762a1bSJed Brown #define CELL_CORNER  1
64c4762a1bSJed Brown #define BC_ANALYTIC  0
65c4762a1bSJed Brown #define BC_NOSTRESS  1
66c4762a1bSJed Brown #define BC_EXPERMNT  2
67c4762a1bSJed Brown #define ADVECT_FV    0
68c4762a1bSJed Brown #define ADVECT_FROMM 1
69c4762a1bSJed Brown #define PLATE_SLAB   0
70c4762a1bSJed Brown #define PLATE_LID    1
71c4762a1bSJed Brown #define EPS_ZERO     0.00000001
72c4762a1bSJed Brown 
73c4762a1bSJed Brown typedef struct { /* holds the variables to be solved for */
74c4762a1bSJed Brown   PetscScalar u, w, p, T;
75c4762a1bSJed Brown } Field;
76c4762a1bSJed Brown 
77c4762a1bSJed Brown typedef struct { /* parameters needed to compute viscosity */
78c4762a1bSJed Brown   PetscReal A, n, Estar, Vstar;
79c4762a1bSJed Brown } ViscParam;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown typedef struct { /* physical and miscelaneous parameters */
82c4762a1bSJed Brown   PetscReal width, depth, scaled_width, scaled_depth, peclet, potentialT;
83c4762a1bSJed Brown   PetscReal slab_dip, slab_age, slab_velocity, kappa, z_scale;
84c4762a1bSJed Brown   PetscReal c, d, sb, cb, skt, visc_cutoff, lid_age, eta0, continuation;
85c4762a1bSJed Brown   PetscReal L, V, lid_depth, fault_depth;
86c4762a1bSJed Brown   ViscParam diffusion, dislocation;
87c4762a1bSJed Brown   PetscInt  ivisc, adv_scheme, ibound, output_ivisc;
88c4762a1bSJed Brown   PetscBool quiet, param_test, output_to_file, pv_analytic;
89c4762a1bSJed Brown   PetscBool interrupted, stop_solve, toggle_kspmon, kspmon;
90c4762a1bSJed Brown   char      filename[PETSC_MAX_PATH_LEN];
91c4762a1bSJed Brown } Parameter;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown typedef struct { /* grid parameters */
94c4762a1bSJed Brown   DMBoundaryType  bx, by;
95c4762a1bSJed Brown   DMDAStencilType stencil;
96c4762a1bSJed Brown   PetscInt        corner, ni, nj, jlid, jfault, inose;
97c4762a1bSJed Brown   PetscInt        dof, stencil_width, mglevels;
98c4762a1bSJed Brown   PetscReal       dx, dz;
99c4762a1bSJed Brown } GridInfo;
100c4762a1bSJed Brown 
101c4762a1bSJed Brown typedef struct { /* application context */
102c4762a1bSJed Brown   Vec        x, Xguess;
103c4762a1bSJed Brown   Parameter *param;
104c4762a1bSJed Brown   GridInfo  *grid;
105c4762a1bSJed Brown } AppCtx;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /* Callback functions (static interface) */
108c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(DMDALocalInfo *, Field **, Field **, void *);
109c4762a1bSJed Brown 
110c4762a1bSJed Brown /* Main routines */
111c4762a1bSJed Brown extern PetscErrorCode SetParams(Parameter *, GridInfo *);
112c4762a1bSJed Brown extern PetscErrorCode ReportParams(Parameter *, GridInfo *);
113c4762a1bSJed Brown extern PetscErrorCode Initialize(DM);
114c4762a1bSJed Brown extern PetscErrorCode UpdateSolution(SNES, AppCtx *, PetscInt *);
115c4762a1bSJed Brown extern PetscErrorCode DoOutput(SNES, PetscInt);
116c4762a1bSJed Brown 
117c4762a1bSJed Brown /* Post-processing & misc */
118c4762a1bSJed Brown extern PetscErrorCode ViscosityField(DM, Vec, Vec);
119c4762a1bSJed Brown extern PetscErrorCode StressField(DM);
120c4762a1bSJed Brown extern PetscErrorCode SNESConverged_Interactive(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
121c4762a1bSJed Brown extern PetscErrorCode InteractiveHandler(int, void *);
122c4762a1bSJed Brown 
123c4762a1bSJed Brown /*-----------------------------------------------------------------------*/
124c4762a1bSJed Brown int main(int argc, char **argv)
125c4762a1bSJed Brown /*-----------------------------------------------------------------------*/
126c4762a1bSJed Brown {
127c4762a1bSJed Brown   SNES      snes;
128c4762a1bSJed Brown   AppCtx   *user; /* user-defined work context */
129c4762a1bSJed Brown   Parameter param;
130c4762a1bSJed Brown   GridInfo  grid;
131c4762a1bSJed Brown   PetscInt  nits;
132c4762a1bSJed Brown   MPI_Comm  comm;
133c4762a1bSJed Brown   DM        da;
134c4762a1bSJed Brown 
135327415f7SBarry Smith   PetscFunctionBeginUser;
1369566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
137c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-file", "ex30_output");
138c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-snes_monitor_short", NULL);
139c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-snes_max_it", "20");
140c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-ksp_max_it", "1500");
141c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-ksp_gmres_restart", "300");
142c4762a1bSJed Brown   PetscOptionsInsert(NULL, &argc, &argv, NULL);
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147c4762a1bSJed Brown      Set up the problem parameters.
148c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1499566063dSJacob Faibussowitsch   PetscCall(SetParams(&param, &grid));
1509566063dSJacob Faibussowitsch   PetscCall(ReportParams(&param, &grid));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1549566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &snes));
1559566063dSJacob Faibussowitsch   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));
1569566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1579566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1589566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
1599566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "x-velocity"));
1609566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "y-velocity"));
1619566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 2, "pressure"));
1629566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 3, "temperature"));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165c4762a1bSJed Brown      Create user context, set problem data, create vector data structures.
166c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1679566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
168c4762a1bSJed Brown   user->param = &param;
169c4762a1bSJed Brown   user->grid  = &grid;
1709566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, user));
1719566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &(user->Xguess)));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174c4762a1bSJed Brown      Set up the SNES solver with callback functions.
175c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1769566063dSJacob Faibussowitsch   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, (void *)user));
1779566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch   PetscCall(SNESSetConvergenceTest(snes, SNESConverged_Interactive, (void *)user, NULL));
1809566063dSJacob Faibussowitsch   PetscCall(PetscPushSignalHandler(InteractiveHandler, (void *)user));
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183c4762a1bSJed Brown      Initialize and solve the nonlinear system
184c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1859566063dSJacob Faibussowitsch   PetscCall(Initialize(da));
1869566063dSJacob Faibussowitsch   PetscCall(UpdateSolution(snes, user, &nits));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189c4762a1bSJed Brown      Output variables.
190c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1919566063dSJacob Faibussowitsch   PetscCall(DoOutput(snes, nits));
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194c4762a1bSJed Brown      Free work space.
195c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Xguess));
1979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->x));
1989566063dSJacob Faibussowitsch   PetscCall(PetscFree(user));
1999566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2009566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
2019566063dSJacob Faibussowitsch   PetscCall(PetscPopSignalHandler());
2029566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
203b122ec5aSJacob Faibussowitsch   return 0;
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown /*=====================================================================
207c4762a1bSJed Brown   PETSc INTERACTION FUNCTIONS (initialize & call SNESSolve)
208c4762a1bSJed Brown   =====================================================================*/
209c4762a1bSJed Brown 
210c4762a1bSJed Brown /*---------------------------------------------------------------------*/
211c4762a1bSJed Brown /*  manages solve: adaptive continuation method  */
2129371c9d4SSatish Balay PetscErrorCode UpdateSolution(SNES snes, AppCtx *user, PetscInt *nits) {
213c4762a1bSJed Brown   KSP                 ksp;
214c4762a1bSJed Brown   PC                  pc;
215c4762a1bSJed Brown   SNESConvergedReason reason    = SNES_CONVERGED_ITERATING;
216c4762a1bSJed Brown   Parameter          *param     = user->param;
217c4762a1bSJed Brown   PetscReal           cont_incr = 0.3;
218c4762a1bSJed Brown   PetscInt            its;
219c4762a1bSJed Brown   PetscBool           q = PETSC_FALSE;
220c4762a1bSJed Brown   DM                  dm;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   PetscFunctionBeginUser;
2239566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
2249566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &user->x));
2259566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
2269566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
2279566063dSJacob Faibussowitsch   PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   *nits = 0;
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   /* Isoviscous solve */
232c4762a1bSJed Brown   if (param->ivisc == VISC_CONST && !param->stop_solve) {
233c4762a1bSJed Brown     param->ivisc = VISC_CONST;
234c4762a1bSJed Brown 
2359566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, 0, user->x));
2369566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes, &reason));
2379566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(snes, &its));
238c4762a1bSJed Brown     *nits += its;
2399566063dSJacob Faibussowitsch     PetscCall(VecCopy(user->x, user->Xguess));
240c4762a1bSJed Brown     if (param->stop_solve) goto done;
241c4762a1bSJed Brown   }
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* Olivine diffusion creep */
244c4762a1bSJed Brown   if (param->ivisc >= VISC_DIFN && !param->stop_solve) {
2459566063dSJacob Faibussowitsch     if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Computing Variable Viscosity Solution\n"));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown     /* continuation method on viscosity cutoff */
248c4762a1bSJed Brown     for (param->continuation = 0.0;; param->continuation += cont_incr) {
2499566063dSJacob Faibussowitsch       if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Continuation parameter = %g\n", (double)param->continuation));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown       /* solve the non-linear system */
2529566063dSJacob Faibussowitsch       PetscCall(VecCopy(user->Xguess, user->x));
2539566063dSJacob Faibussowitsch       PetscCall(SNESSolve(snes, 0, user->x));
2549566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes, &reason));
2559566063dSJacob Faibussowitsch       PetscCall(SNESGetIterationNumber(snes, &its));
256c4762a1bSJed Brown       *nits += its;
25763a3b9bcSJacob Faibussowitsch       if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SNES iterations: %" PetscInt_FMT ", Cumulative: %" PetscInt_FMT "\n", its, *nits));
258c4762a1bSJed Brown       if (param->stop_solve) goto done;
259c4762a1bSJed Brown 
260c4762a1bSJed Brown       if (reason < 0) {
261c4762a1bSJed Brown         /* NOT converged */
262c4762a1bSJed Brown         cont_incr = -PetscAbsReal(cont_incr) / 2.0;
263c4762a1bSJed Brown         if (PetscAbsReal(cont_incr) < 0.01) goto done;
264c4762a1bSJed Brown 
265c4762a1bSJed Brown       } else {
266c4762a1bSJed Brown         /* converged */
2679566063dSJacob Faibussowitsch         PetscCall(VecCopy(user->x, user->Xguess));
268c4762a1bSJed Brown         if (param->continuation >= 1.0) goto done;
269c4762a1bSJed Brown         if (its <= 3) cont_incr = 0.30001;
270c4762a1bSJed Brown         else if (its <= 8) cont_incr = 0.15001;
271c4762a1bSJed Brown         else cont_incr = 0.10001;
272c4762a1bSJed Brown 
273c4762a1bSJed Brown         if (param->continuation + cont_incr > 1.0) cont_incr = 1.0 - param->continuation;
274c4762a1bSJed Brown       } /* endif reason<0 */
275c4762a1bSJed Brown     }
276c4762a1bSJed Brown   }
277c4762a1bSJed Brown done:
2789566063dSJacob Faibussowitsch   if (param->stop_solve && !q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: stopping solve.\n"));
2799566063dSJacob Faibussowitsch   if (reason < 0 && !q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "FAILED TO CONVERGE: stopping solve.\n"));
280c4762a1bSJed Brown   PetscFunctionReturn(0);
281c4762a1bSJed Brown }
282c4762a1bSJed Brown 
283c4762a1bSJed Brown /*=====================================================================
284c4762a1bSJed Brown   PHYSICS FUNCTIONS (compute the discrete residual)
285c4762a1bSJed Brown   =====================================================================*/
286c4762a1bSJed Brown 
287c4762a1bSJed Brown /*---------------------------------------------------------------------*/
2889fbee547SJacob Faibussowitsch static inline PetscScalar UInterp(Field **x, PetscInt i, PetscInt j)
289c4762a1bSJed Brown /*---------------------------------------------------------------------*/
290c4762a1bSJed Brown {
291c4762a1bSJed Brown   return 0.25 * (x[j][i].u + x[j + 1][i].u + x[j][i + 1].u + x[j + 1][i + 1].u);
292c4762a1bSJed Brown }
293c4762a1bSJed Brown 
294c4762a1bSJed Brown /*---------------------------------------------------------------------*/
2959fbee547SJacob Faibussowitsch static inline PetscScalar WInterp(Field **x, PetscInt i, PetscInt j)
296c4762a1bSJed Brown /*---------------------------------------------------------------------*/
297c4762a1bSJed Brown {
298c4762a1bSJed Brown   return 0.25 * (x[j][i].w + x[j + 1][i].w + x[j][i + 1].w + x[j + 1][i + 1].w);
299c4762a1bSJed Brown }
300c4762a1bSJed Brown 
301c4762a1bSJed Brown /*---------------------------------------------------------------------*/
3029fbee547SJacob Faibussowitsch static inline PetscScalar PInterp(Field **x, PetscInt i, PetscInt j)
303c4762a1bSJed Brown /*---------------------------------------------------------------------*/
304c4762a1bSJed Brown {
305c4762a1bSJed Brown   return 0.25 * (x[j][i].p + x[j + 1][i].p + x[j][i + 1].p + x[j + 1][i + 1].p);
306c4762a1bSJed Brown }
307c4762a1bSJed Brown 
308c4762a1bSJed Brown /*---------------------------------------------------------------------*/
3099fbee547SJacob Faibussowitsch static inline PetscScalar TInterp(Field **x, PetscInt i, PetscInt j)
310c4762a1bSJed Brown /*---------------------------------------------------------------------*/
311c4762a1bSJed Brown {
312c4762a1bSJed Brown   return 0.25 * (x[j][i].T + x[j + 1][i].T + x[j][i + 1].T + x[j + 1][i + 1].T);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown /*---------------------------------------------------------------------*/
316c4762a1bSJed Brown /*  isoviscous analytic solution for IC */
3179fbee547SJacob Faibussowitsch static inline PetscScalar HorizVelocity(PetscInt i, PetscInt j, AppCtx *user)
318c4762a1bSJed Brown /*---------------------------------------------------------------------*/
319c4762a1bSJed Brown {
320c4762a1bSJed Brown   Parameter  *param = user->param;
321c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
322c4762a1bSJed Brown   PetscScalar st, ct, th, c = param->c, d = param->d;
323c4762a1bSJed Brown   PetscReal   x, z, r;
324c4762a1bSJed Brown 
3259371c9d4SSatish Balay   x  = (i - grid->jlid) * grid->dx;
3269371c9d4SSatish Balay   z  = (j - grid->jlid - 0.5) * grid->dz;
327c4762a1bSJed Brown   r  = PetscSqrtReal(x * x + z * z);
328c4762a1bSJed Brown   st = z / r;
329c4762a1bSJed Brown   ct = x / r;
330c4762a1bSJed Brown   th = PetscAtanReal(z / x);
331c4762a1bSJed Brown   return ct * (c * th * st + d * (st + th * ct)) + st * (c * (st - th * ct) + d * th * st);
332c4762a1bSJed Brown }
333c4762a1bSJed Brown 
334c4762a1bSJed Brown /*---------------------------------------------------------------------*/
335c4762a1bSJed Brown /*  isoviscous analytic solution for IC */
3369fbee547SJacob Faibussowitsch static inline PetscScalar VertVelocity(PetscInt i, PetscInt j, AppCtx *user)
337c4762a1bSJed Brown /*---------------------------------------------------------------------*/
338c4762a1bSJed Brown {
339c4762a1bSJed Brown   Parameter  *param = user->param;
340c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
341c4762a1bSJed Brown   PetscScalar st, ct, th, c = param->c, d = param->d;
342c4762a1bSJed Brown   PetscReal   x, z, r;
343c4762a1bSJed Brown 
3449371c9d4SSatish Balay   x  = (i - grid->jlid - 0.5) * grid->dx;
3459371c9d4SSatish Balay   z  = (j - grid->jlid) * grid->dz;
3469371c9d4SSatish Balay   r  = PetscSqrtReal(x * x + z * z);
3479371c9d4SSatish Balay   st = z / r;
3489371c9d4SSatish Balay   ct = x / r;
3499371c9d4SSatish Balay   th = PetscAtanReal(z / x);
350c4762a1bSJed Brown   return st * (c * th * st + d * (st + th * ct)) - ct * (c * (st - th * ct) + d * th * st);
351c4762a1bSJed Brown }
352c4762a1bSJed Brown 
353c4762a1bSJed Brown /*---------------------------------------------------------------------*/
354c4762a1bSJed Brown /*  isoviscous analytic solution for IC */
3559fbee547SJacob Faibussowitsch static inline PetscScalar Pressure(PetscInt i, PetscInt j, AppCtx *user)
356c4762a1bSJed Brown /*---------------------------------------------------------------------*/
357c4762a1bSJed Brown {
358c4762a1bSJed Brown   Parameter  *param = user->param;
359c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
360c4762a1bSJed Brown   PetscScalar x, z, r, st, ct, c = param->c, d = param->d;
361c4762a1bSJed Brown 
3629371c9d4SSatish Balay   x  = (i - grid->jlid - 0.5) * grid->dx;
3639371c9d4SSatish Balay   z  = (j - grid->jlid - 0.5) * grid->dz;
3649371c9d4SSatish Balay   r  = PetscSqrtReal(x * x + z * z);
3659371c9d4SSatish Balay   st = z / r;
3669371c9d4SSatish Balay   ct = x / r;
367c4762a1bSJed Brown   return (-2.0 * (c * ct - d * st) / r);
368c4762a1bSJed Brown }
369c4762a1bSJed Brown 
370c4762a1bSJed Brown /*  computes the second invariant of the strain rate tensor */
3719fbee547SJacob Faibussowitsch static inline PetscScalar CalcSecInv(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
372c4762a1bSJed Brown /*---------------------------------------------------------------------*/
373c4762a1bSJed Brown {
374c4762a1bSJed Brown   Parameter  *param = user->param;
375c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
376c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1;
377c4762a1bSJed Brown   PetscScalar uN, uS, uE, uW, wN, wS, wE, wW;
378c4762a1bSJed Brown   PetscScalar eps11, eps12, eps22;
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   if (i < j) return EPS_ZERO;
381c4762a1bSJed Brown   if (i == ilim) i--;
382c4762a1bSJed Brown   if (j == jlim) j--;
383c4762a1bSJed Brown 
384c4762a1bSJed Brown   if (ipos == CELL_CENTER) { /* on cell center */
385c4762a1bSJed Brown     if (j <= grid->jlid) return EPS_ZERO;
386c4762a1bSJed Brown 
3879371c9d4SSatish Balay     uE = x[j][i].u;
3889371c9d4SSatish Balay     uW = x[j][i - 1].u;
3899371c9d4SSatish Balay     wN = x[j][i].w;
3909371c9d4SSatish Balay     wS = x[j - 1][i].w;
391c4762a1bSJed Brown     wE = WInterp(x, i, j - 1);
392c4762a1bSJed Brown     if (i == j) {
3939371c9d4SSatish Balay       uN = param->cb;
3949371c9d4SSatish Balay       wW = param->sb;
395c4762a1bSJed Brown     } else {
3969371c9d4SSatish Balay       uN = UInterp(x, i - 1, j);
3979371c9d4SSatish Balay       wW = WInterp(x, i - 1, j - 1);
398c4762a1bSJed Brown     }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown     if (j == grid->jlid + 1) uS = 0.0;
401c4762a1bSJed Brown     else uS = UInterp(x, i - 1, j - 1);
402c4762a1bSJed Brown 
403c4762a1bSJed Brown   } else { /* on CELL_CORNER */
404c4762a1bSJed Brown     if (j < grid->jlid) return EPS_ZERO;
405c4762a1bSJed Brown 
4069371c9d4SSatish Balay     uN = x[j + 1][i].u;
4079371c9d4SSatish Balay     uS = x[j][i].u;
4089371c9d4SSatish Balay     wE = x[j][i + 1].w;
4099371c9d4SSatish Balay     wW = x[j][i].w;
410c4762a1bSJed Brown     if (i == j) {
411c4762a1bSJed Brown       wN = param->sb;
412c4762a1bSJed Brown       uW = param->cb;
413c4762a1bSJed Brown     } else {
414c4762a1bSJed Brown       wN = WInterp(x, i, j);
415c4762a1bSJed Brown       uW = UInterp(x, i - 1, j);
416c4762a1bSJed Brown     }
417c4762a1bSJed Brown 
418c4762a1bSJed Brown     if (j == grid->jlid) {
4199371c9d4SSatish Balay       uE = 0.0;
4209371c9d4SSatish Balay       uW = 0.0;
421c4762a1bSJed Brown       uS = -uN;
422c4762a1bSJed Brown       wS = -wN;
423c4762a1bSJed Brown     } else {
424c4762a1bSJed Brown       uE = UInterp(x, i, j);
425c4762a1bSJed Brown       wS = WInterp(x, i, j - 1);
426c4762a1bSJed Brown     }
427c4762a1bSJed Brown   }
428c4762a1bSJed Brown 
4299371c9d4SSatish Balay   eps11 = (uE - uW) / grid->dx;
4309371c9d4SSatish Balay   eps22 = (wN - wS) / grid->dz;
431c4762a1bSJed Brown   eps12 = 0.5 * ((uN - uS) / grid->dz + (wE - wW) / grid->dx);
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   return PetscSqrtReal(0.5 * (eps11 * eps11 + 2.0 * eps12 * eps12 + eps22 * eps22));
434c4762a1bSJed Brown }
435c4762a1bSJed Brown 
436c4762a1bSJed Brown /*---------------------------------------------------------------------*/
437c4762a1bSJed Brown /*  computes the shear viscosity */
4389fbee547SJacob Faibussowitsch static inline PetscScalar Viscosity(PetscScalar T, PetscScalar eps, PetscScalar z, Parameter *param)
439c4762a1bSJed Brown /*---------------------------------------------------------------------*/
440c4762a1bSJed Brown {
441c4762a1bSJed Brown   PetscReal   result = 0.0;
442c4762a1bSJed Brown   ViscParam   difn = param->diffusion, disl = param->dislocation;
443c4762a1bSJed Brown   PetscInt    iVisc     = param->ivisc;
444c4762a1bSJed Brown   PetscScalar eps_scale = param->V / (param->L * 1000.0);
445c4762a1bSJed Brown   PetscScalar strain_power, v1, v2, P;
446c4762a1bSJed Brown   PetscScalar rho_g = 32340.0, R = 8.3144;
447c4762a1bSJed Brown 
448c4762a1bSJed Brown   P = rho_g * (z * param->L * 1000.0); /* Pa */
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   if (iVisc == VISC_CONST) {
451c4762a1bSJed Brown     /* constant viscosity */
452c4762a1bSJed Brown     return 1.0;
453c4762a1bSJed Brown   } else if (iVisc == VISC_DIFN) {
454c4762a1bSJed Brown     /* diffusion creep rheology */
455c4762a1bSJed Brown     result = PetscRealPart((difn.A * PetscExpScalar((difn.Estar + P * difn.Vstar) / R / (T + 273.0)) / param->eta0));
456c4762a1bSJed Brown   } else if (iVisc == VISC_DISL) {
457c4762a1bSJed Brown     /* dislocation creep rheology */
458c4762a1bSJed Brown     strain_power = PetscPowScalar(eps * eps_scale, (1.0 - disl.n) / disl.n);
459c4762a1bSJed Brown 
460c4762a1bSJed Brown     result = PetscRealPart(disl.A * PetscExpScalar((disl.Estar + P * disl.Vstar) / disl.n / R / (T + 273.0)) * strain_power / param->eta0);
461c4762a1bSJed Brown   } else if (iVisc == VISC_FULL) {
462c4762a1bSJed Brown     /* dislocation/diffusion creep rheology */
463c4762a1bSJed Brown     strain_power = PetscPowScalar(eps * eps_scale, (1.0 - disl.n) / disl.n);
464c4762a1bSJed Brown 
465c4762a1bSJed Brown     v1 = difn.A * PetscExpScalar((difn.Estar + P * difn.Vstar) / R / (T + 273.0)) / param->eta0;
466c4762a1bSJed Brown     v2 = disl.A * PetscExpScalar((disl.Estar + P * disl.Vstar) / disl.n / R / (T + 273.0)) * strain_power / param->eta0;
467c4762a1bSJed Brown 
468c4762a1bSJed Brown     result = PetscRealPart(1.0 / (1.0 / v1 + 1.0 / v2));
469c4762a1bSJed Brown   }
470c4762a1bSJed Brown 
471c4762a1bSJed Brown   /* max viscosity is param->eta0 */
472c4762a1bSJed Brown   result = PetscMin(result, 1.0);
473c4762a1bSJed Brown   /* min viscosity is param->visc_cutoff */
474c4762a1bSJed Brown   result = PetscMax(result, param->visc_cutoff);
475c4762a1bSJed Brown   /* continuation method */
476c4762a1bSJed Brown   result = PetscPowReal(result, param->continuation);
477c4762a1bSJed Brown   return result;
478c4762a1bSJed Brown }
479c4762a1bSJed Brown 
480c4762a1bSJed Brown /*---------------------------------------------------------------------*/
481c4762a1bSJed Brown /*  computes the residual of the x-component of eqn (1) above */
4829fbee547SJacob Faibussowitsch static inline PetscScalar XMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
483c4762a1bSJed Brown /*---------------------------------------------------------------------*/
484c4762a1bSJed Brown {
485c4762a1bSJed Brown   Parameter  *param = user->param;
486c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
487c4762a1bSJed Brown   PetscScalar dx = grid->dx, dz = grid->dz;
488c4762a1bSJed Brown   PetscScalar etaN, etaS, etaE, etaW, epsN = 0.0, epsS = 0.0, epsE = 0.0, epsW = 0.0;
489c4762a1bSJed Brown   PetscScalar TE = 0.0, TN = 0.0, TS = 0.0, TW = 0.0, dPdx, residual, z_scale;
490c4762a1bSJed Brown   PetscScalar dudxW, dudxE, dudzN, dudzS, dwdxN, dwdxS;
491c4762a1bSJed Brown   PetscInt    jlim = grid->nj - 1;
492c4762a1bSJed Brown 
493c4762a1bSJed Brown   z_scale = param->z_scale;
494c4762a1bSJed Brown 
495c4762a1bSJed Brown   if (param->ivisc == VISC_DIFN || param->ivisc >= VISC_DISL) { /* viscosity is T-dependent */
496c4762a1bSJed Brown     TS = param->potentialT * TInterp(x, i, j - 1) * PetscExpScalar((j - 1.0) * dz * z_scale);
497c4762a1bSJed Brown     if (j == jlim) TN = TS;
498c4762a1bSJed Brown     else TN = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
499c4762a1bSJed Brown     TW = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
500c4762a1bSJed Brown     TE = param->potentialT * x[j][i + 1].T * PetscExpScalar((j - 0.5) * dz * z_scale);
501c4762a1bSJed Brown     if (param->ivisc >= VISC_DISL) { /* olivine dislocation creep */
502c4762a1bSJed Brown       epsN = CalcSecInv(x, i, j, CELL_CORNER, user);
503c4762a1bSJed Brown       epsS = CalcSecInv(x, i, j - 1, CELL_CORNER, user);
504c4762a1bSJed Brown       epsE = CalcSecInv(x, i + 1, j, CELL_CENTER, user);
505c4762a1bSJed Brown       epsW = CalcSecInv(x, i, j, CELL_CENTER, user);
506c4762a1bSJed Brown     }
507c4762a1bSJed Brown   }
508c4762a1bSJed Brown   etaN = Viscosity(TN, epsN, dz * (j + 0.5), param);
509c4762a1bSJed Brown   etaS = Viscosity(TS, epsS, dz * (j - 0.5), param);
510c4762a1bSJed Brown   etaW = Viscosity(TW, epsW, dz * j, param);
511c4762a1bSJed Brown   etaE = Viscosity(TE, epsE, dz * j, param);
512c4762a1bSJed Brown 
513c4762a1bSJed Brown   dPdx = (x[j][i + 1].p - x[j][i].p) / dx;
514c4762a1bSJed Brown   if (j == jlim) dudzN = etaN * (x[j][i].w - x[j][i + 1].w) / dx;
515c4762a1bSJed Brown   else dudzN = etaN * (x[j + 1][i].u - x[j][i].u) / dz;
516c4762a1bSJed Brown   dudzS = etaS * (x[j][i].u - x[j - 1][i].u) / dz;
517c4762a1bSJed Brown   dudxE = etaE * (x[j][i + 1].u - x[j][i].u) / dx;
518c4762a1bSJed Brown   dudxW = etaW * (x[j][i].u - x[j][i - 1].u) / dx;
519c4762a1bSJed Brown 
520c4762a1bSJed Brown   residual = -dPdx /* X-MOMENTUM EQUATION*/
5219371c9d4SSatish Balay            + (dudxE - dudxW) / dx + (dudzN - dudzS) / dz;
522c4762a1bSJed Brown 
523c4762a1bSJed Brown   if (param->ivisc != VISC_CONST) {
524c4762a1bSJed Brown     dwdxN = etaN * (x[j][i + 1].w - x[j][i].w) / dx;
525c4762a1bSJed Brown     dwdxS = etaS * (x[j - 1][i + 1].w - x[j - 1][i].w) / dx;
526c4762a1bSJed Brown 
527c4762a1bSJed Brown     residual += (dudxE - dudxW) / dx + (dwdxN - dwdxS) / dz;
528c4762a1bSJed Brown   }
529c4762a1bSJed Brown 
530c4762a1bSJed Brown   return residual;
531c4762a1bSJed Brown }
532c4762a1bSJed Brown 
533c4762a1bSJed Brown /*---------------------------------------------------------------------*/
534c4762a1bSJed Brown /*  computes the residual of the z-component of eqn (1) above */
5359fbee547SJacob Faibussowitsch static inline PetscScalar ZMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
536c4762a1bSJed Brown /*---------------------------------------------------------------------*/
537c4762a1bSJed Brown {
538c4762a1bSJed Brown   Parameter  *param = user->param;
539c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
540c4762a1bSJed Brown   PetscScalar dx = grid->dx, dz = grid->dz;
541c4762a1bSJed Brown   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;
542c4762a1bSJed Brown   PetscScalar TE = 0.0, TN = 0.0, TS = 0.0, TW = 0.0, dPdz, residual, z_scale;
543c4762a1bSJed Brown   PetscScalar dudzE, dudzW, dwdxW, dwdxE, dwdzN, dwdzS;
544c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1;
545c4762a1bSJed Brown 
546c4762a1bSJed Brown   /* geometric and other parameters */
547c4762a1bSJed Brown   z_scale = param->z_scale;
548c4762a1bSJed Brown 
549c4762a1bSJed Brown   /* viscosity */
550c4762a1bSJed Brown   if (param->ivisc == VISC_DIFN || param->ivisc >= VISC_DISL) { /* viscosity is T-dependent */
551c4762a1bSJed Brown     TN = param->potentialT * x[j + 1][i].T * PetscExpScalar((j + 0.5) * dz * z_scale);
552c4762a1bSJed Brown     TS = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
553c4762a1bSJed Brown     TW = param->potentialT * TInterp(x, i - 1, j) * PetscExpScalar(j * dz * z_scale);
554c4762a1bSJed Brown     if (i == ilim) TE = TW;
555c4762a1bSJed Brown     else TE = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
556c4762a1bSJed Brown     if (param->ivisc >= VISC_DISL) { /* olivine dislocation creep */
557c4762a1bSJed Brown       epsN = CalcSecInv(x, i, j + 1, CELL_CENTER, user);
558c4762a1bSJed Brown       epsS = CalcSecInv(x, i, j, CELL_CENTER, user);
559c4762a1bSJed Brown       epsE = CalcSecInv(x, i, j, CELL_CORNER, user);
560c4762a1bSJed Brown       epsW = CalcSecInv(x, i - 1, j, CELL_CORNER, user);
561c4762a1bSJed Brown     }
562c4762a1bSJed Brown   }
563c4762a1bSJed Brown   etaN = Viscosity(TN, epsN, dz * (j + 1.0), param);
564c4762a1bSJed Brown   etaS = Viscosity(TS, epsS, dz * (j + 0.0), param);
565c4762a1bSJed Brown   etaW = Viscosity(TW, epsW, dz * (j + 0.5), param);
566c4762a1bSJed Brown   etaE = Viscosity(TE, epsE, dz * (j + 0.5), param);
567c4762a1bSJed Brown 
568c4762a1bSJed Brown   dPdz  = (x[j + 1][i].p - x[j][i].p) / dz;
569c4762a1bSJed Brown   dwdzN = etaN * (x[j + 1][i].w - x[j][i].w) / dz;
570c4762a1bSJed Brown   dwdzS = etaS * (x[j][i].w - x[j - 1][i].w) / dz;
571c4762a1bSJed Brown   if (i == ilim) dwdxE = etaE * (x[j][i].u - x[j + 1][i].u) / dz;
572c4762a1bSJed Brown   else dwdxE = etaE * (x[j][i + 1].w - x[j][i].w) / dx;
573c4762a1bSJed Brown   dwdxW = 2.0 * etaW * (x[j][i].w - x[j][i - 1].w) / dx;
574c4762a1bSJed Brown 
575c4762a1bSJed Brown   /* Z-MOMENTUM */
576c4762a1bSJed Brown   residual = -dPdz /* constant viscosity terms */
5779371c9d4SSatish Balay            + (dwdzN - dwdzS) / dz + (dwdxE - dwdxW) / dx;
578c4762a1bSJed Brown 
579c4762a1bSJed Brown   if (param->ivisc != VISC_CONST) {
580c4762a1bSJed Brown     dudzE = etaE * (x[j + 1][i].u - x[j][i].u) / dz;
581c4762a1bSJed Brown     dudzW = etaW * (x[j + 1][i - 1].u - x[j][i - 1].u) / dz;
582c4762a1bSJed Brown 
583c4762a1bSJed Brown     residual += (dwdzN - dwdzS) / dz + (dudzE - dudzW) / dx;
584c4762a1bSJed Brown   }
585c4762a1bSJed Brown 
586c4762a1bSJed Brown   return residual;
587c4762a1bSJed Brown }
588c4762a1bSJed Brown 
589c4762a1bSJed Brown /*---------------------------------------------------------------------*/
590c4762a1bSJed Brown /*  computes the residual of eqn (2) above */
5919fbee547SJacob Faibussowitsch static inline PetscScalar ContinuityResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
592c4762a1bSJed Brown /*---------------------------------------------------------------------*/
593c4762a1bSJed Brown {
594c4762a1bSJed Brown   GridInfo   *grid = user->grid;
595c4762a1bSJed Brown   PetscScalar uE, uW, wN, wS, dudx, dwdz;
596c4762a1bSJed Brown 
5979371c9d4SSatish Balay   uW   = x[j][i - 1].u;
5989371c9d4SSatish Balay   uE   = x[j][i].u;
5999371c9d4SSatish Balay   dudx = (uE - uW) / grid->dx;
6009371c9d4SSatish Balay   wS   = x[j - 1][i].w;
6019371c9d4SSatish Balay   wN   = x[j][i].w;
6029371c9d4SSatish Balay   dwdz = (wN - wS) / grid->dz;
603c4762a1bSJed Brown 
604c4762a1bSJed Brown   return dudx + dwdz;
605c4762a1bSJed Brown }
606c4762a1bSJed Brown 
607c4762a1bSJed Brown /*---------------------------------------------------------------------*/
608c4762a1bSJed Brown /*  computes the residual of eqn (3) above */
6099fbee547SJacob Faibussowitsch static inline PetscScalar EnergyResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
610c4762a1bSJed Brown /*---------------------------------------------------------------------*/
611c4762a1bSJed Brown {
612c4762a1bSJed Brown   Parameter  *param = user->param;
613c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
614c4762a1bSJed Brown   PetscScalar dx = grid->dx, dz = grid->dz;
615c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1, jlid = grid->jlid;
616c4762a1bSJed Brown   PetscScalar TE, TN, TS, TW, residual;
617c4762a1bSJed Brown   PetscScalar uE, uW, wN, wS;
618c4762a1bSJed Brown   PetscScalar fN, fS, fE, fW, dTdxW, dTdxE, dTdzN, dTdzS;
619c4762a1bSJed Brown 
620c4762a1bSJed Brown   dTdzN = (x[j + 1][i].T - x[j][i].T) / dz;
621c4762a1bSJed Brown   dTdzS = (x[j][i].T - x[j - 1][i].T) / dz;
622c4762a1bSJed Brown   dTdxE = (x[j][i + 1].T - x[j][i].T) / dx;
623c4762a1bSJed Brown   dTdxW = (x[j][i].T - x[j][i - 1].T) / dx;
624c4762a1bSJed Brown 
625c4762a1bSJed Brown   residual = ((dTdzN - dTdzS) / dz + /* diffusion term */
6269371c9d4SSatish Balay               (dTdxE - dTdxW) / dx) *
6279371c9d4SSatish Balay              dx * dz / param->peclet;
628c4762a1bSJed Brown 
629c4762a1bSJed Brown   if (j <= jlid && i >= j) {
630c4762a1bSJed Brown     /* don't advect in the lid */
631c4762a1bSJed Brown     return residual;
632c4762a1bSJed Brown   } else if (i < j) {
633c4762a1bSJed Brown     /* beneath the slab sfc */
634c4762a1bSJed Brown     uW = uE = param->cb;
635c4762a1bSJed Brown     wS = wN = param->sb;
636c4762a1bSJed Brown   } else {
637c4762a1bSJed Brown     /* advect in the slab and wedge */
6389371c9d4SSatish Balay     uW = x[j][i - 1].u;
6399371c9d4SSatish Balay     uE = x[j][i].u;
6409371c9d4SSatish Balay     wS = x[j - 1][i].w;
6419371c9d4SSatish Balay     wN = x[j][i].w;
642c4762a1bSJed Brown   }
643c4762a1bSJed Brown 
644c4762a1bSJed Brown   if (param->adv_scheme == ADVECT_FV || i == ilim - 1 || j == jlim - 1 || i == 1 || j == 1) {
645c4762a1bSJed Brown     /* finite volume advection */
646c4762a1bSJed Brown     TS = (x[j][i].T + x[j - 1][i].T) / 2.0;
647c4762a1bSJed Brown     TN = (x[j][i].T + x[j + 1][i].T) / 2.0;
648c4762a1bSJed Brown     TE = (x[j][i].T + x[j][i + 1].T) / 2.0;
649c4762a1bSJed Brown     TW = (x[j][i].T + x[j][i - 1].T) / 2.0;
6509371c9d4SSatish Balay     fN = wN * TN * dx;
6519371c9d4SSatish Balay     fS = wS * TS * dx;
6529371c9d4SSatish Balay     fE = uE * TE * dz;
6539371c9d4SSatish Balay     fW = uW * TW * dz;
654c4762a1bSJed Brown 
655c4762a1bSJed Brown   } else {
656c4762a1bSJed Brown     /* Fromm advection scheme */
6579371c9d4SSatish Balay     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;
6589371c9d4SSatish Balay     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;
6599371c9d4SSatish Balay     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;
6609371c9d4SSatish Balay     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;
661c4762a1bSJed Brown   }
662c4762a1bSJed Brown 
663c4762a1bSJed Brown   residual -= (fE - fW + fN - fS);
664c4762a1bSJed Brown 
665c4762a1bSJed Brown   return residual;
666c4762a1bSJed Brown }
667c4762a1bSJed Brown 
668c4762a1bSJed Brown /*---------------------------------------------------------------------*/
669c4762a1bSJed Brown /*  computes the shear stress---used on the boundaries */
6709fbee547SJacob Faibussowitsch static inline PetscScalar ShearStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
671c4762a1bSJed Brown /*---------------------------------------------------------------------*/
672c4762a1bSJed Brown {
673c4762a1bSJed Brown   Parameter  *param = user->param;
674c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
675c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1;
676c4762a1bSJed Brown   PetscScalar uN, uS, wE, wW;
677c4762a1bSJed Brown 
678c4762a1bSJed Brown   if (j <= grid->jlid || i < j || i == ilim || j == jlim) return EPS_ZERO;
679c4762a1bSJed Brown 
680c4762a1bSJed Brown   if (ipos == CELL_CENTER) { /* on cell center */
681c4762a1bSJed Brown 
682c4762a1bSJed Brown     wE = WInterp(x, i, j - 1);
683c4762a1bSJed Brown     if (i == j) {
684c4762a1bSJed Brown       wW = param->sb;
685c4762a1bSJed Brown       uN = param->cb;
686c4762a1bSJed Brown     } else {
687c4762a1bSJed Brown       wW = WInterp(x, i - 1, j - 1);
688c4762a1bSJed Brown       uN = UInterp(x, i - 1, j);
689c4762a1bSJed Brown     }
690c4762a1bSJed Brown     if (j == grid->jlid + 1) uS = 0.0;
691c4762a1bSJed Brown     else uS = UInterp(x, i - 1, j - 1);
692c4762a1bSJed Brown 
693c4762a1bSJed Brown   } else { /* on cell corner */
694c4762a1bSJed Brown 
6959371c9d4SSatish Balay     uN = x[j + 1][i].u;
6969371c9d4SSatish Balay     uS = x[j][i].u;
6979371c9d4SSatish Balay     wW = x[j][i].w;
6989371c9d4SSatish Balay     wE = x[j][i + 1].w;
699c4762a1bSJed Brown   }
700c4762a1bSJed Brown 
701c4762a1bSJed Brown   return (uN - uS) / grid->dz + (wE - wW) / grid->dx;
702c4762a1bSJed Brown }
703c4762a1bSJed Brown 
704c4762a1bSJed Brown /*---------------------------------------------------------------------*/
705c4762a1bSJed Brown /*  computes the normal stress---used on the boundaries */
7069fbee547SJacob Faibussowitsch static inline PetscScalar XNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
707c4762a1bSJed Brown /*---------------------------------------------------------------------*/
708c4762a1bSJed Brown {
709c4762a1bSJed Brown   Parameter  *param = user->param;
710c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
711c4762a1bSJed Brown   PetscScalar dx = grid->dx, dz = grid->dz;
712c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1, ivisc;
713c4762a1bSJed Brown   PetscScalar epsC = 0.0, etaC, TC, uE, uW, pC, z_scale;
714c4762a1bSJed Brown   if (i < j || j <= grid->jlid) return EPS_ZERO;
715c4762a1bSJed Brown 
7169371c9d4SSatish Balay   ivisc   = param->ivisc;
7179371c9d4SSatish Balay   z_scale = param->z_scale;
718c4762a1bSJed Brown 
719c4762a1bSJed Brown   if (ipos == CELL_CENTER) { /* on cell center */
720c4762a1bSJed Brown 
721c4762a1bSJed Brown     TC = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
722c4762a1bSJed Brown     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CENTER, user);
723c4762a1bSJed Brown     etaC = Viscosity(TC, epsC, dz * j, param);
724c4762a1bSJed Brown 
7259371c9d4SSatish Balay     uW = x[j][i - 1].u;
7269371c9d4SSatish Balay     uE = x[j][i].u;
727c4762a1bSJed Brown     pC = x[j][i].p;
728c4762a1bSJed Brown 
729c4762a1bSJed Brown   } else { /* on cell corner */
730c4762a1bSJed Brown     if (i == ilim || j == jlim) return EPS_ZERO;
731c4762a1bSJed Brown 
732c4762a1bSJed Brown     TC = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
733c4762a1bSJed Brown     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CORNER, user);
734c4762a1bSJed Brown     etaC = Viscosity(TC, epsC, dz * (j + 0.5), param);
735c4762a1bSJed Brown 
736c4762a1bSJed Brown     if (i == j) uW = param->sb;
737c4762a1bSJed Brown     else uW = UInterp(x, i - 1, j);
7389371c9d4SSatish Balay     uE = UInterp(x, i, j);
7399371c9d4SSatish Balay     pC = PInterp(x, i, j);
740c4762a1bSJed Brown   }
741c4762a1bSJed Brown 
742c4762a1bSJed Brown   return 2.0 * etaC * (uE - uW) / dx - pC;
743c4762a1bSJed Brown }
744c4762a1bSJed Brown 
745c4762a1bSJed Brown /*---------------------------------------------------------------------*/
746c4762a1bSJed Brown /*  computes the normal stress---used on the boundaries */
7479fbee547SJacob Faibussowitsch static inline PetscScalar ZNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
748c4762a1bSJed Brown /*---------------------------------------------------------------------*/
749c4762a1bSJed Brown {
750c4762a1bSJed Brown   Parameter  *param = user->param;
751c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
752c4762a1bSJed Brown   PetscScalar dz    = grid->dz;
753c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1, ivisc;
754c4762a1bSJed Brown   PetscScalar epsC = 0.0, etaC, TC;
755c4762a1bSJed Brown   PetscScalar pC, wN, wS, z_scale;
756c4762a1bSJed Brown   if (i < j || j <= grid->jlid) return EPS_ZERO;
757c4762a1bSJed Brown 
7589371c9d4SSatish Balay   ivisc   = param->ivisc;
7599371c9d4SSatish Balay   z_scale = param->z_scale;
760c4762a1bSJed Brown 
761c4762a1bSJed Brown   if (ipos == CELL_CENTER) { /* on cell center */
762c4762a1bSJed Brown 
763c4762a1bSJed Brown     TC = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
764c4762a1bSJed Brown     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CENTER, user);
765c4762a1bSJed Brown     etaC = Viscosity(TC, epsC, dz * j, param);
7669371c9d4SSatish Balay     wN   = x[j][i].w;
7679371c9d4SSatish Balay     wS   = x[j - 1][i].w;
7689371c9d4SSatish Balay     pC   = x[j][i].p;
769c4762a1bSJed Brown 
770c4762a1bSJed Brown   } else { /* on cell corner */
771c4762a1bSJed Brown     if ((i == ilim) || (j == jlim)) return EPS_ZERO;
772c4762a1bSJed Brown 
773c4762a1bSJed Brown     TC = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
774c4762a1bSJed Brown     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CORNER, user);
775c4762a1bSJed Brown     etaC = Viscosity(TC, epsC, dz * (j + 0.5), param);
776c4762a1bSJed Brown     if (i == j) wN = param->sb;
777c4762a1bSJed Brown     else wN = WInterp(x, i, j);
7789371c9d4SSatish Balay     wS = WInterp(x, i, j - 1);
7799371c9d4SSatish Balay     pC = PInterp(x, i, j);
780c4762a1bSJed Brown   }
781c4762a1bSJed Brown 
782c4762a1bSJed Brown   return 2.0 * etaC * (wN - wS) / dz - pC;
783c4762a1bSJed Brown }
784c4762a1bSJed Brown 
785c4762a1bSJed Brown /*---------------------------------------------------------------------*/
786c4762a1bSJed Brown 
787c4762a1bSJed Brown /*=====================================================================
788c4762a1bSJed Brown   INITIALIZATION, POST-PROCESSING AND OUTPUT FUNCTIONS
789c4762a1bSJed Brown   =====================================================================*/
790c4762a1bSJed Brown 
791c4762a1bSJed Brown /*---------------------------------------------------------------------*/
792c4762a1bSJed Brown /* initializes the problem parameters and checks for
793c4762a1bSJed Brown    command line changes */
794c4762a1bSJed Brown PetscErrorCode SetParams(Parameter *param, GridInfo *grid)
795c4762a1bSJed Brown /*---------------------------------------------------------------------*/
796c4762a1bSJed Brown {
797c4762a1bSJed Brown   PetscReal SEC_PER_YR                     = 3600.00 * 24.00 * 365.2500;
798c4762a1bSJed Brown   PetscReal alpha_g_on_cp_units_inverse_km = 4.0e-5 * 9.8;
799c4762a1bSJed Brown 
800c4762a1bSJed Brown   /* domain geometry */
801c4762a1bSJed Brown   param->slab_dip    = 45.0;
802c4762a1bSJed Brown   param->width       = 320.0; /* km */
803c4762a1bSJed Brown   param->depth       = 300.0; /* km */
804c4762a1bSJed Brown   param->lid_depth   = 35.0;  /* km */
805c4762a1bSJed Brown   param->fault_depth = 35.0;  /* km */
806c4762a1bSJed Brown 
8079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_dip", &(param->slab_dip), NULL));
8089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-width", &(param->width), NULL));
8099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-depth", &(param->depth), NULL));
8109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-lid_depth", &(param->lid_depth), NULL));
8119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-fault_depth", &(param->fault_depth), NULL));
812c4762a1bSJed Brown 
813c4762a1bSJed Brown   param->slab_dip = param->slab_dip * PETSC_PI / 180.0; /* radians */
814c4762a1bSJed Brown 
815c4762a1bSJed Brown   /* grid information */
8169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-jfault", &(grid->jfault), NULL));
817c4762a1bSJed Brown   grid->ni = 82;
8189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ni", &(grid->ni), NULL));
819c4762a1bSJed Brown 
820c4762a1bSJed Brown   grid->dx     = param->width / ((PetscReal)(grid->ni - 2)); /* km */
821c4762a1bSJed Brown   grid->dz     = grid->dx * PetscTanReal(param->slab_dip);   /* km */
822c4762a1bSJed Brown   grid->nj     = (PetscInt)(param->depth / grid->dz + 3.0);  /* gridpoints*/
823c4762a1bSJed Brown   param->depth = grid->dz * (grid->nj - 2);                  /* km */
824c4762a1bSJed Brown   grid->inose  = 0;                                          /* gridpoints*/
8259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-inose", &(grid->inose), NULL));
826c4762a1bSJed Brown   grid->bx            = DM_BOUNDARY_NONE;
827c4762a1bSJed Brown   grid->by            = DM_BOUNDARY_NONE;
828c4762a1bSJed Brown   grid->stencil       = DMDA_STENCIL_BOX;
829c4762a1bSJed Brown   grid->dof           = 4;
830c4762a1bSJed Brown   grid->stencil_width = 2;
831c4762a1bSJed Brown   grid->mglevels      = 1;
832c4762a1bSJed Brown 
833c4762a1bSJed Brown   /* boundary conditions */
834c4762a1bSJed Brown   param->pv_analytic = PETSC_FALSE;
835c4762a1bSJed Brown   param->ibound      = BC_NOSTRESS;
8369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ibound", &(param->ibound), NULL));
837c4762a1bSJed Brown 
838c4762a1bSJed Brown   /* physical constants */
839c4762a1bSJed Brown   param->slab_velocity = 5.0;       /* cm/yr */
840c4762a1bSJed Brown   param->slab_age      = 50.0;      /* Ma */
841c4762a1bSJed Brown   param->lid_age       = 50.0;      /* Ma */
842c4762a1bSJed Brown   param->kappa         = 0.7272e-6; /* m^2/sec */
843c4762a1bSJed Brown   param->potentialT    = 1300.0;    /* degrees C */
844c4762a1bSJed Brown 
8459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_velocity", &(param->slab_velocity), NULL));
8469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_age", &(param->slab_age), NULL));
8479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-lid_age", &(param->lid_age), NULL));
8489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", &(param->kappa), NULL));
8499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-potentialT", &(param->potentialT), NULL));
850c4762a1bSJed Brown 
851c4762a1bSJed Brown   /* viscosity */
852c4762a1bSJed Brown   param->ivisc        = 3;    /* 0=isovisc, 1=difn creep, 2=disl creep, 3=full */
853c4762a1bSJed Brown   param->eta0         = 1e24; /* Pa-s */
854c4762a1bSJed Brown   param->visc_cutoff  = 0.0;  /* factor of eta_0 */
855c4762a1bSJed Brown   param->continuation = 1.0;
856c4762a1bSJed Brown 
857c4762a1bSJed Brown   /* constants for diffusion creep */
858c4762a1bSJed Brown   param->diffusion.A     = 1.8e7; /* Pa-s */
859c4762a1bSJed Brown   param->diffusion.n     = 1.0;   /* dim'less */
860c4762a1bSJed Brown   param->diffusion.Estar = 375e3; /* J/mol */
861c4762a1bSJed Brown   param->diffusion.Vstar = 5e-6;  /* m^3/mol */
862c4762a1bSJed Brown 
863c4762a1bSJed Brown   /* constants for param->dislocationocation creep */
864c4762a1bSJed Brown   param->dislocation.A     = 2.8969e4; /* Pa-s */
865c4762a1bSJed Brown   param->dislocation.n     = 3.5;      /* dim'less */
866c4762a1bSJed Brown   param->dislocation.Estar = 530e3;    /* J/mol */
867c4762a1bSJed Brown   param->dislocation.Vstar = 14e-6;    /* m^3/mol */
868c4762a1bSJed Brown 
8699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ivisc", &(param->ivisc), NULL));
8709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-visc_cutoff", &(param->visc_cutoff), NULL));
871c4762a1bSJed Brown 
872c4762a1bSJed Brown   param->output_ivisc = param->ivisc;
873c4762a1bSJed Brown 
8749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-output_ivisc", &(param->output_ivisc), NULL));
8759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-vstar", &(param->dislocation.Vstar), NULL));
876c4762a1bSJed Brown 
877c4762a1bSJed Brown   /* output options */
878c4762a1bSJed Brown   param->quiet      = PETSC_FALSE;
879c4762a1bSJed Brown   param->param_test = PETSC_FALSE;
880c4762a1bSJed Brown 
8819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-quiet", &(param->quiet)));
8829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test", &(param->param_test)));
8839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-file", param->filename, sizeof(param->filename), &(param->output_to_file)));
884c4762a1bSJed Brown 
885c4762a1bSJed Brown   /* advection */
886c4762a1bSJed Brown   param->adv_scheme = ADVECT_FROMM; /* advection scheme: 0=finite vol, 1=Fromm */
887c4762a1bSJed Brown 
8889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-adv_scheme", &(param->adv_scheme), NULL));
889c4762a1bSJed Brown 
890c4762a1bSJed Brown   /* misc. flags */
891c4762a1bSJed Brown   param->stop_solve    = PETSC_FALSE;
892c4762a1bSJed Brown   param->interrupted   = PETSC_FALSE;
893c4762a1bSJed Brown   param->kspmon        = PETSC_FALSE;
894c4762a1bSJed Brown   param->toggle_kspmon = PETSC_FALSE;
895c4762a1bSJed Brown 
896c4762a1bSJed Brown   /* derived parameters for slab angle */
897c4762a1bSJed Brown   param->sb = PetscSinReal(param->slab_dip);
898c4762a1bSJed Brown   param->cb = PetscCosReal(param->slab_dip);
899c4762a1bSJed Brown   param->c  = param->slab_dip * param->sb / (param->slab_dip * param->slab_dip - param->sb * param->sb);
900c4762a1bSJed Brown   param->d  = (param->slab_dip * param->cb - param->sb) / (param->slab_dip * param->slab_dip - param->sb * param->sb);
901c4762a1bSJed Brown 
902c4762a1bSJed Brown   /* length, velocity and time scale for non-dimensionalization */
903c4762a1bSJed Brown   param->L = PetscMin(param->width, param->depth);      /* km */
904c4762a1bSJed Brown   param->V = param->slab_velocity / 100.0 / SEC_PER_YR; /* m/sec */
905c4762a1bSJed Brown 
906c4762a1bSJed Brown   /* other unit conversions and derived parameters */
907c4762a1bSJed Brown   param->scaled_width = param->width / param->L;                   /* dim'less */
908c4762a1bSJed Brown   param->scaled_depth = param->depth / param->L;                   /* dim'less */
909c4762a1bSJed Brown   param->lid_depth    = param->lid_depth / param->L;               /* dim'less */
910c4762a1bSJed Brown   param->fault_depth  = param->fault_depth / param->L;             /* dim'less */
911c4762a1bSJed Brown   grid->dx            = grid->dx / param->L;                       /* dim'less */
912c4762a1bSJed Brown   grid->dz            = grid->dz / param->L;                       /* dim'less */
913c4762a1bSJed Brown   grid->jlid          = (PetscInt)(param->lid_depth / grid->dz);   /* gridcells */
914c4762a1bSJed Brown   grid->jfault        = (PetscInt)(param->fault_depth / grid->dz); /* gridcells */
915c4762a1bSJed Brown   param->lid_depth    = grid->jlid * grid->dz;                     /* dim'less */
916c4762a1bSJed Brown   param->fault_depth  = grid->jfault * grid->dz;                   /* dim'less */
917c4762a1bSJed Brown   grid->corner        = grid->jlid + 1;                            /* gridcells */
918c4762a1bSJed Brown   param->peclet       = param->V                                   /* m/sec */
919c4762a1bSJed Brown                 * param->L * 1000.0                                /* m */
920c4762a1bSJed Brown                 / param->kappa;                                    /* m^2/sec */
921c4762a1bSJed Brown   param->z_scale = param->L * alpha_g_on_cp_units_inverse_km;
922c4762a1bSJed Brown   param->skt     = PetscSqrtReal(param->kappa * param->slab_age * SEC_PER_YR);
9239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-peclet", &(param->peclet), NULL));
924c4762a1bSJed Brown 
9255f80ce2aSJacob Faibussowitsch   return 0;
926c4762a1bSJed Brown }
927c4762a1bSJed Brown 
928c4762a1bSJed Brown /*---------------------------------------------------------------------*/
929c4762a1bSJed Brown /*  prints a report of the problem parameters to stdout */
930c4762a1bSJed Brown PetscErrorCode ReportParams(Parameter *param, GridInfo *grid)
931c4762a1bSJed Brown /*---------------------------------------------------------------------*/
932c4762a1bSJed Brown {
933c4762a1bSJed Brown   char date[30];
934c4762a1bSJed Brown 
9359566063dSJacob Faibussowitsch   PetscCall(PetscGetDate(date, 30));
936c4762a1bSJed Brown 
937c4762a1bSJed Brown   if (!(param->quiet)) {
9389566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---------------------BEGIN ex30 PARAM REPORT-------------------\n"));
9399566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Domain: \n"));
9409566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Width = %g km,         Depth = %g km\n", (double)param->width, (double)param->depth));
9419566063dSJacob Faibussowitsch     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));
9429566063dSJacob Faibussowitsch     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)));
943c4762a1bSJed Brown 
9449566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nGrid: \n"));
94563a3b9bcSJacob Faibussowitsch     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)));
94663a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  jlid = %3" PetscInt_FMT "              jfault = %3" PetscInt_FMT " \n", grid->jlid, grid->jfault));
9479566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Pe = %g\n", (double)param->peclet));
948c4762a1bSJed Brown 
9499566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nRheology:"));
950c4762a1bSJed Brown     if (param->ivisc == VISC_CONST) {
9519566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Isoviscous \n"));
952*48a46eb9SPierre Jolivet       if (param->pv_analytic) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Pressure and Velocity prescribed! \n"));
953c4762a1bSJed Brown     } else if (param->ivisc == VISC_DIFN) {
9549566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Diffusion Creep (T-Dependent Newtonian) \n"));
9559566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
956c4762a1bSJed Brown     } else if (param->ivisc == VISC_DISL) {
9579566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Dislocation Creep (T-Dependent Non-Newtonian) \n"));
9589566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
959c4762a1bSJed Brown     } else if (param->ivisc == VISC_FULL) {
9609566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Full Rheology \n"));
9619566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
962c4762a1bSJed Brown     } else {
9639566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Invalid! \n"));
9645f80ce2aSJacob Faibussowitsch       return 1;
965c4762a1bSJed Brown     }
966c4762a1bSJed Brown 
9679566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Boundary condition:"));
968c4762a1bSJed Brown     if (param->ibound == BC_ANALYTIC) {
9699566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Isoviscous Analytic Dirichlet \n"));
970c4762a1bSJed Brown     } else if (param->ibound == BC_NOSTRESS) {
9719566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Stress-Free (normal & shear stress)\n"));
972c4762a1bSJed Brown     } else if (param->ibound == BC_EXPERMNT) {
9739566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Experimental boundary condition \n"));
974c4762a1bSJed Brown     } else {
9759566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Invalid! \n"));
9765f80ce2aSJacob Faibussowitsch       return 1;
977c4762a1bSJed Brown     }
978c4762a1bSJed Brown 
97960d4fc61SSatish Balay     if (param->output_to_file) {
980c4762a1bSJed Brown #if defined(PETSC_HAVE_MATLAB_ENGINE)
9819566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output Destination:       Mat file \"%s\"\n", param->filename));
982c4762a1bSJed Brown #else
9839566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output Destination:       PETSc binary file \"%s\"\n", param->filename));
984c4762a1bSJed Brown #endif
98560d4fc61SSatish Balay     }
986*48a46eb9SPierre Jolivet     if (param->output_ivisc != param->ivisc) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Output viscosity: -ivisc %" PetscInt_FMT "\n", param->output_ivisc));
987c4762a1bSJed Brown 
9889566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---------------------END ex30 PARAM REPORT---------------------\n"));
989c4762a1bSJed Brown   }
990c4762a1bSJed Brown   if (param->param_test) PetscEnd();
9915f80ce2aSJacob Faibussowitsch   return 0;
992c4762a1bSJed Brown }
993c4762a1bSJed Brown 
994c4762a1bSJed Brown /* ------------------------------------------------------------------- */
995a5b23f4aSJose E. Roman /*  generates an initial guess using the analytic solution for isoviscous
996c4762a1bSJed Brown     corner flow */
997c4762a1bSJed Brown PetscErrorCode Initialize(DM da)
998c4762a1bSJed Brown /* ------------------------------------------------------------------- */
999c4762a1bSJed Brown {
1000c4762a1bSJed Brown   AppCtx    *user;
1001c4762a1bSJed Brown   Parameter *param;
1002c4762a1bSJed Brown   GridInfo  *grid;
1003c4762a1bSJed Brown   PetscInt   i, j, is, js, im, jm;
1004c4762a1bSJed Brown   Field    **x;
1005c4762a1bSJed Brown   Vec        Xguess;
1006c4762a1bSJed Brown 
1007c4762a1bSJed Brown   /* Get the fine grid */
10089566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &user));
1009c4762a1bSJed Brown   Xguess = user->Xguess;
1010c4762a1bSJed Brown   param  = user->param;
1011c4762a1bSJed Brown   grid   = user->grid;
10129566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
10139566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xguess, (void **)&x));
1014c4762a1bSJed Brown 
1015c4762a1bSJed Brown   /* Compute initial guess */
1016c4762a1bSJed Brown   for (j = js; j < js + jm; j++) {
1017c4762a1bSJed Brown     for (i = is; i < is + im; i++) {
1018c4762a1bSJed Brown       if (i < j) x[j][i].u = param->cb;
1019c4762a1bSJed Brown       else if (j <= grid->jlid) x[j][i].u = 0.0;
1020c4762a1bSJed Brown       else x[j][i].u = HorizVelocity(i, j, user);
1021c4762a1bSJed Brown 
1022c4762a1bSJed Brown       if (i <= j) x[j][i].w = param->sb;
1023c4762a1bSJed Brown       else if (j <= grid->jlid) x[j][i].w = 0.0;
1024c4762a1bSJed Brown       else x[j][i].w = VertVelocity(i, j, user);
1025c4762a1bSJed Brown 
1026c4762a1bSJed Brown       if (i < j || j <= grid->jlid) x[j][i].p = 0.0;
1027c4762a1bSJed Brown       else x[j][i].p = Pressure(i, j, user);
1028c4762a1bSJed Brown 
1029c4762a1bSJed Brown       x[j][i].T = PetscMin(grid->dz * (j - 0.5), 1.0);
1030c4762a1bSJed Brown     }
1031c4762a1bSJed Brown   }
1032c4762a1bSJed Brown 
1033c4762a1bSJed Brown   /* Restore x to Xguess */
10349566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xguess, (void **)&x));
1035c4762a1bSJed Brown 
1036c4762a1bSJed Brown   return 0;
1037c4762a1bSJed Brown }
1038c4762a1bSJed Brown 
1039c4762a1bSJed Brown /*---------------------------------------------------------------------*/
1040c4762a1bSJed Brown /*  controls output to a file */
1041c4762a1bSJed Brown PetscErrorCode DoOutput(SNES snes, PetscInt its)
1042c4762a1bSJed Brown /*---------------------------------------------------------------------*/
1043c4762a1bSJed Brown {
1044c4762a1bSJed Brown   AppCtx     *user;
1045c4762a1bSJed Brown   Parameter  *param;
1046c4762a1bSJed Brown   GridInfo   *grid;
1047c4762a1bSJed Brown   PetscInt    ivt;
1048c4762a1bSJed Brown   PetscMPIInt rank;
1049c4762a1bSJed Brown   PetscViewer viewer;
1050c4762a1bSJed Brown   Vec         res, pars;
1051c4762a1bSJed Brown   MPI_Comm    comm;
1052c4762a1bSJed Brown   DM          da;
1053c4762a1bSJed Brown 
10549566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
10559566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &user));
1056c4762a1bSJed Brown   param = user->param;
1057c4762a1bSJed Brown   grid  = user->grid;
1058c4762a1bSJed Brown   ivt   = param->ivisc;
1059c4762a1bSJed Brown 
1060c4762a1bSJed Brown   param->ivisc = param->output_ivisc;
1061c4762a1bSJed Brown 
1062c4762a1bSJed Brown   /* compute final residual and final viscosity/strain rate fields */
10639566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &res, NULL, NULL));
10649566063dSJacob Faibussowitsch   PetscCall(ViscosityField(da, user->x, user->Xguess));
1065c4762a1bSJed Brown 
1066c4762a1bSJed Brown   /* get the communicator and the rank of the processor */
10679566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
10689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1069c4762a1bSJed Brown 
1070c4762a1bSJed Brown   if (param->output_to_file) { /* send output to binary file */
10719566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, &pars));
1072dd400576SPatrick Sanan     if (rank == 0) { /* on processor 0 */
10739566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(pars, 20, PETSC_DETERMINE));
10749566063dSJacob Faibussowitsch       PetscCall(VecSetFromOptions(pars));
10759566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 0, (PetscScalar)(grid->ni), INSERT_VALUES));
10769566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 1, (PetscScalar)(grid->nj), INSERT_VALUES));
10779566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 2, (PetscScalar)(grid->dx), INSERT_VALUES));
10789566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 3, (PetscScalar)(grid->dz), INSERT_VALUES));
10799566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 4, (PetscScalar)(param->L), INSERT_VALUES));
10809566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 5, (PetscScalar)(param->V), INSERT_VALUES));
1081c4762a1bSJed Brown       /* skipped 6 intentionally */
10829566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 7, (PetscScalar)(param->slab_dip), INSERT_VALUES));
10839566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 8, (PetscScalar)(grid->jlid), INSERT_VALUES));
10849566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 9, (PetscScalar)(param->lid_depth), INSERT_VALUES));
10859566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 10, (PetscScalar)(grid->jfault), INSERT_VALUES));
10869566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 11, (PetscScalar)(param->fault_depth), INSERT_VALUES));
10879566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 12, (PetscScalar)(param->potentialT), INSERT_VALUES));
10889566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 13, (PetscScalar)(param->ivisc), INSERT_VALUES));
10899566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 14, (PetscScalar)(param->visc_cutoff), INSERT_VALUES));
10909566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 15, (PetscScalar)(param->ibound), INSERT_VALUES));
10919566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 16, (PetscScalar)(its), INSERT_VALUES));
1092c4762a1bSJed Brown     } else { /* on some other processor */
10939566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(pars, 0, PETSC_DETERMINE));
10949566063dSJacob Faibussowitsch       PetscCall(VecSetFromOptions(pars));
1095c4762a1bSJed Brown     }
10969371c9d4SSatish Balay     PetscCall(VecAssemblyBegin(pars));
10979371c9d4SSatish Balay     PetscCall(VecAssemblyEnd(pars));
1098c4762a1bSJed Brown 
1099c4762a1bSJed Brown     /* create viewer */
1100c4762a1bSJed Brown #if defined(PETSC_HAVE_MATLAB_ENGINE)
11019566063dSJacob Faibussowitsch     PetscCall(PetscViewerMatlabOpen(PETSC_COMM_WORLD, param->filename, FILE_MODE_WRITE, &viewer));
1102c4762a1bSJed Brown #else
11039566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, param->filename, FILE_MODE_WRITE, &viewer));
1104c4762a1bSJed Brown #endif
1105c4762a1bSJed Brown 
1106c4762a1bSJed Brown     /* send vectors to viewer */
11079566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)res, "res"));
11089566063dSJacob Faibussowitsch     PetscCall(VecView(res, viewer));
11099566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)user->x, "out"));
11109566063dSJacob Faibussowitsch     PetscCall(VecView(user->x, viewer));
11119566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)(user->Xguess), "aux"));
11129566063dSJacob Faibussowitsch     PetscCall(VecView(user->Xguess, viewer));
11139566063dSJacob Faibussowitsch     PetscCall(StressField(da)); /* compute stress fields */
11149566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)(user->Xguess), "str"));
11159566063dSJacob Faibussowitsch     PetscCall(VecView(user->Xguess, viewer));
11169566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)pars, "par"));
11179566063dSJacob Faibussowitsch     PetscCall(VecView(pars, viewer));
1118c4762a1bSJed Brown 
1119c4762a1bSJed Brown     /* destroy viewer and vector */
11209566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
11219566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pars));
1122c4762a1bSJed Brown   }
1123c4762a1bSJed Brown 
1124c4762a1bSJed Brown   param->ivisc = ivt;
1125c4762a1bSJed Brown   return 0;
1126c4762a1bSJed Brown }
1127c4762a1bSJed Brown 
1128c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1129c4762a1bSJed Brown /* Compute both the second invariant of the strain rate tensor and the viscosity, at both cell centers and cell corners */
1130c4762a1bSJed Brown PetscErrorCode ViscosityField(DM da, Vec X, Vec V)
1131c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1132c4762a1bSJed Brown {
1133c4762a1bSJed Brown   AppCtx    *user;
1134c4762a1bSJed Brown   Parameter *param;
1135c4762a1bSJed Brown   GridInfo  *grid;
1136c4762a1bSJed Brown   Vec        localX;
1137c4762a1bSJed Brown   Field    **v, **x;
1138c4762a1bSJed Brown   PetscReal  eps, /* dx,*/ dz, T, epsC, TC;
1139c4762a1bSJed Brown   PetscInt   i, j, is, js, im, jm, ilim, jlim, ivt;
1140c4762a1bSJed Brown 
1141c4762a1bSJed Brown   PetscFunctionBeginUser;
11429566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &user));
1143c4762a1bSJed Brown   param        = user->param;
1144c4762a1bSJed Brown   grid         = user->grid;
1145c4762a1bSJed Brown   ivt          = param->ivisc;
1146c4762a1bSJed Brown   param->ivisc = param->output_ivisc;
1147c4762a1bSJed Brown 
11489566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
11499566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
11509566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
11519566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, (void **)&x));
11529566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, V, (void **)&v));
1153c4762a1bSJed Brown 
1154c4762a1bSJed Brown   /* Parameters */
1155c4762a1bSJed Brown   /* dx = grid->dx; */ dz = grid->dz;
1156c4762a1bSJed Brown 
11579371c9d4SSatish Balay   ilim = grid->ni - 1;
11589371c9d4SSatish Balay   jlim = grid->nj - 1;
1159c4762a1bSJed Brown 
1160c4762a1bSJed Brown   /* Compute real temperature, strain rate and viscosity */
11619566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
1162c4762a1bSJed Brown   for (j = js; j < js + jm; j++) {
1163c4762a1bSJed Brown     for (i = is; i < is + im; i++) {
1164c4762a1bSJed Brown       T = PetscRealPart(param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * param->z_scale));
1165c4762a1bSJed Brown       if (i < ilim && j < jlim) {
1166c4762a1bSJed Brown         TC = PetscRealPart(param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * param->z_scale));
1167c4762a1bSJed Brown       } else {
1168c4762a1bSJed Brown         TC = T;
1169c4762a1bSJed Brown       }
1170c4762a1bSJed Brown       eps  = PetscRealPart((CalcSecInv(x, i, j, CELL_CENTER, user)));
1171c4762a1bSJed Brown       epsC = PetscRealPart(CalcSecInv(x, i, j, CELL_CORNER, user));
1172c4762a1bSJed Brown 
1173c4762a1bSJed Brown       v[j][i].u = eps;
1174c4762a1bSJed Brown       v[j][i].w = epsC;
1175c4762a1bSJed Brown       v[j][i].p = Viscosity(T, eps, dz * (j - 0.5), param);
1176c4762a1bSJed Brown       v[j][i].T = Viscosity(TC, epsC, dz * j, param);
1177c4762a1bSJed Brown     }
1178c4762a1bSJed Brown   }
11799566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, V, (void **)&v));
11809566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, (void **)&x));
11819566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
1182c4762a1bSJed Brown 
1183c4762a1bSJed Brown   param->ivisc = ivt;
1184c4762a1bSJed Brown   PetscFunctionReturn(0);
1185c4762a1bSJed Brown }
1186c4762a1bSJed Brown 
1187c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1188c4762a1bSJed Brown /* post-processing: compute stress everywhere */
1189c4762a1bSJed Brown PetscErrorCode StressField(DM da)
1190c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1191c4762a1bSJed Brown {
1192c4762a1bSJed Brown   AppCtx  *user;
1193c4762a1bSJed Brown   PetscInt i, j, is, js, im, jm;
1194c4762a1bSJed Brown   Vec      locVec;
1195c4762a1bSJed Brown   Field  **x, **y;
1196c4762a1bSJed Brown 
11979566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &user));
1198c4762a1bSJed Brown 
1199c4762a1bSJed Brown   /* Get the fine grid of Xguess and X */
12009566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
12019566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, user->Xguess, (void **)&x));
1202c4762a1bSJed Brown 
12039566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &locVec));
12049566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, user->x, INSERT_VALUES, locVec));
12059566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, user->x, INSERT_VALUES, locVec));
12069566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, locVec, (void **)&y));
1207c4762a1bSJed Brown 
1208c4762a1bSJed Brown   /* Compute stress on the corner points */
1209c4762a1bSJed Brown   for (j = js; j < js + jm; j++) {
1210c4762a1bSJed Brown     for (i = is; i < is + im; i++) {
1211c4762a1bSJed Brown       x[j][i].u = ShearStress(y, i, j, CELL_CENTER, user);
1212c4762a1bSJed Brown       x[j][i].w = ShearStress(y, i, j, CELL_CORNER, user);
1213c4762a1bSJed Brown       x[j][i].p = XNormalStress(y, i, j, CELL_CENTER, user);
1214c4762a1bSJed Brown       x[j][i].T = ZNormalStress(y, i, j, CELL_CENTER, user);
1215c4762a1bSJed Brown     }
1216c4762a1bSJed Brown   }
1217c4762a1bSJed Brown 
1218c4762a1bSJed Brown   /* Restore the fine grid of Xguess and X */
12199566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, user->Xguess, (void **)&x));
12209566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, locVec, (void **)&y));
12219566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &locVec));
1222c4762a1bSJed Brown   return 0;
1223c4762a1bSJed Brown }
1224c4762a1bSJed Brown 
1225c4762a1bSJed Brown /*=====================================================================
1226c4762a1bSJed Brown   UTILITY FUNCTIONS
1227c4762a1bSJed Brown   =====================================================================*/
1228c4762a1bSJed Brown 
1229c4762a1bSJed Brown /*---------------------------------------------------------------------*/
1230c4762a1bSJed Brown /* returns the velocity of the subducting slab and handles fault nodes
1231c4762a1bSJed Brown    for BC */
12329fbee547SJacob Faibussowitsch static inline PetscScalar SlabVel(char c, PetscInt i, PetscInt j, AppCtx *user)
1233c4762a1bSJed Brown /*---------------------------------------------------------------------*/
1234c4762a1bSJed Brown {
1235c4762a1bSJed Brown   Parameter *param = user->param;
1236c4762a1bSJed Brown   GridInfo  *grid  = user->grid;
1237c4762a1bSJed Brown 
1238c4762a1bSJed Brown   if (c == 'U' || c == 'u') {
1239c4762a1bSJed Brown     if (i < j - 1) return param->cb;
1240c4762a1bSJed Brown     else if (j <= grid->jfault) return 0.0;
1241c4762a1bSJed Brown     else return param->cb;
1242c4762a1bSJed Brown 
1243c4762a1bSJed Brown   } else {
1244c4762a1bSJed Brown     if (i < j) return param->sb;
1245c4762a1bSJed Brown     else if (j <= grid->jfault) return 0.0;
1246c4762a1bSJed Brown     else return param->sb;
1247c4762a1bSJed Brown   }
1248c4762a1bSJed Brown }
1249c4762a1bSJed Brown 
1250c4762a1bSJed Brown /*---------------------------------------------------------------------*/
1251c4762a1bSJed Brown /*  solution to diffusive half-space cooling model for BC */
12529fbee547SJacob Faibussowitsch static inline PetscScalar PlateModel(PetscInt j, PetscInt plate, AppCtx *user)
1253c4762a1bSJed Brown /*---------------------------------------------------------------------*/
1254c4762a1bSJed Brown {
1255c4762a1bSJed Brown   Parameter  *param = user->param;
1256c4762a1bSJed Brown   PetscScalar z;
1257c4762a1bSJed Brown   if (plate == PLATE_LID) z = (j - 0.5) * user->grid->dz;
1258c4762a1bSJed Brown   else z = (j - 0.5) * user->grid->dz * param->cb; /* PLATE_SLAB */
1259c4762a1bSJed Brown #if defined(PETSC_HAVE_ERF)
1260c4762a1bSJed Brown   return (PetscReal)(erf((double)PetscRealPart(z * param->L / 2.0 / param->skt)));
1261c4762a1bSJed Brown #else
1262c4762a1bSJed Brown   (*PetscErrorPrintf)("erf() not available on this machine\n");
1263c4762a1bSJed Brown   MPI_Abort(PETSC_COMM_SELF, 1);
1264c4762a1bSJed Brown #endif
1265c4762a1bSJed Brown }
1266c4762a1bSJed Brown 
1267c4762a1bSJed Brown /*=====================================================================
1268c4762a1bSJed Brown   INTERACTIVE SIGNAL HANDLING
1269c4762a1bSJed Brown   =====================================================================*/
1270c4762a1bSJed Brown 
1271c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1272c4762a1bSJed Brown PetscErrorCode SNESConverged_Interactive(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *ctx)
1273c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1274c4762a1bSJed Brown {
1275c4762a1bSJed Brown   AppCtx    *user  = (AppCtx *)ctx;
1276c4762a1bSJed Brown   Parameter *param = user->param;
1277c4762a1bSJed Brown   KSP        ksp;
1278c4762a1bSJed Brown 
1279c4762a1bSJed Brown   PetscFunctionBeginUser;
1280c4762a1bSJed Brown   if (param->interrupted) {
1281c4762a1bSJed Brown     param->interrupted = PETSC_FALSE;
12829566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: exiting SNES solve. \n"));
1283c4762a1bSJed Brown     *reason = SNES_CONVERGED_FNORM_ABS;
1284c4762a1bSJed Brown     PetscFunctionReturn(0);
1285c4762a1bSJed Brown   } else if (param->toggle_kspmon) {
1286c4762a1bSJed Brown     param->toggle_kspmon = PETSC_FALSE;
1287c4762a1bSJed Brown 
12889566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
1289c4762a1bSJed Brown 
1290c4762a1bSJed Brown     if (param->kspmon) {
12919566063dSJacob Faibussowitsch       PetscCall(KSPMonitorCancel(ksp));
1292c4762a1bSJed Brown 
1293c4762a1bSJed Brown       param->kspmon = PETSC_FALSE;
12949566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: deactivating ksp singular value monitor. \n"));
1295c4762a1bSJed Brown     } else {
1296c4762a1bSJed Brown       PetscViewerAndFormat *vf;
12979566063dSJacob Faibussowitsch       PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf));
12989566063dSJacob Faibussowitsch       PetscCall(KSPMonitorSet(ksp, (PetscErrorCode(*)(KSP, PetscInt, PetscReal, void *))KSPMonitorSingularValue, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
1299c4762a1bSJed Brown 
1300c4762a1bSJed Brown       param->kspmon = PETSC_TRUE;
13019566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: activating ksp singular value monitor. \n"));
1302c4762a1bSJed Brown     }
1303c4762a1bSJed Brown   }
1304c4762a1bSJed Brown   PetscFunctionReturn(SNESConvergedDefault(snes, it, xnorm, snorm, fnorm, reason, ctx));
1305c4762a1bSJed Brown }
1306c4762a1bSJed Brown 
1307c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1308c4762a1bSJed Brown #include <signal.h>
1309c4762a1bSJed Brown PetscErrorCode InteractiveHandler(int signum, void *ctx)
1310c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1311c4762a1bSJed Brown {
1312c4762a1bSJed Brown   AppCtx    *user  = (AppCtx *)ctx;
1313c4762a1bSJed Brown   Parameter *param = user->param;
1314c4762a1bSJed Brown 
1315c4762a1bSJed Brown   if (signum == SIGILL) {
1316c4762a1bSJed Brown     param->toggle_kspmon = PETSC_TRUE;
1317c4762a1bSJed Brown #if !defined(PETSC_MISSING_SIGCONT)
1318c4762a1bSJed Brown   } else if (signum == SIGCONT) {
1319c4762a1bSJed Brown     param->interrupted = PETSC_TRUE;
1320c4762a1bSJed Brown #endif
1321c4762a1bSJed Brown #if !defined(PETSC_MISSING_SIGURG)
1322c4762a1bSJed Brown   } else if (signum == SIGURG) {
1323c4762a1bSJed Brown     param->stop_solve = PETSC_TRUE;
1324c4762a1bSJed Brown #endif
1325c4762a1bSJed Brown   }
1326c4762a1bSJed Brown   return 0;
1327c4762a1bSJed Brown }
1328c4762a1bSJed Brown 
1329c4762a1bSJed Brown /*---------------------------------------------------------------------*/
1330c4762a1bSJed Brown /*  main call-back function that computes the processor-local piece
1331c4762a1bSJed Brown     of the residual */
1332c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field **x, Field **f, void *ptr)
1333c4762a1bSJed Brown /*---------------------------------------------------------------------*/
1334c4762a1bSJed Brown {
1335c4762a1bSJed Brown   AppCtx     *user  = (AppCtx *)ptr;
1336c4762a1bSJed Brown   Parameter  *param = user->param;
1337c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
1338c4762a1bSJed Brown   PetscScalar mag_w, mag_u;
1339c4762a1bSJed Brown   PetscInt    i, j, mx, mz, ilim, jlim;
1340c4762a1bSJed Brown   PetscInt    is, ie, js, je, ibound; /* ,ivisc */
1341c4762a1bSJed Brown 
1342c4762a1bSJed Brown   PetscFunctionBeginUser;
1343c4762a1bSJed Brown   /* Define global and local grid parameters */
13449371c9d4SSatish Balay   mx   = info->mx;
13459371c9d4SSatish Balay   mz   = info->my;
13469371c9d4SSatish Balay   ilim = mx - 1;
13479371c9d4SSatish Balay   jlim = mz - 1;
13489371c9d4SSatish Balay   is   = info->xs;
13499371c9d4SSatish Balay   ie   = info->xs + info->xm;
13509371c9d4SSatish Balay   js   = info->ys;
13519371c9d4SSatish Balay   je   = info->ys + info->ym;
1352c4762a1bSJed Brown 
1353c4762a1bSJed Brown   /* Define geometric and numeric parameters */
1354c4762a1bSJed Brown   /* ivisc = param->ivisc; */ ibound = param->ibound;
1355c4762a1bSJed Brown 
1356c4762a1bSJed Brown   for (j = js; j < je; j++) {
1357c4762a1bSJed Brown     for (i = is; i < ie; i++) {
1358c4762a1bSJed Brown       /************* X-MOMENTUM/VELOCITY *************/
1359c4762a1bSJed Brown       if (i < j) f[j][i].u = x[j][i].u - SlabVel('U', i, j, user);
1360c4762a1bSJed Brown       else if (j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1361c4762a1bSJed Brown         /* in the lithospheric lid */
1362c4762a1bSJed Brown         f[j][i].u = x[j][i].u - 0.0;
1363c4762a1bSJed Brown       } else if (i == ilim) {
1364c4762a1bSJed Brown         /* on the right side boundary */
1365c4762a1bSJed Brown         if (ibound == BC_ANALYTIC) {
1366c4762a1bSJed Brown           f[j][i].u = x[j][i].u - HorizVelocity(i, j, user);
1367c4762a1bSJed Brown         } else {
1368c4762a1bSJed Brown           f[j][i].u = XNormalStress(x, i, j, CELL_CENTER, user) - EPS_ZERO;
1369c4762a1bSJed Brown         }
1370c4762a1bSJed Brown 
1371c4762a1bSJed Brown       } else if (j == jlim) {
1372c4762a1bSJed Brown         /* on the bottom boundary */
1373c4762a1bSJed Brown         if (ibound == BC_ANALYTIC) {
1374c4762a1bSJed Brown           f[j][i].u = x[j][i].u - HorizVelocity(i, j, user);
1375c4762a1bSJed Brown         } else if (ibound == BC_NOSTRESS) {
1376c4762a1bSJed Brown           f[j][i].u = XMomentumResidual(x, i, j, user);
1377c4762a1bSJed Brown         } else {
1378c4762a1bSJed Brown           /* experimental boundary condition */
1379c4762a1bSJed Brown         }
1380c4762a1bSJed Brown 
1381c4762a1bSJed Brown       } else {
1382c4762a1bSJed Brown         /* in the mantle wedge */
1383c4762a1bSJed Brown         f[j][i].u = XMomentumResidual(x, i, j, user);
1384c4762a1bSJed Brown       }
1385c4762a1bSJed Brown 
1386c4762a1bSJed Brown       /************* Z-MOMENTUM/VELOCITY *************/
1387c4762a1bSJed Brown       if (i <= j) {
1388c4762a1bSJed Brown         f[j][i].w = x[j][i].w - SlabVel('W', i, j, user);
1389c4762a1bSJed Brown 
1390c4762a1bSJed Brown       } else if (j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1391c4762a1bSJed Brown         /* in the lithospheric lid */
1392c4762a1bSJed Brown         f[j][i].w = x[j][i].w - 0.0;
1393c4762a1bSJed Brown 
1394c4762a1bSJed Brown       } else if (j == jlim) {
1395c4762a1bSJed Brown         /* on the bottom boundary */
1396c4762a1bSJed Brown         if (ibound == BC_ANALYTIC) {
1397c4762a1bSJed Brown           f[j][i].w = x[j][i].w - VertVelocity(i, j, user);
1398c4762a1bSJed Brown         } else {
1399c4762a1bSJed Brown           f[j][i].w = ZNormalStress(x, i, j, CELL_CENTER, user) - EPS_ZERO;
1400c4762a1bSJed Brown         }
1401c4762a1bSJed Brown 
1402c4762a1bSJed Brown       } else if (i == ilim) {
1403c4762a1bSJed Brown         /* on the right side boundary */
1404c4762a1bSJed Brown         if (ibound == BC_ANALYTIC) {
1405c4762a1bSJed Brown           f[j][i].w = x[j][i].w - VertVelocity(i, j, user);
1406c4762a1bSJed Brown         } else if (ibound == BC_NOSTRESS) {
1407c4762a1bSJed Brown           f[j][i].w = ZMomentumResidual(x, i, j, user);
1408c4762a1bSJed Brown         } else {
1409c4762a1bSJed Brown           /* experimental boundary condition */
1410c4762a1bSJed Brown         }
1411c4762a1bSJed Brown 
1412c4762a1bSJed Brown       } else {
1413c4762a1bSJed Brown         /* in the mantle wedge */
1414c4762a1bSJed Brown         f[j][i].w = ZMomentumResidual(x, i, j, user);
1415c4762a1bSJed Brown       }
1416c4762a1bSJed Brown 
1417c4762a1bSJed Brown       /************* CONTINUITY/PRESSURE *************/
1418c4762a1bSJed Brown       if (i < j || j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1419c4762a1bSJed Brown         /* in the lid or slab */
1420c4762a1bSJed Brown         f[j][i].p = x[j][i].p;
1421c4762a1bSJed Brown 
1422c4762a1bSJed Brown       } else if ((i == ilim || j == jlim) && ibound == BC_ANALYTIC) {
1423c4762a1bSJed Brown         /* on an analytic boundary */
1424c4762a1bSJed Brown         f[j][i].p = x[j][i].p - Pressure(i, j, user);
1425c4762a1bSJed Brown 
1426c4762a1bSJed Brown       } else {
1427c4762a1bSJed Brown         /* in the mantle wedge */
1428c4762a1bSJed Brown         f[j][i].p = ContinuityResidual(x, i, j, user);
1429c4762a1bSJed Brown       }
1430c4762a1bSJed Brown 
1431c4762a1bSJed Brown       /************* TEMPERATURE *************/
1432c4762a1bSJed Brown       if (j == 0) {
1433c4762a1bSJed Brown         /* on the surface */
1434c4762a1bSJed Brown         f[j][i].T = x[j][i].T + x[j + 1][i].T + PetscMax(PetscRealPart(x[j][i].T), 0.0);
1435c4762a1bSJed Brown 
1436c4762a1bSJed Brown       } else if (i == 0) {
1437c4762a1bSJed Brown         /* slab inflow boundary */
1438c4762a1bSJed Brown         f[j][i].T = x[j][i].T - PlateModel(j, PLATE_SLAB, user);
1439c4762a1bSJed Brown 
1440c4762a1bSJed Brown       } else if (i == ilim) {
1441c4762a1bSJed Brown         /* right side boundary */
1442c4762a1bSJed Brown         mag_u     = 1.0 - PetscPowRealInt((1.0 - PetscMax(PetscMin(PetscRealPart(x[j][i - 1].u) / param->cb, 1.0), 0.0)), 5);
1443c4762a1bSJed Brown         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);
1444c4762a1bSJed Brown 
1445c4762a1bSJed Brown       } else if (j == jlim) {
1446c4762a1bSJed Brown         /* bottom boundary */
1447c4762a1bSJed Brown         mag_w     = 1.0 - PetscPowRealInt((1.0 - PetscMax(PetscMin(PetscRealPart(x[j - 1][i].w) / param->sb, 1.0), 0.0)), 5);
1448c4762a1bSJed Brown         f[j][i].T = x[j][i].T - mag_w * x[j - 1][i - 1].T - (1.0 - mag_w);
1449c4762a1bSJed Brown 
1450c4762a1bSJed Brown       } else {
1451c4762a1bSJed Brown         /* in the mantle wedge */
1452c4762a1bSJed Brown         f[j][i].T = EnergyResidual(x, i, j, user);
1453c4762a1bSJed Brown       }
1454c4762a1bSJed Brown     }
1455c4762a1bSJed Brown   }
1456c4762a1bSJed Brown   PetscFunctionReturn(0);
1457c4762a1bSJed Brown }
1458c4762a1bSJed Brown 
1459c4762a1bSJed Brown /*TEST
1460c4762a1bSJed Brown 
1461c4762a1bSJed Brown    build:
1462c4762a1bSJed Brown       requires: !complex erf
1463c4762a1bSJed Brown 
1464c4762a1bSJed Brown    test:
1465c4762a1bSJed Brown       args: -ni 18
1466c4762a1bSJed Brown       filter: grep -v Destination
1467c4762a1bSJed Brown       requires: !single
1468c4762a1bSJed Brown 
1469c4762a1bSJed Brown TEST*/
1470