xref: /petsc/src/snes/tutorials/ex30.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 
123*f6dfbefdSBarry Smith int main(int argc, char **argv) {
124c4762a1bSJed Brown   SNES      snes;
125c4762a1bSJed Brown   AppCtx   *user; /* user-defined work context */
126c4762a1bSJed Brown   Parameter param;
127c4762a1bSJed Brown   GridInfo  grid;
128c4762a1bSJed Brown   PetscInt  nits;
129c4762a1bSJed Brown   MPI_Comm  comm;
130c4762a1bSJed Brown   DM        da;
131c4762a1bSJed Brown 
132327415f7SBarry Smith   PetscFunctionBeginUser;
1339566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
134c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-file", "ex30_output");
135c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-snes_monitor_short", NULL);
136c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-snes_max_it", "20");
137c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-ksp_max_it", "1500");
138c4762a1bSJed Brown   PetscOptionsSetValue(NULL, "-ksp_gmres_restart", "300");
139c4762a1bSJed Brown   PetscOptionsInsert(NULL, &argc, &argv, NULL);
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown      Set up the problem parameters.
145c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1469566063dSJacob Faibussowitsch   PetscCall(SetParams(&param, &grid));
1479566063dSJacob Faibussowitsch   PetscCall(ReportParams(&param, &grid));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1519566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &snes));
1529566063dSJacob 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));
1539566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1549566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1559566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
1569566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "x-velocity"));
1579566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "y-velocity"));
1589566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 2, "pressure"));
1599566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 3, "temperature"));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162c4762a1bSJed Brown      Create user context, set problem data, create vector data structures.
163c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1649566063dSJacob Faibussowitsch   PetscCall(PetscNew(&user));
165c4762a1bSJed Brown   user->param = &param;
166c4762a1bSJed Brown   user->grid  = &grid;
1679566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, user));
1689566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &(user->Xguess)));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171c4762a1bSJed Brown      Set up the SNES solver with callback functions.
172c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1739566063dSJacob Faibussowitsch   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, (void *)user));
1749566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch   PetscCall(SNESSetConvergenceTest(snes, SNESConverged_Interactive, (void *)user, NULL));
1779566063dSJacob Faibussowitsch   PetscCall(PetscPushSignalHandler(InteractiveHandler, (void *)user));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown      Initialize and solve the nonlinear system
181c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1829566063dSJacob Faibussowitsch   PetscCall(Initialize(da));
1839566063dSJacob Faibussowitsch   PetscCall(UpdateSolution(snes, user, &nits));
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186c4762a1bSJed Brown      Output variables.
187c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1889566063dSJacob Faibussowitsch   PetscCall(DoOutput(snes, nits));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191c4762a1bSJed Brown      Free work space.
192c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Xguess));
1949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->x));
1959566063dSJacob Faibussowitsch   PetscCall(PetscFree(user));
1969566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1979566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1989566063dSJacob Faibussowitsch   PetscCall(PetscPopSignalHandler());
1999566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
200b122ec5aSJacob Faibussowitsch   return 0;
201c4762a1bSJed Brown }
202c4762a1bSJed Brown 
203c4762a1bSJed Brown /*=====================================================================
204c4762a1bSJed Brown   PETSc INTERACTION FUNCTIONS (initialize & call SNESSolve)
205c4762a1bSJed Brown   =====================================================================*/
206c4762a1bSJed Brown 
207c4762a1bSJed Brown /*  manages solve: adaptive continuation method  */
2089371c9d4SSatish Balay PetscErrorCode UpdateSolution(SNES snes, AppCtx *user, PetscInt *nits) {
209c4762a1bSJed Brown   KSP                 ksp;
210c4762a1bSJed Brown   PC                  pc;
211c4762a1bSJed Brown   SNESConvergedReason reason    = SNES_CONVERGED_ITERATING;
212c4762a1bSJed Brown   Parameter          *param     = user->param;
213c4762a1bSJed Brown   PetscReal           cont_incr = 0.3;
214c4762a1bSJed Brown   PetscInt            its;
215c4762a1bSJed Brown   PetscBool           q = PETSC_FALSE;
216c4762a1bSJed Brown   DM                  dm;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   PetscFunctionBeginUser;
2199566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
2209566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &user->x));
2219566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
2229566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
2239566063dSJacob Faibussowitsch   PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   *nits = 0;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /* Isoviscous solve */
228c4762a1bSJed Brown   if (param->ivisc == VISC_CONST && !param->stop_solve) {
229c4762a1bSJed Brown     param->ivisc = VISC_CONST;
230c4762a1bSJed Brown 
2319566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, 0, user->x));
2329566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes, &reason));
2339566063dSJacob Faibussowitsch     PetscCall(SNESGetIterationNumber(snes, &its));
234c4762a1bSJed Brown     *nits += its;
2359566063dSJacob Faibussowitsch     PetscCall(VecCopy(user->x, user->Xguess));
236c4762a1bSJed Brown     if (param->stop_solve) goto done;
237c4762a1bSJed Brown   }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   /* Olivine diffusion creep */
240c4762a1bSJed Brown   if (param->ivisc >= VISC_DIFN && !param->stop_solve) {
2419566063dSJacob Faibussowitsch     if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Computing Variable Viscosity Solution\n"));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown     /* continuation method on viscosity cutoff */
244c4762a1bSJed Brown     for (param->continuation = 0.0;; param->continuation += cont_incr) {
2459566063dSJacob Faibussowitsch       if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Continuation parameter = %g\n", (double)param->continuation));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown       /* solve the non-linear system */
2489566063dSJacob Faibussowitsch       PetscCall(VecCopy(user->Xguess, user->x));
2499566063dSJacob Faibussowitsch       PetscCall(SNESSolve(snes, 0, user->x));
2509566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes, &reason));
2519566063dSJacob Faibussowitsch       PetscCall(SNESGetIterationNumber(snes, &its));
252c4762a1bSJed Brown       *nits += its;
25363a3b9bcSJacob Faibussowitsch       if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SNES iterations: %" PetscInt_FMT ", Cumulative: %" PetscInt_FMT "\n", its, *nits));
254c4762a1bSJed Brown       if (param->stop_solve) goto done;
255c4762a1bSJed Brown 
256c4762a1bSJed Brown       if (reason < 0) {
257c4762a1bSJed Brown         /* NOT converged */
258c4762a1bSJed Brown         cont_incr = -PetscAbsReal(cont_incr) / 2.0;
259c4762a1bSJed Brown         if (PetscAbsReal(cont_incr) < 0.01) goto done;
260c4762a1bSJed Brown 
261c4762a1bSJed Brown       } else {
262c4762a1bSJed Brown         /* converged */
2639566063dSJacob Faibussowitsch         PetscCall(VecCopy(user->x, user->Xguess));
264c4762a1bSJed Brown         if (param->continuation >= 1.0) goto done;
265c4762a1bSJed Brown         if (its <= 3) cont_incr = 0.30001;
266c4762a1bSJed Brown         else if (its <= 8) cont_incr = 0.15001;
267c4762a1bSJed Brown         else cont_incr = 0.10001;
268c4762a1bSJed Brown 
269c4762a1bSJed Brown         if (param->continuation + cont_incr > 1.0) cont_incr = 1.0 - param->continuation;
270c4762a1bSJed Brown       } /* endif reason<0 */
271c4762a1bSJed Brown     }
272c4762a1bSJed Brown   }
273c4762a1bSJed Brown done:
2749566063dSJacob Faibussowitsch   if (param->stop_solve && !q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: stopping solve.\n"));
2759566063dSJacob Faibussowitsch   if (reason < 0 && !q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "FAILED TO CONVERGE: stopping solve.\n"));
276c4762a1bSJed Brown   PetscFunctionReturn(0);
277c4762a1bSJed Brown }
278c4762a1bSJed Brown 
279c4762a1bSJed Brown /*=====================================================================
280c4762a1bSJed Brown   PHYSICS FUNCTIONS (compute the discrete residual)
281c4762a1bSJed Brown   =====================================================================*/
282c4762a1bSJed Brown 
283*f6dfbefdSBarry Smith static inline PetscScalar UInterp(Field **x, PetscInt i, PetscInt j) {
284c4762a1bSJed Brown   return 0.25 * (x[j][i].u + x[j + 1][i].u + x[j][i + 1].u + x[j + 1][i + 1].u);
285c4762a1bSJed Brown }
286c4762a1bSJed Brown 
287*f6dfbefdSBarry Smith static inline PetscScalar WInterp(Field **x, PetscInt i, PetscInt j) {
288c4762a1bSJed Brown   return 0.25 * (x[j][i].w + x[j + 1][i].w + x[j][i + 1].w + x[j + 1][i + 1].w);
289c4762a1bSJed Brown }
290c4762a1bSJed Brown 
291*f6dfbefdSBarry Smith static inline PetscScalar PInterp(Field **x, PetscInt i, PetscInt j) {
292c4762a1bSJed Brown   return 0.25 * (x[j][i].p + x[j + 1][i].p + x[j][i + 1].p + x[j + 1][i + 1].p);
293c4762a1bSJed Brown }
294c4762a1bSJed Brown 
295*f6dfbefdSBarry Smith static inline PetscScalar TInterp(Field **x, PetscInt i, PetscInt j) {
296c4762a1bSJed Brown   return 0.25 * (x[j][i].T + x[j + 1][i].T + x[j][i + 1].T + x[j + 1][i + 1].T);
297c4762a1bSJed Brown }
298c4762a1bSJed Brown 
299c4762a1bSJed Brown /*  isoviscous analytic solution for IC */
300*f6dfbefdSBarry Smith static inline PetscScalar HorizVelocity(PetscInt i, PetscInt j, AppCtx *user) {
301c4762a1bSJed Brown   Parameter  *param = user->param;
302c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
303c4762a1bSJed Brown   PetscScalar st, ct, th, c = param->c, d = param->d;
304c4762a1bSJed Brown   PetscReal   x, z, r;
305c4762a1bSJed Brown 
3069371c9d4SSatish Balay   x  = (i - grid->jlid) * grid->dx;
3079371c9d4SSatish Balay   z  = (j - grid->jlid - 0.5) * grid->dz;
308c4762a1bSJed Brown   r  = PetscSqrtReal(x * x + z * z);
309c4762a1bSJed Brown   st = z / r;
310c4762a1bSJed Brown   ct = x / r;
311c4762a1bSJed Brown   th = PetscAtanReal(z / x);
312c4762a1bSJed Brown   return ct * (c * th * st + d * (st + th * ct)) + st * (c * (st - th * ct) + d * th * st);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown /*  isoviscous analytic solution for IC */
3169fbee547SJacob Faibussowitsch static inline PetscScalar VertVelocity(PetscInt i, PetscInt j, AppCtx *user)
317*f6dfbefdSBarry Smith 
318c4762a1bSJed Brown {
319c4762a1bSJed Brown   Parameter  *param = user->param;
320c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
321c4762a1bSJed Brown   PetscScalar st, ct, th, c = param->c, d = param->d;
322c4762a1bSJed Brown   PetscReal   x, z, r;
323c4762a1bSJed Brown 
3249371c9d4SSatish Balay   x  = (i - grid->jlid - 0.5) * grid->dx;
3259371c9d4SSatish Balay   z  = (j - grid->jlid) * grid->dz;
3269371c9d4SSatish Balay   r  = PetscSqrtReal(x * x + z * z);
3279371c9d4SSatish Balay   st = z / r;
3289371c9d4SSatish Balay   ct = x / r;
3299371c9d4SSatish Balay   th = PetscAtanReal(z / x);
330c4762a1bSJed Brown   return st * (c * th * st + d * (st + th * ct)) - ct * (c * (st - th * ct) + d * th * st);
331c4762a1bSJed Brown }
332c4762a1bSJed Brown 
333c4762a1bSJed Brown /*  isoviscous analytic solution for IC */
334*f6dfbefdSBarry Smith static inline PetscScalar Pressure(PetscInt i, PetscInt j, AppCtx *user) {
335c4762a1bSJed Brown   Parameter  *param = user->param;
336c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
337c4762a1bSJed Brown   PetscScalar x, z, r, st, ct, c = param->c, d = param->d;
338c4762a1bSJed Brown 
3399371c9d4SSatish Balay   x  = (i - grid->jlid - 0.5) * grid->dx;
3409371c9d4SSatish Balay   z  = (j - grid->jlid - 0.5) * grid->dz;
3419371c9d4SSatish Balay   r  = PetscSqrtReal(x * x + z * z);
3429371c9d4SSatish Balay   st = z / r;
3439371c9d4SSatish Balay   ct = x / r;
344c4762a1bSJed Brown   return (-2.0 * (c * ct - d * st) / r);
345c4762a1bSJed Brown }
346c4762a1bSJed Brown 
347c4762a1bSJed Brown /*  computes the second invariant of the strain rate tensor */
348*f6dfbefdSBarry Smith static inline PetscScalar CalcSecInv(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user) {
349c4762a1bSJed Brown   Parameter  *param = user->param;
350c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
351c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1;
352c4762a1bSJed Brown   PetscScalar uN, uS, uE, uW, wN, wS, wE, wW;
353c4762a1bSJed Brown   PetscScalar eps11, eps12, eps22;
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   if (i < j) return EPS_ZERO;
356c4762a1bSJed Brown   if (i == ilim) i--;
357c4762a1bSJed Brown   if (j == jlim) j--;
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   if (ipos == CELL_CENTER) { /* on cell center */
360c4762a1bSJed Brown     if (j <= grid->jlid) return EPS_ZERO;
361c4762a1bSJed Brown 
3629371c9d4SSatish Balay     uE = x[j][i].u;
3639371c9d4SSatish Balay     uW = x[j][i - 1].u;
3649371c9d4SSatish Balay     wN = x[j][i].w;
3659371c9d4SSatish Balay     wS = x[j - 1][i].w;
366c4762a1bSJed Brown     wE = WInterp(x, i, j - 1);
367c4762a1bSJed Brown     if (i == j) {
3689371c9d4SSatish Balay       uN = param->cb;
3699371c9d4SSatish Balay       wW = param->sb;
370c4762a1bSJed Brown     } else {
3719371c9d4SSatish Balay       uN = UInterp(x, i - 1, j);
3729371c9d4SSatish Balay       wW = WInterp(x, i - 1, j - 1);
373c4762a1bSJed Brown     }
374c4762a1bSJed Brown 
375c4762a1bSJed Brown     if (j == grid->jlid + 1) uS = 0.0;
376c4762a1bSJed Brown     else uS = UInterp(x, i - 1, j - 1);
377c4762a1bSJed Brown 
378c4762a1bSJed Brown   } else { /* on CELL_CORNER */
379c4762a1bSJed Brown     if (j < grid->jlid) return EPS_ZERO;
380c4762a1bSJed Brown 
3819371c9d4SSatish Balay     uN = x[j + 1][i].u;
3829371c9d4SSatish Balay     uS = x[j][i].u;
3839371c9d4SSatish Balay     wE = x[j][i + 1].w;
3849371c9d4SSatish Balay     wW = x[j][i].w;
385c4762a1bSJed Brown     if (i == j) {
386c4762a1bSJed Brown       wN = param->sb;
387c4762a1bSJed Brown       uW = param->cb;
388c4762a1bSJed Brown     } else {
389c4762a1bSJed Brown       wN = WInterp(x, i, j);
390c4762a1bSJed Brown       uW = UInterp(x, i - 1, j);
391c4762a1bSJed Brown     }
392c4762a1bSJed Brown 
393c4762a1bSJed Brown     if (j == grid->jlid) {
3949371c9d4SSatish Balay       uE = 0.0;
3959371c9d4SSatish Balay       uW = 0.0;
396c4762a1bSJed Brown       uS = -uN;
397c4762a1bSJed Brown       wS = -wN;
398c4762a1bSJed Brown     } else {
399c4762a1bSJed Brown       uE = UInterp(x, i, j);
400c4762a1bSJed Brown       wS = WInterp(x, i, j - 1);
401c4762a1bSJed Brown     }
402c4762a1bSJed Brown   }
403c4762a1bSJed Brown 
4049371c9d4SSatish Balay   eps11 = (uE - uW) / grid->dx;
4059371c9d4SSatish Balay   eps22 = (wN - wS) / grid->dz;
406c4762a1bSJed Brown   eps12 = 0.5 * ((uN - uS) / grid->dz + (wE - wW) / grid->dx);
407c4762a1bSJed Brown 
408c4762a1bSJed Brown   return PetscSqrtReal(0.5 * (eps11 * eps11 + 2.0 * eps12 * eps12 + eps22 * eps22));
409c4762a1bSJed Brown }
410c4762a1bSJed Brown 
411c4762a1bSJed Brown /*  computes the shear viscosity */
412*f6dfbefdSBarry Smith static inline PetscScalar Viscosity(PetscScalar T, PetscScalar eps, PetscScalar z, Parameter *param) {
413c4762a1bSJed Brown   PetscReal   result = 0.0;
414c4762a1bSJed Brown   ViscParam   difn = param->diffusion, disl = param->dislocation;
415c4762a1bSJed Brown   PetscInt    iVisc     = param->ivisc;
416c4762a1bSJed Brown   PetscScalar eps_scale = param->V / (param->L * 1000.0);
417c4762a1bSJed Brown   PetscScalar strain_power, v1, v2, P;
418c4762a1bSJed Brown   PetscScalar rho_g = 32340.0, R = 8.3144;
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   P = rho_g * (z * param->L * 1000.0); /* Pa */
421c4762a1bSJed Brown 
422c4762a1bSJed Brown   if (iVisc == VISC_CONST) {
423c4762a1bSJed Brown     /* constant viscosity */
424c4762a1bSJed Brown     return 1.0;
425c4762a1bSJed Brown   } else if (iVisc == VISC_DIFN) {
426c4762a1bSJed Brown     /* diffusion creep rheology */
427c4762a1bSJed Brown     result = PetscRealPart((difn.A * PetscExpScalar((difn.Estar + P * difn.Vstar) / R / (T + 273.0)) / param->eta0));
428c4762a1bSJed Brown   } else if (iVisc == VISC_DISL) {
429c4762a1bSJed Brown     /* dislocation creep rheology */
430c4762a1bSJed Brown     strain_power = PetscPowScalar(eps * eps_scale, (1.0 - disl.n) / disl.n);
431c4762a1bSJed Brown 
432c4762a1bSJed Brown     result = PetscRealPart(disl.A * PetscExpScalar((disl.Estar + P * disl.Vstar) / disl.n / R / (T + 273.0)) * strain_power / param->eta0);
433c4762a1bSJed Brown   } else if (iVisc == VISC_FULL) {
434c4762a1bSJed Brown     /* dislocation/diffusion creep rheology */
435c4762a1bSJed Brown     strain_power = PetscPowScalar(eps * eps_scale, (1.0 - disl.n) / disl.n);
436c4762a1bSJed Brown 
437c4762a1bSJed Brown     v1 = difn.A * PetscExpScalar((difn.Estar + P * difn.Vstar) / R / (T + 273.0)) / param->eta0;
438c4762a1bSJed Brown     v2 = disl.A * PetscExpScalar((disl.Estar + P * disl.Vstar) / disl.n / R / (T + 273.0)) * strain_power / param->eta0;
439c4762a1bSJed Brown 
440c4762a1bSJed Brown     result = PetscRealPart(1.0 / (1.0 / v1 + 1.0 / v2));
441c4762a1bSJed Brown   }
442c4762a1bSJed Brown 
443c4762a1bSJed Brown   /* max viscosity is param->eta0 */
444c4762a1bSJed Brown   result = PetscMin(result, 1.0);
445c4762a1bSJed Brown   /* min viscosity is param->visc_cutoff */
446c4762a1bSJed Brown   result = PetscMax(result, param->visc_cutoff);
447c4762a1bSJed Brown   /* continuation method */
448c4762a1bSJed Brown   result = PetscPowReal(result, param->continuation);
449c4762a1bSJed Brown   return result;
450c4762a1bSJed Brown }
451c4762a1bSJed Brown 
452c4762a1bSJed Brown /*  computes the residual of the x-component of eqn (1) above */
453*f6dfbefdSBarry Smith static inline PetscScalar XMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user) {
454c4762a1bSJed Brown   Parameter  *param = user->param;
455c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
456c4762a1bSJed Brown   PetscScalar dx = grid->dx, dz = grid->dz;
457c4762a1bSJed Brown   PetscScalar etaN, etaS, etaE, etaW, epsN = 0.0, epsS = 0.0, epsE = 0.0, epsW = 0.0;
458c4762a1bSJed Brown   PetscScalar TE = 0.0, TN = 0.0, TS = 0.0, TW = 0.0, dPdx, residual, z_scale;
459c4762a1bSJed Brown   PetscScalar dudxW, dudxE, dudzN, dudzS, dwdxN, dwdxS;
460c4762a1bSJed Brown   PetscInt    jlim = grid->nj - 1;
461c4762a1bSJed Brown 
462c4762a1bSJed Brown   z_scale = param->z_scale;
463c4762a1bSJed Brown 
464c4762a1bSJed Brown   if (param->ivisc == VISC_DIFN || param->ivisc >= VISC_DISL) { /* viscosity is T-dependent */
465c4762a1bSJed Brown     TS = param->potentialT * TInterp(x, i, j - 1) * PetscExpScalar((j - 1.0) * dz * z_scale);
466c4762a1bSJed Brown     if (j == jlim) TN = TS;
467c4762a1bSJed Brown     else TN = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
468c4762a1bSJed Brown     TW = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
469c4762a1bSJed Brown     TE = param->potentialT * x[j][i + 1].T * PetscExpScalar((j - 0.5) * dz * z_scale);
470c4762a1bSJed Brown     if (param->ivisc >= VISC_DISL) { /* olivine dislocation creep */
471c4762a1bSJed Brown       epsN = CalcSecInv(x, i, j, CELL_CORNER, user);
472c4762a1bSJed Brown       epsS = CalcSecInv(x, i, j - 1, CELL_CORNER, user);
473c4762a1bSJed Brown       epsE = CalcSecInv(x, i + 1, j, CELL_CENTER, user);
474c4762a1bSJed Brown       epsW = CalcSecInv(x, i, j, CELL_CENTER, user);
475c4762a1bSJed Brown     }
476c4762a1bSJed Brown   }
477c4762a1bSJed Brown   etaN = Viscosity(TN, epsN, dz * (j + 0.5), param);
478c4762a1bSJed Brown   etaS = Viscosity(TS, epsS, dz * (j - 0.5), param);
479c4762a1bSJed Brown   etaW = Viscosity(TW, epsW, dz * j, param);
480c4762a1bSJed Brown   etaE = Viscosity(TE, epsE, dz * j, param);
481c4762a1bSJed Brown 
482c4762a1bSJed Brown   dPdx = (x[j][i + 1].p - x[j][i].p) / dx;
483c4762a1bSJed Brown   if (j == jlim) dudzN = etaN * (x[j][i].w - x[j][i + 1].w) / dx;
484c4762a1bSJed Brown   else dudzN = etaN * (x[j + 1][i].u - x[j][i].u) / dz;
485c4762a1bSJed Brown   dudzS = etaS * (x[j][i].u - x[j - 1][i].u) / dz;
486c4762a1bSJed Brown   dudxE = etaE * (x[j][i + 1].u - x[j][i].u) / dx;
487c4762a1bSJed Brown   dudxW = etaW * (x[j][i].u - x[j][i - 1].u) / dx;
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   residual = -dPdx /* X-MOMENTUM EQUATION*/
4909371c9d4SSatish Balay            + (dudxE - dudxW) / dx + (dudzN - dudzS) / dz;
491c4762a1bSJed Brown 
492c4762a1bSJed Brown   if (param->ivisc != VISC_CONST) {
493c4762a1bSJed Brown     dwdxN = etaN * (x[j][i + 1].w - x[j][i].w) / dx;
494c4762a1bSJed Brown     dwdxS = etaS * (x[j - 1][i + 1].w - x[j - 1][i].w) / dx;
495c4762a1bSJed Brown 
496c4762a1bSJed Brown     residual += (dudxE - dudxW) / dx + (dwdxN - dwdxS) / dz;
497c4762a1bSJed Brown   }
498c4762a1bSJed Brown 
499c4762a1bSJed Brown   return residual;
500c4762a1bSJed Brown }
501c4762a1bSJed Brown 
502c4762a1bSJed Brown /*  computes the residual of the z-component of eqn (1) above */
5039fbee547SJacob Faibussowitsch static inline PetscScalar ZMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
504*f6dfbefdSBarry Smith 
505c4762a1bSJed Brown {
506c4762a1bSJed Brown   Parameter  *param = user->param;
507c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
508c4762a1bSJed Brown   PetscScalar dx = grid->dx, dz = grid->dz;
509c4762a1bSJed 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;
510c4762a1bSJed Brown   PetscScalar TE = 0.0, TN = 0.0, TS = 0.0, TW = 0.0, dPdz, residual, z_scale;
511c4762a1bSJed Brown   PetscScalar dudzE, dudzW, dwdxW, dwdxE, dwdzN, dwdzS;
512c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1;
513c4762a1bSJed Brown 
514c4762a1bSJed Brown   /* geometric and other parameters */
515c4762a1bSJed Brown   z_scale = param->z_scale;
516c4762a1bSJed Brown 
517c4762a1bSJed Brown   /* viscosity */
518c4762a1bSJed Brown   if (param->ivisc == VISC_DIFN || param->ivisc >= VISC_DISL) { /* viscosity is T-dependent */
519c4762a1bSJed Brown     TN = param->potentialT * x[j + 1][i].T * PetscExpScalar((j + 0.5) * dz * z_scale);
520c4762a1bSJed Brown     TS = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
521c4762a1bSJed Brown     TW = param->potentialT * TInterp(x, i - 1, j) * PetscExpScalar(j * dz * z_scale);
522c4762a1bSJed Brown     if (i == ilim) TE = TW;
523c4762a1bSJed Brown     else TE = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
524c4762a1bSJed Brown     if (param->ivisc >= VISC_DISL) { /* olivine dislocation creep */
525c4762a1bSJed Brown       epsN = CalcSecInv(x, i, j + 1, CELL_CENTER, user);
526c4762a1bSJed Brown       epsS = CalcSecInv(x, i, j, CELL_CENTER, user);
527c4762a1bSJed Brown       epsE = CalcSecInv(x, i, j, CELL_CORNER, user);
528c4762a1bSJed Brown       epsW = CalcSecInv(x, i - 1, j, CELL_CORNER, user);
529c4762a1bSJed Brown     }
530c4762a1bSJed Brown   }
531c4762a1bSJed Brown   etaN = Viscosity(TN, epsN, dz * (j + 1.0), param);
532c4762a1bSJed Brown   etaS = Viscosity(TS, epsS, dz * (j + 0.0), param);
533c4762a1bSJed Brown   etaW = Viscosity(TW, epsW, dz * (j + 0.5), param);
534c4762a1bSJed Brown   etaE = Viscosity(TE, epsE, dz * (j + 0.5), param);
535c4762a1bSJed Brown 
536c4762a1bSJed Brown   dPdz  = (x[j + 1][i].p - x[j][i].p) / dz;
537c4762a1bSJed Brown   dwdzN = etaN * (x[j + 1][i].w - x[j][i].w) / dz;
538c4762a1bSJed Brown   dwdzS = etaS * (x[j][i].w - x[j - 1][i].w) / dz;
539c4762a1bSJed Brown   if (i == ilim) dwdxE = etaE * (x[j][i].u - x[j + 1][i].u) / dz;
540c4762a1bSJed Brown   else dwdxE = etaE * (x[j][i + 1].w - x[j][i].w) / dx;
541c4762a1bSJed Brown   dwdxW = 2.0 * etaW * (x[j][i].w - x[j][i - 1].w) / dx;
542c4762a1bSJed Brown 
543c4762a1bSJed Brown   /* Z-MOMENTUM */
544c4762a1bSJed Brown   residual = -dPdz /* constant viscosity terms */
5459371c9d4SSatish Balay            + (dwdzN - dwdzS) / dz + (dwdxE - dwdxW) / dx;
546c4762a1bSJed Brown 
547c4762a1bSJed Brown   if (param->ivisc != VISC_CONST) {
548c4762a1bSJed Brown     dudzE = etaE * (x[j + 1][i].u - x[j][i].u) / dz;
549c4762a1bSJed Brown     dudzW = etaW * (x[j + 1][i - 1].u - x[j][i - 1].u) / dz;
550c4762a1bSJed Brown 
551c4762a1bSJed Brown     residual += (dwdzN - dwdzS) / dz + (dudzE - dudzW) / dx;
552c4762a1bSJed Brown   }
553c4762a1bSJed Brown 
554c4762a1bSJed Brown   return residual;
555c4762a1bSJed Brown }
556c4762a1bSJed Brown 
557c4762a1bSJed Brown /*  computes the residual of eqn (2) above */
558*f6dfbefdSBarry Smith static inline PetscScalar ContinuityResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user) {
559c4762a1bSJed Brown   GridInfo   *grid = user->grid;
560c4762a1bSJed Brown   PetscScalar uE, uW, wN, wS, dudx, dwdz;
561c4762a1bSJed Brown 
5629371c9d4SSatish Balay   uW   = x[j][i - 1].u;
5639371c9d4SSatish Balay   uE   = x[j][i].u;
5649371c9d4SSatish Balay   dudx = (uE - uW) / grid->dx;
5659371c9d4SSatish Balay   wS   = x[j - 1][i].w;
5669371c9d4SSatish Balay   wN   = x[j][i].w;
5679371c9d4SSatish Balay   dwdz = (wN - wS) / grid->dz;
568c4762a1bSJed Brown 
569c4762a1bSJed Brown   return dudx + dwdz;
570c4762a1bSJed Brown }
571c4762a1bSJed Brown 
572c4762a1bSJed Brown /*  computes the residual of eqn (3) above */
573*f6dfbefdSBarry Smith static inline PetscScalar EnergyResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user) {
574c4762a1bSJed Brown   Parameter  *param = user->param;
575c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
576c4762a1bSJed Brown   PetscScalar dx = grid->dx, dz = grid->dz;
577c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1, jlid = grid->jlid;
578c4762a1bSJed Brown   PetscScalar TE, TN, TS, TW, residual;
579c4762a1bSJed Brown   PetscScalar uE, uW, wN, wS;
580c4762a1bSJed Brown   PetscScalar fN, fS, fE, fW, dTdxW, dTdxE, dTdzN, dTdzS;
581c4762a1bSJed Brown 
582c4762a1bSJed Brown   dTdzN = (x[j + 1][i].T - x[j][i].T) / dz;
583c4762a1bSJed Brown   dTdzS = (x[j][i].T - x[j - 1][i].T) / dz;
584c4762a1bSJed Brown   dTdxE = (x[j][i + 1].T - x[j][i].T) / dx;
585c4762a1bSJed Brown   dTdxW = (x[j][i].T - x[j][i - 1].T) / dx;
586c4762a1bSJed Brown 
587c4762a1bSJed Brown   residual = ((dTdzN - dTdzS) / dz + /* diffusion term */
5889371c9d4SSatish Balay               (dTdxE - dTdxW) / dx) *
5899371c9d4SSatish Balay              dx * dz / param->peclet;
590c4762a1bSJed Brown 
591c4762a1bSJed Brown   if (j <= jlid && i >= j) {
592c4762a1bSJed Brown     /* don't advect in the lid */
593c4762a1bSJed Brown     return residual;
594c4762a1bSJed Brown   } else if (i < j) {
595c4762a1bSJed Brown     /* beneath the slab sfc */
596c4762a1bSJed Brown     uW = uE = param->cb;
597c4762a1bSJed Brown     wS = wN = param->sb;
598c4762a1bSJed Brown   } else {
599c4762a1bSJed Brown     /* advect in the slab and wedge */
6009371c9d4SSatish Balay     uW = x[j][i - 1].u;
6019371c9d4SSatish Balay     uE = x[j][i].u;
6029371c9d4SSatish Balay     wS = x[j - 1][i].w;
6039371c9d4SSatish Balay     wN = x[j][i].w;
604c4762a1bSJed Brown   }
605c4762a1bSJed Brown 
606c4762a1bSJed Brown   if (param->adv_scheme == ADVECT_FV || i == ilim - 1 || j == jlim - 1 || i == 1 || j == 1) {
607c4762a1bSJed Brown     /* finite volume advection */
608c4762a1bSJed Brown     TS = (x[j][i].T + x[j - 1][i].T) / 2.0;
609c4762a1bSJed Brown     TN = (x[j][i].T + x[j + 1][i].T) / 2.0;
610c4762a1bSJed Brown     TE = (x[j][i].T + x[j][i + 1].T) / 2.0;
611c4762a1bSJed Brown     TW = (x[j][i].T + x[j][i - 1].T) / 2.0;
6129371c9d4SSatish Balay     fN = wN * TN * dx;
6139371c9d4SSatish Balay     fS = wS * TS * dx;
6149371c9d4SSatish Balay     fE = uE * TE * dz;
6159371c9d4SSatish Balay     fW = uW * TW * dz;
616c4762a1bSJed Brown 
617c4762a1bSJed Brown   } else {
618c4762a1bSJed Brown     /* Fromm advection scheme */
6199371c9d4SSatish 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;
6209371c9d4SSatish 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;
6219371c9d4SSatish 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;
6229371c9d4SSatish 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;
623c4762a1bSJed Brown   }
624c4762a1bSJed Brown 
625c4762a1bSJed Brown   residual -= (fE - fW + fN - fS);
626c4762a1bSJed Brown 
627c4762a1bSJed Brown   return residual;
628c4762a1bSJed Brown }
629c4762a1bSJed Brown 
630c4762a1bSJed Brown /*  computes the shear stress---used on the boundaries */
631*f6dfbefdSBarry Smith static inline PetscScalar ShearStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user) {
632c4762a1bSJed Brown   Parameter  *param = user->param;
633c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
634c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1;
635c4762a1bSJed Brown   PetscScalar uN, uS, wE, wW;
636c4762a1bSJed Brown 
637c4762a1bSJed Brown   if (j <= grid->jlid || i < j || i == ilim || j == jlim) return EPS_ZERO;
638c4762a1bSJed Brown 
639c4762a1bSJed Brown   if (ipos == CELL_CENTER) { /* on cell center */
640c4762a1bSJed Brown 
641c4762a1bSJed Brown     wE = WInterp(x, i, j - 1);
642c4762a1bSJed Brown     if (i == j) {
643c4762a1bSJed Brown       wW = param->sb;
644c4762a1bSJed Brown       uN = param->cb;
645c4762a1bSJed Brown     } else {
646c4762a1bSJed Brown       wW = WInterp(x, i - 1, j - 1);
647c4762a1bSJed Brown       uN = UInterp(x, i - 1, j);
648c4762a1bSJed Brown     }
649c4762a1bSJed Brown     if (j == grid->jlid + 1) uS = 0.0;
650c4762a1bSJed Brown     else uS = UInterp(x, i - 1, j - 1);
651c4762a1bSJed Brown 
652c4762a1bSJed Brown   } else { /* on cell corner */
653c4762a1bSJed Brown 
6549371c9d4SSatish Balay     uN = x[j + 1][i].u;
6559371c9d4SSatish Balay     uS = x[j][i].u;
6569371c9d4SSatish Balay     wW = x[j][i].w;
6579371c9d4SSatish Balay     wE = x[j][i + 1].w;
658c4762a1bSJed Brown   }
659c4762a1bSJed Brown 
660c4762a1bSJed Brown   return (uN - uS) / grid->dz + (wE - wW) / grid->dx;
661c4762a1bSJed Brown }
662c4762a1bSJed Brown 
663c4762a1bSJed Brown /*  computes the normal stress---used on the boundaries */
664*f6dfbefdSBarry Smith static inline PetscScalar XNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user) {
665c4762a1bSJed Brown   Parameter  *param = user->param;
666c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
667c4762a1bSJed Brown   PetscScalar dx = grid->dx, dz = grid->dz;
668c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1, ivisc;
669c4762a1bSJed Brown   PetscScalar epsC = 0.0, etaC, TC, uE, uW, pC, z_scale;
670c4762a1bSJed Brown   if (i < j || j <= grid->jlid) return EPS_ZERO;
671c4762a1bSJed Brown 
6729371c9d4SSatish Balay   ivisc   = param->ivisc;
6739371c9d4SSatish Balay   z_scale = param->z_scale;
674c4762a1bSJed Brown 
675c4762a1bSJed Brown   if (ipos == CELL_CENTER) { /* on cell center */
676c4762a1bSJed Brown 
677c4762a1bSJed Brown     TC = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
678c4762a1bSJed Brown     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CENTER, user);
679c4762a1bSJed Brown     etaC = Viscosity(TC, epsC, dz * j, param);
680c4762a1bSJed Brown 
6819371c9d4SSatish Balay     uW = x[j][i - 1].u;
6829371c9d4SSatish Balay     uE = x[j][i].u;
683c4762a1bSJed Brown     pC = x[j][i].p;
684c4762a1bSJed Brown 
685c4762a1bSJed Brown   } else { /* on cell corner */
686c4762a1bSJed Brown     if (i == ilim || j == jlim) return EPS_ZERO;
687c4762a1bSJed Brown 
688c4762a1bSJed Brown     TC = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
689c4762a1bSJed Brown     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CORNER, user);
690c4762a1bSJed Brown     etaC = Viscosity(TC, epsC, dz * (j + 0.5), param);
691c4762a1bSJed Brown 
692c4762a1bSJed Brown     if (i == j) uW = param->sb;
693c4762a1bSJed Brown     else uW = UInterp(x, i - 1, j);
6949371c9d4SSatish Balay     uE = UInterp(x, i, j);
6959371c9d4SSatish Balay     pC = PInterp(x, i, j);
696c4762a1bSJed Brown   }
697c4762a1bSJed Brown 
698c4762a1bSJed Brown   return 2.0 * etaC * (uE - uW) / dx - pC;
699c4762a1bSJed Brown }
700c4762a1bSJed Brown 
701c4762a1bSJed Brown /*  computes the normal stress---used on the boundaries */
702*f6dfbefdSBarry Smith static inline PetscScalar ZNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user) {
703c4762a1bSJed Brown   Parameter  *param = user->param;
704c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
705c4762a1bSJed Brown   PetscScalar dz    = grid->dz;
706c4762a1bSJed Brown   PetscInt    ilim = grid->ni - 1, jlim = grid->nj - 1, ivisc;
707c4762a1bSJed Brown   PetscScalar epsC = 0.0, etaC, TC;
708c4762a1bSJed Brown   PetscScalar pC, wN, wS, z_scale;
709c4762a1bSJed Brown   if (i < j || j <= grid->jlid) return EPS_ZERO;
710c4762a1bSJed Brown 
7119371c9d4SSatish Balay   ivisc   = param->ivisc;
7129371c9d4SSatish Balay   z_scale = param->z_scale;
713c4762a1bSJed Brown 
714c4762a1bSJed Brown   if (ipos == CELL_CENTER) { /* on cell center */
715c4762a1bSJed Brown 
716c4762a1bSJed Brown     TC = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
717c4762a1bSJed Brown     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CENTER, user);
718c4762a1bSJed Brown     etaC = Viscosity(TC, epsC, dz * j, param);
7199371c9d4SSatish Balay     wN   = x[j][i].w;
7209371c9d4SSatish Balay     wS   = x[j - 1][i].w;
7219371c9d4SSatish Balay     pC   = x[j][i].p;
722c4762a1bSJed Brown 
723c4762a1bSJed Brown   } else { /* on cell corner */
724c4762a1bSJed Brown     if ((i == ilim) || (j == jlim)) return EPS_ZERO;
725c4762a1bSJed Brown 
726c4762a1bSJed Brown     TC = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
727c4762a1bSJed Brown     if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CORNER, user);
728c4762a1bSJed Brown     etaC = Viscosity(TC, epsC, dz * (j + 0.5), param);
729c4762a1bSJed Brown     if (i == j) wN = param->sb;
730c4762a1bSJed Brown     else wN = WInterp(x, i, j);
7319371c9d4SSatish Balay     wS = WInterp(x, i, j - 1);
7329371c9d4SSatish Balay     pC = PInterp(x, i, j);
733c4762a1bSJed Brown   }
734c4762a1bSJed Brown 
735c4762a1bSJed Brown   return 2.0 * etaC * (wN - wS) / dz - pC;
736c4762a1bSJed Brown }
737c4762a1bSJed Brown 
738c4762a1bSJed Brown /*=====================================================================
739c4762a1bSJed Brown   INITIALIZATION, POST-PROCESSING AND OUTPUT FUNCTIONS
740c4762a1bSJed Brown   =====================================================================*/
741c4762a1bSJed Brown 
742c4762a1bSJed Brown /* initializes the problem parameters and checks for
743c4762a1bSJed Brown    command line changes */
744*f6dfbefdSBarry Smith PetscErrorCode SetParams(Parameter *param, GridInfo *grid) {
745c4762a1bSJed Brown   PetscReal SEC_PER_YR                     = 3600.00 * 24.00 * 365.2500;
746c4762a1bSJed Brown   PetscReal alpha_g_on_cp_units_inverse_km = 4.0e-5 * 9.8;
747c4762a1bSJed Brown 
748c4762a1bSJed Brown   /* domain geometry */
749c4762a1bSJed Brown   param->slab_dip    = 45.0;
750c4762a1bSJed Brown   param->width       = 320.0; /* km */
751c4762a1bSJed Brown   param->depth       = 300.0; /* km */
752c4762a1bSJed Brown   param->lid_depth   = 35.0;  /* km */
753c4762a1bSJed Brown   param->fault_depth = 35.0;  /* km */
754c4762a1bSJed Brown 
7559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_dip", &(param->slab_dip), NULL));
7569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-width", &(param->width), NULL));
7579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-depth", &(param->depth), NULL));
7589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-lid_depth", &(param->lid_depth), NULL));
7599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-fault_depth", &(param->fault_depth), NULL));
760c4762a1bSJed Brown 
761c4762a1bSJed Brown   param->slab_dip = param->slab_dip * PETSC_PI / 180.0; /* radians */
762c4762a1bSJed Brown 
763c4762a1bSJed Brown   /* grid information */
7649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-jfault", &(grid->jfault), NULL));
765c4762a1bSJed Brown   grid->ni = 82;
7669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ni", &(grid->ni), NULL));
767c4762a1bSJed Brown 
768c4762a1bSJed Brown   grid->dx     = param->width / ((PetscReal)(grid->ni - 2)); /* km */
769c4762a1bSJed Brown   grid->dz     = grid->dx * PetscTanReal(param->slab_dip);   /* km */
770c4762a1bSJed Brown   grid->nj     = (PetscInt)(param->depth / grid->dz + 3.0);  /* gridpoints*/
771c4762a1bSJed Brown   param->depth = grid->dz * (grid->nj - 2);                  /* km */
772c4762a1bSJed Brown   grid->inose  = 0;                                          /* gridpoints*/
7739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-inose", &(grid->inose), NULL));
774c4762a1bSJed Brown   grid->bx            = DM_BOUNDARY_NONE;
775c4762a1bSJed Brown   grid->by            = DM_BOUNDARY_NONE;
776c4762a1bSJed Brown   grid->stencil       = DMDA_STENCIL_BOX;
777c4762a1bSJed Brown   grid->dof           = 4;
778c4762a1bSJed Brown   grid->stencil_width = 2;
779c4762a1bSJed Brown   grid->mglevels      = 1;
780c4762a1bSJed Brown 
781c4762a1bSJed Brown   /* boundary conditions */
782c4762a1bSJed Brown   param->pv_analytic = PETSC_FALSE;
783c4762a1bSJed Brown   param->ibound      = BC_NOSTRESS;
7849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ibound", &(param->ibound), NULL));
785c4762a1bSJed Brown 
786c4762a1bSJed Brown   /* physical constants */
787c4762a1bSJed Brown   param->slab_velocity = 5.0;       /* cm/yr */
788c4762a1bSJed Brown   param->slab_age      = 50.0;      /* Ma */
789c4762a1bSJed Brown   param->lid_age       = 50.0;      /* Ma */
790c4762a1bSJed Brown   param->kappa         = 0.7272e-6; /* m^2/sec */
791c4762a1bSJed Brown   param->potentialT    = 1300.0;    /* degrees C */
792c4762a1bSJed Brown 
7939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_velocity", &(param->slab_velocity), NULL));
7949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_age", &(param->slab_age), NULL));
7959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-lid_age", &(param->lid_age), NULL));
7969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", &(param->kappa), NULL));
7979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-potentialT", &(param->potentialT), NULL));
798c4762a1bSJed Brown 
799c4762a1bSJed Brown   /* viscosity */
800c4762a1bSJed Brown   param->ivisc        = 3;    /* 0=isovisc, 1=difn creep, 2=disl creep, 3=full */
801c4762a1bSJed Brown   param->eta0         = 1e24; /* Pa-s */
802c4762a1bSJed Brown   param->visc_cutoff  = 0.0;  /* factor of eta_0 */
803c4762a1bSJed Brown   param->continuation = 1.0;
804c4762a1bSJed Brown 
805c4762a1bSJed Brown   /* constants for diffusion creep */
806c4762a1bSJed Brown   param->diffusion.A     = 1.8e7; /* Pa-s */
807c4762a1bSJed Brown   param->diffusion.n     = 1.0;   /* dim'less */
808c4762a1bSJed Brown   param->diffusion.Estar = 375e3; /* J/mol */
809c4762a1bSJed Brown   param->diffusion.Vstar = 5e-6;  /* m^3/mol */
810c4762a1bSJed Brown 
811c4762a1bSJed Brown   /* constants for param->dislocationocation creep */
812c4762a1bSJed Brown   param->dislocation.A     = 2.8969e4; /* Pa-s */
813c4762a1bSJed Brown   param->dislocation.n     = 3.5;      /* dim'less */
814c4762a1bSJed Brown   param->dislocation.Estar = 530e3;    /* J/mol */
815c4762a1bSJed Brown   param->dislocation.Vstar = 14e-6;    /* m^3/mol */
816c4762a1bSJed Brown 
8179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ivisc", &(param->ivisc), NULL));
8189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-visc_cutoff", &(param->visc_cutoff), NULL));
819c4762a1bSJed Brown 
820c4762a1bSJed Brown   param->output_ivisc = param->ivisc;
821c4762a1bSJed Brown 
8229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-output_ivisc", &(param->output_ivisc), NULL));
8239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-vstar", &(param->dislocation.Vstar), NULL));
824c4762a1bSJed Brown 
825c4762a1bSJed Brown   /* output options */
826c4762a1bSJed Brown   param->quiet      = PETSC_FALSE;
827c4762a1bSJed Brown   param->param_test = PETSC_FALSE;
828c4762a1bSJed Brown 
8299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-quiet", &(param->quiet)));
8309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test", &(param->param_test)));
8319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-file", param->filename, sizeof(param->filename), &(param->output_to_file)));
832c4762a1bSJed Brown 
833c4762a1bSJed Brown   /* advection */
834c4762a1bSJed Brown   param->adv_scheme = ADVECT_FROMM; /* advection scheme: 0=finite vol, 1=Fromm */
835c4762a1bSJed Brown 
8369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-adv_scheme", &(param->adv_scheme), NULL));
837c4762a1bSJed Brown 
838c4762a1bSJed Brown   /* misc. flags */
839c4762a1bSJed Brown   param->stop_solve    = PETSC_FALSE;
840c4762a1bSJed Brown   param->interrupted   = PETSC_FALSE;
841c4762a1bSJed Brown   param->kspmon        = PETSC_FALSE;
842c4762a1bSJed Brown   param->toggle_kspmon = PETSC_FALSE;
843c4762a1bSJed Brown 
844c4762a1bSJed Brown   /* derived parameters for slab angle */
845c4762a1bSJed Brown   param->sb = PetscSinReal(param->slab_dip);
846c4762a1bSJed Brown   param->cb = PetscCosReal(param->slab_dip);
847c4762a1bSJed Brown   param->c  = param->slab_dip * param->sb / (param->slab_dip * param->slab_dip - param->sb * param->sb);
848c4762a1bSJed Brown   param->d  = (param->slab_dip * param->cb - param->sb) / (param->slab_dip * param->slab_dip - param->sb * param->sb);
849c4762a1bSJed Brown 
850c4762a1bSJed Brown   /* length, velocity and time scale for non-dimensionalization */
851c4762a1bSJed Brown   param->L = PetscMin(param->width, param->depth);      /* km */
852c4762a1bSJed Brown   param->V = param->slab_velocity / 100.0 / SEC_PER_YR; /* m/sec */
853c4762a1bSJed Brown 
854c4762a1bSJed Brown   /* other unit conversions and derived parameters */
855c4762a1bSJed Brown   param->scaled_width = param->width / param->L;                   /* dim'less */
856c4762a1bSJed Brown   param->scaled_depth = param->depth / param->L;                   /* dim'less */
857c4762a1bSJed Brown   param->lid_depth    = param->lid_depth / param->L;               /* dim'less */
858c4762a1bSJed Brown   param->fault_depth  = param->fault_depth / param->L;             /* dim'less */
859c4762a1bSJed Brown   grid->dx            = grid->dx / param->L;                       /* dim'less */
860c4762a1bSJed Brown   grid->dz            = grid->dz / param->L;                       /* dim'less */
861c4762a1bSJed Brown   grid->jlid          = (PetscInt)(param->lid_depth / grid->dz);   /* gridcells */
862c4762a1bSJed Brown   grid->jfault        = (PetscInt)(param->fault_depth / grid->dz); /* gridcells */
863c4762a1bSJed Brown   param->lid_depth    = grid->jlid * grid->dz;                     /* dim'less */
864c4762a1bSJed Brown   param->fault_depth  = grid->jfault * grid->dz;                   /* dim'less */
865c4762a1bSJed Brown   grid->corner        = grid->jlid + 1;                            /* gridcells */
866c4762a1bSJed Brown   param->peclet       = param->V                                   /* m/sec */
867c4762a1bSJed Brown                 * param->L * 1000.0                                /* m */
868c4762a1bSJed Brown                 / param->kappa;                                    /* m^2/sec */
869c4762a1bSJed Brown   param->z_scale = param->L * alpha_g_on_cp_units_inverse_km;
870c4762a1bSJed Brown   param->skt     = PetscSqrtReal(param->kappa * param->slab_age * SEC_PER_YR);
8719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-peclet", &(param->peclet), NULL));
872c4762a1bSJed Brown 
8735f80ce2aSJacob Faibussowitsch   return 0;
874c4762a1bSJed Brown }
875c4762a1bSJed Brown 
876c4762a1bSJed Brown /*  prints a report of the problem parameters to stdout */
877*f6dfbefdSBarry Smith PetscErrorCode ReportParams(Parameter *param, GridInfo *grid) {
878c4762a1bSJed Brown   char date[30];
879c4762a1bSJed Brown 
8809566063dSJacob Faibussowitsch   PetscCall(PetscGetDate(date, 30));
881c4762a1bSJed Brown 
882c4762a1bSJed Brown   if (!(param->quiet)) {
8839566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---------------------BEGIN ex30 PARAM REPORT-------------------\n"));
8849566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Domain: \n"));
8859566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Width = %g km,         Depth = %g km\n", (double)param->width, (double)param->depth));
8869566063dSJacob 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));
8879566063dSJacob 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)));
888c4762a1bSJed Brown 
8899566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nGrid: \n"));
89063a3b9bcSJacob 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)));
89163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  jlid = %3" PetscInt_FMT "              jfault = %3" PetscInt_FMT " \n", grid->jlid, grid->jfault));
8929566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  Pe = %g\n", (double)param->peclet));
893c4762a1bSJed Brown 
8949566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nRheology:"));
895c4762a1bSJed Brown     if (param->ivisc == VISC_CONST) {
8969566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Isoviscous \n"));
89748a46eb9SPierre Jolivet       if (param->pv_analytic) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Pressure and Velocity prescribed! \n"));
898c4762a1bSJed Brown     } else if (param->ivisc == VISC_DIFN) {
8999566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Diffusion Creep (T-Dependent Newtonian) \n"));
9009566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
901c4762a1bSJed Brown     } else if (param->ivisc == VISC_DISL) {
9029566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Dislocation Creep (T-Dependent Non-Newtonian) \n"));
9039566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
904c4762a1bSJed Brown     } else if (param->ivisc == VISC_FULL) {
9059566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Full Rheology \n"));
9069566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
907c4762a1bSJed Brown     } else {
9089566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                 Invalid! \n"));
9095f80ce2aSJacob Faibussowitsch       return 1;
910c4762a1bSJed Brown     }
911c4762a1bSJed Brown 
9129566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Boundary condition:"));
913c4762a1bSJed Brown     if (param->ibound == BC_ANALYTIC) {
9149566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Isoviscous Analytic Dirichlet \n"));
915c4762a1bSJed Brown     } else if (param->ibound == BC_NOSTRESS) {
9169566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Stress-Free (normal & shear stress)\n"));
917c4762a1bSJed Brown     } else if (param->ibound == BC_EXPERMNT) {
9189566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Experimental boundary condition \n"));
919c4762a1bSJed Brown     } else {
9209566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "       Invalid! \n"));
9215f80ce2aSJacob Faibussowitsch       return 1;
922c4762a1bSJed Brown     }
923c4762a1bSJed Brown 
92460d4fc61SSatish Balay     if (param->output_to_file) {
925c4762a1bSJed Brown #if defined(PETSC_HAVE_MATLAB_ENGINE)
9269566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output Destination:       Mat file \"%s\"\n", param->filename));
927c4762a1bSJed Brown #else
9289566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output Destination:       PETSc binary file \"%s\"\n", param->filename));
929c4762a1bSJed Brown #endif
93060d4fc61SSatish Balay     }
93148a46eb9SPierre Jolivet     if (param->output_ivisc != param->ivisc) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "                          Output viscosity: -ivisc %" PetscInt_FMT "\n", param->output_ivisc));
932c4762a1bSJed Brown 
9339566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---------------------END ex30 PARAM REPORT---------------------\n"));
934c4762a1bSJed Brown   }
935c4762a1bSJed Brown   if (param->param_test) PetscEnd();
9365f80ce2aSJacob Faibussowitsch   return 0;
937c4762a1bSJed Brown }
938c4762a1bSJed Brown 
939c4762a1bSJed Brown /* ------------------------------------------------------------------- */
940a5b23f4aSJose E. Roman /*  generates an initial guess using the analytic solution for isoviscous
941c4762a1bSJed Brown     corner flow */
942c4762a1bSJed Brown PetscErrorCode Initialize(DM da)
943c4762a1bSJed Brown /* ------------------------------------------------------------------- */
944c4762a1bSJed Brown {
945c4762a1bSJed Brown   AppCtx    *user;
946c4762a1bSJed Brown   Parameter *param;
947c4762a1bSJed Brown   GridInfo  *grid;
948c4762a1bSJed Brown   PetscInt   i, j, is, js, im, jm;
949c4762a1bSJed Brown   Field    **x;
950c4762a1bSJed Brown   Vec        Xguess;
951c4762a1bSJed Brown 
952c4762a1bSJed Brown   /* Get the fine grid */
9539566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &user));
954c4762a1bSJed Brown   Xguess = user->Xguess;
955c4762a1bSJed Brown   param  = user->param;
956c4762a1bSJed Brown   grid   = user->grid;
9579566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
9589566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xguess, (void **)&x));
959c4762a1bSJed Brown 
960c4762a1bSJed Brown   /* Compute initial guess */
961c4762a1bSJed Brown   for (j = js; j < js + jm; j++) {
962c4762a1bSJed Brown     for (i = is; i < is + im; i++) {
963c4762a1bSJed Brown       if (i < j) x[j][i].u = param->cb;
964c4762a1bSJed Brown       else if (j <= grid->jlid) x[j][i].u = 0.0;
965c4762a1bSJed Brown       else x[j][i].u = HorizVelocity(i, j, user);
966c4762a1bSJed Brown 
967c4762a1bSJed Brown       if (i <= j) x[j][i].w = param->sb;
968c4762a1bSJed Brown       else if (j <= grid->jlid) x[j][i].w = 0.0;
969c4762a1bSJed Brown       else x[j][i].w = VertVelocity(i, j, user);
970c4762a1bSJed Brown 
971c4762a1bSJed Brown       if (i < j || j <= grid->jlid) x[j][i].p = 0.0;
972c4762a1bSJed Brown       else x[j][i].p = Pressure(i, j, user);
973c4762a1bSJed Brown 
974c4762a1bSJed Brown       x[j][i].T = PetscMin(grid->dz * (j - 0.5), 1.0);
975c4762a1bSJed Brown     }
976c4762a1bSJed Brown   }
977c4762a1bSJed Brown 
978c4762a1bSJed Brown   /* Restore x to Xguess */
9799566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xguess, (void **)&x));
980c4762a1bSJed Brown 
981c4762a1bSJed Brown   return 0;
982c4762a1bSJed Brown }
983c4762a1bSJed Brown 
984c4762a1bSJed Brown /*  controls output to a file */
985*f6dfbefdSBarry Smith PetscErrorCode DoOutput(SNES snes, PetscInt its) {
986c4762a1bSJed Brown   AppCtx     *user;
987c4762a1bSJed Brown   Parameter  *param;
988c4762a1bSJed Brown   GridInfo   *grid;
989c4762a1bSJed Brown   PetscInt    ivt;
990c4762a1bSJed Brown   PetscMPIInt rank;
991c4762a1bSJed Brown   PetscViewer viewer;
992c4762a1bSJed Brown   Vec         res, pars;
993c4762a1bSJed Brown   MPI_Comm    comm;
994c4762a1bSJed Brown   DM          da;
995c4762a1bSJed Brown 
9969566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
9979566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &user));
998c4762a1bSJed Brown   param = user->param;
999c4762a1bSJed Brown   grid  = user->grid;
1000c4762a1bSJed Brown   ivt   = param->ivisc;
1001c4762a1bSJed Brown 
1002c4762a1bSJed Brown   param->ivisc = param->output_ivisc;
1003c4762a1bSJed Brown 
1004c4762a1bSJed Brown   /* compute final residual and final viscosity/strain rate fields */
10059566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &res, NULL, NULL));
10069566063dSJacob Faibussowitsch   PetscCall(ViscosityField(da, user->x, user->Xguess));
1007c4762a1bSJed Brown 
1008c4762a1bSJed Brown   /* get the communicator and the rank of the processor */
10099566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
10109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1011c4762a1bSJed Brown 
1012c4762a1bSJed Brown   if (param->output_to_file) { /* send output to binary file */
10139566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, &pars));
1014dd400576SPatrick Sanan     if (rank == 0) { /* on processor 0 */
10159566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(pars, 20, PETSC_DETERMINE));
10169566063dSJacob Faibussowitsch       PetscCall(VecSetFromOptions(pars));
10179566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 0, (PetscScalar)(grid->ni), INSERT_VALUES));
10189566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 1, (PetscScalar)(grid->nj), INSERT_VALUES));
10199566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 2, (PetscScalar)(grid->dx), INSERT_VALUES));
10209566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 3, (PetscScalar)(grid->dz), INSERT_VALUES));
10219566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 4, (PetscScalar)(param->L), INSERT_VALUES));
10229566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 5, (PetscScalar)(param->V), INSERT_VALUES));
1023c4762a1bSJed Brown       /* skipped 6 intentionally */
10249566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 7, (PetscScalar)(param->slab_dip), INSERT_VALUES));
10259566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 8, (PetscScalar)(grid->jlid), INSERT_VALUES));
10269566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 9, (PetscScalar)(param->lid_depth), INSERT_VALUES));
10279566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 10, (PetscScalar)(grid->jfault), INSERT_VALUES));
10289566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 11, (PetscScalar)(param->fault_depth), INSERT_VALUES));
10299566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 12, (PetscScalar)(param->potentialT), INSERT_VALUES));
10309566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 13, (PetscScalar)(param->ivisc), INSERT_VALUES));
10319566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 14, (PetscScalar)(param->visc_cutoff), INSERT_VALUES));
10329566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 15, (PetscScalar)(param->ibound), INSERT_VALUES));
10339566063dSJacob Faibussowitsch       PetscCall(VecSetValue(pars, 16, (PetscScalar)(its), INSERT_VALUES));
1034c4762a1bSJed Brown     } else { /* on some other processor */
10359566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(pars, 0, PETSC_DETERMINE));
10369566063dSJacob Faibussowitsch       PetscCall(VecSetFromOptions(pars));
1037c4762a1bSJed Brown     }
10389371c9d4SSatish Balay     PetscCall(VecAssemblyBegin(pars));
10399371c9d4SSatish Balay     PetscCall(VecAssemblyEnd(pars));
1040c4762a1bSJed Brown 
1041c4762a1bSJed Brown     /* create viewer */
1042c4762a1bSJed Brown #if defined(PETSC_HAVE_MATLAB_ENGINE)
10439566063dSJacob Faibussowitsch     PetscCall(PetscViewerMatlabOpen(PETSC_COMM_WORLD, param->filename, FILE_MODE_WRITE, &viewer));
1044c4762a1bSJed Brown #else
10459566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, param->filename, FILE_MODE_WRITE, &viewer));
1046c4762a1bSJed Brown #endif
1047c4762a1bSJed Brown 
1048c4762a1bSJed Brown     /* send vectors to viewer */
10499566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)res, "res"));
10509566063dSJacob Faibussowitsch     PetscCall(VecView(res, viewer));
10519566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)user->x, "out"));
10529566063dSJacob Faibussowitsch     PetscCall(VecView(user->x, viewer));
10539566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)(user->Xguess), "aux"));
10549566063dSJacob Faibussowitsch     PetscCall(VecView(user->Xguess, viewer));
10559566063dSJacob Faibussowitsch     PetscCall(StressField(da)); /* compute stress fields */
10569566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)(user->Xguess), "str"));
10579566063dSJacob Faibussowitsch     PetscCall(VecView(user->Xguess, viewer));
10589566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)pars, "par"));
10599566063dSJacob Faibussowitsch     PetscCall(VecView(pars, viewer));
1060c4762a1bSJed Brown 
1061c4762a1bSJed Brown     /* destroy viewer and vector */
10629566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
10639566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&pars));
1064c4762a1bSJed Brown   }
1065c4762a1bSJed Brown 
1066c4762a1bSJed Brown   param->ivisc = ivt;
1067c4762a1bSJed Brown   return 0;
1068c4762a1bSJed Brown }
1069c4762a1bSJed Brown 
1070c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1071c4762a1bSJed Brown /* Compute both the second invariant of the strain rate tensor and the viscosity, at both cell centers and cell corners */
1072c4762a1bSJed Brown PetscErrorCode ViscosityField(DM da, Vec X, Vec V)
1073c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1074c4762a1bSJed Brown {
1075c4762a1bSJed Brown   AppCtx    *user;
1076c4762a1bSJed Brown   Parameter *param;
1077c4762a1bSJed Brown   GridInfo  *grid;
1078c4762a1bSJed Brown   Vec        localX;
1079c4762a1bSJed Brown   Field    **v, **x;
1080c4762a1bSJed Brown   PetscReal  eps, /* dx,*/ dz, T, epsC, TC;
1081c4762a1bSJed Brown   PetscInt   i, j, is, js, im, jm, ilim, jlim, ivt;
1082c4762a1bSJed Brown 
1083c4762a1bSJed Brown   PetscFunctionBeginUser;
10849566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &user));
1085c4762a1bSJed Brown   param        = user->param;
1086c4762a1bSJed Brown   grid         = user->grid;
1087c4762a1bSJed Brown   ivt          = param->ivisc;
1088c4762a1bSJed Brown   param->ivisc = param->output_ivisc;
1089c4762a1bSJed Brown 
10909566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
10919566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
10929566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
10939566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, (void **)&x));
10949566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, V, (void **)&v));
1095c4762a1bSJed Brown 
1096c4762a1bSJed Brown   /* Parameters */
1097c4762a1bSJed Brown   /* dx = grid->dx; */ dz = grid->dz;
1098c4762a1bSJed Brown 
10999371c9d4SSatish Balay   ilim = grid->ni - 1;
11009371c9d4SSatish Balay   jlim = grid->nj - 1;
1101c4762a1bSJed Brown 
1102c4762a1bSJed Brown   /* Compute real temperature, strain rate and viscosity */
11039566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
1104c4762a1bSJed Brown   for (j = js; j < js + jm; j++) {
1105c4762a1bSJed Brown     for (i = is; i < is + im; i++) {
1106c4762a1bSJed Brown       T = PetscRealPart(param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * param->z_scale));
1107c4762a1bSJed Brown       if (i < ilim && j < jlim) {
1108c4762a1bSJed Brown         TC = PetscRealPart(param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * param->z_scale));
1109c4762a1bSJed Brown       } else {
1110c4762a1bSJed Brown         TC = T;
1111c4762a1bSJed Brown       }
1112c4762a1bSJed Brown       eps  = PetscRealPart((CalcSecInv(x, i, j, CELL_CENTER, user)));
1113c4762a1bSJed Brown       epsC = PetscRealPart(CalcSecInv(x, i, j, CELL_CORNER, user));
1114c4762a1bSJed Brown 
1115c4762a1bSJed Brown       v[j][i].u = eps;
1116c4762a1bSJed Brown       v[j][i].w = epsC;
1117c4762a1bSJed Brown       v[j][i].p = Viscosity(T, eps, dz * (j - 0.5), param);
1118c4762a1bSJed Brown       v[j][i].T = Viscosity(TC, epsC, dz * j, param);
1119c4762a1bSJed Brown     }
1120c4762a1bSJed Brown   }
11219566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, V, (void **)&v));
11229566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, (void **)&x));
11239566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
1124c4762a1bSJed Brown 
1125c4762a1bSJed Brown   param->ivisc = ivt;
1126c4762a1bSJed Brown   PetscFunctionReturn(0);
1127c4762a1bSJed Brown }
1128c4762a1bSJed Brown 
1129c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1130c4762a1bSJed Brown /* post-processing: compute stress everywhere */
1131c4762a1bSJed Brown PetscErrorCode StressField(DM da)
1132c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1133c4762a1bSJed Brown {
1134c4762a1bSJed Brown   AppCtx  *user;
1135c4762a1bSJed Brown   PetscInt i, j, is, js, im, jm;
1136c4762a1bSJed Brown   Vec      locVec;
1137c4762a1bSJed Brown   Field  **x, **y;
1138c4762a1bSJed Brown 
11399566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &user));
1140c4762a1bSJed Brown 
1141c4762a1bSJed Brown   /* Get the fine grid of Xguess and X */
11429566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
11439566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, user->Xguess, (void **)&x));
1144c4762a1bSJed Brown 
11459566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &locVec));
11469566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, user->x, INSERT_VALUES, locVec));
11479566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, user->x, INSERT_VALUES, locVec));
11489566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, locVec, (void **)&y));
1149c4762a1bSJed Brown 
1150c4762a1bSJed Brown   /* Compute stress on the corner points */
1151c4762a1bSJed Brown   for (j = js; j < js + jm; j++) {
1152c4762a1bSJed Brown     for (i = is; i < is + im; i++) {
1153c4762a1bSJed Brown       x[j][i].u = ShearStress(y, i, j, CELL_CENTER, user);
1154c4762a1bSJed Brown       x[j][i].w = ShearStress(y, i, j, CELL_CORNER, user);
1155c4762a1bSJed Brown       x[j][i].p = XNormalStress(y, i, j, CELL_CENTER, user);
1156c4762a1bSJed Brown       x[j][i].T = ZNormalStress(y, i, j, CELL_CENTER, user);
1157c4762a1bSJed Brown     }
1158c4762a1bSJed Brown   }
1159c4762a1bSJed Brown 
1160c4762a1bSJed Brown   /* Restore the fine grid of Xguess and X */
11619566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, user->Xguess, (void **)&x));
11629566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, locVec, (void **)&y));
11639566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &locVec));
1164c4762a1bSJed Brown   return 0;
1165c4762a1bSJed Brown }
1166c4762a1bSJed Brown 
1167c4762a1bSJed Brown /*=====================================================================
1168c4762a1bSJed Brown   UTILITY FUNCTIONS
1169c4762a1bSJed Brown   =====================================================================*/
1170c4762a1bSJed Brown 
1171*f6dfbefdSBarry Smith /* returns the velocity of the subducting slab and handles fault nodes for BC */
1172*f6dfbefdSBarry Smith static inline PetscScalar SlabVel(char c, PetscInt i, PetscInt j, AppCtx *user) {
1173c4762a1bSJed Brown   Parameter *param = user->param;
1174c4762a1bSJed Brown   GridInfo  *grid  = user->grid;
1175c4762a1bSJed Brown 
1176c4762a1bSJed Brown   if (c == 'U' || c == 'u') {
1177c4762a1bSJed Brown     if (i < j - 1) return param->cb;
1178c4762a1bSJed Brown     else if (j <= grid->jfault) return 0.0;
1179c4762a1bSJed Brown     else return param->cb;
1180c4762a1bSJed Brown 
1181c4762a1bSJed Brown   } else {
1182c4762a1bSJed Brown     if (i < j) return param->sb;
1183c4762a1bSJed Brown     else if (j <= grid->jfault) return 0.0;
1184c4762a1bSJed Brown     else return param->sb;
1185c4762a1bSJed Brown   }
1186c4762a1bSJed Brown }
1187c4762a1bSJed Brown 
1188c4762a1bSJed Brown /*  solution to diffusive half-space cooling model for BC */
1189*f6dfbefdSBarry Smith static inline PetscScalar PlateModel(PetscInt j, PetscInt plate, AppCtx *user) {
1190c4762a1bSJed Brown   Parameter  *param = user->param;
1191c4762a1bSJed Brown   PetscScalar z;
1192c4762a1bSJed Brown   if (plate == PLATE_LID) z = (j - 0.5) * user->grid->dz;
1193c4762a1bSJed Brown   else z = (j - 0.5) * user->grid->dz * param->cb; /* PLATE_SLAB */
1194c4762a1bSJed Brown #if defined(PETSC_HAVE_ERF)
1195c4762a1bSJed Brown   return (PetscReal)(erf((double)PetscRealPart(z * param->L / 2.0 / param->skt)));
1196c4762a1bSJed Brown #else
1197c4762a1bSJed Brown   (*PetscErrorPrintf)("erf() not available on this machine\n");
1198c4762a1bSJed Brown   MPI_Abort(PETSC_COMM_SELF, 1);
1199c4762a1bSJed Brown #endif
1200c4762a1bSJed Brown }
1201c4762a1bSJed Brown 
1202c4762a1bSJed Brown /*=====================================================================
1203c4762a1bSJed Brown   INTERACTIVE SIGNAL HANDLING
1204c4762a1bSJed Brown   =====================================================================*/
1205c4762a1bSJed Brown 
1206c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1207c4762a1bSJed Brown PetscErrorCode SNESConverged_Interactive(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *ctx)
1208c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1209c4762a1bSJed Brown {
1210c4762a1bSJed Brown   AppCtx    *user  = (AppCtx *)ctx;
1211c4762a1bSJed Brown   Parameter *param = user->param;
1212c4762a1bSJed Brown   KSP        ksp;
1213c4762a1bSJed Brown 
1214c4762a1bSJed Brown   PetscFunctionBeginUser;
1215c4762a1bSJed Brown   if (param->interrupted) {
1216c4762a1bSJed Brown     param->interrupted = PETSC_FALSE;
12179566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: exiting SNES solve. \n"));
1218c4762a1bSJed Brown     *reason = SNES_CONVERGED_FNORM_ABS;
1219c4762a1bSJed Brown     PetscFunctionReturn(0);
1220c4762a1bSJed Brown   } else if (param->toggle_kspmon) {
1221c4762a1bSJed Brown     param->toggle_kspmon = PETSC_FALSE;
1222c4762a1bSJed Brown 
12239566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
1224c4762a1bSJed Brown 
1225c4762a1bSJed Brown     if (param->kspmon) {
12269566063dSJacob Faibussowitsch       PetscCall(KSPMonitorCancel(ksp));
1227c4762a1bSJed Brown 
1228c4762a1bSJed Brown       param->kspmon = PETSC_FALSE;
12299566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: deactivating ksp singular value monitor. \n"));
1230c4762a1bSJed Brown     } else {
1231c4762a1bSJed Brown       PetscViewerAndFormat *vf;
12329566063dSJacob Faibussowitsch       PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf));
12339566063dSJacob Faibussowitsch       PetscCall(KSPMonitorSet(ksp, (PetscErrorCode(*)(KSP, PetscInt, PetscReal, void *))KSPMonitorSingularValue, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
1234c4762a1bSJed Brown 
1235c4762a1bSJed Brown       param->kspmon = PETSC_TRUE;
12369566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: activating ksp singular value monitor. \n"));
1237c4762a1bSJed Brown     }
1238c4762a1bSJed Brown   }
123911cc89d2SBarry Smith   PetscCall(SNESConvergedDefault(snes, it, xnorm, snorm, fnorm, reason, ctx));
124011cc89d2SBarry Smith   PetscFunctionReturn(0);
1241c4762a1bSJed Brown }
1242c4762a1bSJed Brown 
1243c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1244c4762a1bSJed Brown #include <signal.h>
1245c4762a1bSJed Brown PetscErrorCode InteractiveHandler(int signum, void *ctx)
1246c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1247c4762a1bSJed Brown {
1248c4762a1bSJed Brown   AppCtx    *user  = (AppCtx *)ctx;
1249c4762a1bSJed Brown   Parameter *param = user->param;
1250c4762a1bSJed Brown 
1251c4762a1bSJed Brown   if (signum == SIGILL) {
1252c4762a1bSJed Brown     param->toggle_kspmon = PETSC_TRUE;
1253c4762a1bSJed Brown #if !defined(PETSC_MISSING_SIGCONT)
1254c4762a1bSJed Brown   } else if (signum == SIGCONT) {
1255c4762a1bSJed Brown     param->interrupted = PETSC_TRUE;
1256c4762a1bSJed Brown #endif
1257c4762a1bSJed Brown #if !defined(PETSC_MISSING_SIGURG)
1258c4762a1bSJed Brown   } else if (signum == SIGURG) {
1259c4762a1bSJed Brown     param->stop_solve = PETSC_TRUE;
1260c4762a1bSJed Brown #endif
1261c4762a1bSJed Brown   }
1262c4762a1bSJed Brown   return 0;
1263c4762a1bSJed Brown }
1264c4762a1bSJed Brown 
1265*f6dfbefdSBarry Smith /*  main call-back function that computes the processor-local piece of the residual */
1266*f6dfbefdSBarry Smith PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field **x, Field **f, void *ptr) {
1267c4762a1bSJed Brown   AppCtx     *user  = (AppCtx *)ptr;
1268c4762a1bSJed Brown   Parameter  *param = user->param;
1269c4762a1bSJed Brown   GridInfo   *grid  = user->grid;
1270c4762a1bSJed Brown   PetscScalar mag_w, mag_u;
1271c4762a1bSJed Brown   PetscInt    i, j, mx, mz, ilim, jlim;
1272c4762a1bSJed Brown   PetscInt    is, ie, js, je, ibound; /* ,ivisc */
1273c4762a1bSJed Brown 
1274c4762a1bSJed Brown   PetscFunctionBeginUser;
1275c4762a1bSJed Brown   /* Define global and local grid parameters */
12769371c9d4SSatish Balay   mx   = info->mx;
12779371c9d4SSatish Balay   mz   = info->my;
12789371c9d4SSatish Balay   ilim = mx - 1;
12799371c9d4SSatish Balay   jlim = mz - 1;
12809371c9d4SSatish Balay   is   = info->xs;
12819371c9d4SSatish Balay   ie   = info->xs + info->xm;
12829371c9d4SSatish Balay   js   = info->ys;
12839371c9d4SSatish Balay   je   = info->ys + info->ym;
1284c4762a1bSJed Brown 
1285c4762a1bSJed Brown   /* Define geometric and numeric parameters */
1286c4762a1bSJed Brown   /* ivisc = param->ivisc; */ ibound = param->ibound;
1287c4762a1bSJed Brown 
1288c4762a1bSJed Brown   for (j = js; j < je; j++) {
1289c4762a1bSJed Brown     for (i = is; i < ie; i++) {
1290c4762a1bSJed Brown       /************* X-MOMENTUM/VELOCITY *************/
1291c4762a1bSJed Brown       if (i < j) f[j][i].u = x[j][i].u - SlabVel('U', i, j, user);
1292c4762a1bSJed Brown       else if (j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1293c4762a1bSJed Brown         /* in the lithospheric lid */
1294c4762a1bSJed Brown         f[j][i].u = x[j][i].u - 0.0;
1295c4762a1bSJed Brown       } else if (i == ilim) {
1296c4762a1bSJed Brown         /* on the right side boundary */
1297c4762a1bSJed Brown         if (ibound == BC_ANALYTIC) {
1298c4762a1bSJed Brown           f[j][i].u = x[j][i].u - HorizVelocity(i, j, user);
1299c4762a1bSJed Brown         } else {
1300c4762a1bSJed Brown           f[j][i].u = XNormalStress(x, i, j, CELL_CENTER, user) - EPS_ZERO;
1301c4762a1bSJed Brown         }
1302c4762a1bSJed Brown 
1303c4762a1bSJed Brown       } else if (j == jlim) {
1304c4762a1bSJed Brown         /* on the bottom boundary */
1305c4762a1bSJed Brown         if (ibound == BC_ANALYTIC) {
1306c4762a1bSJed Brown           f[j][i].u = x[j][i].u - HorizVelocity(i, j, user);
1307c4762a1bSJed Brown         } else if (ibound == BC_NOSTRESS) {
1308c4762a1bSJed Brown           f[j][i].u = XMomentumResidual(x, i, j, user);
1309c4762a1bSJed Brown         } else {
1310c4762a1bSJed Brown           /* experimental boundary condition */
1311c4762a1bSJed Brown         }
1312c4762a1bSJed Brown 
1313c4762a1bSJed Brown       } else {
1314c4762a1bSJed Brown         /* in the mantle wedge */
1315c4762a1bSJed Brown         f[j][i].u = XMomentumResidual(x, i, j, user);
1316c4762a1bSJed Brown       }
1317c4762a1bSJed Brown 
1318c4762a1bSJed Brown       /************* Z-MOMENTUM/VELOCITY *************/
1319c4762a1bSJed Brown       if (i <= j) {
1320c4762a1bSJed Brown         f[j][i].w = x[j][i].w - SlabVel('W', i, j, user);
1321c4762a1bSJed Brown 
1322c4762a1bSJed Brown       } else if (j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1323c4762a1bSJed Brown         /* in the lithospheric lid */
1324c4762a1bSJed Brown         f[j][i].w = x[j][i].w - 0.0;
1325c4762a1bSJed Brown 
1326c4762a1bSJed Brown       } else if (j == jlim) {
1327c4762a1bSJed Brown         /* on the bottom boundary */
1328c4762a1bSJed Brown         if (ibound == BC_ANALYTIC) {
1329c4762a1bSJed Brown           f[j][i].w = x[j][i].w - VertVelocity(i, j, user);
1330c4762a1bSJed Brown         } else {
1331c4762a1bSJed Brown           f[j][i].w = ZNormalStress(x, i, j, CELL_CENTER, user) - EPS_ZERO;
1332c4762a1bSJed Brown         }
1333c4762a1bSJed Brown 
1334c4762a1bSJed Brown       } else if (i == ilim) {
1335c4762a1bSJed Brown         /* on the right side boundary */
1336c4762a1bSJed Brown         if (ibound == BC_ANALYTIC) {
1337c4762a1bSJed Brown           f[j][i].w = x[j][i].w - VertVelocity(i, j, user);
1338c4762a1bSJed Brown         } else if (ibound == BC_NOSTRESS) {
1339c4762a1bSJed Brown           f[j][i].w = ZMomentumResidual(x, i, j, user);
1340c4762a1bSJed Brown         } else {
1341c4762a1bSJed Brown           /* experimental boundary condition */
1342c4762a1bSJed Brown         }
1343c4762a1bSJed Brown 
1344c4762a1bSJed Brown       } else {
1345c4762a1bSJed Brown         /* in the mantle wedge */
1346c4762a1bSJed Brown         f[j][i].w = ZMomentumResidual(x, i, j, user);
1347c4762a1bSJed Brown       }
1348c4762a1bSJed Brown 
1349c4762a1bSJed Brown       /************* CONTINUITY/PRESSURE *************/
1350c4762a1bSJed Brown       if (i < j || j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1351c4762a1bSJed Brown         /* in the lid or slab */
1352c4762a1bSJed Brown         f[j][i].p = x[j][i].p;
1353c4762a1bSJed Brown 
1354c4762a1bSJed Brown       } else if ((i == ilim || j == jlim) && ibound == BC_ANALYTIC) {
1355c4762a1bSJed Brown         /* on an analytic boundary */
1356c4762a1bSJed Brown         f[j][i].p = x[j][i].p - Pressure(i, j, user);
1357c4762a1bSJed Brown 
1358c4762a1bSJed Brown       } else {
1359c4762a1bSJed Brown         /* in the mantle wedge */
1360c4762a1bSJed Brown         f[j][i].p = ContinuityResidual(x, i, j, user);
1361c4762a1bSJed Brown       }
1362c4762a1bSJed Brown 
1363c4762a1bSJed Brown       /************* TEMPERATURE *************/
1364c4762a1bSJed Brown       if (j == 0) {
1365c4762a1bSJed Brown         /* on the surface */
1366c4762a1bSJed Brown         f[j][i].T = x[j][i].T + x[j + 1][i].T + PetscMax(PetscRealPart(x[j][i].T), 0.0);
1367c4762a1bSJed Brown 
1368c4762a1bSJed Brown       } else if (i == 0) {
1369c4762a1bSJed Brown         /* slab inflow boundary */
1370c4762a1bSJed Brown         f[j][i].T = x[j][i].T - PlateModel(j, PLATE_SLAB, user);
1371c4762a1bSJed Brown 
1372c4762a1bSJed Brown       } else if (i == ilim) {
1373c4762a1bSJed Brown         /* right side boundary */
1374c4762a1bSJed Brown         mag_u     = 1.0 - PetscPowRealInt((1.0 - PetscMax(PetscMin(PetscRealPart(x[j][i - 1].u) / param->cb, 1.0), 0.0)), 5);
1375c4762a1bSJed 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);
1376c4762a1bSJed Brown 
1377c4762a1bSJed Brown       } else if (j == jlim) {
1378c4762a1bSJed Brown         /* bottom boundary */
1379c4762a1bSJed Brown         mag_w     = 1.0 - PetscPowRealInt((1.0 - PetscMax(PetscMin(PetscRealPart(x[j - 1][i].w) / param->sb, 1.0), 0.0)), 5);
1380c4762a1bSJed Brown         f[j][i].T = x[j][i].T - mag_w * x[j - 1][i - 1].T - (1.0 - mag_w);
1381c4762a1bSJed Brown 
1382c4762a1bSJed Brown       } else {
1383c4762a1bSJed Brown         /* in the mantle wedge */
1384c4762a1bSJed Brown         f[j][i].T = EnergyResidual(x, i, j, user);
1385c4762a1bSJed Brown       }
1386c4762a1bSJed Brown     }
1387c4762a1bSJed Brown   }
1388c4762a1bSJed Brown   PetscFunctionReturn(0);
1389c4762a1bSJed Brown }
1390c4762a1bSJed Brown 
1391c4762a1bSJed Brown /*TEST
1392c4762a1bSJed Brown 
1393c4762a1bSJed Brown    build:
1394c4762a1bSJed Brown       requires: !complex erf
1395c4762a1bSJed Brown 
1396c4762a1bSJed Brown    test:
1397c4762a1bSJed Brown       args: -ni 18
1398c4762a1bSJed Brown       filter: grep -v Destination
1399c4762a1bSJed Brown       requires: !single
1400c4762a1bSJed Brown 
1401c4762a1bSJed Brown TEST*/
1402