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
81da81f932SPierre Jolivet typedef struct { /* physical and miscellaneous 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
main(int argc,char ** argv)123d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
124d71ae5a4SJacob Faibussowitsch {
125c4762a1bSJed Brown SNES snes;
126c4762a1bSJed Brown AppCtx *user; /* user-defined work context */
127c4762a1bSJed Brown Parameter param;
128c4762a1bSJed Brown GridInfo grid;
129c4762a1bSJed Brown PetscInt nits;
130c4762a1bSJed Brown MPI_Comm comm;
131c4762a1bSJed Brown DM da;
132c4762a1bSJed Brown
133327415f7SBarry Smith PetscFunctionBeginUser;
134c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1353ba16761SJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-file", "ex30_output"));
1363ba16761SJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-snes_monitor_short", NULL));
1373ba16761SJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-snes_max_it", "20"));
1383ba16761SJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-ksp_max_it", "1500"));
1393ba16761SJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-ksp_gmres_restart", "300"));
1403ba16761SJacob Faibussowitsch PetscCall(PetscOptionsInsert(NULL, &argc, &argv, NULL));
141c4762a1bSJed Brown
142c4762a1bSJed Brown comm = PETSC_COMM_WORLD;
143c4762a1bSJed Brown
144c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown Set up the problem parameters.
146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1479566063dSJacob Faibussowitsch PetscCall(SetParams(¶m, &grid));
1489566063dSJacob Faibussowitsch PetscCall(ReportParams(¶m, &grid));
149c4762a1bSJed Brown
150c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1529566063dSJacob Faibussowitsch PetscCall(SNESCreate(comm, &snes));
1539566063dSJacob 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));
1549566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
1559566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
1569566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, da));
1579566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "x-velocity"));
1589566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "y-velocity"));
1599566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 2, "pressure"));
1609566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 3, "temperature"));
161c4762a1bSJed Brown
162c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163c4762a1bSJed Brown Create user context, set problem data, create vector data structures.
164c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1659566063dSJacob Faibussowitsch PetscCall(PetscNew(&user));
166c4762a1bSJed Brown user->param = ¶m;
167c4762a1bSJed Brown user->grid = &grid;
1689566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da, user));
169f4f49eeaSPierre Jolivet PetscCall(DMCreateGlobalVector(da, &user->Xguess));
170c4762a1bSJed Brown
171c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172c4762a1bSJed Brown Set up the SNES solver with callback functions.
173c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1749566063dSJacob Faibussowitsch PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode (*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, (void *)user));
1759566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
176c4762a1bSJed Brown
1779566063dSJacob Faibussowitsch PetscCall(SNESSetConvergenceTest(snes, SNESConverged_Interactive, (void *)user, NULL));
1789566063dSJacob Faibussowitsch PetscCall(PetscPushSignalHandler(InteractiveHandler, (void *)user));
179c4762a1bSJed Brown
180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181c4762a1bSJed Brown Initialize and solve the nonlinear system
182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1839566063dSJacob Faibussowitsch PetscCall(Initialize(da));
1849566063dSJacob Faibussowitsch PetscCall(UpdateSolution(snes, user, &nits));
185c4762a1bSJed Brown
186c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown Output variables.
188c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1899566063dSJacob Faibussowitsch PetscCall(DoOutput(snes, nits));
190c4762a1bSJed Brown
191c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192c4762a1bSJed Brown Free work space.
193c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->Xguess));
1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->x));
1969566063dSJacob Faibussowitsch PetscCall(PetscFree(user));
1979566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
1989566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
1999566063dSJacob Faibussowitsch PetscCall(PetscPopSignalHandler());
2009566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
201b122ec5aSJacob Faibussowitsch return 0;
202c4762a1bSJed Brown }
203c4762a1bSJed Brown
204c4762a1bSJed Brown /*=====================================================================
205c4762a1bSJed Brown PETSc INTERACTION FUNCTIONS (initialize & call SNESSolve)
206c4762a1bSJed Brown =====================================================================*/
207c4762a1bSJed Brown
208c4762a1bSJed Brown /* manages solve: adaptive continuation method */
UpdateSolution(SNES snes,AppCtx * user,PetscInt * nits)209d71ae5a4SJacob Faibussowitsch PetscErrorCode UpdateSolution(SNES snes, AppCtx *user, PetscInt *nits)
210d71ae5a4SJacob Faibussowitsch {
211c4762a1bSJed Brown KSP ksp;
212c4762a1bSJed Brown PC pc;
213c4762a1bSJed Brown SNESConvergedReason reason = SNES_CONVERGED_ITERATING;
214c4762a1bSJed Brown Parameter *param = user->param;
215c4762a1bSJed Brown PetscReal cont_incr = 0.3;
216c4762a1bSJed Brown PetscInt its;
217c4762a1bSJed Brown PetscBool q = PETSC_FALSE;
218c4762a1bSJed Brown DM dm;
219c4762a1bSJed Brown
220c4762a1bSJed Brown PetscFunctionBeginUser;
2219566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm));
2229566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &user->x));
2239566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp));
2249566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc));
2259566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE));
226c4762a1bSJed Brown
227c4762a1bSJed Brown *nits = 0;
228c4762a1bSJed Brown
229c4762a1bSJed Brown /* Isoviscous solve */
230c4762a1bSJed Brown if (param->ivisc == VISC_CONST && !param->stop_solve) {
231c4762a1bSJed Brown param->ivisc = VISC_CONST;
232c4762a1bSJed Brown
2339566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, 0, user->x));
2349566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes, &reason));
2359566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its));
236c4762a1bSJed Brown *nits += its;
2379566063dSJacob Faibussowitsch PetscCall(VecCopy(user->x, user->Xguess));
238c4762a1bSJed Brown if (param->stop_solve) goto done;
239c4762a1bSJed Brown }
240c4762a1bSJed Brown
241c4762a1bSJed Brown /* Olivine diffusion creep */
242c4762a1bSJed Brown if (param->ivisc >= VISC_DIFN && !param->stop_solve) {
2439566063dSJacob Faibussowitsch if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Computing Variable Viscosity Solution\n"));
244c4762a1bSJed Brown
245c4762a1bSJed Brown /* continuation method on viscosity cutoff */
246c4762a1bSJed Brown for (param->continuation = 0.0;; param->continuation += cont_incr) {
2479566063dSJacob Faibussowitsch if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Continuation parameter = %g\n", (double)param->continuation));
248c4762a1bSJed Brown
249c4762a1bSJed Brown /* solve the non-linear system */
2509566063dSJacob Faibussowitsch PetscCall(VecCopy(user->Xguess, user->x));
2519566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, 0, user->x));
2529566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes, &reason));
2539566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its));
254c4762a1bSJed Brown *nits += its;
25563a3b9bcSJacob Faibussowitsch if (!q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " SNES iterations: %" PetscInt_FMT ", Cumulative: %" PetscInt_FMT "\n", its, *nits));
256c4762a1bSJed Brown if (param->stop_solve) goto done;
257c4762a1bSJed Brown
258c4762a1bSJed Brown if (reason < 0) {
259c4762a1bSJed Brown /* NOT converged */
260c4762a1bSJed Brown cont_incr = -PetscAbsReal(cont_incr) / 2.0;
261c4762a1bSJed Brown if (PetscAbsReal(cont_incr) < 0.01) goto done;
262c4762a1bSJed Brown
263c4762a1bSJed Brown } else {
264c4762a1bSJed Brown /* converged */
2659566063dSJacob Faibussowitsch PetscCall(VecCopy(user->x, user->Xguess));
266c4762a1bSJed Brown if (param->continuation >= 1.0) goto done;
267c4762a1bSJed Brown if (its <= 3) cont_incr = 0.30001;
268c4762a1bSJed Brown else if (its <= 8) cont_incr = 0.15001;
269c4762a1bSJed Brown else cont_incr = 0.10001;
270c4762a1bSJed Brown
271c4762a1bSJed Brown if (param->continuation + cont_incr > 1.0) cont_incr = 1.0 - param->continuation;
272c4762a1bSJed Brown } /* endif reason<0 */
273c4762a1bSJed Brown }
274c4762a1bSJed Brown }
275c4762a1bSJed Brown done:
2769566063dSJacob Faibussowitsch if (param->stop_solve && !q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: stopping solve.\n"));
2779566063dSJacob Faibussowitsch if (reason < 0 && !q) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "FAILED TO CONVERGE: stopping solve.\n"));
2783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
279c4762a1bSJed Brown }
280c4762a1bSJed Brown
281c4762a1bSJed Brown /*=====================================================================
282c4762a1bSJed Brown PHYSICS FUNCTIONS (compute the discrete residual)
283c4762a1bSJed Brown =====================================================================*/
284c4762a1bSJed Brown
UInterp(Field ** x,PetscInt i,PetscInt j)285d71ae5a4SJacob Faibussowitsch static inline PetscScalar UInterp(Field **x, PetscInt i, PetscInt j)
286d71ae5a4SJacob Faibussowitsch {
287c4762a1bSJed Brown return 0.25 * (x[j][i].u + x[j + 1][i].u + x[j][i + 1].u + x[j + 1][i + 1].u);
288c4762a1bSJed Brown }
289c4762a1bSJed Brown
WInterp(Field ** x,PetscInt i,PetscInt j)290d71ae5a4SJacob Faibussowitsch static inline PetscScalar WInterp(Field **x, PetscInt i, PetscInt j)
291d71ae5a4SJacob Faibussowitsch {
292c4762a1bSJed Brown return 0.25 * (x[j][i].w + x[j + 1][i].w + x[j][i + 1].w + x[j + 1][i + 1].w);
293c4762a1bSJed Brown }
294c4762a1bSJed Brown
PInterp(Field ** x,PetscInt i,PetscInt j)295d71ae5a4SJacob Faibussowitsch static inline PetscScalar PInterp(Field **x, PetscInt i, PetscInt j)
296d71ae5a4SJacob Faibussowitsch {
297c4762a1bSJed Brown return 0.25 * (x[j][i].p + x[j + 1][i].p + x[j][i + 1].p + x[j + 1][i + 1].p);
298c4762a1bSJed Brown }
299c4762a1bSJed Brown
TInterp(Field ** x,PetscInt i,PetscInt j)300d71ae5a4SJacob Faibussowitsch static inline PetscScalar TInterp(Field **x, PetscInt i, PetscInt j)
301d71ae5a4SJacob Faibussowitsch {
302c4762a1bSJed Brown return 0.25 * (x[j][i].T + x[j + 1][i].T + x[j][i + 1].T + x[j + 1][i + 1].T);
303c4762a1bSJed Brown }
304c4762a1bSJed Brown
305c4762a1bSJed Brown /* isoviscous analytic solution for IC */
HorizVelocity(PetscInt i,PetscInt j,AppCtx * user)306d71ae5a4SJacob Faibussowitsch static inline PetscScalar HorizVelocity(PetscInt i, PetscInt j, AppCtx *user)
307d71ae5a4SJacob Faibussowitsch {
308c4762a1bSJed Brown Parameter *param = user->param;
309c4762a1bSJed Brown GridInfo *grid = user->grid;
310c4762a1bSJed Brown PetscScalar st, ct, th, c = param->c, d = param->d;
311c4762a1bSJed Brown PetscReal x, z, r;
312c4762a1bSJed Brown
3139371c9d4SSatish Balay x = (i - grid->jlid) * grid->dx;
3149371c9d4SSatish Balay z = (j - grid->jlid - 0.5) * grid->dz;
315c4762a1bSJed Brown r = PetscSqrtReal(x * x + z * z);
316c4762a1bSJed Brown st = z / r;
317c4762a1bSJed Brown ct = x / r;
318c4762a1bSJed Brown th = PetscAtanReal(z / x);
319c4762a1bSJed Brown return ct * (c * th * st + d * (st + th * ct)) + st * (c * (st - th * ct) + d * th * st);
320c4762a1bSJed Brown }
321c4762a1bSJed Brown
322c4762a1bSJed Brown /* isoviscous analytic solution for IC */
VertVelocity(PetscInt i,PetscInt j,AppCtx * user)3239fbee547SJacob Faibussowitsch static inline PetscScalar VertVelocity(PetscInt i, PetscInt j, AppCtx *user)
324c4762a1bSJed Brown {
325c4762a1bSJed Brown Parameter *param = user->param;
326c4762a1bSJed Brown GridInfo *grid = user->grid;
327c4762a1bSJed Brown PetscScalar st, ct, th, c = param->c, d = param->d;
328c4762a1bSJed Brown PetscReal x, z, r;
329c4762a1bSJed Brown
3309371c9d4SSatish Balay x = (i - grid->jlid - 0.5) * grid->dx;
3319371c9d4SSatish Balay z = (j - grid->jlid) * grid->dz;
3329371c9d4SSatish Balay r = PetscSqrtReal(x * x + z * z);
3339371c9d4SSatish Balay st = z / r;
3349371c9d4SSatish Balay ct = x / r;
3359371c9d4SSatish Balay th = PetscAtanReal(z / x);
336c4762a1bSJed Brown return st * (c * th * st + d * (st + th * ct)) - ct * (c * (st - th * ct) + d * th * st);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown
339c4762a1bSJed Brown /* isoviscous analytic solution for IC */
Pressure(PetscInt i,PetscInt j,AppCtx * user)340d71ae5a4SJacob Faibussowitsch static inline PetscScalar Pressure(PetscInt i, PetscInt j, AppCtx *user)
341d71ae5a4SJacob Faibussowitsch {
342c4762a1bSJed Brown Parameter *param = user->param;
343c4762a1bSJed Brown GridInfo *grid = user->grid;
344c4762a1bSJed Brown PetscScalar x, z, r, st, ct, c = param->c, d = param->d;
345c4762a1bSJed Brown
3469371c9d4SSatish Balay x = (i - grid->jlid - 0.5) * grid->dx;
3479371c9d4SSatish Balay z = (j - grid->jlid - 0.5) * grid->dz;
3489371c9d4SSatish Balay r = PetscSqrtReal(x * x + z * z);
3499371c9d4SSatish Balay st = z / r;
3509371c9d4SSatish Balay ct = x / r;
3514ad8454bSPierre Jolivet return -2.0 * (c * ct - d * st) / r;
352c4762a1bSJed Brown }
353c4762a1bSJed Brown
354c4762a1bSJed Brown /* computes the second invariant of the strain rate tensor */
CalcSecInv(Field ** x,PetscInt i,PetscInt j,PetscInt ipos,AppCtx * user)355d71ae5a4SJacob Faibussowitsch static inline PetscScalar CalcSecInv(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
356d71ae5a4SJacob Faibussowitsch {
357c4762a1bSJed Brown Parameter *param = user->param;
358c4762a1bSJed Brown GridInfo *grid = user->grid;
359c4762a1bSJed Brown PetscInt ilim = grid->ni - 1, jlim = grid->nj - 1;
360c4762a1bSJed Brown PetscScalar uN, uS, uE, uW, wN, wS, wE, wW;
361c4762a1bSJed Brown PetscScalar eps11, eps12, eps22;
362c4762a1bSJed Brown
363c4762a1bSJed Brown if (i < j) return EPS_ZERO;
364c4762a1bSJed Brown if (i == ilim) i--;
365c4762a1bSJed Brown if (j == jlim) j--;
366c4762a1bSJed Brown
367c4762a1bSJed Brown if (ipos == CELL_CENTER) { /* on cell center */
368c4762a1bSJed Brown if (j <= grid->jlid) return EPS_ZERO;
369c4762a1bSJed Brown
3709371c9d4SSatish Balay uE = x[j][i].u;
3719371c9d4SSatish Balay uW = x[j][i - 1].u;
3729371c9d4SSatish Balay wN = x[j][i].w;
3739371c9d4SSatish Balay wS = x[j - 1][i].w;
374c4762a1bSJed Brown wE = WInterp(x, i, j - 1);
375c4762a1bSJed Brown if (i == j) {
3769371c9d4SSatish Balay uN = param->cb;
3779371c9d4SSatish Balay wW = param->sb;
378c4762a1bSJed Brown } else {
3799371c9d4SSatish Balay uN = UInterp(x, i - 1, j);
3809371c9d4SSatish Balay wW = WInterp(x, i - 1, j - 1);
381c4762a1bSJed Brown }
382c4762a1bSJed Brown
383c4762a1bSJed Brown if (j == grid->jlid + 1) uS = 0.0;
384c4762a1bSJed Brown else uS = UInterp(x, i - 1, j - 1);
385c4762a1bSJed Brown
386c4762a1bSJed Brown } else { /* on CELL_CORNER */
387c4762a1bSJed Brown if (j < grid->jlid) return EPS_ZERO;
388c4762a1bSJed Brown
3899371c9d4SSatish Balay uN = x[j + 1][i].u;
3909371c9d4SSatish Balay uS = x[j][i].u;
3919371c9d4SSatish Balay wE = x[j][i + 1].w;
3929371c9d4SSatish Balay wW = x[j][i].w;
393c4762a1bSJed Brown if (i == j) {
394c4762a1bSJed Brown wN = param->sb;
395c4762a1bSJed Brown uW = param->cb;
396c4762a1bSJed Brown } else {
397c4762a1bSJed Brown wN = WInterp(x, i, j);
398c4762a1bSJed Brown uW = UInterp(x, i - 1, j);
399c4762a1bSJed Brown }
400c4762a1bSJed Brown
401c4762a1bSJed Brown if (j == grid->jlid) {
4029371c9d4SSatish Balay uE = 0.0;
4039371c9d4SSatish Balay uW = 0.0;
404c4762a1bSJed Brown uS = -uN;
405c4762a1bSJed Brown wS = -wN;
406c4762a1bSJed Brown } else {
407c4762a1bSJed Brown uE = UInterp(x, i, j);
408c4762a1bSJed Brown wS = WInterp(x, i, j - 1);
409c4762a1bSJed Brown }
410c4762a1bSJed Brown }
411c4762a1bSJed Brown
4129371c9d4SSatish Balay eps11 = (uE - uW) / grid->dx;
4139371c9d4SSatish Balay eps22 = (wN - wS) / grid->dz;
414c4762a1bSJed Brown eps12 = 0.5 * ((uN - uS) / grid->dz + (wE - wW) / grid->dx);
415c4762a1bSJed Brown
416c4762a1bSJed Brown return PetscSqrtReal(0.5 * (eps11 * eps11 + 2.0 * eps12 * eps12 + eps22 * eps22));
417c4762a1bSJed Brown }
418c4762a1bSJed Brown
419c4762a1bSJed Brown /* computes the shear viscosity */
Viscosity(PetscScalar T,PetscScalar eps,PetscScalar z,Parameter * param)420d71ae5a4SJacob Faibussowitsch static inline PetscScalar Viscosity(PetscScalar T, PetscScalar eps, PetscScalar z, Parameter *param)
421d71ae5a4SJacob Faibussowitsch {
422c4762a1bSJed Brown PetscReal result = 0.0;
423c4762a1bSJed Brown ViscParam difn = param->diffusion, disl = param->dislocation;
424c4762a1bSJed Brown PetscInt iVisc = param->ivisc;
425c4762a1bSJed Brown PetscScalar eps_scale = param->V / (param->L * 1000.0);
426c4762a1bSJed Brown PetscScalar strain_power, v1, v2, P;
427c4762a1bSJed Brown PetscScalar rho_g = 32340.0, R = 8.3144;
428c4762a1bSJed Brown
429c4762a1bSJed Brown P = rho_g * (z * param->L * 1000.0); /* Pa */
430c4762a1bSJed Brown
431c4762a1bSJed Brown if (iVisc == VISC_CONST) {
432c4762a1bSJed Brown /* constant viscosity */
433c4762a1bSJed Brown return 1.0;
434c4762a1bSJed Brown } else if (iVisc == VISC_DIFN) {
435c4762a1bSJed Brown /* diffusion creep rheology */
436f4f49eeaSPierre Jolivet result = PetscRealPart(difn.A * PetscExpScalar((difn.Estar + P * difn.Vstar) / R / (T + 273.0)) / param->eta0);
437c4762a1bSJed Brown } else if (iVisc == VISC_DISL) {
438c4762a1bSJed Brown /* dislocation creep rheology */
439c4762a1bSJed Brown strain_power = PetscPowScalar(eps * eps_scale, (1.0 - disl.n) / disl.n);
440c4762a1bSJed Brown
441c4762a1bSJed Brown result = PetscRealPart(disl.A * PetscExpScalar((disl.Estar + P * disl.Vstar) / disl.n / R / (T + 273.0)) * strain_power / param->eta0);
442c4762a1bSJed Brown } else if (iVisc == VISC_FULL) {
443c4762a1bSJed Brown /* dislocation/diffusion creep rheology */
444c4762a1bSJed Brown strain_power = PetscPowScalar(eps * eps_scale, (1.0 - disl.n) / disl.n);
445c4762a1bSJed Brown
446c4762a1bSJed Brown v1 = difn.A * PetscExpScalar((difn.Estar + P * difn.Vstar) / R / (T + 273.0)) / param->eta0;
447c4762a1bSJed Brown v2 = disl.A * PetscExpScalar((disl.Estar + P * disl.Vstar) / disl.n / R / (T + 273.0)) * strain_power / param->eta0;
448c4762a1bSJed Brown
449c4762a1bSJed Brown result = PetscRealPart(1.0 / (1.0 / v1 + 1.0 / v2));
450c4762a1bSJed Brown }
451c4762a1bSJed Brown
452c4762a1bSJed Brown /* max viscosity is param->eta0 */
453c4762a1bSJed Brown result = PetscMin(result, 1.0);
454c4762a1bSJed Brown /* min viscosity is param->visc_cutoff */
455c4762a1bSJed Brown result = PetscMax(result, param->visc_cutoff);
456c4762a1bSJed Brown /* continuation method */
457c4762a1bSJed Brown result = PetscPowReal(result, param->continuation);
458c4762a1bSJed Brown return result;
459c4762a1bSJed Brown }
460c4762a1bSJed Brown
461c4762a1bSJed Brown /* computes the residual of the x-component of eqn (1) above */
XMomentumResidual(Field ** x,PetscInt i,PetscInt j,AppCtx * user)462d71ae5a4SJacob Faibussowitsch static inline PetscScalar XMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
463d71ae5a4SJacob Faibussowitsch {
464c4762a1bSJed Brown Parameter *param = user->param;
465c4762a1bSJed Brown GridInfo *grid = user->grid;
466c4762a1bSJed Brown PetscScalar dx = grid->dx, dz = grid->dz;
467c4762a1bSJed Brown PetscScalar etaN, etaS, etaE, etaW, epsN = 0.0, epsS = 0.0, epsE = 0.0, epsW = 0.0;
468c4762a1bSJed Brown PetscScalar TE = 0.0, TN = 0.0, TS = 0.0, TW = 0.0, dPdx, residual, z_scale;
469c4762a1bSJed Brown PetscScalar dudxW, dudxE, dudzN, dudzS, dwdxN, dwdxS;
470c4762a1bSJed Brown PetscInt jlim = grid->nj - 1;
471c4762a1bSJed Brown
472c4762a1bSJed Brown z_scale = param->z_scale;
473c4762a1bSJed Brown
474c4762a1bSJed Brown if (param->ivisc == VISC_DIFN || param->ivisc >= VISC_DISL) { /* viscosity is T-dependent */
475c4762a1bSJed Brown TS = param->potentialT * TInterp(x, i, j - 1) * PetscExpScalar((j - 1.0) * dz * z_scale);
476c4762a1bSJed Brown if (j == jlim) TN = TS;
477c4762a1bSJed Brown else TN = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
478c4762a1bSJed Brown TW = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
479c4762a1bSJed Brown TE = param->potentialT * x[j][i + 1].T * PetscExpScalar((j - 0.5) * dz * z_scale);
480c4762a1bSJed Brown if (param->ivisc >= VISC_DISL) { /* olivine dislocation creep */
481c4762a1bSJed Brown epsN = CalcSecInv(x, i, j, CELL_CORNER, user);
482c4762a1bSJed Brown epsS = CalcSecInv(x, i, j - 1, CELL_CORNER, user);
483c4762a1bSJed Brown epsE = CalcSecInv(x, i + 1, j, CELL_CENTER, user);
484c4762a1bSJed Brown epsW = CalcSecInv(x, i, j, CELL_CENTER, user);
485c4762a1bSJed Brown }
486c4762a1bSJed Brown }
487c4762a1bSJed Brown etaN = Viscosity(TN, epsN, dz * (j + 0.5), param);
488c4762a1bSJed Brown etaS = Viscosity(TS, epsS, dz * (j - 0.5), param);
489c4762a1bSJed Brown etaW = Viscosity(TW, epsW, dz * j, param);
490c4762a1bSJed Brown etaE = Viscosity(TE, epsE, dz * j, param);
491c4762a1bSJed Brown
492c4762a1bSJed Brown dPdx = (x[j][i + 1].p - x[j][i].p) / dx;
493c4762a1bSJed Brown if (j == jlim) dudzN = etaN * (x[j][i].w - x[j][i + 1].w) / dx;
494c4762a1bSJed Brown else dudzN = etaN * (x[j + 1][i].u - x[j][i].u) / dz;
495c4762a1bSJed Brown dudzS = etaS * (x[j][i].u - x[j - 1][i].u) / dz;
496c4762a1bSJed Brown dudxE = etaE * (x[j][i + 1].u - x[j][i].u) / dx;
497c4762a1bSJed Brown dudxW = etaW * (x[j][i].u - x[j][i - 1].u) / dx;
498c4762a1bSJed Brown
499c4762a1bSJed Brown residual = -dPdx /* X-MOMENTUM EQUATION*/
5009371c9d4SSatish Balay + (dudxE - dudxW) / dx + (dudzN - dudzS) / dz;
501c4762a1bSJed Brown
502c4762a1bSJed Brown if (param->ivisc != VISC_CONST) {
503c4762a1bSJed Brown dwdxN = etaN * (x[j][i + 1].w - x[j][i].w) / dx;
504c4762a1bSJed Brown dwdxS = etaS * (x[j - 1][i + 1].w - x[j - 1][i].w) / dx;
505c4762a1bSJed Brown
506c4762a1bSJed Brown residual += (dudxE - dudxW) / dx + (dwdxN - dwdxS) / dz;
507c4762a1bSJed Brown }
508c4762a1bSJed Brown
509c4762a1bSJed Brown return residual;
510c4762a1bSJed Brown }
511c4762a1bSJed Brown
512c4762a1bSJed Brown /* computes the residual of the z-component of eqn (1) above */
ZMomentumResidual(Field ** x,PetscInt i,PetscInt j,AppCtx * user)5139fbee547SJacob Faibussowitsch static inline PetscScalar ZMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
514c4762a1bSJed Brown {
515c4762a1bSJed Brown Parameter *param = user->param;
516c4762a1bSJed Brown GridInfo *grid = user->grid;
517c4762a1bSJed Brown PetscScalar dx = grid->dx, dz = grid->dz;
518c4762a1bSJed 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;
519c4762a1bSJed Brown PetscScalar TE = 0.0, TN = 0.0, TS = 0.0, TW = 0.0, dPdz, residual, z_scale;
520c4762a1bSJed Brown PetscScalar dudzE, dudzW, dwdxW, dwdxE, dwdzN, dwdzS;
521c4762a1bSJed Brown PetscInt ilim = grid->ni - 1;
522c4762a1bSJed Brown
523c4762a1bSJed Brown /* geometric and other parameters */
524c4762a1bSJed Brown z_scale = param->z_scale;
525c4762a1bSJed Brown
526c4762a1bSJed Brown /* viscosity */
527c4762a1bSJed Brown if (param->ivisc == VISC_DIFN || param->ivisc >= VISC_DISL) { /* viscosity is T-dependent */
528c4762a1bSJed Brown TN = param->potentialT * x[j + 1][i].T * PetscExpScalar((j + 0.5) * dz * z_scale);
529c4762a1bSJed Brown TS = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
530c4762a1bSJed Brown TW = param->potentialT * TInterp(x, i - 1, j) * PetscExpScalar(j * dz * z_scale);
531c4762a1bSJed Brown if (i == ilim) TE = TW;
532c4762a1bSJed Brown else TE = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
533c4762a1bSJed Brown if (param->ivisc >= VISC_DISL) { /* olivine dislocation creep */
534c4762a1bSJed Brown epsN = CalcSecInv(x, i, j + 1, CELL_CENTER, user);
535c4762a1bSJed Brown epsS = CalcSecInv(x, i, j, CELL_CENTER, user);
536c4762a1bSJed Brown epsE = CalcSecInv(x, i, j, CELL_CORNER, user);
537c4762a1bSJed Brown epsW = CalcSecInv(x, i - 1, j, CELL_CORNER, user);
538c4762a1bSJed Brown }
539c4762a1bSJed Brown }
540c4762a1bSJed Brown etaN = Viscosity(TN, epsN, dz * (j + 1.0), param);
541c4762a1bSJed Brown etaS = Viscosity(TS, epsS, dz * (j + 0.0), param);
542c4762a1bSJed Brown etaW = Viscosity(TW, epsW, dz * (j + 0.5), param);
543c4762a1bSJed Brown etaE = Viscosity(TE, epsE, dz * (j + 0.5), param);
544c4762a1bSJed Brown
545c4762a1bSJed Brown dPdz = (x[j + 1][i].p - x[j][i].p) / dz;
546c4762a1bSJed Brown dwdzN = etaN * (x[j + 1][i].w - x[j][i].w) / dz;
547c4762a1bSJed Brown dwdzS = etaS * (x[j][i].w - x[j - 1][i].w) / dz;
548c4762a1bSJed Brown if (i == ilim) dwdxE = etaE * (x[j][i].u - x[j + 1][i].u) / dz;
549c4762a1bSJed Brown else dwdxE = etaE * (x[j][i + 1].w - x[j][i].w) / dx;
550c4762a1bSJed Brown dwdxW = 2.0 * etaW * (x[j][i].w - x[j][i - 1].w) / dx;
551c4762a1bSJed Brown
552c4762a1bSJed Brown /* Z-MOMENTUM */
553c4762a1bSJed Brown residual = -dPdz /* constant viscosity terms */
5549371c9d4SSatish Balay + (dwdzN - dwdzS) / dz + (dwdxE - dwdxW) / dx;
555c4762a1bSJed Brown
556c4762a1bSJed Brown if (param->ivisc != VISC_CONST) {
557c4762a1bSJed Brown dudzE = etaE * (x[j + 1][i].u - x[j][i].u) / dz;
558c4762a1bSJed Brown dudzW = etaW * (x[j + 1][i - 1].u - x[j][i - 1].u) / dz;
559c4762a1bSJed Brown
560c4762a1bSJed Brown residual += (dwdzN - dwdzS) / dz + (dudzE - dudzW) / dx;
561c4762a1bSJed Brown }
562c4762a1bSJed Brown
563c4762a1bSJed Brown return residual;
564c4762a1bSJed Brown }
565c4762a1bSJed Brown
566c4762a1bSJed Brown /* computes the residual of eqn (2) above */
ContinuityResidual(Field ** x,PetscInt i,PetscInt j,AppCtx * user)567d71ae5a4SJacob Faibussowitsch static inline PetscScalar ContinuityResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
568d71ae5a4SJacob Faibussowitsch {
569c4762a1bSJed Brown GridInfo *grid = user->grid;
570c4762a1bSJed Brown PetscScalar uE, uW, wN, wS, dudx, dwdz;
571c4762a1bSJed Brown
5729371c9d4SSatish Balay uW = x[j][i - 1].u;
5739371c9d4SSatish Balay uE = x[j][i].u;
5749371c9d4SSatish Balay dudx = (uE - uW) / grid->dx;
5759371c9d4SSatish Balay wS = x[j - 1][i].w;
5769371c9d4SSatish Balay wN = x[j][i].w;
5779371c9d4SSatish Balay dwdz = (wN - wS) / grid->dz;
578c4762a1bSJed Brown
579c4762a1bSJed Brown return dudx + dwdz;
580c4762a1bSJed Brown }
581c4762a1bSJed Brown
582c4762a1bSJed Brown /* computes the residual of eqn (3) above */
EnergyResidual(Field ** x,PetscInt i,PetscInt j,AppCtx * user)583d71ae5a4SJacob Faibussowitsch static inline PetscScalar EnergyResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
584d71ae5a4SJacob Faibussowitsch {
585c4762a1bSJed Brown Parameter *param = user->param;
586c4762a1bSJed Brown GridInfo *grid = user->grid;
587c4762a1bSJed Brown PetscScalar dx = grid->dx, dz = grid->dz;
588c4762a1bSJed Brown PetscInt ilim = grid->ni - 1, jlim = grid->nj - 1, jlid = grid->jlid;
589c4762a1bSJed Brown PetscScalar TE, TN, TS, TW, residual;
590c4762a1bSJed Brown PetscScalar uE, uW, wN, wS;
591c4762a1bSJed Brown PetscScalar fN, fS, fE, fW, dTdxW, dTdxE, dTdzN, dTdzS;
592c4762a1bSJed Brown
593c4762a1bSJed Brown dTdzN = (x[j + 1][i].T - x[j][i].T) / dz;
594c4762a1bSJed Brown dTdzS = (x[j][i].T - x[j - 1][i].T) / dz;
595c4762a1bSJed Brown dTdxE = (x[j][i + 1].T - x[j][i].T) / dx;
596c4762a1bSJed Brown dTdxW = (x[j][i].T - x[j][i - 1].T) / dx;
597c4762a1bSJed Brown
598c4762a1bSJed Brown residual = ((dTdzN - dTdzS) / dz + /* diffusion term */
5999371c9d4SSatish Balay (dTdxE - dTdxW) / dx) *
6009371c9d4SSatish Balay dx * dz / param->peclet;
601c4762a1bSJed Brown
602c4762a1bSJed Brown if (j <= jlid && i >= j) {
603c4762a1bSJed Brown /* don't advect in the lid */
604c4762a1bSJed Brown return residual;
605c4762a1bSJed Brown } else if (i < j) {
606c4762a1bSJed Brown /* beneath the slab sfc */
607c4762a1bSJed Brown uW = uE = param->cb;
608c4762a1bSJed Brown wS = wN = param->sb;
609c4762a1bSJed Brown } else {
610c4762a1bSJed Brown /* advect in the slab and wedge */
6119371c9d4SSatish Balay uW = x[j][i - 1].u;
6129371c9d4SSatish Balay uE = x[j][i].u;
6139371c9d4SSatish Balay wS = x[j - 1][i].w;
6149371c9d4SSatish Balay wN = x[j][i].w;
615c4762a1bSJed Brown }
616c4762a1bSJed Brown
617c4762a1bSJed Brown if (param->adv_scheme == ADVECT_FV || i == ilim - 1 || j == jlim - 1 || i == 1 || j == 1) {
618c4762a1bSJed Brown /* finite volume advection */
619c4762a1bSJed Brown TS = (x[j][i].T + x[j - 1][i].T) / 2.0;
620c4762a1bSJed Brown TN = (x[j][i].T + x[j + 1][i].T) / 2.0;
621c4762a1bSJed Brown TE = (x[j][i].T + x[j][i + 1].T) / 2.0;
622c4762a1bSJed Brown TW = (x[j][i].T + x[j][i - 1].T) / 2.0;
6239371c9d4SSatish Balay fN = wN * TN * dx;
6249371c9d4SSatish Balay fS = wS * TS * dx;
6259371c9d4SSatish Balay fE = uE * TE * dz;
6269371c9d4SSatish Balay fW = uW * TW * dz;
627c4762a1bSJed Brown
628c4762a1bSJed Brown } else {
629c4762a1bSJed Brown /* Fromm advection scheme */
6309371c9d4SSatish 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;
6319371c9d4SSatish 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;
6329371c9d4SSatish 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;
6339371c9d4SSatish 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;
634c4762a1bSJed Brown }
635c4762a1bSJed Brown
636c4762a1bSJed Brown residual -= (fE - fW + fN - fS);
637c4762a1bSJed Brown
638c4762a1bSJed Brown return residual;
639c4762a1bSJed Brown }
640c4762a1bSJed Brown
641c4762a1bSJed Brown /* computes the shear stress---used on the boundaries */
ShearStress(Field ** x,PetscInt i,PetscInt j,PetscInt ipos,AppCtx * user)642d71ae5a4SJacob Faibussowitsch static inline PetscScalar ShearStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
643d71ae5a4SJacob Faibussowitsch {
644c4762a1bSJed Brown Parameter *param = user->param;
645c4762a1bSJed Brown GridInfo *grid = user->grid;
646c4762a1bSJed Brown PetscInt ilim = grid->ni - 1, jlim = grid->nj - 1;
647c4762a1bSJed Brown PetscScalar uN, uS, wE, wW;
648c4762a1bSJed Brown
649c4762a1bSJed Brown if (j <= grid->jlid || i < j || i == ilim || j == jlim) return EPS_ZERO;
650c4762a1bSJed Brown
651c4762a1bSJed Brown if (ipos == CELL_CENTER) { /* on cell center */
652c4762a1bSJed Brown
653c4762a1bSJed Brown wE = WInterp(x, i, j - 1);
654c4762a1bSJed Brown if (i == j) {
655c4762a1bSJed Brown wW = param->sb;
656c4762a1bSJed Brown uN = param->cb;
657c4762a1bSJed Brown } else {
658c4762a1bSJed Brown wW = WInterp(x, i - 1, j - 1);
659c4762a1bSJed Brown uN = UInterp(x, i - 1, j);
660c4762a1bSJed Brown }
661c4762a1bSJed Brown if (j == grid->jlid + 1) uS = 0.0;
662c4762a1bSJed Brown else uS = UInterp(x, i - 1, j - 1);
663c4762a1bSJed Brown
664c4762a1bSJed Brown } else { /* on cell corner */
665c4762a1bSJed Brown
6669371c9d4SSatish Balay uN = x[j + 1][i].u;
6679371c9d4SSatish Balay uS = x[j][i].u;
6689371c9d4SSatish Balay wW = x[j][i].w;
6699371c9d4SSatish Balay wE = x[j][i + 1].w;
670c4762a1bSJed Brown }
671c4762a1bSJed Brown
672c4762a1bSJed Brown return (uN - uS) / grid->dz + (wE - wW) / grid->dx;
673c4762a1bSJed Brown }
674c4762a1bSJed Brown
675c4762a1bSJed Brown /* computes the normal stress---used on the boundaries */
XNormalStress(Field ** x,PetscInt i,PetscInt j,PetscInt ipos,AppCtx * user)676d71ae5a4SJacob Faibussowitsch static inline PetscScalar XNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
677d71ae5a4SJacob Faibussowitsch {
678c4762a1bSJed Brown Parameter *param = user->param;
679c4762a1bSJed Brown GridInfo *grid = user->grid;
680c4762a1bSJed Brown PetscScalar dx = grid->dx, dz = grid->dz;
681c4762a1bSJed Brown PetscInt ilim = grid->ni - 1, jlim = grid->nj - 1, ivisc;
682c4762a1bSJed Brown PetscScalar epsC = 0.0, etaC, TC, uE, uW, pC, z_scale;
683c4762a1bSJed Brown if (i < j || j <= grid->jlid) return EPS_ZERO;
684c4762a1bSJed Brown
6859371c9d4SSatish Balay ivisc = param->ivisc;
6869371c9d4SSatish Balay z_scale = param->z_scale;
687c4762a1bSJed Brown
688c4762a1bSJed Brown if (ipos == CELL_CENTER) { /* on cell center */
689c4762a1bSJed Brown
690c4762a1bSJed Brown TC = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
691c4762a1bSJed Brown if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CENTER, user);
692c4762a1bSJed Brown etaC = Viscosity(TC, epsC, dz * j, param);
693c4762a1bSJed Brown
6949371c9d4SSatish Balay uW = x[j][i - 1].u;
6959371c9d4SSatish Balay uE = x[j][i].u;
696c4762a1bSJed Brown pC = x[j][i].p;
697c4762a1bSJed Brown
698c4762a1bSJed Brown } else { /* on cell corner */
699c4762a1bSJed Brown if (i == ilim || j == jlim) return EPS_ZERO;
700c4762a1bSJed Brown
701c4762a1bSJed Brown TC = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
702c4762a1bSJed Brown if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CORNER, user);
703c4762a1bSJed Brown etaC = Viscosity(TC, epsC, dz * (j + 0.5), param);
704c4762a1bSJed Brown
705c4762a1bSJed Brown if (i == j) uW = param->sb;
706c4762a1bSJed Brown else uW = UInterp(x, i - 1, j);
7079371c9d4SSatish Balay uE = UInterp(x, i, j);
7089371c9d4SSatish Balay pC = PInterp(x, i, j);
709c4762a1bSJed Brown }
710c4762a1bSJed Brown
711c4762a1bSJed Brown return 2.0 * etaC * (uE - uW) / dx - pC;
712c4762a1bSJed Brown }
713c4762a1bSJed Brown
714c4762a1bSJed Brown /* computes the normal stress---used on the boundaries */
ZNormalStress(Field ** x,PetscInt i,PetscInt j,PetscInt ipos,AppCtx * user)715d71ae5a4SJacob Faibussowitsch static inline PetscScalar ZNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
716d71ae5a4SJacob Faibussowitsch {
717c4762a1bSJed Brown Parameter *param = user->param;
718c4762a1bSJed Brown GridInfo *grid = user->grid;
719c4762a1bSJed Brown PetscScalar dz = grid->dz;
720c4762a1bSJed Brown PetscInt ilim = grid->ni - 1, jlim = grid->nj - 1, ivisc;
721c4762a1bSJed Brown PetscScalar epsC = 0.0, etaC, TC;
722c4762a1bSJed Brown PetscScalar pC, wN, wS, z_scale;
723c4762a1bSJed Brown if (i < j || j <= grid->jlid) return EPS_ZERO;
724c4762a1bSJed Brown
7259371c9d4SSatish Balay ivisc = param->ivisc;
7269371c9d4SSatish Balay z_scale = param->z_scale;
727c4762a1bSJed Brown
728c4762a1bSJed Brown if (ipos == CELL_CENTER) { /* on cell center */
729c4762a1bSJed Brown
730c4762a1bSJed Brown TC = param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * z_scale);
731c4762a1bSJed Brown if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CENTER, user);
732c4762a1bSJed Brown etaC = Viscosity(TC, epsC, dz * j, param);
7339371c9d4SSatish Balay wN = x[j][i].w;
7349371c9d4SSatish Balay wS = x[j - 1][i].w;
7359371c9d4SSatish Balay pC = x[j][i].p;
736c4762a1bSJed Brown
737c4762a1bSJed Brown } else { /* on cell corner */
738c4762a1bSJed Brown if ((i == ilim) || (j == jlim)) return EPS_ZERO;
739c4762a1bSJed Brown
740c4762a1bSJed Brown TC = param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * z_scale);
741c4762a1bSJed Brown if (ivisc >= VISC_DISL) epsC = CalcSecInv(x, i, j, CELL_CORNER, user);
742c4762a1bSJed Brown etaC = Viscosity(TC, epsC, dz * (j + 0.5), param);
743c4762a1bSJed Brown if (i == j) wN = param->sb;
744c4762a1bSJed Brown else wN = WInterp(x, i, j);
7459371c9d4SSatish Balay wS = WInterp(x, i, j - 1);
7469371c9d4SSatish Balay pC = PInterp(x, i, j);
747c4762a1bSJed Brown }
748c4762a1bSJed Brown
749c4762a1bSJed Brown return 2.0 * etaC * (wN - wS) / dz - pC;
750c4762a1bSJed Brown }
751c4762a1bSJed Brown
752c4762a1bSJed Brown /*=====================================================================
753c4762a1bSJed Brown INITIALIZATION, POST-PROCESSING AND OUTPUT FUNCTIONS
754c4762a1bSJed Brown =====================================================================*/
755c4762a1bSJed Brown
756c4762a1bSJed Brown /* initializes the problem parameters and checks for
757c4762a1bSJed Brown command line changes */
SetParams(Parameter * param,GridInfo * grid)758d71ae5a4SJacob Faibussowitsch PetscErrorCode SetParams(Parameter *param, GridInfo *grid)
759d71ae5a4SJacob Faibussowitsch {
760c4762a1bSJed Brown PetscReal SEC_PER_YR = 3600.00 * 24.00 * 365.2500;
761c4762a1bSJed Brown PetscReal alpha_g_on_cp_units_inverse_km = 4.0e-5 * 9.8;
762c4762a1bSJed Brown
7633ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
764c4762a1bSJed Brown /* domain geometry */
765c4762a1bSJed Brown param->slab_dip = 45.0;
766c4762a1bSJed Brown param->width = 320.0; /* km */
767c4762a1bSJed Brown param->depth = 300.0; /* km */
768c4762a1bSJed Brown param->lid_depth = 35.0; /* km */
769c4762a1bSJed Brown param->fault_depth = 35.0; /* km */
770c4762a1bSJed Brown
771f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_dip", ¶m->slab_dip, NULL));
772f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-width", ¶m->width, NULL));
773f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-depth", ¶m->depth, NULL));
774f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-lid_depth", ¶m->lid_depth, NULL));
775f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-fault_depth", ¶m->fault_depth, NULL));
776c4762a1bSJed Brown
777c4762a1bSJed Brown param->slab_dip = param->slab_dip * PETSC_PI / 180.0; /* radians */
778c4762a1bSJed Brown
779c4762a1bSJed Brown /* grid information */
780f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetInt(NULL, NULL, "-jfault", &grid->jfault, NULL));
781c4762a1bSJed Brown grid->ni = 82;
782f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetInt(NULL, NULL, "-ni", &grid->ni, NULL));
783c4762a1bSJed Brown
784c4762a1bSJed Brown grid->dx = param->width / ((PetscReal)(grid->ni - 2)); /* km */
785c4762a1bSJed Brown grid->dz = grid->dx * PetscTanReal(param->slab_dip); /* km */
786c4762a1bSJed Brown grid->nj = (PetscInt)(param->depth / grid->dz + 3.0); /* gridpoints*/
787c4762a1bSJed Brown param->depth = grid->dz * (grid->nj - 2); /* km */
788c4762a1bSJed Brown grid->inose = 0; /* gridpoints*/
789f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetInt(NULL, NULL, "-inose", &grid->inose, NULL));
790c4762a1bSJed Brown grid->bx = DM_BOUNDARY_NONE;
791c4762a1bSJed Brown grid->by = DM_BOUNDARY_NONE;
792c4762a1bSJed Brown grid->stencil = DMDA_STENCIL_BOX;
793c4762a1bSJed Brown grid->dof = 4;
794c4762a1bSJed Brown grid->stencil_width = 2;
795c4762a1bSJed Brown grid->mglevels = 1;
796c4762a1bSJed Brown
797c4762a1bSJed Brown /* boundary conditions */
798c4762a1bSJed Brown param->pv_analytic = PETSC_FALSE;
799c4762a1bSJed Brown param->ibound = BC_NOSTRESS;
800f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetInt(NULL, NULL, "-ibound", ¶m->ibound, NULL));
801c4762a1bSJed Brown
802c4762a1bSJed Brown /* physical constants */
803c4762a1bSJed Brown param->slab_velocity = 5.0; /* cm/yr */
804c4762a1bSJed Brown param->slab_age = 50.0; /* Ma */
805c4762a1bSJed Brown param->lid_age = 50.0; /* Ma */
806c4762a1bSJed Brown param->kappa = 0.7272e-6; /* m^2/sec */
807c4762a1bSJed Brown param->potentialT = 1300.0; /* degrees C */
808c4762a1bSJed Brown
809f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_velocity", ¶m->slab_velocity, NULL));
810f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-slab_age", ¶m->slab_age, NULL));
811f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-lid_age", ¶m->lid_age, NULL));
812f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", ¶m->kappa, NULL));
813f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-potentialT", ¶m->potentialT, NULL));
814c4762a1bSJed Brown
815c4762a1bSJed Brown /* viscosity */
816c4762a1bSJed Brown param->ivisc = 3; /* 0=isovisc, 1=difn creep, 2=disl creep, 3=full */
817c4762a1bSJed Brown param->eta0 = 1e24; /* Pa-s */
818c4762a1bSJed Brown param->visc_cutoff = 0.0; /* factor of eta_0 */
819c4762a1bSJed Brown param->continuation = 1.0;
820c4762a1bSJed Brown
821c4762a1bSJed Brown /* constants for diffusion creep */
822c4762a1bSJed Brown param->diffusion.A = 1.8e7; /* Pa-s */
823c4762a1bSJed Brown param->diffusion.n = 1.0; /* dim'less */
824c4762a1bSJed Brown param->diffusion.Estar = 375e3; /* J/mol */
825c4762a1bSJed Brown param->diffusion.Vstar = 5e-6; /* m^3/mol */
826c4762a1bSJed Brown
827c4762a1bSJed Brown /* constants for param->dislocationocation creep */
828c4762a1bSJed Brown param->dislocation.A = 2.8969e4; /* Pa-s */
829c4762a1bSJed Brown param->dislocation.n = 3.5; /* dim'less */
830c4762a1bSJed Brown param->dislocation.Estar = 530e3; /* J/mol */
831c4762a1bSJed Brown param->dislocation.Vstar = 14e-6; /* m^3/mol */
832c4762a1bSJed Brown
833f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetInt(NULL, NULL, "-ivisc", ¶m->ivisc, NULL));
834f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-visc_cutoff", ¶m->visc_cutoff, NULL));
835c4762a1bSJed Brown
836c4762a1bSJed Brown param->output_ivisc = param->ivisc;
837c4762a1bSJed Brown
838f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetInt(NULL, NULL, "-output_ivisc", ¶m->output_ivisc, NULL));
839f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-vstar", ¶m->dislocation.Vstar, NULL));
840c4762a1bSJed Brown
841c4762a1bSJed Brown /* output options */
842c4762a1bSJed Brown param->quiet = PETSC_FALSE;
843c4762a1bSJed Brown param->param_test = PETSC_FALSE;
844c4762a1bSJed Brown
845f4f49eeaSPierre Jolivet PetscCall(PetscOptionsHasName(NULL, NULL, "-quiet", ¶m->quiet));
846f4f49eeaSPierre Jolivet PetscCall(PetscOptionsHasName(NULL, NULL, "-test", ¶m->param_test));
847f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetString(NULL, NULL, "-file", param->filename, sizeof(param->filename), ¶m->output_to_file));
848c4762a1bSJed Brown
849c4762a1bSJed Brown /* advection */
850c4762a1bSJed Brown param->adv_scheme = ADVECT_FROMM; /* advection scheme: 0=finite vol, 1=Fromm */
851c4762a1bSJed Brown
852f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetInt(NULL, NULL, "-adv_scheme", ¶m->adv_scheme, NULL));
853c4762a1bSJed Brown
854c4762a1bSJed Brown /* misc. flags */
855c4762a1bSJed Brown param->stop_solve = PETSC_FALSE;
856c4762a1bSJed Brown param->interrupted = PETSC_FALSE;
857c4762a1bSJed Brown param->kspmon = PETSC_FALSE;
858c4762a1bSJed Brown param->toggle_kspmon = PETSC_FALSE;
859c4762a1bSJed Brown
860c4762a1bSJed Brown /* derived parameters for slab angle */
861c4762a1bSJed Brown param->sb = PetscSinReal(param->slab_dip);
862c4762a1bSJed Brown param->cb = PetscCosReal(param->slab_dip);
863c4762a1bSJed Brown param->c = param->slab_dip * param->sb / (param->slab_dip * param->slab_dip - param->sb * param->sb);
864c4762a1bSJed Brown param->d = (param->slab_dip * param->cb - param->sb) / (param->slab_dip * param->slab_dip - param->sb * param->sb);
865c4762a1bSJed Brown
866c4762a1bSJed Brown /* length, velocity and time scale for non-dimensionalization */
867c4762a1bSJed Brown param->L = PetscMin(param->width, param->depth); /* km */
868c4762a1bSJed Brown param->V = param->slab_velocity / 100.0 / SEC_PER_YR; /* m/sec */
869c4762a1bSJed Brown
870c4762a1bSJed Brown /* other unit conversions and derived parameters */
871c4762a1bSJed Brown param->scaled_width = param->width / param->L; /* dim'less */
872c4762a1bSJed Brown param->scaled_depth = param->depth / param->L; /* dim'less */
873c4762a1bSJed Brown param->lid_depth = param->lid_depth / param->L; /* dim'less */
874c4762a1bSJed Brown param->fault_depth = param->fault_depth / param->L; /* dim'less */
875c4762a1bSJed Brown grid->dx = grid->dx / param->L; /* dim'less */
876c4762a1bSJed Brown grid->dz = grid->dz / param->L; /* dim'less */
877c4762a1bSJed Brown grid->jlid = (PetscInt)(param->lid_depth / grid->dz); /* gridcells */
878c4762a1bSJed Brown grid->jfault = (PetscInt)(param->fault_depth / grid->dz); /* gridcells */
879c4762a1bSJed Brown param->lid_depth = grid->jlid * grid->dz; /* dim'less */
880c4762a1bSJed Brown param->fault_depth = grid->jfault * grid->dz; /* dim'less */
881c4762a1bSJed Brown grid->corner = grid->jlid + 1; /* gridcells */
882c4762a1bSJed Brown param->peclet = param->V /* m/sec */
883c4762a1bSJed Brown * param->L * 1000.0 /* m */
884c4762a1bSJed Brown / param->kappa; /* m^2/sec */
885c4762a1bSJed Brown param->z_scale = param->L * alpha_g_on_cp_units_inverse_km;
886c4762a1bSJed Brown param->skt = PetscSqrtReal(param->kappa * param->slab_age * SEC_PER_YR);
887f4f49eeaSPierre Jolivet PetscCall(PetscOptionsGetReal(NULL, NULL, "-peclet", ¶m->peclet, NULL));
8883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
889c4762a1bSJed Brown }
890c4762a1bSJed Brown
891c4762a1bSJed Brown /* prints a report of the problem parameters to stdout */
ReportParams(Parameter * param,GridInfo * grid)892d71ae5a4SJacob Faibussowitsch PetscErrorCode ReportParams(Parameter *param, GridInfo *grid)
893d71ae5a4SJacob Faibussowitsch {
894c4762a1bSJed Brown char date[30];
895c4762a1bSJed Brown
8963ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
8979566063dSJacob Faibussowitsch PetscCall(PetscGetDate(date, 30));
898c4762a1bSJed Brown
899f4f49eeaSPierre Jolivet if (!param->quiet) {
9009566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---------------------BEGIN ex30 PARAM REPORT-------------------\n"));
9019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Domain: \n"));
9029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Width = %g km, Depth = %g km\n", (double)param->width, (double)param->depth));
9039566063dSJacob 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));
9049566063dSJacob 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)));
905c4762a1bSJed Brown
9069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nGrid: \n"));
90763a3b9bcSJacob 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)));
90863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " jlid = %3" PetscInt_FMT " jfault = %3" PetscInt_FMT " \n", grid->jlid, grid->jfault));
9099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Pe = %g\n", (double)param->peclet));
910c4762a1bSJed Brown
9119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nRheology:"));
912c4762a1bSJed Brown if (param->ivisc == VISC_CONST) {
9139566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Isoviscous \n"));
91448a46eb9SPierre Jolivet if (param->pv_analytic) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Pressure and Velocity prescribed! \n"));
915c4762a1bSJed Brown } else if (param->ivisc == VISC_DIFN) {
9169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Diffusion Creep (T-Dependent Newtonian) \n"));
9179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
918c4762a1bSJed Brown } else if (param->ivisc == VISC_DISL) {
9199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Dislocation Creep (T-Dependent Non-Newtonian) \n"));
9209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
921c4762a1bSJed Brown } else if (param->ivisc == VISC_FULL) {
9229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Full Rheology \n"));
9239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Viscosity range: %g--%g Pa-sec \n", (double)param->eta0, (double)(param->visc_cutoff * param->eta0)));
924c4762a1bSJed Brown } else {
9259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Invalid! \n"));
9263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_ERR_ARG_WRONG);
927c4762a1bSJed Brown }
928c4762a1bSJed Brown
9299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Boundary condition:"));
930c4762a1bSJed Brown if (param->ibound == BC_ANALYTIC) {
9319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Isoviscous Analytic Dirichlet \n"));
932c4762a1bSJed Brown } else if (param->ibound == BC_NOSTRESS) {
9339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Stress-Free (normal & shear stress)\n"));
934c4762a1bSJed Brown } else if (param->ibound == BC_EXPERMNT) {
9359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Experimental boundary condition \n"));
936c4762a1bSJed Brown } else {
9379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Invalid! \n"));
9383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_ERR_ARG_WRONG);
939c4762a1bSJed Brown }
940c4762a1bSJed Brown
94160d4fc61SSatish Balay if (param->output_to_file) {
942d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
9439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output Destination: Mat file \"%s\"\n", param->filename));
944c4762a1bSJed Brown #else
9459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Output Destination: PETSc binary file \"%s\"\n", param->filename));
946c4762a1bSJed Brown #endif
94760d4fc61SSatish Balay }
94848a46eb9SPierre Jolivet if (param->output_ivisc != param->ivisc) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Output viscosity: -ivisc %" PetscInt_FMT "\n", param->output_ivisc));
949c4762a1bSJed Brown
9509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---------------------END ex30 PARAM REPORT---------------------\n"));
951c4762a1bSJed Brown }
9523ba16761SJacob Faibussowitsch if (param->param_test) PetscCall(PetscEnd());
9533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
954c4762a1bSJed Brown }
955c4762a1bSJed Brown
956c4762a1bSJed Brown /* ------------------------------------------------------------------- */
957a5b23f4aSJose E. Roman /* generates an initial guess using the analytic solution for isoviscous
958c4762a1bSJed Brown corner flow */
Initialize(DM da)959c4762a1bSJed Brown PetscErrorCode Initialize(DM da)
960c4762a1bSJed Brown /* ------------------------------------------------------------------- */
961c4762a1bSJed Brown {
962c4762a1bSJed Brown AppCtx *user;
963c4762a1bSJed Brown Parameter *param;
964c4762a1bSJed Brown GridInfo *grid;
965c4762a1bSJed Brown PetscInt i, j, is, js, im, jm;
966c4762a1bSJed Brown Field **x;
967c4762a1bSJed Brown Vec Xguess;
968c4762a1bSJed Brown
9693ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
970c4762a1bSJed Brown /* Get the fine grid */
9719566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(da, &user));
972c4762a1bSJed Brown Xguess = user->Xguess;
973c4762a1bSJed Brown param = user->param;
974c4762a1bSJed Brown grid = user->grid;
9759566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
9769566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xguess, (void **)&x));
977c4762a1bSJed Brown
978c4762a1bSJed Brown /* Compute initial guess */
979c4762a1bSJed Brown for (j = js; j < js + jm; j++) {
980c4762a1bSJed Brown for (i = is; i < is + im; i++) {
981c4762a1bSJed Brown if (i < j) x[j][i].u = param->cb;
982c4762a1bSJed Brown else if (j <= grid->jlid) x[j][i].u = 0.0;
983c4762a1bSJed Brown else x[j][i].u = HorizVelocity(i, j, user);
984c4762a1bSJed Brown
985c4762a1bSJed Brown if (i <= j) x[j][i].w = param->sb;
986c4762a1bSJed Brown else if (j <= grid->jlid) x[j][i].w = 0.0;
987c4762a1bSJed Brown else x[j][i].w = VertVelocity(i, j, user);
988c4762a1bSJed Brown
989c4762a1bSJed Brown if (i < j || j <= grid->jlid) x[j][i].p = 0.0;
990c4762a1bSJed Brown else x[j][i].p = Pressure(i, j, user);
991c4762a1bSJed Brown
992c4762a1bSJed Brown x[j][i].T = PetscMin(grid->dz * (j - 0.5), 1.0);
993c4762a1bSJed Brown }
994c4762a1bSJed Brown }
995c4762a1bSJed Brown
996c4762a1bSJed Brown /* Restore x to Xguess */
9979566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xguess, (void **)&x));
9983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
999c4762a1bSJed Brown }
1000c4762a1bSJed Brown
1001c4762a1bSJed Brown /* controls output to a file */
DoOutput(SNES snes,PetscInt its)1002d71ae5a4SJacob Faibussowitsch PetscErrorCode DoOutput(SNES snes, PetscInt its)
1003d71ae5a4SJacob Faibussowitsch {
1004c4762a1bSJed Brown AppCtx *user;
1005c4762a1bSJed Brown Parameter *param;
1006c4762a1bSJed Brown GridInfo *grid;
1007c4762a1bSJed Brown PetscInt ivt;
1008c4762a1bSJed Brown PetscMPIInt rank;
1009c4762a1bSJed Brown PetscViewer viewer;
1010c4762a1bSJed Brown Vec res, pars;
1011c4762a1bSJed Brown MPI_Comm comm;
1012c4762a1bSJed Brown DM da;
1013c4762a1bSJed Brown
10143ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
10159566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &da));
10169566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(da, &user));
1017c4762a1bSJed Brown param = user->param;
1018c4762a1bSJed Brown grid = user->grid;
1019c4762a1bSJed Brown ivt = param->ivisc;
1020c4762a1bSJed Brown
1021c4762a1bSJed Brown param->ivisc = param->output_ivisc;
1022c4762a1bSJed Brown
1023c4762a1bSJed Brown /* compute final residual and final viscosity/strain rate fields */
10249566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &res, NULL, NULL));
10259566063dSJacob Faibussowitsch PetscCall(ViscosityField(da, user->x, user->Xguess));
1026c4762a1bSJed Brown
1027c4762a1bSJed Brown /* get the communicator and the rank of the processor */
10289566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
10299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
1030c4762a1bSJed Brown
1031c4762a1bSJed Brown if (param->output_to_file) { /* send output to binary file */
10329566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &pars));
1033dd400576SPatrick Sanan if (rank == 0) { /* on processor 0 */
10349566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pars, 20, PETSC_DETERMINE));
10359566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pars));
1036f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 0, (PetscScalar)grid->ni, INSERT_VALUES));
1037f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 1, (PetscScalar)grid->nj, INSERT_VALUES));
1038f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 2, (PetscScalar)grid->dx, INSERT_VALUES));
1039f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 3, (PetscScalar)grid->dz, INSERT_VALUES));
1040f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 4, (PetscScalar)param->L, INSERT_VALUES));
1041f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 5, (PetscScalar)param->V, INSERT_VALUES));
1042c4762a1bSJed Brown /* skipped 6 intentionally */
1043f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 7, (PetscScalar)param->slab_dip, INSERT_VALUES));
1044f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 8, (PetscScalar)grid->jlid, INSERT_VALUES));
1045f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 9, (PetscScalar)param->lid_depth, INSERT_VALUES));
1046f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 10, (PetscScalar)grid->jfault, INSERT_VALUES));
1047f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 11, (PetscScalar)param->fault_depth, INSERT_VALUES));
1048f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 12, (PetscScalar)param->potentialT, INSERT_VALUES));
1049f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 13, (PetscScalar)param->ivisc, INSERT_VALUES));
1050f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 14, (PetscScalar)param->visc_cutoff, INSERT_VALUES));
1051f4f49eeaSPierre Jolivet PetscCall(VecSetValue(pars, 15, (PetscScalar)param->ibound, INSERT_VALUES));
105257508eceSPierre Jolivet PetscCall(VecSetValue(pars, 16, (PetscScalar)its, INSERT_VALUES));
1053c4762a1bSJed Brown } else { /* on some other processor */
10549566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pars, 0, PETSC_DETERMINE));
10559566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pars));
1056c4762a1bSJed Brown }
10579371c9d4SSatish Balay PetscCall(VecAssemblyBegin(pars));
10589371c9d4SSatish Balay PetscCall(VecAssemblyEnd(pars));
1059c4762a1bSJed Brown
1060c4762a1bSJed Brown /* create viewer */
1061d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
10629566063dSJacob Faibussowitsch PetscCall(PetscViewerMatlabOpen(PETSC_COMM_WORLD, param->filename, FILE_MODE_WRITE, &viewer));
1063c4762a1bSJed Brown #else
10649566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, param->filename, FILE_MODE_WRITE, &viewer));
1065c4762a1bSJed Brown #endif
1066c4762a1bSJed Brown
1067c4762a1bSJed Brown /* send vectors to viewer */
10689566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)res, "res"));
10699566063dSJacob Faibussowitsch PetscCall(VecView(res, viewer));
10709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)user->x, "out"));
10719566063dSJacob Faibussowitsch PetscCall(VecView(user->x, viewer));
1072f4f49eeaSPierre Jolivet PetscCall(PetscObjectSetName((PetscObject)user->Xguess, "aux"));
10739566063dSJacob Faibussowitsch PetscCall(VecView(user->Xguess, viewer));
10749566063dSJacob Faibussowitsch PetscCall(StressField(da)); /* compute stress fields */
1075f4f49eeaSPierre Jolivet PetscCall(PetscObjectSetName((PetscObject)user->Xguess, "str"));
10769566063dSJacob Faibussowitsch PetscCall(VecView(user->Xguess, viewer));
10779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)pars, "par"));
10789566063dSJacob Faibussowitsch PetscCall(VecView(pars, viewer));
1079c4762a1bSJed Brown
1080c4762a1bSJed Brown /* destroy viewer and vector */
10819566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
10829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pars));
1083c4762a1bSJed Brown }
1084c4762a1bSJed Brown
1085c4762a1bSJed Brown param->ivisc = ivt;
10863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1087c4762a1bSJed Brown }
1088c4762a1bSJed Brown
1089c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1090c4762a1bSJed Brown /* Compute both the second invariant of the strain rate tensor and the viscosity, at both cell centers and cell corners */
ViscosityField(DM da,Vec X,Vec V)1091c4762a1bSJed Brown PetscErrorCode ViscosityField(DM da, Vec X, Vec V)
1092c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1093c4762a1bSJed Brown {
1094c4762a1bSJed Brown AppCtx *user;
1095c4762a1bSJed Brown Parameter *param;
1096c4762a1bSJed Brown GridInfo *grid;
1097c4762a1bSJed Brown Vec localX;
1098c4762a1bSJed Brown Field **v, **x;
1099c4762a1bSJed Brown PetscReal eps, /* dx,*/ dz, T, epsC, TC;
1100c4762a1bSJed Brown PetscInt i, j, is, js, im, jm, ilim, jlim, ivt;
1101c4762a1bSJed Brown
1102c4762a1bSJed Brown PetscFunctionBeginUser;
11039566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(da, &user));
1104c4762a1bSJed Brown param = user->param;
1105c4762a1bSJed Brown grid = user->grid;
1106c4762a1bSJed Brown ivt = param->ivisc;
1107c4762a1bSJed Brown param->ivisc = param->output_ivisc;
1108c4762a1bSJed Brown
11099566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX));
11109566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
11119566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
11129566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localX, (void **)&x));
11139566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, V, (void **)&v));
1114c4762a1bSJed Brown
1115c4762a1bSJed Brown /* Parameters */
1116c4762a1bSJed Brown /* dx = grid->dx; */ dz = grid->dz;
1117c4762a1bSJed Brown
11189371c9d4SSatish Balay ilim = grid->ni - 1;
11199371c9d4SSatish Balay jlim = grid->nj - 1;
1120c4762a1bSJed Brown
1121c4762a1bSJed Brown /* Compute real temperature, strain rate and viscosity */
11229566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
1123c4762a1bSJed Brown for (j = js; j < js + jm; j++) {
1124c4762a1bSJed Brown for (i = is; i < is + im; i++) {
1125c4762a1bSJed Brown T = PetscRealPart(param->potentialT * x[j][i].T * PetscExpScalar((j - 0.5) * dz * param->z_scale));
1126c4762a1bSJed Brown if (i < ilim && j < jlim) {
1127c4762a1bSJed Brown TC = PetscRealPart(param->potentialT * TInterp(x, i, j) * PetscExpScalar(j * dz * param->z_scale));
1128c4762a1bSJed Brown } else {
1129c4762a1bSJed Brown TC = T;
1130c4762a1bSJed Brown }
1131f4f49eeaSPierre Jolivet eps = PetscRealPart(CalcSecInv(x, i, j, CELL_CENTER, user));
1132c4762a1bSJed Brown epsC = PetscRealPart(CalcSecInv(x, i, j, CELL_CORNER, user));
1133c4762a1bSJed Brown
1134c4762a1bSJed Brown v[j][i].u = eps;
1135c4762a1bSJed Brown v[j][i].w = epsC;
1136c4762a1bSJed Brown v[j][i].p = Viscosity(T, eps, dz * (j - 0.5), param);
1137c4762a1bSJed Brown v[j][i].T = Viscosity(TC, epsC, dz * j, param);
1138c4762a1bSJed Brown }
1139c4762a1bSJed Brown }
11409566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, V, (void **)&v));
11419566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localX, (void **)&x));
11429566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX));
1143c4762a1bSJed Brown
1144c4762a1bSJed Brown param->ivisc = ivt;
11453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1146c4762a1bSJed Brown }
1147c4762a1bSJed Brown
1148c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1149c4762a1bSJed Brown /* post-processing: compute stress everywhere */
StressField(DM da)1150c4762a1bSJed Brown PetscErrorCode StressField(DM da)
1151c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1152c4762a1bSJed Brown {
1153c4762a1bSJed Brown AppCtx *user;
1154c4762a1bSJed Brown PetscInt i, j, is, js, im, jm;
1155c4762a1bSJed Brown Vec locVec;
1156c4762a1bSJed Brown Field **x, **y;
1157c4762a1bSJed Brown
11583ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
11599566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(da, &user));
1160c4762a1bSJed Brown
1161c4762a1bSJed Brown /* Get the fine grid of Xguess and X */
11629566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &is, &js, NULL, &im, &jm, NULL));
11639566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, user->Xguess, (void **)&x));
1164c4762a1bSJed Brown
11659566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &locVec));
11669566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, user->x, INSERT_VALUES, locVec));
11679566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, user->x, INSERT_VALUES, locVec));
11689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, locVec, (void **)&y));
1169c4762a1bSJed Brown
1170c4762a1bSJed Brown /* Compute stress on the corner points */
1171c4762a1bSJed Brown for (j = js; j < js + jm; j++) {
1172c4762a1bSJed Brown for (i = is; i < is + im; i++) {
1173c4762a1bSJed Brown x[j][i].u = ShearStress(y, i, j, CELL_CENTER, user);
1174c4762a1bSJed Brown x[j][i].w = ShearStress(y, i, j, CELL_CORNER, user);
1175c4762a1bSJed Brown x[j][i].p = XNormalStress(y, i, j, CELL_CENTER, user);
1176c4762a1bSJed Brown x[j][i].T = ZNormalStress(y, i, j, CELL_CENTER, user);
1177c4762a1bSJed Brown }
1178c4762a1bSJed Brown }
1179c4762a1bSJed Brown
1180c4762a1bSJed Brown /* Restore the fine grid of Xguess and X */
11819566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, user->Xguess, (void **)&x));
11829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, locVec, (void **)&y));
11839566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &locVec));
11843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1185c4762a1bSJed Brown }
1186c4762a1bSJed Brown
1187c4762a1bSJed Brown /*=====================================================================
1188c4762a1bSJed Brown UTILITY FUNCTIONS
1189c4762a1bSJed Brown =====================================================================*/
1190c4762a1bSJed Brown
1191f6dfbefdSBarry Smith /* returns the velocity of the subducting slab and handles fault nodes for BC */
SlabVel(char c,PetscInt i,PetscInt j,AppCtx * user)1192d71ae5a4SJacob Faibussowitsch static inline PetscScalar SlabVel(char c, PetscInt i, PetscInt j, AppCtx *user)
1193d71ae5a4SJacob Faibussowitsch {
1194c4762a1bSJed Brown Parameter *param = user->param;
1195c4762a1bSJed Brown GridInfo *grid = user->grid;
1196c4762a1bSJed Brown
1197c4762a1bSJed Brown if (c == 'U' || c == 'u') {
1198c4762a1bSJed Brown if (i < j - 1) return param->cb;
1199c4762a1bSJed Brown else if (j <= grid->jfault) return 0.0;
1200c4762a1bSJed Brown else return param->cb;
1201c4762a1bSJed Brown
1202c4762a1bSJed Brown } else {
1203c4762a1bSJed Brown if (i < j) return param->sb;
1204c4762a1bSJed Brown else if (j <= grid->jfault) return 0.0;
1205c4762a1bSJed Brown else return param->sb;
1206c4762a1bSJed Brown }
1207c4762a1bSJed Brown }
1208c4762a1bSJed Brown
1209c4762a1bSJed Brown /* solution to diffusive half-space cooling model for BC */
PlateModel(PetscInt j,PetscInt plate,AppCtx * user)1210d71ae5a4SJacob Faibussowitsch static inline PetscScalar PlateModel(PetscInt j, PetscInt plate, AppCtx *user)
1211d71ae5a4SJacob Faibussowitsch {
1212c4762a1bSJed Brown Parameter *param = user->param;
1213c4762a1bSJed Brown PetscScalar z;
1214c4762a1bSJed Brown if (plate == PLATE_LID) z = (j - 0.5) * user->grid->dz;
1215c4762a1bSJed Brown else z = (j - 0.5) * user->grid->dz * param->cb; /* PLATE_SLAB */
1216c4762a1bSJed Brown #if defined(PETSC_HAVE_ERF)
1217c4762a1bSJed Brown return (PetscReal)(erf((double)PetscRealPart(z * param->L / 2.0 / param->skt)));
1218c4762a1bSJed Brown #else
1219c4762a1bSJed Brown (*PetscErrorPrintf)("erf() not available on this machine\n");
1220c4762a1bSJed Brown MPI_Abort(PETSC_COMM_SELF, 1);
1221c4762a1bSJed Brown #endif
1222c4762a1bSJed Brown }
1223c4762a1bSJed Brown
1224c4762a1bSJed Brown /*=====================================================================
1225c4762a1bSJed Brown INTERACTIVE SIGNAL HANDLING
1226c4762a1bSJed Brown =====================================================================*/
1227c4762a1bSJed Brown
1228c4762a1bSJed Brown /* ------------------------------------------------------------------- */
SNESConverged_Interactive(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason * reason,PetscCtx ctx)1229*2a8381b2SBarry Smith PetscErrorCode SNESConverged_Interactive(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, PetscCtx ctx)
1230c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1231c4762a1bSJed Brown {
1232c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx;
1233c4762a1bSJed Brown Parameter *param = user->param;
1234c4762a1bSJed Brown KSP ksp;
1235c4762a1bSJed Brown
1236c4762a1bSJed Brown PetscFunctionBeginUser;
1237c4762a1bSJed Brown if (param->interrupted) {
1238c4762a1bSJed Brown param->interrupted = PETSC_FALSE;
12399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: exiting SNES solve. \n"));
1240c4762a1bSJed Brown *reason = SNES_CONVERGED_FNORM_ABS;
12413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1242c4762a1bSJed Brown } else if (param->toggle_kspmon) {
1243c4762a1bSJed Brown param->toggle_kspmon = PETSC_FALSE;
1244c4762a1bSJed Brown
12459566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp));
1246c4762a1bSJed Brown
1247c4762a1bSJed Brown if (param->kspmon) {
12489566063dSJacob Faibussowitsch PetscCall(KSPMonitorCancel(ksp));
1249c4762a1bSJed Brown
1250c4762a1bSJed Brown param->kspmon = PETSC_FALSE;
12519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: deactivating ksp singular value monitor. \n"));
1252c4762a1bSJed Brown } else {
1253c4762a1bSJed Brown PetscViewerAndFormat *vf;
12549566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_DEFAULT, &vf));
125549abdd8aSBarry Smith PetscCall(KSPMonitorSet(ksp, (PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *))KSPMonitorSingularValue, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
1256c4762a1bSJed Brown
1257c4762a1bSJed Brown param->kspmon = PETSC_TRUE;
12589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "USER SIGNAL: activating ksp singular value monitor. \n"));
1259c4762a1bSJed Brown }
1260c4762a1bSJed Brown }
126111cc89d2SBarry Smith PetscCall(SNESConvergedDefault(snes, it, xnorm, snorm, fnorm, reason, ctx));
12623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1263c4762a1bSJed Brown }
1264c4762a1bSJed Brown
1265c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1266c4762a1bSJed Brown #include <signal.h>
InteractiveHandler(int signum,PetscCtx ctx)1267*2a8381b2SBarry Smith PetscErrorCode InteractiveHandler(int signum, PetscCtx ctx)
1268c4762a1bSJed Brown /* ------------------------------------------------------------------- */
1269c4762a1bSJed Brown {
1270c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx;
1271c4762a1bSJed Brown Parameter *param = user->param;
1272c4762a1bSJed Brown
1273c4762a1bSJed Brown if (signum == SIGILL) {
1274c4762a1bSJed Brown param->toggle_kspmon = PETSC_TRUE;
1275c4762a1bSJed Brown #if !defined(PETSC_MISSING_SIGCONT)
1276c4762a1bSJed Brown } else if (signum == SIGCONT) {
1277c4762a1bSJed Brown param->interrupted = PETSC_TRUE;
1278c4762a1bSJed Brown #endif
1279c4762a1bSJed Brown #if !defined(PETSC_MISSING_SIGURG)
1280c4762a1bSJed Brown } else if (signum == SIGURG) {
1281c4762a1bSJed Brown param->stop_solve = PETSC_TRUE;
1282c4762a1bSJed Brown #endif
1283c4762a1bSJed Brown }
12843ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
1285c4762a1bSJed Brown }
1286c4762a1bSJed Brown
1287f6dfbefdSBarry Smith /* main call-back function that computes the processor-local piece of the residual */
FormFunctionLocal(DMDALocalInfo * info,Field ** x,Field ** f,void * ptr)1288d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field **x, Field **f, void *ptr)
1289d71ae5a4SJacob Faibussowitsch {
1290c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
1291c4762a1bSJed Brown Parameter *param = user->param;
1292c4762a1bSJed Brown GridInfo *grid = user->grid;
1293c4762a1bSJed Brown PetscScalar mag_w, mag_u;
1294c4762a1bSJed Brown PetscInt i, j, mx, mz, ilim, jlim;
1295c4762a1bSJed Brown PetscInt is, ie, js, je, ibound; /* ,ivisc */
1296c4762a1bSJed Brown
1297c4762a1bSJed Brown PetscFunctionBeginUser;
1298c4762a1bSJed Brown /* Define global and local grid parameters */
12999371c9d4SSatish Balay mx = info->mx;
13009371c9d4SSatish Balay mz = info->my;
13019371c9d4SSatish Balay ilim = mx - 1;
13029371c9d4SSatish Balay jlim = mz - 1;
13039371c9d4SSatish Balay is = info->xs;
13049371c9d4SSatish Balay ie = info->xs + info->xm;
13059371c9d4SSatish Balay js = info->ys;
13069371c9d4SSatish Balay je = info->ys + info->ym;
1307c4762a1bSJed Brown
1308c4762a1bSJed Brown /* Define geometric and numeric parameters */
1309c4762a1bSJed Brown /* ivisc = param->ivisc; */ ibound = param->ibound;
1310c4762a1bSJed Brown
1311c4762a1bSJed Brown for (j = js; j < je; j++) {
1312c4762a1bSJed Brown for (i = is; i < ie; i++) {
1313c4762a1bSJed Brown /************* X-MOMENTUM/VELOCITY *************/
1314c4762a1bSJed Brown if (i < j) f[j][i].u = x[j][i].u - SlabVel('U', i, j, user);
1315c4762a1bSJed Brown else if (j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1316c4762a1bSJed Brown /* in the lithospheric lid */
1317c4762a1bSJed Brown f[j][i].u = x[j][i].u - 0.0;
1318c4762a1bSJed Brown } else if (i == ilim) {
1319c4762a1bSJed Brown /* on the right side boundary */
1320c4762a1bSJed Brown if (ibound == BC_ANALYTIC) {
1321c4762a1bSJed Brown f[j][i].u = x[j][i].u - HorizVelocity(i, j, user);
1322c4762a1bSJed Brown } else {
1323c4762a1bSJed Brown f[j][i].u = XNormalStress(x, i, j, CELL_CENTER, user) - EPS_ZERO;
1324c4762a1bSJed Brown }
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].u = x[j][i].u - HorizVelocity(i, j, user);
1330c4762a1bSJed Brown } else if (ibound == BC_NOSTRESS) {
1331c4762a1bSJed Brown f[j][i].u = XMomentumResidual(x, i, j, user);
1332c4762a1bSJed Brown } else {
1333c4762a1bSJed Brown /* experimental boundary condition */
1334c4762a1bSJed Brown }
1335c4762a1bSJed Brown
1336c4762a1bSJed Brown } else {
1337c4762a1bSJed Brown /* in the mantle wedge */
1338c4762a1bSJed Brown f[j][i].u = XMomentumResidual(x, i, j, user);
1339c4762a1bSJed Brown }
1340c4762a1bSJed Brown
1341c4762a1bSJed Brown /************* Z-MOMENTUM/VELOCITY *************/
1342c4762a1bSJed Brown if (i <= j) {
1343c4762a1bSJed Brown f[j][i].w = x[j][i].w - SlabVel('W', i, j, user);
1344c4762a1bSJed Brown
1345c4762a1bSJed Brown } else if (j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1346c4762a1bSJed Brown /* in the lithospheric lid */
1347c4762a1bSJed Brown f[j][i].w = x[j][i].w - 0.0;
1348c4762a1bSJed Brown
1349c4762a1bSJed Brown } else if (j == jlim) {
1350c4762a1bSJed Brown /* on the bottom boundary */
1351c4762a1bSJed Brown if (ibound == BC_ANALYTIC) {
1352c4762a1bSJed Brown f[j][i].w = x[j][i].w - VertVelocity(i, j, user);
1353c4762a1bSJed Brown } else {
1354c4762a1bSJed Brown f[j][i].w = ZNormalStress(x, i, j, CELL_CENTER, user) - EPS_ZERO;
1355c4762a1bSJed Brown }
1356c4762a1bSJed Brown
1357c4762a1bSJed Brown } else if (i == ilim) {
1358c4762a1bSJed Brown /* on the right side boundary */
1359c4762a1bSJed Brown if (ibound == BC_ANALYTIC) {
1360c4762a1bSJed Brown f[j][i].w = x[j][i].w - VertVelocity(i, j, user);
1361c4762a1bSJed Brown } else if (ibound == BC_NOSTRESS) {
1362c4762a1bSJed Brown f[j][i].w = ZMomentumResidual(x, i, j, user);
1363c4762a1bSJed Brown } else {
1364c4762a1bSJed Brown /* experimental boundary condition */
1365c4762a1bSJed Brown }
1366c4762a1bSJed Brown
1367c4762a1bSJed Brown } else {
1368c4762a1bSJed Brown /* in the mantle wedge */
1369c4762a1bSJed Brown f[j][i].w = ZMomentumResidual(x, i, j, user);
1370c4762a1bSJed Brown }
1371c4762a1bSJed Brown
1372c4762a1bSJed Brown /************* CONTINUITY/PRESSURE *************/
1373c4762a1bSJed Brown if (i < j || j <= grid->jlid || (j < grid->corner + grid->inose && i < grid->corner + grid->inose)) {
1374c4762a1bSJed Brown /* in the lid or slab */
1375c4762a1bSJed Brown f[j][i].p = x[j][i].p;
1376c4762a1bSJed Brown
1377c4762a1bSJed Brown } else if ((i == ilim || j == jlim) && ibound == BC_ANALYTIC) {
1378c4762a1bSJed Brown /* on an analytic boundary */
1379c4762a1bSJed Brown f[j][i].p = x[j][i].p - Pressure(i, j, user);
1380c4762a1bSJed Brown
1381c4762a1bSJed Brown } else {
1382c4762a1bSJed Brown /* in the mantle wedge */
1383c4762a1bSJed Brown f[j][i].p = ContinuityResidual(x, i, j, user);
1384c4762a1bSJed Brown }
1385c4762a1bSJed Brown
1386c4762a1bSJed Brown /************* TEMPERATURE *************/
1387c4762a1bSJed Brown if (j == 0) {
1388c4762a1bSJed Brown /* on the surface */
1389c4762a1bSJed Brown f[j][i].T = x[j][i].T + x[j + 1][i].T + PetscMax(PetscRealPart(x[j][i].T), 0.0);
1390c4762a1bSJed Brown
1391c4762a1bSJed Brown } else if (i == 0) {
1392c4762a1bSJed Brown /* slab inflow boundary */
1393c4762a1bSJed Brown f[j][i].T = x[j][i].T - PlateModel(j, PLATE_SLAB, user);
1394c4762a1bSJed Brown
1395c4762a1bSJed Brown } else if (i == ilim) {
1396c4762a1bSJed Brown /* right side boundary */
139757508eceSPierre Jolivet mag_u = 1.0 - PetscPowRealInt(1.0 - PetscMax(PetscMin(PetscRealPart(x[j][i - 1].u) / param->cb, 1.0), 0.0), 5);
1398c4762a1bSJed 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);
1399c4762a1bSJed Brown
1400c4762a1bSJed Brown } else if (j == jlim) {
1401c4762a1bSJed Brown /* bottom boundary */
140257508eceSPierre Jolivet mag_w = 1.0 - PetscPowRealInt(1.0 - PetscMax(PetscMin(PetscRealPart(x[j - 1][i].w) / param->sb, 1.0), 0.0), 5);
1403c4762a1bSJed Brown f[j][i].T = x[j][i].T - mag_w * x[j - 1][i - 1].T - (1.0 - mag_w);
1404c4762a1bSJed Brown
1405c4762a1bSJed Brown } else {
1406c4762a1bSJed Brown /* in the mantle wedge */
1407c4762a1bSJed Brown f[j][i].T = EnergyResidual(x, i, j, user);
1408c4762a1bSJed Brown }
1409c4762a1bSJed Brown }
1410c4762a1bSJed Brown }
14113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1412c4762a1bSJed Brown }
1413c4762a1bSJed Brown
1414c4762a1bSJed Brown /*TEST
1415c4762a1bSJed Brown
1416c4762a1bSJed Brown build:
1417c4762a1bSJed Brown requires: !complex erf
1418c4762a1bSJed Brown
1419c4762a1bSJed Brown test:
1420308041e4SStefano Zampini args: -ni 18 -fp_trap 0
1421c4762a1bSJed Brown filter: grep -v Destination
1422c4762a1bSJed Brown requires: !single
1423c4762a1bSJed Brown
1424c4762a1bSJed Brown TEST*/
1425