Actual source code: ex30.c
1: static char help[] =
2: "ex30: Steady-state 2D subduction flow, pressure and temperature solver.\n\
3: The flow is driven by the subducting slab.\n\
4: ---------------------------------ex30 help---------------------------------\n\
5: -OPTION <DEFAULT> = (UNITS) DESCRIPTION.\n\n\
6: -width <320> = (km) width of domain.\n\
7: -depth <300> = (km) depth of domain.\n\
8: -slab_dip <45> = (degrees) dip angle of the slab (determines the grid aspect ratio).\n\
9: -lid_depth <35> = (km) depth of the static conductive lid.\n\
10: -fault_depth <35> = (km) depth of slab-wedge mechanical coupling\n\
11: ( fault dept >= lid depth ).\n\
12: \n\
13: -ni <82> = grid cells in x-direction. (nj adjusts to accommodate\n\
14: the slab dip & depth). DO NOT USE -da_grid_x option!!!\n\
15: -ivisc <3> = rheology option.\n\
16: 0 --- constant viscosity.\n\
17: 1 --- olivine diffusion creep rheology (T&P-dependent, newtonian).\n\
18: 2 --- olivine dislocation creep rheology (T&P-dependent, non-newtonian).\n\
19: 3 --- Full mantle rheology, combination of 1 & 2.\n\
20: \n\
21: -slab_velocity <5> = (cm/year) convergence rate of slab into subduction zone.\n\
22: -slab_age <50> = (million yrs) age of slab for thermal profile boundary condition.\n\
23: -lid_age <50> = (million yrs) age of lid for thermal profile boundary condition.\n\
24: \n\
25: FOR OTHER PARAMETER OPTIONS AND THEIR DEFAULT VALUES, see SetParams() in ex30.c.\n\
26: ---------------------------------ex30 help---------------------------------\n";
29: /* ------------------------------------------------------------------------
30:
31: This PETSc 2.2.0 example by Richard F. Katz
32: http://www.ldeo.columbia.edu/~katz/
34: The problem is modeled by the partial differential equation system
35:
36: (1) -Grad(P) + Div[Eta (Grad(v) + Grad(v)^T)] = 0
37: (2) Div(U,W) = 0
38: (3) dT/dt + Div(vT) - 1/Pe Del^2(T) = 0
39: (4) Eta(T,Eps_dot) = constant if ivisc==0
40: = diffusion creep (T,P-dependent) if ivisc==1
41: = dislocation creep (T,P,v-dependent) if ivisc==2
42: = mantle viscosity (difn & disl) if ivisc==3
44: which is uniformly discretized on a staggered mesh:
45: -------w_ij------
46: | |
47: u_i-1j P,T_ij u_ij
48: | |
49: ------w_ij-1-----
51: ------------------------------------------------------------------------- */
53: #include petscsnes.h
54: #include petscda.h
56: #define VISC_CONST 0
57: #define VISC_DIFN 1
58: #define VISC_DISL 2
59: #define VISC_FULL 3
60: #define CELL_CENTER 0
61: #define CELL_CORNER 1
62: #define BC_ANALYTIC 0
63: #define BC_NOSTRESS 1
64: #define BC_EXPERMNT 2
65: #define ADVECT_FV 0
66: #define ADVECT_FROMM 1
67: #define PLATE_SLAB 0
68: #define PLATE_LID 1
69: #define EPS_ZERO 0.00000001
71: typedef struct { /* holds the variables to be solved for */
72: PetscScalar u,w,p,T;
73: } Field;
75: typedef struct { /* parameters needed to compute viscosity */
76: PetscReal A,n,Estar,Vstar;
77: } ViscParam;
79: typedef struct { /* physical and miscelaneous parameters */
80: PetscReal width, depth, scaled_width, scaled_depth, peclet, potentialT;
81: PetscReal slab_dip, slab_age, slab_velocity, kappa, z_scale;
82: PetscReal c, d, sb, cb, skt, visc_cutoff, lid_age, eta0, continuation;
83: PetscReal L, V, lid_depth, fault_depth;
84: ViscParam diffusion, dislocation;
85: PetscInt ivisc, adv_scheme, ibound, output_ivisc;
86: PetscTruth quiet, param_test, output_to_file, pv_analytic;
87: PetscTruth interrupted, stop_solve, toggle_kspmon, kspmon;
88: char filename[PETSC_MAX_PATH_LEN];
89: } Parameter;
91: typedef struct { /* grid parameters */
92: DAPeriodicType periodic;
93: DAStencilType stencil;
94: PetscInt corner,ni,nj,jlid,jfault,inose;
95: PetscInt dof,stencil_width,mglevels;
96: PassiveScalar dx,dz;
97: } GridInfo;
99: typedef struct { /* application context */
100: Vec Xguess;
101: Parameter *param;
102: GridInfo *grid;
103: } AppCtx;
105: /* Callback functions (static interface) */
109: /* Main routines */
116: /* Physics subroutines */
127: /* Utilities for interpolation, ICs and BCs */
138: /* Post-processing & misc */
145: /*-----------------------------------------------------------------------*/
148: int main(int argc,char **argv)
149: /*-----------------------------------------------------------------------*/
150: {
151: DMMG *dmmg; /* multilevel grid structure */
152: AppCtx *user; /* user-defined work context */
153: Parameter param;
154: GridInfo grid;
155: PetscInt nits;
157: MPI_Comm comm;
158: DA da;
160: PetscInitialize(&argc,&argv,(char *)0,help);
161: PetscOptionsSetValue("-file","ex30_output");
162: PetscOptionsSetValue("-snes_monitor",PETSC_NULL);
163: PetscOptionsSetValue("-snes_max_it","20");
164: PetscOptionsSetValue("-ksp_max_it","1500");
165: PetscOptionsSetValue("-ksp_gmres_restart","300");
166: PetscOptionsInsert(&argc,&argv,PETSC_NULL);
168: comm = PETSC_COMM_WORLD;
170: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171: Set up the problem parameters.
172: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173: SetParams(¶m,&grid);
174: ReportParams(¶m,&grid);
176: #if 0
177: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178: Create user context, set problem data, create vector data structures.
179: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180: PetscMalloc(sizeof(AppCtx),&user);
181: user->param = ¶m;
182: user->grid = &grid;
184: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185: Create distributed array multigrid object (DMMG) to manage parallel grid and vectors
186: for principal unknowns (x) and governing residuals (f)
187: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
188: DMMGCreate(comm,grid.mglevels,user,&dmmg);
189: DACreate2d(comm,grid.periodic,grid.stencil,grid.ni,grid.nj,PETSC_DECIDE,PETSC_DECIDE,grid.dof,grid.stencil_width,0,0,&da);
190: DMMGSetDM(dmmg,(DM)da);
191: DADestroy(da);
192: DASetFieldName(DMMGGetDA(dmmg),0,"x-velocity");
193: DASetFieldName(DMMGGetDA(dmmg),1,"y-velocity");
194: DASetFieldName(DMMGGetDA(dmmg),2,"pressure");
195: DASetFieldName(DMMGGetDA(dmmg),3,"temperature");
196: VecDuplicate(dmmg[0]->x, &(user->Xguess));
197: #else
198: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199: Create distributed array multigrid object (DMMG) to manage parallel grid and vectors
200: for principal unknowns (x) and governing residuals (f)
201: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202: DMMGCreate(comm,grid.mglevels,&user,&dmmg);
203: DACreate2d(comm,grid.periodic,grid.stencil,grid.ni,grid.nj,PETSC_DECIDE,PETSC_DECIDE,grid.dof,grid.stencil_width,0,0,&da);
204: DMMGSetDM(dmmg,(DM)da);
205: DADestroy(da);
206: DASetFieldName(DMMGGetDA(dmmg),0,"x-velocity");
207: DASetFieldName(DMMGGetDA(dmmg),1,"y-velocity");
208: DASetFieldName(DMMGGetDA(dmmg),2,"pressure");
209: DASetFieldName(DMMGGetDA(dmmg),3,"temperature");
211: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212: Create user context, set problem data, create vector data structures.
213: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
214: PetscMalloc(sizeof(AppCtx),&user);
215: user->param = ¶m;
216: user->grid = &grid;
217: dmmg[0]->user = user;
218: VecDuplicate(dmmg[0]->x, &(user->Xguess));
219: #endif
221: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222: Set up the SNES solver with callback functions.
223: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224: DMMGSetSNESLocal(dmmg,FormFunctionLocal,0,0,0);
225: DMMGSetInitialGuess(dmmg,FormInitialGuess);
226: SNESSetConvergenceTest(DMMGGetSNES(dmmg),SNESConverged_Interactive,(void*)user);
227: PetscPushSignalHandler(InteractiveHandler,(void*)user);
228:
229: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230: Initialize and solve the nonlinear system
231: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
232: Initialize(dmmg);
233: UpdateSolution(dmmg,user,&nits);
234:
235: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236: Output variables.
237: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238: DoOutput(dmmg,nits);
239:
240: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241: Free work space.
242: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243: VecDestroy(user->Xguess);
244: PetscFree(user);
245: DMMGDestroy(dmmg);
246:
247: PetscFinalize();
248: return 0;
249: }
251: /*=====================================================================
252: PETSc INTERACTION FUNCTIONS (initialize & call SNESSolve)
253: =====================================================================*/
255: /*---------------------------------------------------------------------*/
258: /* manages solve: adaptive continuation method */
259: PetscErrorCode UpdateSolution(DMMG *dmmg, AppCtx *user, PetscInt *nits)
260: {
261: SNES snes;
262: KSP ksp;
263: PC pc;
264: SNESConvergedReason reason;
265: Parameter *param = user->param;
266: PassiveScalar cont_incr=0.3;
267: PetscInt its;
268: PetscErrorCode ierr;
269: PetscTruth q = PETSC_FALSE;
272: snes = DMMGGetSNES(dmmg);
273: SNESGetKSP(snes,&ksp);
274: KSPGetPC(ksp, &pc);
275: KSPSetComputeSingularValues(ksp, PETSC_TRUE);
277: *nits=0;
279: /* Isoviscous solve */
280: if (param->ivisc == VISC_CONST && !param->stop_solve) {
281: param->ivisc = VISC_CONST;
282: DMMGSolve(dmmg);
283: VecCopy(DMMGGetx(dmmg),user->Xguess);
284: SNESGetIterationNumber(snes, &its);
285: *nits +=its;
286: if (param->stop_solve) goto done;
287: }
289: /* Olivine diffusion creep */
290: if (param->ivisc >= VISC_DIFN && !param->stop_solve) {
291: if (!q) PetscPrintf(PETSC_COMM_WORLD,"Computing Variable Viscosity Solution\n");
293: /* continuation method on viscosity cutoff */
294: for (param->continuation = 0.0; param->continuation<=1.0;) {
295: if (!q) PetscPrintf(PETSC_COMM_WORLD," Continuation parameter = %g\n", param->continuation);
297: /* solve the non-linear system */
298: DMMGSolve(dmmg);
299: SNESGetConvergedReason(snes,&reason);
300: SNESGetIterationNumber(snes,&its);
301: *nits += its;
302: if (!q) PetscPrintf(PETSC_COMM_WORLD," Newton iterations: %D, Cumulative: %D\n", its, *nits);
303: if (param->stop_solve || param->continuation==1.0) goto done;
305: if (reason<0) {
306: /* NOT converged */
307: cont_incr = cont_incr/2.0;
308: param->continuation -= cont_incr;
309: if (cont_incr<0.01) goto done;
311: } else {
312: /* converged */
313: VecCopy(DMMGGetx(dmmg),user->Xguess);
314: if (its<=3) {
315: cont_incr = 0.30001;
316: } else if (its<=8) {
317: cont_incr = 0.15001;
318: } else {
319: cont_incr = 0.10001;
320: }
321: param->continuation = PetscMin(param->continuation+cont_incr,1.0);
322: } /* endif reason<0 */
323: }
324: }
325: done:
326: if (param->stop_solve && !q) PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: stopping solve.\n");
327: if (reason<0 && !q) PetscPrintf(PETSC_COMM_WORLD,"FAILED TO CONVERGE: stopping solve.\n");
328: return(0);
329: }
331: /* ------------------------------------------------------------------- */
334: /* used by SNESSolve to get an initial guess for the solution X */
335: PetscErrorCode FormInitialGuess(DMMG dmmg,Vec X)
336: /* ------------------------------------------------------------------- */
337: {
338: AppCtx *user = (AppCtx*)dmmg->user;
341: VecCopy(user->Xguess, X);
342: return 0;
343: }
345: /*=====================================================================
346: PHYSICS FUNCTIONS (compute the discrete residual)
347: =====================================================================*/
349: /*---------------------------------------------------------------------*/
352: /* main call-back function that computes the processor-local piece
353: of the residual */
354: PetscErrorCode FormFunctionLocal(DALocalInfo *info,Field **x,Field **f,void *ptr)
355: /*---------------------------------------------------------------------*/
356: {
357: AppCtx *user = (AppCtx*)ptr;
358: Parameter *param = user->param;
359: GridInfo *grid = user->grid;
360: PetscScalar mag_w, mag_u;
361: PetscInt i,j,mx,mz,ilim,jlim;
362: PetscInt is,ie,js,je,ivisc,ibound;
366: /* Define global and local grid parameters */
367: mx = info->mx; mz = info->my;
368: ilim = mx-1; jlim = mz-1;
369: is = info->xs; ie = info->xs+info->xm;
370: js = info->ys; je = info->ys+info->ym;
372: /* Define geometric and numeric parameters */
373: ivisc = param->ivisc; ibound = param->ibound;
375: for (j=js; j<je; j++) {
376: for (i=is; i<ie; i++) {
378: /************* X-MOMENTUM/VELOCITY *************/
379: if (i<j) {
380: f[j][i].u = x[j][i].u - SlabVel('U',i,j,user);
382: } else if (j<=grid->jlid || (j<grid->corner+grid->inose && i<grid->corner+grid->inose)) {
383: /* in the lithospheric lid */
384: f[j][i].u = x[j][i].u - 0.0;
386: } else if (i==ilim) {
387: /* on the right side boundary */
388: if (ibound==BC_ANALYTIC) {
389: f[j][i].u = x[j][i].u - HorizVelocity(i,j,user);
390: } else {
391: f[j][i].u = XNormalStress(x,i,j,CELL_CENTER,user) - EPS_ZERO;
392: }
394: } else if (j==jlim) {
395: /* on the bottom boundary */
396: if (ibound==BC_ANALYTIC) {
397: f[j][i].u = x[j][i].u - HorizVelocity(i,j,user);
398: } else if (ibound==BC_NOSTRESS) {
399: f[j][i].u = XMomentumResidual(x,i,j,user);
400: } else {
401: /* experimental boundary condition */
402: }
404: } else {
405: /* in the mantle wedge */
406: f[j][i].u = XMomentumResidual(x,i,j,user);
407: }
408:
409: /************* Z-MOMENTUM/VELOCITY *************/
410: if (i<=j) {
411: f[j][i].w = x[j][i].w - SlabVel('W',i,j,user);
413: } else if (j<=grid->jlid || (j<grid->corner+grid->inose && i<grid->corner+grid->inose)) {
414: /* in the lithospheric lid */
415: f[j][i].w = x[j][i].w - 0.0;
417: } else if (j==jlim) {
418: /* on the bottom boundary */
419: if (ibound==BC_ANALYTIC) {
420: f[j][i].w = x[j][i].w - VertVelocity(i,j,user);
421: } else {
422: f[j][i].w = ZNormalStress(x,i,j,CELL_CENTER,user) - EPS_ZERO;
423: }
425: } else if (i==ilim) {
426: /* on the right side boundary */
427: if (ibound==BC_ANALYTIC) {
428: f[j][i].w = x[j][i].w - VertVelocity(i,j,user);
429: } else if (ibound==BC_NOSTRESS) {
430: f[j][i].w = ZMomentumResidual(x,i,j,user);
431: } else {
432: /* experimental boundary condition */
433: }
435: } else {
436: /* in the mantle wedge */
437: f[j][i].w = ZMomentumResidual(x,i,j,user);
438: }
440: /************* CONTINUITY/PRESSURE *************/
441: if (i<j || j<=grid->jlid || (j<grid->corner+grid->inose && i<grid->corner+grid->inose)) {
442: /* in the lid or slab */
443: f[j][i].p = x[j][i].p;
444:
445: } else if ((i==ilim || j==jlim) && ibound==BC_ANALYTIC) {
446: /* on an analytic boundary */
447: f[j][i].p = x[j][i].p - Pressure(i,j,user);
449: } else {
450: /* in the mantle wedge */
451: f[j][i].p = ContinuityResidual(x,i,j,user);
452: }
454: /************* TEMPERATURE *************/
455: if (j==0) {
456: /* on the surface */
457: f[j][i].T = x[j][i].T + x[j+1][i].T + PetscMax(x[j][i].T,0.0);
459: } else if (i==0) {
460: /* slab inflow boundary */
461: f[j][i].T = x[j][i].T - PlateModel(j,PLATE_SLAB,user);
463: } else if (i==ilim) {
464: /* right side boundary */
465: mag_u = 1.0 - pow( (1.0-PetscMax(PetscMin(x[j][i-1].u/param->cb,1.0),0.0)), 5.0 );
466: 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);
468: } else if (j==jlim) {
469: /* bottom boundary */
470: mag_w = 1.0 - pow( (1.0-PetscMax(PetscMin(x[j-1][i].w/param->sb,1.0),0.0)), 5.0 );
471: f[j][i].T = x[j][i].T - mag_w*x[j-1][i-1].T - (1.0-mag_w);
473: } else {
474: /* in the mantle wedge */
475: f[j][i].T = EnergyResidual(x,i,j,user);
476: }
477: }
478: }
479: return(0);
480: }
482: /*---------------------------------------------------------------------*/
485: /* computes the residual of the x-component of eqn (1) above */
486: PetscScalar XMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
487: /*---------------------------------------------------------------------*/
488: {
489: Parameter *param=user->param;
490: GridInfo *grid =user->grid;
491: PetscScalar dx = grid->dx, dz=grid->dz;
492: PetscScalar etaN,etaS,etaE,etaW,epsN=0.0,epsS=0.0,epsE=0.0,epsW=0.0;
493: PetscScalar TE=0.0,TN=0.0,TS=0.0,TW=0.0, dPdx, residual, z_scale;
494: PetscScalar dudxW,dudxE,dudzN,dudzS,dwdxN,dwdxS;
495: PetscInt jlim = grid->nj-1;
497: z_scale = param->z_scale;
499: if ( param->ivisc==VISC_DIFN || param->ivisc>=VISC_DISL ) { /* viscosity is T-dependent */
500: TS = param->potentialT * TInterp(x,i,j-1) * exp( (j-1.0)*dz*z_scale );
501: if (j==jlim) TN = TS;
502: else TN = param->potentialT * TInterp(x,i,j) * exp( j *dz*z_scale );
503: TW = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*z_scale );
504: TE = param->potentialT * x[j][i+1].T * exp( (j-0.5)*dz*z_scale );
505: if (param->ivisc>=VISC_DISL) { /* olivine dislocation creep */
506: epsN = CalcSecInv(x,i,j, CELL_CORNER,user);
507: epsS = CalcSecInv(x,i,j-1,CELL_CORNER,user);
508: epsE = CalcSecInv(x,i+1,j,CELL_CENTER,user);
509: epsW = CalcSecInv(x,i,j, CELL_CENTER,user);
510: }
511: }
512: etaN = Viscosity(TN,epsN,dz*(j+0.5),param);
513: etaS = Viscosity(TS,epsS,dz*(j-0.5),param);
514: etaW = Viscosity(TW,epsW,dz*j,param);
515: etaE = Viscosity(TE,epsE,dz*j,param);
517: dPdx = ( x[j][i+1].p - x[j][i].p )/dx;
518: if (j==jlim) dudzN = etaN * ( x[j][i].w - x[j][i+1].w )/dx;
519: else dudzN = etaN * ( x[j+1][i].u - x[j][i].u )/dz;
520: dudzS = etaS * ( x[j][i].u - x[j-1][i].u )/dz;
521: dudxE = etaE * ( x[j][i+1].u - x[j][i].u )/dx;
522: dudxW = etaW * ( x[j][i].u - x[j][i-1].u )/dx;
524: residual = -dPdx /* X-MOMENTUM EQUATION*/
525: +( dudxE - dudxW )/dx
526: +( dudzN - dudzS )/dz;
528: if ( param->ivisc!=VISC_CONST ) {
529: dwdxN = etaN * ( x[j][i+1].w - x[j][i].w )/dx;
530: dwdxS = etaS * ( x[j-1][i+1].w - x[j-1][i].w )/dx;
532: residual += ( dudxE - dudxW )/dx + ( dwdxN - dwdxS )/dz;
533: }
535: return residual;
536: }
538: /*---------------------------------------------------------------------*/
541: /* computes the residual of the z-component of eqn (1) above */
542: PetscScalar ZMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
543: /*---------------------------------------------------------------------*/
544: {
545: Parameter *param=user->param;
546: GridInfo *grid =user->grid;
547: PetscScalar dx = grid->dx, dz=grid->dz;
548: 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;
549: PetscScalar TE=0.0,TN=0.0,TS=0.0,TW=0.0, dPdz, residual,z_scale;
550: PetscScalar dudzE,dudzW,dwdxW,dwdxE,dwdzN,dwdzS;
551: PetscInt ilim = grid->ni-1;
553: /* geometric and other parameters */
554: z_scale = param->z_scale;
555:
556: /* viscosity */
557: if ( param->ivisc==VISC_DIFN || param->ivisc>=VISC_DISL ) { /* viscosity is T-dependent */
558: TN = param->potentialT * x[j+1][i].T * exp( (j+0.5)*dz*z_scale );
559: TS = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*z_scale );
560: TW = param->potentialT * TInterp(x,i-1,j) * exp( j *dz*z_scale );
561: if (i==ilim) TE = TW;
562: else TE = param->potentialT * TInterp(x,i,j) * exp( j*dz*z_scale );
563: if (param->ivisc>=VISC_DISL) { /* olivine dislocation creep */
564: epsN = CalcSecInv(x,i,j+1,CELL_CENTER,user);
565: epsS = CalcSecInv(x,i,j, CELL_CENTER,user);
566: epsE = CalcSecInv(x,i,j, CELL_CORNER,user);
567: epsW = CalcSecInv(x,i-1,j,CELL_CORNER,user);
568: }
569: }
570: etaN = Viscosity(TN,epsN,dz*(j+1),param);
571: etaS = Viscosity(TS,epsS,dz*j,param);
572: etaW = Viscosity(TW,epsW,dz*(j+0.5),param);
573: etaE = Viscosity(TE,epsE,dz*(j+0.5),param);
575: dPdz = ( x[j+1][i].p - x[j][i].p )/dz;
576: dwdzN = etaN * ( x[j+1][i].w - x[j][i].w )/dz;
577: dwdzS = etaS * ( x[j][i].w - x[j-1][i].w )/dz;
578: if (i==ilim) dwdxE = etaE * ( x[j][i].u - x[j+1][i].u )/dz;
579: else dwdxE = etaE * ( x[j][i+1].w - x[j][i].w )/dx;
580: dwdxW = 2.0*etaW * ( x[j][i].w - x[j][i-1].w )/dx;
581:
582: /* Z-MOMENTUM */
583: residual = -dPdz /* constant viscosity terms */
584: +( dwdzN - dwdzS )/dz
585: +( dwdxE - dwdxW )/dx;
587: if ( param->ivisc!=VISC_CONST ) {
588: dudzE = etaE * ( x[j+1][i].u - x[j][i].u )/dz;
589: dudzW = etaW * ( x[j+1][i-1].u - x[j][i-1].u )/dz;
591: residual += ( dwdzN - dwdzS )/dz + ( dudzE - dudzW )/dx;
592: }
594: return residual;
595: }
597: /*---------------------------------------------------------------------*/
600: /* computes the residual of eqn (2) above */
601: PetscScalar ContinuityResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
602: /*---------------------------------------------------------------------*/
603: {
604: GridInfo *grid =user->grid;
605: PetscScalar uE,uW,wN,wS,dudx,dwdz;
607: uW = x[j][i-1].u; uE = x[j][i].u; dudx = ( uE - uW )/grid->dx;
608: wS = x[j-1][i].w; wN = x[j][i].w; dwdz = ( wN - wS )/grid->dz;
610: return dudx + dwdz;
611: }
613: /*---------------------------------------------------------------------*/
616: /* computes the residual of eqn (3) above */
617: PetscScalar EnergyResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
618: /*---------------------------------------------------------------------*/
619: {
620: Parameter *param=user->param;
621: GridInfo *grid =user->grid;
622: PetscScalar dx = grid->dx, dz=grid->dz;
623: PetscInt ilim=grid->ni-1, jlim=grid->nj-1, jlid=grid->jlid;
624: PetscScalar TE, TN, TS, TW, residual;
625: PetscScalar uE,uW,wN,wS;
626: PetscScalar fN,fS,fE,fW,dTdxW,dTdxE,dTdzN,dTdzS;
628: dTdzN = ( x[j+1][i].T - x[j][i].T )/dz;
629: dTdzS = ( x[j][i].T - x[j-1][i].T )/dz;
630: dTdxE = ( x[j][i+1].T - x[j][i].T )/dx;
631: dTdxW = ( x[j][i].T - x[j][i-1].T )/dx;
632:
633: residual = ( ( dTdzN - dTdzS )/dz + /* diffusion term */
634: ( dTdxE - dTdxW )/dx )*dx*dz/param->peclet;
636: if (j<=jlid && i>=j) {
637: /* don't advect in the lid */
638: return residual;
640: } else if (i<j) {
641: /* beneath the slab sfc */
642: uW = uE = param->cb;
643: wS = wN = param->sb;
645: } else {
646: /* advect in the slab and wedge */
647: uW = x[j][i-1].u; uE = x[j][i].u;
648: wS = x[j-1][i].w; wN = x[j][i].w;
649: }
651: if ( param->adv_scheme==ADVECT_FV || i==ilim-1 || j==jlim-1 || i==1 || j==1 ) {
652: /* finite volume advection */
653: TS = ( x[j][i].T + x[j-1][i].T )/2.0;
654: TN = ( x[j][i].T + x[j+1][i].T )/2.0;
655: TE = ( x[j][i].T + x[j][i+1].T )/2.0;
656: TW = ( x[j][i].T + x[j][i-1].T )/2.0;
657: fN = wN*TN*dx; fS = wS*TS*dx;
658: fE = uE*TE*dz; fW = uW*TW*dz;
659:
660: } else {
661: /* Fromm advection scheme */
662: 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
663: - fabs(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;
664: 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
665: - fabs(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;
666: 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
667: - fabs(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;
668: 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
669: - fabs(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;
670: }
671:
672: residual -= ( fE - fW + fN - fS );
674: return residual;
675: }
677: /*---------------------------------------------------------------------*/
680: /* computes the shear stress---used on the boundaries */
681: PetscScalar ShearStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
682: /*---------------------------------------------------------------------*/
683: {
684: Parameter *param=user->param;
685: GridInfo *grid=user->grid;
686: PetscInt ilim=grid->ni-1, jlim=grid->nj-1;
687: PetscScalar uN, uS, wE, wW;
689: if ( j<=grid->jlid || i<j || i==ilim || j==jlim ) return EPS_ZERO;
691: if (ipos==CELL_CENTER) { /* on cell center */
693: wE = WInterp(x,i,j-1);
694: if (i==j) { wW = param->sb; uN = param->cb;}
695: else { wW = WInterp(x,i-1,j-1); uN = UInterp(x,i-1,j); }
696: if (j==grid->jlid+1) uS = 0.0;
697: else uS = UInterp(x,i-1,j-1);
699: } else { /* on cell corner */
701: uN = x[j+1][i].u; uS = x[j][i].u;
702: wW = x[j][i].w; wE = x[j][i+1].w;
704: }
706: return (uN-uS)/grid->dz + (wE-wW)/grid->dx;
707: }
709: /*---------------------------------------------------------------------*/
712: /* computes the normal stress---used on the boundaries */
713: PetscScalar XNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
714: /*---------------------------------------------------------------------*/
715: {
716: Parameter *param=user->param;
717: GridInfo *grid =user->grid;
718: PetscScalar dx = grid->dx, dz=grid->dz;
719: PetscInt ilim=grid->ni-1, jlim=grid->nj-1, ivisc;
720: PetscScalar epsC=0.0, etaC, TC, uE, uW, pC, z_scale;
721: if (i<j || j<=grid->jlid) return EPS_ZERO;
723: ivisc=param->ivisc; z_scale = param->z_scale;
725: if (ipos==CELL_CENTER) { /* on cell center */
727: TC = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*z_scale );
728: if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CENTER,user);
729: etaC = Viscosity(TC,epsC,dz*j,param);
731: uW = x[j][i-1].u; uE = x[j][i].u;
732: pC = x[j][i].p;
734: } else { /* on cell corner */
735: if ( i==ilim || j==jlim ) return EPS_ZERO;
737: TC = param->potentialT * TInterp(x,i,j) * exp( j*dz*z_scale );
738: if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CORNER,user);
739: etaC = Viscosity(TC,epsC,dz*(j+0.5),param);
741: if (i==j) uW = param->sb;
742: else uW = UInterp(x,i-1,j);
743: uE = UInterp(x,i,j); pC = PInterp(x,i,j);
744: }
745:
746: return 2.0*etaC*(uE-uW)/dx - pC;
747: }
749: /*---------------------------------------------------------------------*/
752: /* computes the normal stress---used on the boundaries */
753: PetscScalar ZNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
754: /*---------------------------------------------------------------------*/
755: {
756: Parameter *param=user->param;
757: GridInfo *grid =user->grid;
758: PetscScalar dz=grid->dz;
759: PetscInt ilim=grid->ni-1, jlim=grid->nj-1, ivisc;
760: PetscScalar epsC=0.0, etaC, TC;
761: PetscScalar pC, wN, wS, z_scale;
762: if (i<j || j<=grid->jlid) return EPS_ZERO;
764: ivisc=param->ivisc; z_scale = param->z_scale;
766: if (ipos==CELL_CENTER) { /* on cell center */
768: TC = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*z_scale );
769: if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CENTER,user);
770: etaC = Viscosity(TC,epsC,dz*j,param);
771: wN = x[j][i].w; wS = x[j-1][i].w; pC = x[j][i].p;
773: } else { /* on cell corner */
774: if ( (i==ilim) || (j==jlim) ) return EPS_ZERO;
776: TC = param->potentialT * TInterp(x,i,j) * exp( j*dz*z_scale );
777: if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CORNER,user);
778: etaC = Viscosity(TC,epsC,dz*(j+0.5),param);
779: if (i==j) wN = param->sb;
780: else wN = WInterp(x,i,j);
781: wS = WInterp(x,i,j-1); pC = PInterp(x,i,j);
782: }
784: return 2.0*etaC*(wN-wS)/dz - pC;
785: }
787: /*---------------------------------------------------------------------*/
790: /* computes the second invariant of the strain rate tensor */
791: PetscScalar CalcSecInv(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
792: /*---------------------------------------------------------------------*/
793: {
794: Parameter *param = user->param;
795: GridInfo *grid = user->grid;
796: PetscInt ilim=grid->ni-1, jlim=grid->nj-1;
797: PetscScalar uN,uS,uE,uW,wN,wS,wE,wW;
798: PetscScalar eps11, eps12, eps22;
800: if (i<j) return EPS_ZERO;
801: if (i==ilim) i--; if (j==jlim) j--;
803: if (ipos==CELL_CENTER) { /* on cell center */
804: if (j<=grid->jlid) return EPS_ZERO;
806: uE = x[j][i].u; uW = x[j][i-1].u;
807: wN = x[j][i].w; wS = x[j-1][i].w;
808: wE = WInterp(x,i,j-1);
809: if (i==j) { uN = param->cb; wW = param->sb; }
810: else { uN = UInterp(x,i-1,j); wW = WInterp(x,i-1,j-1); }
812: if (j==grid->jlid+1) uS = 0.0;
813: else uS = UInterp(x,i-1,j-1);
815: } else { /* on CELL_CORNER */
816: if (j<grid->jlid) return EPS_ZERO;
818: uN = x[j+1][i].u; uS = x[j][i].u;
819: wE = x[j][i+1].w; wW = x[j][i].w;
820: if (i==j) { wN = param->sb; uW = param->cb; }
821: else { wN = WInterp(x,i,j); uW = UInterp(x,i-1,j); }
823: if (j==grid->jlid) {
824: uE = 0.0; uW = 0.0;
825: uS = -uN;
826: wS = -wN;
827: } else {
828: uE = UInterp(x,i,j);
829: wS = WInterp(x,i,j-1);
830: }
831: }
833: eps11 = (uE-uW)/grid->dx; eps22 = (wN-wS)/grid->dz;
834: eps12 = 0.5*((uN-uS)/grid->dz + (wE-wW)/grid->dx);
836: return sqrt( 0.5*( eps11*eps11 + 2.0*eps12*eps12 + eps22*eps22 ) );
837: }
839: /*---------------------------------------------------------------------*/
842: /* computes the shear viscosity */
843: PetscScalar Viscosity(PetscScalar T, PetscScalar eps, PassiveScalar z,
844: Parameter *param)
845: /*---------------------------------------------------------------------*/
846: {
847: PetscScalar result=0.0;
848: ViscParam difn=param->diffusion, disl=param->dislocation;
849: PetscInt iVisc=param->ivisc;
850: double eps_scale=param->V/(param->L*1000.0);
851: double strain_power, v1, v2, P;
852: double rho_g = 32340.0, R=8.3144;
854: P = rho_g*(z*param->L*1000.0); /* Pa */
856: if (iVisc==VISC_CONST) {
857: /* constant viscosity */
858: return 1.0;
860: } else if (iVisc==VISC_DIFN) {
861: /* diffusion creep rheology */
862: result = difn.A*PetscExpScalar((difn.Estar + P*difn.Vstar)/R/(T+273.0))/param->eta0;
864: } else if (iVisc==VISC_DISL) {
865: /* dislocation creep rheology */
866: strain_power = pow( eps*eps_scale, (1.0-disl.n)/disl.n );
867: result = disl.A*PetscExpScalar((disl.Estar + P*disl.Vstar)/disl.n/R/(T+273.0))*strain_power/param->eta0;
869: } else if (iVisc==VISC_FULL) {
870: /* dislocation/diffusion creep rheology */
871: strain_power = pow( eps*eps_scale, (1.0-disl.n)/disl.n );
872: v1 = difn.A*PetscExpScalar((difn.Estar + P*difn.Vstar)/R/(T+273.0))/param->eta0;
873: v2 = disl.A*PetscExpScalar((disl.Estar + P*disl.Vstar)/disl.n/R/(T+273.0))*strain_power/param->eta0;
874: result = 1.0/(1.0/v1 + 1.0/v2);
875: }
877: /* max viscosity is param->eta0 */
878: result = PetscMin( result, 1.0 );
879: /* min viscosity is param->visc_cutoff */
880: result = PetscMax( result, param->visc_cutoff );
881: /* continuation method */
882: result = pow(result,param->continuation);
883: return result;
884: }
886: /*=====================================================================
887: INITIALIZATION, POST-PROCESSING AND OUTPUT FUNCTIONS
888: =====================================================================*/
890: /*---------------------------------------------------------------------*/
893: /* initializes the problem parameters and checks for
894: command line changes */
895: PetscErrorCode SetParams(Parameter *param, GridInfo *grid)
896: /*---------------------------------------------------------------------*/
897: {
898: PetscErrorCode ierr, ierr_out=0;
899: PetscReal SEC_PER_YR = 3600.00*24.00*365.2500;
900: PetscReal PI = 3.14159265358979323846;
901: PetscReal alpha_g_on_cp_units_inverse_km=4.0e-5*9.8;
902:
903: /* domain geometry */
904: param->slab_dip = 45.0;
905: param->width = 320.0; /* km */
906: param->depth = 300.0; /* km */
907: param->lid_depth = 35.0; /* km */
908: param->fault_depth = 35.0; /* km */
909: PetscOptionsGetReal(PETSC_NULL,"-slab_dip",&(param->slab_dip),PETSC_NULL);
910: PetscOptionsGetReal(PETSC_NULL,"-width",&(param->width),PETSC_NULL);
911: PetscOptionsGetReal(PETSC_NULL,"-depth",&(param->depth),PETSC_NULL);
912: PetscOptionsGetReal(PETSC_NULL,"-lid_depth",&(param->lid_depth),PETSC_NULL);
913: PetscOptionsGetReal(PETSC_NULL,"-fault_depth",&(param->fault_depth),PETSC_NULL);
914: param->slab_dip = param->slab_dip*PI/180.0; /* radians */
916: /* grid information */
917: PetscOptionsGetInt(PETSC_NULL, "-jfault",&(grid->jfault),PETSC_NULL);
918: grid->ni = 82;
919: PetscOptionsGetInt(PETSC_NULL, "-ni",&(grid->ni),PETSC_NULL);
920: grid->dx = param->width/((double)(grid->ni-2)); /* km */
921: grid->dz = grid->dx*tan(param->slab_dip); /* km */
922: grid->nj = (PetscInt)(param->depth/grid->dz + 3.0); /* gridpoints*/
923: param->depth = grid->dz*(grid->nj-2); /* km */
924: grid->inose = 0; /* gridpoints*/
925: PetscOptionsGetInt(PETSC_NULL,"-inose",&(grid->inose),PETSC_NULL);
926: grid->periodic = DA_NONPERIODIC;
927: grid->stencil = DA_STENCIL_BOX;
928: grid->dof = 4;
929: grid->stencil_width = 2;
930: grid->mglevels = 1;
932: /* boundary conditions */
933: param->pv_analytic = PETSC_FALSE;
934: param->ibound = BC_NOSTRESS;
935: PetscOptionsGetInt(PETSC_NULL,"-ibound",&(param->ibound),PETSC_NULL);
937: /* physical constants */
938: param->slab_velocity = 5.0; /* cm/yr */
939: param->slab_age = 50.0; /* Ma */
940: param->lid_age = 50.0; /* Ma */
941: param->kappa = 0.7272e-6; /* m^2/sec */
942: param->potentialT = 1300.0; /* degrees C */
943: PetscOptionsGetReal(PETSC_NULL,"-slab_velocity",&(param->slab_velocity),PETSC_NULL);
944: PetscOptionsGetReal(PETSC_NULL,"-slab_age",&(param->slab_age),PETSC_NULL);
945: PetscOptionsGetReal(PETSC_NULL,"-lid_age",&(param->lid_age),PETSC_NULL);
946: PetscOptionsGetReal(PETSC_NULL,"-kappa",&(param->kappa),PETSC_NULL);
947: PetscOptionsGetReal(PETSC_NULL,"-potentialT",&(param->potentialT),PETSC_NULL);
949: /* viscosity */
950: param->ivisc = 3; /* 0=isovisc, 1=difn creep, 2=disl creep, 3=full */
951: param->eta0 = 1e24; /* Pa-s */
952: param->visc_cutoff = 0.0; /* factor of eta_0 */
953: param->continuation = 1.0;
954: /* constants for diffusion creep */
955: param->diffusion.A = 1.8e7; /* Pa-s */
956: param->diffusion.n = 1.0; /* dim'less */
957: param->diffusion.Estar = 375e3; /* J/mol */
958: param->diffusion.Vstar = 5e-6; /* m^3/mol */
959: /* constants for param->dislocationocation creep */
960: param->dislocation.A = 2.8969e4; /* Pa-s */
961: param->dislocation.n = 3.5; /* dim'less */
962: param->dislocation.Estar = 530e3; /* J/mol */
963: param->dislocation.Vstar = 14e-6; /* m^3/mol */
964: PetscOptionsGetInt(PETSC_NULL, "-ivisc",&(param->ivisc),PETSC_NULL);
965: PetscOptionsGetReal(PETSC_NULL,"-visc_cutoff",&(param->visc_cutoff),PETSC_NULL);
966: param->output_ivisc = param->ivisc;
967: PetscOptionsGetInt(PETSC_NULL,"-output_ivisc",&(param->output_ivisc),PETSC_NULL);
968: PetscOptionsGetReal(PETSC_NULL,"-vstar",&(param->dislocation.Vstar),PETSC_NULL);
970: /* output options */
971: param->quiet = PETSC_FALSE;
972: param->param_test = PETSC_FALSE;
973: PetscOptionsHasName(PETSC_NULL,"-quiet",&(param->quiet));
974: PetscOptionsHasName(PETSC_NULL,"-test",&(param->param_test));
975: PetscOptionsGetString(PETSC_NULL,"-file",param->filename,PETSC_MAX_PATH_LEN,&(param->output_to_file));
977: /* advection */
978: param->adv_scheme = ADVECT_FROMM; /* advection scheme: 0=finite vol, 1=Fromm */
979: PetscOptionsGetInt(PETSC_NULL,"-adv_scheme",&(param->adv_scheme),PETSC_NULL);
981: /* misc. flags */
982: param->stop_solve = PETSC_FALSE;
983: param->interrupted = PETSC_FALSE;
984: param->kspmon = PETSC_FALSE;
985: param->toggle_kspmon = PETSC_FALSE;
987: /* derived parameters for slab angle */
988: param->sb = sin(param->slab_dip);
989: param->cb = cos(param->slab_dip);
990: param->c = param->slab_dip*param->sb/(param->slab_dip*param->slab_dip-param->sb*param->sb);
991: param->d = (param->slab_dip*param->cb-param->sb)/(param->slab_dip*param->slab_dip-param->sb*param->sb);
993: /* length, velocity and time scale for non-dimensionalization */
994: param->L = PetscMin(param->width,param->depth); /* km */
995: param->V = param->slab_velocity/100.0/SEC_PER_YR; /* m/sec */
997: /* other unit conversions and derived parameters */
998: param->scaled_width = param->width/param->L; /* dim'less */
999: param->scaled_depth = param->depth/param->L; /* dim'less */
1000: param->lid_depth = param->lid_depth/param->L; /* dim'less */
1001: param->fault_depth = param->fault_depth/param->L; /* dim'less */
1002: grid->dx = grid->dx/param->L; /* dim'less */
1003: grid->dz = grid->dz/param->L; /* dim'less */
1004: grid->jlid = (PetscInt)(param->lid_depth/grid->dz); /* gridcells */
1005: grid->jfault = (PetscInt)(param->fault_depth/grid->dz); /* gridcells */
1006: param->lid_depth = grid->jlid*grid->dz; /* dim'less */
1007: param->fault_depth = grid->jfault*grid->dz; /* dim'less */
1008: grid->corner = grid->jlid+1; /* gridcells */
1009: param->peclet = param->V /* m/sec */
1010: * param->L*1000.0 /* m */
1011: / param->kappa; /* m^2/sec */
1012: param->z_scale = param->L * alpha_g_on_cp_units_inverse_km;
1013: param->skt = sqrt(param->kappa*param->slab_age*SEC_PER_YR);
1014: PetscOptionsGetReal(PETSC_NULL,"-peclet",&(param->peclet),PETSC_NULL);
1015:
1016: return ierr_out;
1017: }
1019: /*---------------------------------------------------------------------*/
1022: /* prints a report of the problem parameters to stdout */
1023: PetscErrorCode ReportParams(Parameter *param, GridInfo *grid)
1024: /*---------------------------------------------------------------------*/
1025: {
1026: PetscErrorCode ierr, ierr_out=0;
1027: char date[30];
1028: PetscReal PI = 3.14159265358979323846;
1030: PetscGetDate(date,30);
1032: if ( !(param->quiet) ) {
1033: PetscPrintf(PETSC_COMM_WORLD,"---------------------BEGIN ex30 PARAM REPORT-------------------\n");
1034: PetscPrintf(PETSC_COMM_WORLD," %s\n",&(date[0]));
1036: PetscPrintf(PETSC_COMM_WORLD,"Domain: \n");
1037: PetscPrintf(PETSC_COMM_WORLD," Width = %g km, Depth = %g km\n",param->width,param->depth);
1038: PetscPrintf(PETSC_COMM_WORLD," Slab dip = %g degrees, Slab velocity = %g cm/yr\n",param->slab_dip*180.0/PI,param->slab_velocity);
1039: PetscPrintf(PETSC_COMM_WORLD," Lid depth = %5.2f km, Fault depth = %5.2f km\n",param->lid_depth*param->L,param->fault_depth*param->L);
1041: PetscPrintf(PETSC_COMM_WORLD,"\nGrid: \n");
1042: PetscPrintf(PETSC_COMM_WORLD," [ni,nj] = %D, %D [dx,dz] = %g, %g km\n",grid->ni,grid->nj,grid->dx*param->L,grid->dz*param->L);
1043: PetscPrintf(PETSC_COMM_WORLD," jlid = %3D jfault = %3D \n",grid->jlid,grid->jfault);
1044: PetscPrintf(PETSC_COMM_WORLD," Pe = %g\n",param->peclet);
1046: PetscPrintf(PETSC_COMM_WORLD,"\nRheology:");
1047: if (param->ivisc==VISC_CONST) {
1048: PetscPrintf(PETSC_COMM_WORLD," Isoviscous \n");
1049: if (param->pv_analytic)
1050: PetscPrintf(PETSC_COMM_WORLD," Pressure and Velocity prescribed! \n");
1051: } else if (param->ivisc==VISC_DIFN) {
1052: PetscPrintf(PETSC_COMM_WORLD," Diffusion Creep (T-Dependent Newtonian) \n");
1053: PetscPrintf(PETSC_COMM_WORLD," Viscosity range: %g--%g Pa-sec \n",param->eta0,param->visc_cutoff*param->eta0);
1054: } else if (param->ivisc==VISC_DISL ) {
1055: PetscPrintf(PETSC_COMM_WORLD," Dislocation Creep (T-Dependent Non-Newtonian) \n");
1056: PetscPrintf(PETSC_COMM_WORLD," Viscosity range: %g--%g Pa-sec \n",param->eta0,param->visc_cutoff*param->eta0);
1057: } else if (param->ivisc==VISC_FULL ) {
1058: PetscPrintf(PETSC_COMM_WORLD," Full Rheology \n");
1059: PetscPrintf(PETSC_COMM_WORLD," Viscosity range: %g--%g Pa-sec \n",param->eta0,param->visc_cutoff*param->eta0);
1060: } else {
1061: PetscPrintf(PETSC_COMM_WORLD," Invalid! \n");
1062: ierr_out=1;
1063: }
1065: PetscPrintf(PETSC_COMM_WORLD,"Boundary condition:");
1066: if ( param->ibound==BC_ANALYTIC ) {
1067: PetscPrintf(PETSC_COMM_WORLD," Isoviscous Analytic Dirichlet \n");
1068: } else if ( param->ibound==BC_NOSTRESS ) {
1069: PetscPrintf(PETSC_COMM_WORLD," Stress-Free (normal & shear stress)\n");
1070: } else if ( param->ibound==BC_EXPERMNT ) {
1071: PetscPrintf(PETSC_COMM_WORLD," Experimental boundary condition \n");
1072: } else {
1073: PetscPrintf(PETSC_COMM_WORLD," Invalid! \n");
1074: ierr_out=1;
1075: }
1077: if (param->output_to_file)
1078: #if defined(PETSC_HAVE_MATLAB) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
1079: PetscPrintf(PETSC_COMM_WORLD,"Output Destination: Mat file \"%s\"\n",param->filename);
1080: #else
1081: PetscPrintf(PETSC_COMM_WORLD,"Output Destination: PETSc binary file \"%s\"\n",param->filename);
1082: #endif
1083: if ( param->output_ivisc != param->ivisc )
1084: PetscPrintf(PETSC_COMM_WORLD," Output viscosity: -ivisc %D\n",param->output_ivisc);
1086: PetscPrintf(PETSC_COMM_WORLD,"---------------------END ex30 PARAM REPORT---------------------\n");
1087: }
1088: if ( param->param_test ) PetscEnd();
1089: return ierr_out;
1090: }
1092: /* ------------------------------------------------------------------- */
1095: /* generates an inital guess using the analytic solution for isoviscous
1096: corner flow */
1097: PetscErrorCode Initialize(DMMG *dmmg)
1098: /* ------------------------------------------------------------------- */
1099: {
1100: AppCtx *user = (AppCtx*)dmmg[0]->user;
1101: Parameter *param = user->param;
1102: GridInfo *grid = user->grid;
1103: DA da;
1104: PetscInt i,j,is,js,im,jm;
1106: Field **x;
1108: /* Get the fine grid */
1109: da = (DA)(dmmg[0]->dm);
1110: DAGetCorners(da,&is,&js,PETSC_NULL,&im,&jm,PETSC_NULL);
1111: DAVecGetArray(da,((AppCtx*)dmmg[0]->user)->Xguess,(void**)&x);
1113: /* Compute initial guess */
1114: for (j=js; j<js+jm; j++) {
1115: for (i=is; i<is+im; i++) {
1116: if (i<j) {
1117: x[j][i].u = param->cb;
1118: } else if (j<=grid->jlid) {
1119: x[j][i].u = 0.0;
1120: } else {
1121: x[j][i].u = HorizVelocity(i,j,user);
1122: }
1123: if (i<=j) {
1124: x[j][i].w = param->sb;
1125: } else if (j<=grid->jlid) {
1126: x[j][i].w = 0.0;
1127: } else {
1128: x[j][i].w = VertVelocity(i,j,user);
1129: }
1130: if (i<j || j<=grid->jlid) {
1131: x[j][i].p = 0.0;
1132: } else {
1133: x[j][i].p = Pressure(i,j,user);
1134: }
1135: x[j][i].T = PetscMin(grid->dz*(j-0.5),1.0);
1136: }
1137: }
1139: /* Restore x to Xguess */
1140: DAVecRestoreArray(da,((AppCtx*)dmmg[0]->user)->Xguess,(void**)&x);
1142: return 0;
1143: }
1145: /*---------------------------------------------------------------------*/
1148: /* controls output to a file */
1149: PetscErrorCode DoOutput(DMMG *dmmg, PetscInt its)
1150: /*---------------------------------------------------------------------*/
1151: {
1152: AppCtx *user = (AppCtx*)dmmg[0]->user;
1153: Parameter *param = user->param;
1154: GridInfo *grid = user->grid;
1156: PetscMPIInt rank;
1157: PetscInt ivt=param->ivisc;
1158: PetscViewer viewer;
1159: Vec res, pars;
1160: MPI_Comm comm;
1162: param->ivisc = param->output_ivisc;
1164: /* compute final residual and final viscosity/strain rate fields */
1165: SNESGetFunction(DMMGGetSNES(dmmg), &res, PETSC_NULL, PETSC_NULL);
1166: ViscosityField(DMMGGetDMMG(dmmg), DMMGGetx(dmmg), ((AppCtx *)dmmg[0]->user)->Xguess);
1168: /* get the communicator and the rank of the processor */
1169: PetscObjectGetComm((PetscObject)DMMGGetSNES(dmmg), &comm);
1170: MPI_Comm_rank(comm, &rank);
1172: if (param->output_to_file) { /* send output to binary file */
1173: VecCreate(comm, &pars);
1174: if (rank == 0) { /* on processor 0 */
1175: VecSetSizes(pars, 20, PETSC_DETERMINE);
1176: VecSetFromOptions(pars);
1177: VecSetValue(pars,0, (PetscScalar)(grid->ni),INSERT_VALUES);
1178: VecSetValue(pars,1, (PetscScalar)(grid->nj),INSERT_VALUES);
1179: VecSetValue(pars,2, (PetscScalar)(grid->dx),INSERT_VALUES);
1180: VecSetValue(pars,3, (PetscScalar)(grid->dz),INSERT_VALUES);
1181: VecSetValue(pars,4, (PetscScalar)(param->L),INSERT_VALUES);
1182: VecSetValue(pars,5, (PetscScalar)(param->V),INSERT_VALUES);
1183: /* skipped 6 intentionally */
1184: VecSetValue(pars,7, (PetscScalar)(param->slab_dip),INSERT_VALUES);
1185: VecSetValue(pars,8, (PetscScalar)(grid->jlid),INSERT_VALUES);
1186: VecSetValue(pars,9, (PetscScalar)(param->lid_depth),INSERT_VALUES);
1187: VecSetValue(pars,10,(PetscScalar)(grid->jfault),INSERT_VALUES);
1188: VecSetValue(pars,11,(PetscScalar)(param->fault_depth),INSERT_VALUES);
1189: VecSetValue(pars,12,(PetscScalar)(param->potentialT),INSERT_VALUES);
1190: VecSetValue(pars,13,(PetscScalar)(param->ivisc),INSERT_VALUES);
1191: VecSetValue(pars,14,(PetscScalar)(param->visc_cutoff),INSERT_VALUES);
1192: VecSetValue(pars,15,(PetscScalar)(param->ibound),INSERT_VALUES);
1193: VecSetValue(pars,16,(PetscScalar)(its),INSERT_VALUES);
1194: } else { /* on some other processor */
1195: VecSetSizes(pars, 0, PETSC_DETERMINE);
1196: VecSetFromOptions(pars);
1197: }
1198: VecAssemblyBegin(pars); VecAssemblyEnd(pars);
1200: /* create viewer */
1201: #if defined(PETSC_HAVE_MATLAB) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
1202: PetscViewerMatlabOpen(PETSC_COMM_WORLD,param->filename,PETSC_FILE_CREATE,&viewer);
1203: #else
1204: PetscViewerBinaryOpen(PETSC_COMM_WORLD,param->filename,PETSC_FILE_CREATE,&viewer);
1205: #endif
1207: /* send vectors to viewer */
1208: PetscObjectSetName((PetscObject)res,"res");
1209: VecView(res,viewer);
1210: PetscObjectSetName((PetscObject)DMMGGetx(dmmg),"out");
1211: VecView(DMMGGetx(dmmg), viewer);
1212: PetscObjectSetName((PetscObject)(user->Xguess),"aux");
1213: VecView(user->Xguess, viewer);
1214: StressField(dmmg); /* compute stress fields */
1215: PetscObjectSetName((PetscObject)(user->Xguess),"str");
1216: VecView(user->Xguess, viewer);
1217: PetscObjectSetName((PetscObject)pars,"par");
1218: VecView(pars, viewer);
1219:
1220: /* destroy viewer and vector */
1221: PetscViewerDestroy(viewer);
1222: VecDestroy(pars);
1223: }
1224:
1225: param->ivisc = ivt;
1226: return 0;
1227: }
1229: /* ------------------------------------------------------------------- */
1232: /* Compute both the second invariant of the strain rate tensor and the viscosity, at both cell centers and cell corners */
1233: PetscErrorCode ViscosityField(DMMG dmmg, Vec X, Vec V)
1234: /* ------------------------------------------------------------------- */
1235: {
1236: DA da = (DA) dmmg->dm;
1237: AppCtx *user = (AppCtx *) dmmg->user;
1238: Parameter *param = user->param;
1239: GridInfo *grid = user->grid;
1240: Vec localX;
1241: Field **v, **x;
1242: PassiveReal eps, dx, dz, T, epsC, TC;
1243: PetscInt i,j,is,js,im,jm,ilim,jlim,ivt;
1247: ivt = param->ivisc;
1248: param->ivisc = param->output_ivisc;
1250: DACreateLocalVector(da, &localX);
1251: DAGlobalToLocalBegin(da, X, INSERT_VALUES, localX);
1252: DAGlobalToLocalEnd(da, X, INSERT_VALUES, localX);
1253: DAVecGetArray(da,localX,(void**)&x);
1254: DAVecGetArray(da,V,(void**)&v);
1256: /* Parameters */
1257: dx = grid->dx; dz = grid->dz;
1258: ilim = grid->ni-1; jlim = grid->nj-1;
1260: /* Compute real temperature, strain rate and viscosity */
1261: DAGetCorners(da,&is,&js,PETSC_NULL,&im,&jm,PETSC_NULL);
1262: for (j=js; j<js+jm; j++) {
1263: for (i=is; i<is+im; i++) {
1264: T = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*param->z_scale );
1265: if (i<ilim && j<jlim) {
1266: TC = param->potentialT * TInterp(x,i,j) * exp( j*dz*param->z_scale );
1267: } else {
1268: TC = T;
1269: }
1270: eps = CalcSecInv(x,i,j,CELL_CENTER,user);
1271: epsC = CalcSecInv(x,i,j,CELL_CORNER,user);
1272: v[j][i].u = eps;
1273: v[j][i].w = epsC;
1274: v[j][i].p = Viscosity(T,eps,dz*(j-0.5),param);
1275: v[j][i].T = Viscosity(TC,epsC,dz*j,param);
1276: }
1277: }
1278: DAVecRestoreArray(da,V,(void**)&v);
1279: DAVecRestoreArray(da,localX,(void**)&x);
1280: param->ivisc = ivt;
1281: return(0);
1282: }
1284: /* ------------------------------------------------------------------- */
1287: /* post-processing: compute stress everywhere */
1288: PetscErrorCode StressField(DMMG *dmmg)
1289: /* ------------------------------------------------------------------- */
1290: {
1291: AppCtx *user = (AppCtx*)dmmg[0]->user;
1292: PetscInt i,j,is,js,im,jm;
1294: DA da;
1295: Vec locVec;
1296: Field **x, **y;
1298: /* Get the fine grid of Xguess and X */
1299: da = (DA)(dmmg[0]->dm);
1300: DAGetCorners(da,&is,&js,PETSC_NULL,&im,&jm,PETSC_NULL);
1301: DAVecGetArray(da,((AppCtx*)dmmg[0]->user)->Xguess,(void**)&x);
1303: DACreateLocalVector(da, &locVec);
1304: DAGlobalToLocalBegin(da, DMMGGetx(dmmg), INSERT_VALUES, locVec);
1305: DAGlobalToLocalEnd(da, DMMGGetx(dmmg), INSERT_VALUES, locVec);
1306: DAVecGetArray(da,locVec,(void**)&y);
1308: /* Compute stress on the corner points */
1309: for (j=js; j<js+jm; j++) {
1310: for (i=is; i<is+im; i++) {
1311:
1312: x[j][i].u = ShearStress(y,i,j,CELL_CENTER,user);
1313: x[j][i].w = ShearStress(y,i,j,CELL_CORNER,user);
1314: x[j][i].p = XNormalStress(y,i,j,CELL_CENTER,user);
1315: x[j][i].T = ZNormalStress(y,i,j,CELL_CENTER,user);
1316: }
1317: }
1319: /* Restore the fine grid of Xguess and X */
1320: DAVecRestoreArray(da,((AppCtx*)dmmg[0]->user)->Xguess,(void**)&x);
1321: DAVecRestoreArray(da,locVec,(void**)&y);
1323: return 0;
1324: }
1326: /*=====================================================================
1327: UTILITY FUNCTIONS
1328: =====================================================================*/
1330: /*---------------------------------------------------------------------*/
1333: /* returns the velocity of the subducting slab and handles fault nodes
1334: for BC */
1335: PassiveScalar SlabVel(char c, PetscInt i, PetscInt j, AppCtx *user)
1336: /*---------------------------------------------------------------------*/
1337: {
1338: Parameter *param = user->param;
1339: GridInfo *grid = user->grid;
1341: if (c=='U' || c=='u') {
1342: if (i<j-1) {
1343: return param->cb;
1344: } else if (j<=grid->jfault) {
1345: return 0.0;
1346: } else
1347: return param->cb;
1349: } else {
1350: if (i<j) {
1351: return param->sb;
1352: } else if (j<=grid->jfault) {
1353: return 0.0;
1354: } else
1355: return param->sb;
1356: }
1357: }
1359: /*---------------------------------------------------------------------*/
1362: /* solution to diffusive half-space cooling model for BC */
1363: PassiveScalar PlateModel(PetscInt j, PetscInt plate, AppCtx *user)
1364: /*---------------------------------------------------------------------*/
1365: {
1366: Parameter *param = user->param;
1367: PassiveScalar z;
1368: if (plate==PLATE_LID)
1369: z = (j-0.5)*user->grid->dz;
1370: else /* PLATE_SLAB */
1371: z = (j-0.5)*user->grid->dz*param->cb;
1373: return erf(z*param->L/2.0/param->skt);
1374: }
1376: /*---------------------------------------------------------------------*/
1379: PetscScalar UInterp(Field **x, PetscInt i, PetscInt j)
1380: /*---------------------------------------------------------------------*/
1381: {
1382: return 0.25*(x[j][i].u+x[j+1][i].u+x[j][i+1].u+x[j+1][i+1].u);
1383: }
1385: /*---------------------------------------------------------------------*/
1388: PetscScalar WInterp(Field **x, PetscInt i, PetscInt j)
1389: /*---------------------------------------------------------------------*/
1390: {
1391: return 0.25*(x[j][i].w+x[j+1][i].w+x[j][i+1].w+x[j+1][i+1].w);
1392: }
1394: /*---------------------------------------------------------------------*/
1397: PetscScalar PInterp(Field **x, PetscInt i, PetscInt j)
1398: /*---------------------------------------------------------------------*/
1399: {
1400: return 0.25*(x[j][i].p+x[j+1][i].p+x[j][i+1].p+x[j+1][i+1].p);
1401: }
1403: /*---------------------------------------------------------------------*/
1406: PetscScalar TInterp(Field **x, PetscInt i, PetscInt j)
1407: /*---------------------------------------------------------------------*/
1408: {
1409: return 0.25*(x[j][i].T+x[j+1][i].T+x[j][i+1].T+x[j+1][i+1].T);
1410: }
1412: /*---------------------------------------------------------------------*/
1415: /* isoviscous analytic solution for IC */
1416: PassiveScalar HorizVelocity(PetscInt i, PetscInt j, AppCtx *user)
1417: /*---------------------------------------------------------------------*/
1418: {
1419: Parameter *param = user->param;
1420: GridInfo *grid = user->grid;
1421: PetscScalar x, z, r, st, ct, th, c=param->c, d=param->d;
1422:
1423: x = (i - grid->jlid)*grid->dx; z = (j - grid->jlid - 0.5)*grid->dz;
1424: r = sqrt(x*x+z*z); st = z/r; ct = x/r; th = atan(z/x);
1425: return ct*(c*th*st+d*(st+th*ct)) + st*(c*(st-th*ct)+d*th*st);
1426: }
1428: /*---------------------------------------------------------------------*/
1431: /* isoviscous analytic solution for IC */
1432: PetscScalar VertVelocity(PetscInt i, PetscInt j, AppCtx *user)
1433: /*---------------------------------------------------------------------*/
1434: {
1435: Parameter *param = user->param;
1436: GridInfo *grid = user->grid;
1437: PetscScalar x, z, r, st, ct, th, c=param->c, d=param->d;
1438:
1439: x = (i - grid->jlid - 0.5)*grid->dx; z = (j - grid->jlid)*grid->dz;
1440: r = sqrt(x*x+z*z); st = z/r; ct = x/r; th = atan(z/x);
1441: return st*(c*th*st+d*(st+th*ct)) - ct*(c*(st-th*ct)+d*th*st);
1442: }
1444: /*---------------------------------------------------------------------*/
1447: /* isoviscous analytic solution for IC */
1448: PetscScalar Pressure(PetscInt i, PetscInt j, AppCtx *user)
1449: /*---------------------------------------------------------------------*/
1450: {
1451: Parameter *param = user->param;
1452: GridInfo *grid = user->grid;
1453: PetscScalar x, z, r, st, ct, c=param->c, d=param->d;
1455: x = (i - grid->jlid - 0.5)*grid->dx; z = (j - grid->jlid - 0.5)*grid->dz;
1456: r = sqrt(x*x+z*z); st = z/r; ct = x/r;
1457: return (-2.0*(c*ct-d*st)/r);
1458: }
1460: /* ------------------------------------------------------------------- */
1463: /* utility function */
1464: PetscTruth OptionsHasName(const char pre[],const char name[])
1465: /* ------------------------------------------------------------------- */
1466: {
1467: PetscTruth retval;
1469: PetscOptionsHasName(pre,name,&retval);//
1470: return retval;
1471: }
1473: /*=====================================================================
1474: INTERACTIVE SIGNAL HANDLING
1475: =====================================================================*/
1477: /* ------------------------------------------------------------------- */
1480: PetscErrorCode SNESConverged_Interactive(SNES snes, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *ctx)
1481: /* ------------------------------------------------------------------- */
1482: {
1483: AppCtx *user = (AppCtx *) ctx;
1484: Parameter *param = user->param;
1485: KSP ksp;
1489: if (param->interrupted) {
1490: param->interrupted = PETSC_FALSE;
1491: PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: exiting SNES solve. \n");
1492: *reason = SNES_CONVERGED_FNORM_ABS;
1493: return(0);
1494: } else if (param->toggle_kspmon) {
1495: param->toggle_kspmon = PETSC_FALSE;
1496: SNESGetKSP(snes, &ksp);
1497: if (param->kspmon) {
1498: KSPClearMonitor(ksp);
1499: param->kspmon = PETSC_FALSE;
1500: PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: deactivating ksp singular value monitor. \n");
1501: } else {
1502: KSPSetMonitor(ksp,KSPSingularValueMonitor,PETSC_NULL,PETSC_NULL);
1503: param->kspmon = PETSC_TRUE;
1504: PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: activating ksp singular value monitor. \n");
1505: }
1506: }
1507: PetscFunctionReturn(SNESConverged_LS(snes, xnorm, pnorm, fnorm, reason, ctx));
1508: }
1510: /* ------------------------------------------------------------------- */
1511: #include <signal.h>
1514: PetscErrorCode InteractiveHandler(int signum, void *ctx)
1515: /* ------------------------------------------------------------------- */
1516: {
1517: AppCtx *user = (AppCtx *) ctx;
1518: Parameter *param = user->param;
1520: if (signum == SIGILL) {
1521: param->toggle_kspmon = PETSC_TRUE;
1522: } else if (signum == SIGCONT) {
1523: param->interrupted = PETSC_TRUE;
1524: } else if (signum == SIGURG) {
1525: param->stop_solve = PETSC_TRUE;
1526: }
1527: return 0;
1528: }