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