Actual source code: ex30.c
petsc-3.3-p2 2012-07-13
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";
28: /*F-----------------------------------------------------------------------
This PETSc 2.2.0 example by Richard F. Katz
http://www.ldeo.columbia.edu/~katz/
The problem is modeled by the partial differential equation system
\begin{eqnarray}
-\nabla P + \nabla \cdot [\eta (\nabla v + \nabla v^T)] & = & 0 \\
\nabla \cdot v & = & 0 \\
dT/dt + \nabla \cdot (vT) - 1/Pe \triangle^2(T) & = & 0 \\
\end{eqnarray}
\begin{eqnarray}
\eta(T,Eps\_dot) & = & \hbox{constant } \hbox{if ivisc} ==0 \\
& = & \hbox{diffusion creep (T,P-dependent) } \hbox{if ivisc} ==1 \\
& = & \hbox{dislocation creep (T,P,v-dependent)} \hbox{if ivisc} ==2 \\
& = & \hbox{mantle viscosity (difn and disl) } \hbox{if ivisc} ==3
\end{eqnarray}
which is uniformly discretized on a staggered mesh:
-------$w_{ij}$------
$u_{i-1j}$ $P,T_{ij}$ $u_{ij}$
------$w_{ij-1}$-----
53: ------------------------------------------------------------------------F*/
55: #include <petscsnes.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 miscelaneous 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: DMDABoundaryType 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 *);
122: extern PetscBool OptionsHasName(const char pre[],const char name[]);
124: /*-----------------------------------------------------------------------*/
127: int main(int argc,char **argv)
128: /*-----------------------------------------------------------------------*/
129: {
130: SNES snes;
131: AppCtx *user; /* user-defined work context */
132: Parameter param;
133: GridInfo grid;
134: PetscInt nits;
136: MPI_Comm comm;
137: DM da;
139: PetscInitialize(&argc,&argv,(char *)0,help);
140: PetscOptionsSetValue("-file","ex30_output");
141: PetscOptionsSetValue("-snes_monitor_short",PETSC_NULL);
142: PetscOptionsSetValue("-snes_max_it","20");
143: PetscOptionsSetValue("-ksp_max_it","1500");
144: PetscOptionsSetValue("-ksp_gmres_restart","300");
145: PetscOptionsInsert(&argc,&argv,PETSC_NULL);
147: comm = PETSC_COMM_WORLD;
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Set up the problem parameters.
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152: SetParams(¶m,&grid);
153: ReportParams(¶m,&grid);
155: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157: SNESCreate(comm,&snes);
158: DMDACreate2d(comm,grid.bx,grid.by,grid.stencil,grid.ni,grid.nj,PETSC_DECIDE,PETSC_DECIDE,grid.dof,grid.stencil_width,0,0,&da);
159: SNESSetDM(snes,da);
160: DMDASetFieldName(da,0,"x-velocity");
161: DMDASetFieldName(da,1,"y-velocity");
162: DMDASetFieldName(da,2,"pressure");
163: DMDASetFieldName(da,3,"temperature");
166: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: Create user context, set problem data, create vector data structures.
168: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169: PetscMalloc(sizeof(AppCtx),&user);
170: user->param = ¶m;
171: user->grid = &grid;
172: DMSetApplicationContext(da,user);
173: DMCreateGlobalVector(da,&(user->Xguess));
176: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177: Set up the SNES solver with callback functions.
178: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
179: DMDASetLocalFunction(da,(DMDALocalFunction1)FormFunctionLocal);
180: SNESSetFromOptions(snes);
183: SNESSetConvergenceTest(snes,SNESConverged_Interactive,(void*)user,PETSC_NULL);
184: PetscPushSignalHandler(InteractiveHandler,(void*)user);
185:
186: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187: Initialize and solve the nonlinear system
188: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189: Initialize(da);
190: UpdateSolution(snes,user,&nits);
191:
192: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193: Output variables.
194: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195: DoOutput(snes,nits);
196:
197: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198: Free work space.
199: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200: VecDestroy(&user->Xguess);
201: VecDestroy(&user->x);
202: PetscFree(user);
203: SNESDestroy(&snes);
204: DMDestroy(&da);
205: PetscPopSignalHandler();
206: PetscFinalize();
207: return 0;
208: }
210: /*=====================================================================
211: PETSc INTERACTION FUNCTIONS (initialize & call SNESSolve)
212: =====================================================================*/
214: /*---------------------------------------------------------------------*/
217: /* manages solve: adaptive continuation method */
218: PetscErrorCode UpdateSolution(SNES snes, AppCtx *user, PetscInt *nits)
219: {
220: KSP ksp;
221: PC pc;
222: SNESConvergedReason reason;
223: Parameter *param = user->param;
224: PetscReal cont_incr=0.3;
225: PetscInt its;
226: PetscErrorCode ierr;
227: PetscBool q = PETSC_FALSE;
228: DM dm;
231: SNESGetDM(snes,&dm);
232: DMCreateGlobalVector(dm,&user->x);
233: SNESGetKSP(snes,&ksp);
234: KSPGetPC(ksp, &pc);
235: KSPSetComputeSingularValues(ksp, PETSC_TRUE);
237: *nits=0;
239: /* Isoviscous solve */
240: if (param->ivisc == VISC_CONST && !param->stop_solve) {
241: param->ivisc = VISC_CONST;
242: SNESSolve(snes,0,user->x);
243: VecCopy(user->x,user->Xguess);
244: SNESGetIterationNumber(snes, &its);
245: *nits +=its;
246: if (param->stop_solve) goto done;
247: }
249: /* Olivine diffusion creep */
250: if (param->ivisc >= VISC_DIFN && !param->stop_solve) {
251: if (!q) PetscPrintf(PETSC_COMM_WORLD,"Computing Variable Viscosity Solution\n");
253: /* continuation method on viscosity cutoff */
254: for (param->continuation=0.0; ; param->continuation+=cont_incr) {
255: if (!q) PetscPrintf(PETSC_COMM_WORLD," Continuation parameter = %G\n", param->continuation);
257: /* solve the non-linear system */
258: VecCopy(user->Xguess,user->x);
259: SNESSolve(snes,0,user->x);
260: SNESGetConvergedReason(snes,&reason);
261: SNESGetIterationNumber(snes,&its);
262: *nits += its;
263: if (!q) PetscPrintf(PETSC_COMM_WORLD," SNES iterations: %D, Cumulative: %D\n", its, *nits);
264: if (param->stop_solve) goto done;
266: if (reason<0) {
267: /* NOT converged */
268: cont_incr = -fabs(cont_incr)/2.0;
269: if (fabs(cont_incr)<0.01) goto done;
271: } else {
272: /* converged */
273: VecCopy(user->x,user->Xguess);
274: if (param->continuation >= 1.0) goto done;
275: if (its<=3) {
276: cont_incr = 0.30001;
277: } else if (its<=8) {
278: cont_incr = 0.15001;
279: } else {
280: cont_incr = 0.10001;
281: }
282: if (param->continuation+cont_incr > 1.0) {
283: cont_incr = 1.0 - param->continuation;
284: }
285: } /* endif reason<0 */
286: }
287: }
288: done:
289: if (param->stop_solve && !q) PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: stopping solve.\n");
290: if (reason<0 && !q) PetscPrintf(PETSC_COMM_WORLD,"FAILED TO CONVERGE: stopping solve.\n");
291: return(0);
292: }
295: /*=====================================================================
296: PHYSICS FUNCTIONS (compute the discrete residual)
297: =====================================================================*/
300: /*---------------------------------------------------------------------*/
303: PETSC_STATIC_INLINE PetscScalar UInterp(Field **x, PetscInt i, PetscInt j)
304: /*---------------------------------------------------------------------*/
305: {
306: return 0.25*(x[j][i].u+x[j+1][i].u+x[j][i+1].u+x[j+1][i+1].u);
307: }
309: /*---------------------------------------------------------------------*/
312: PETSC_STATIC_INLINE PetscScalar WInterp(Field **x, PetscInt i, PetscInt j)
313: /*---------------------------------------------------------------------*/
314: {
315: return 0.25*(x[j][i].w+x[j+1][i].w+x[j][i+1].w+x[j+1][i+1].w);
316: }
318: /*---------------------------------------------------------------------*/
321: PETSC_STATIC_INLINE PetscScalar PInterp(Field **x, PetscInt i, PetscInt j)
322: /*---------------------------------------------------------------------*/
323: {
324: return 0.25*(x[j][i].p+x[j+1][i].p+x[j][i+1].p+x[j+1][i+1].p);
325: }
327: /*---------------------------------------------------------------------*/
330: PETSC_STATIC_INLINE PetscScalar TInterp(Field **x, PetscInt i, PetscInt j)
331: /*---------------------------------------------------------------------*/
332: {
333: return 0.25*(x[j][i].T+x[j+1][i].T+x[j][i+1].T+x[j+1][i+1].T);
334: }
336: /*---------------------------------------------------------------------*/
339: /* isoviscous analytic solution for IC */
340: PETSC_STATIC_INLINE PassiveScalar HorizVelocity(PetscInt i, PetscInt j, AppCtx *user)
341: /*---------------------------------------------------------------------*/
342: {
343: Parameter *param = user->param;
344: GridInfo *grid = user->grid;
345: PetscScalar x, z, r, st, ct, th, c=param->c, d=param->d;
346:
347: x = (i - grid->jlid)*grid->dx; z = (j - grid->jlid - 0.5)*grid->dz;
348: r = sqrt(x*x+z*z); st = z/r; ct = x/r; th = atan(z/x);
349: return ct*(c*th*st+d*(st+th*ct)) + st*(c*(st-th*ct)+d*th*st);
350: }
352: /*---------------------------------------------------------------------*/
355: /* isoviscous analytic solution for IC */
356: PETSC_STATIC_INLINE PetscScalar VertVelocity(PetscInt i, PetscInt j, AppCtx *user)
357: /*---------------------------------------------------------------------*/
358: {
359: Parameter *param = user->param;
360: GridInfo *grid = user->grid;
361: PetscScalar x, z, r, st, ct, th, c=param->c, d=param->d;
362:
363: x = (i - grid->jlid - 0.5)*grid->dx; z = (j - grid->jlid)*grid->dz;
364: r = sqrt(x*x+z*z); st = z/r; ct = x/r; th = atan(z/x);
365: return st*(c*th*st+d*(st+th*ct)) - ct*(c*(st-th*ct)+d*th*st);
366: }
368: /*---------------------------------------------------------------------*/
371: /* isoviscous analytic solution for IC */
372: PETSC_STATIC_INLINE PetscScalar Pressure(PetscInt i, PetscInt j, AppCtx *user)
373: /*---------------------------------------------------------------------*/
374: {
375: Parameter *param = user->param;
376: GridInfo *grid = user->grid;
377: PetscScalar x, z, r, st, ct, c=param->c, d=param->d;
379: x = (i - grid->jlid - 0.5)*grid->dx; z = (j - grid->jlid - 0.5)*grid->dz;
380: r = sqrt(x*x+z*z); st = z/r; ct = x/r;
381: return (-2.0*(c*ct-d*st)/r);
382: }
386: /* computes the second invariant of the strain rate tensor */
387: PETSC_STATIC_INLINE PetscScalar CalcSecInv(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
388: /*---------------------------------------------------------------------*/
389: {
390: Parameter *param = user->param;
391: GridInfo *grid = user->grid;
392: PetscInt ilim=grid->ni-1, jlim=grid->nj-1;
393: PetscScalar uN,uS,uE,uW,wN,wS,wE,wW;
394: PetscScalar eps11, eps12, eps22;
396: if (i<j) return EPS_ZERO;
397: if (i==ilim) i--; if (j==jlim) j--;
399: if (ipos==CELL_CENTER) { /* on cell center */
400: if (j<=grid->jlid) return EPS_ZERO;
402: uE = x[j][i].u; uW = x[j][i-1].u;
403: wN = x[j][i].w; wS = x[j-1][i].w;
404: wE = WInterp(x,i,j-1);
405: if (i==j) { uN = param->cb; wW = param->sb; }
406: else { uN = UInterp(x,i-1,j); wW = WInterp(x,i-1,j-1); }
408: if (j==grid->jlid+1) uS = 0.0;
409: else uS = UInterp(x,i-1,j-1);
411: } else { /* on CELL_CORNER */
412: if (j<grid->jlid) return EPS_ZERO;
414: uN = x[j+1][i].u; uS = x[j][i].u;
415: wE = x[j][i+1].w; wW = x[j][i].w;
416: if (i==j) { wN = param->sb; uW = param->cb; }
417: else { wN = WInterp(x,i,j); uW = UInterp(x,i-1,j); }
419: if (j==grid->jlid) {
420: uE = 0.0; uW = 0.0;
421: uS = -uN;
422: wS = -wN;
423: } else {
424: uE = UInterp(x,i,j);
425: wS = WInterp(x,i,j-1);
426: }
427: }
429: eps11 = (uE-uW)/grid->dx; eps22 = (wN-wS)/grid->dz;
430: eps12 = 0.5*((uN-uS)/grid->dz + (wE-wW)/grid->dx);
432: return sqrt( 0.5*( eps11*eps11 + 2.0*eps12*eps12 + eps22*eps22 ) );
433: }
435: /*---------------------------------------------------------------------*/
438: /* computes the shear viscosity */
439: PETSC_STATIC_INLINE PetscScalar Viscosity(PetscScalar T, PetscScalar eps, PassiveScalar z, Parameter *param)
440: /*---------------------------------------------------------------------*/
441: {
442: PetscReal result=0.0;
443: ViscParam difn=param->diffusion, disl=param->dislocation;
444: PetscInt iVisc=param->ivisc;
445: double eps_scale=param->V/(param->L*1000.0);
446: double strain_power, v1, v2, P;
447: double rho_g = 32340.0, R=8.3144;
449: P = rho_g*(z*param->L*1000.0); /* Pa */
451: if (iVisc==VISC_CONST) {
452: /* constant viscosity */
453: return 1.0;
455: } else if (iVisc==VISC_DIFN) {
456: /* diffusion creep rheology */
457: result = difn.A*PetscExpScalar((difn.Estar + P*difn.Vstar)/R/(T+273.0))/param->eta0;
459: } else if (iVisc==VISC_DISL) {
460: /* dislocation creep rheology */
461: strain_power = pow( eps*eps_scale, (1.0-disl.n)/disl.n );
462: result = disl.A*PetscExpScalar((disl.Estar + P*disl.Vstar)/disl.n/R/(T+273.0))*strain_power/param->eta0;
464: } else if (iVisc==VISC_FULL) {
465: /* dislocation/diffusion creep rheology */
466: strain_power = pow( eps*eps_scale, (1.0-disl.n)/disl.n );
467: v1 = difn.A*PetscExpScalar((difn.Estar + P*difn.Vstar)/R/(T+273.0))/param->eta0;
468: v2 = disl.A*PetscExpScalar((disl.Estar + P*disl.Vstar)/disl.n/R/(T+273.0))*strain_power/param->eta0;
469: result = 1.0/(1.0/v1 + 1.0/v2);
470: }
472: /* max viscosity is param->eta0 */
473: result = PetscMin( result, 1.0 );
474: /* min viscosity is param->visc_cutoff */
475: result = PetscMax( result, param->visc_cutoff );
476: /* continuation method */
477: result = pow(result,param->continuation);
478: return result;
479: }
481: /*---------------------------------------------------------------------*/
484: /* computes the residual of the x-component of eqn (1) above */
485: PETSC_STATIC_INLINE PetscScalar XMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
486: /*---------------------------------------------------------------------*/
487: {
488: Parameter *param=user->param;
489: GridInfo *grid =user->grid;
490: PetscScalar dx = grid->dx, dz=grid->dz;
491: PetscScalar etaN,etaS,etaE,etaW,epsN=0.0,epsS=0.0,epsE=0.0,epsW=0.0;
492: PetscScalar TE=0.0,TN=0.0,TS=0.0,TW=0.0, dPdx, residual, z_scale;
493: PetscScalar dudxW,dudxE,dudzN,dudzS,dwdxN,dwdxS;
494: PetscInt jlim = grid->nj-1;
496: z_scale = param->z_scale;
498: if ( param->ivisc==VISC_DIFN || param->ivisc>=VISC_DISL ) { /* viscosity is T-dependent */
499: TS = param->potentialT * TInterp(x,i,j-1) * exp( (j-1.0)*dz*z_scale );
500: if (j==jlim) TN = TS;
501: else TN = param->potentialT * TInterp(x,i,j) * exp( j *dz*z_scale );
502: TW = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*z_scale );
503: TE = param->potentialT * x[j][i+1].T * exp( (j-0.5)*dz*z_scale );
504: if (param->ivisc>=VISC_DISL) { /* olivine dislocation creep */
505: epsN = CalcSecInv(x,i,j, CELL_CORNER,user);
506: epsS = CalcSecInv(x,i,j-1,CELL_CORNER,user);
507: epsE = CalcSecInv(x,i+1,j,CELL_CENTER,user);
508: epsW = CalcSecInv(x,i,j, CELL_CENTER,user);
509: }
510: }
511: etaN = Viscosity(TN,epsN,dz*(j+0.5),param);
512: etaS = Viscosity(TS,epsS,dz*(j-0.5),param);
513: etaW = Viscosity(TW,epsW,dz*j,param);
514: etaE = Viscosity(TE,epsE,dz*j,param);
516: dPdx = ( x[j][i+1].p - x[j][i].p )/dx;
517: if (j==jlim) dudzN = etaN * ( x[j][i].w - x[j][i+1].w )/dx;
518: else dudzN = etaN * ( x[j+1][i].u - x[j][i].u )/dz;
519: dudzS = etaS * ( x[j][i].u - x[j-1][i].u )/dz;
520: dudxE = etaE * ( x[j][i+1].u - x[j][i].u )/dx;
521: dudxW = etaW * ( x[j][i].u - x[j][i-1].u )/dx;
523: residual = -dPdx /* X-MOMENTUM EQUATION*/
524: +( dudxE - dudxW )/dx
525: +( dudzN - dudzS )/dz;
527: if ( param->ivisc!=VISC_CONST ) {
528: dwdxN = etaN * ( x[j][i+1].w - x[j][i].w )/dx;
529: dwdxS = etaS * ( x[j-1][i+1].w - x[j-1][i].w )/dx;
531: residual += ( dudxE - dudxW )/dx + ( dwdxN - dwdxS )/dz;
532: }
534: return residual;
535: }
537: /*---------------------------------------------------------------------*/
540: /* computes the residual of the z-component of eqn (1) above */
541: PETSC_STATIC_INLINE PetscScalar ZMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
542: /*---------------------------------------------------------------------*/
543: {
544: Parameter *param=user->param;
545: GridInfo *grid =user->grid;
546: PetscScalar dx = grid->dx, dz=grid->dz;
547: 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;
548: PetscScalar TE=0.0,TN=0.0,TS=0.0,TW=0.0, dPdz, residual,z_scale;
549: PetscScalar dudzE,dudzW,dwdxW,dwdxE,dwdzN,dwdzS;
550: PetscInt ilim = grid->ni-1;
552: /* geometric and other parameters */
553: z_scale = param->z_scale;
554:
555: /* viscosity */
556: if ( param->ivisc==VISC_DIFN || param->ivisc>=VISC_DISL ) { /* viscosity is T-dependent */
557: TN = param->potentialT * x[j+1][i].T * exp( (j+0.5)*dz*z_scale );
558: TS = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*z_scale );
559: TW = param->potentialT * TInterp(x,i-1,j) * exp( j *dz*z_scale );
560: if (i==ilim) TE = TW;
561: else TE = param->potentialT * TInterp(x,i,j) * exp( j*dz*z_scale );
562: if (param->ivisc>=VISC_DISL) { /* olivine dislocation creep */
563: epsN = CalcSecInv(x,i,j+1,CELL_CENTER,user);
564: epsS = CalcSecInv(x,i,j, CELL_CENTER,user);
565: epsE = CalcSecInv(x,i,j, CELL_CORNER,user);
566: epsW = CalcSecInv(x,i-1,j,CELL_CORNER,user);
567: }
568: }
569: etaN = Viscosity(TN,epsN,dz*(j+1),param);
570: etaS = Viscosity(TS,epsS,dz*j,param);
571: etaW = Viscosity(TW,epsW,dz*(j+0.5),param);
572: etaE = Viscosity(TE,epsE,dz*(j+0.5),param);
574: dPdz = ( x[j+1][i].p - x[j][i].p )/dz;
575: dwdzN = etaN * ( x[j+1][i].w - x[j][i].w )/dz;
576: dwdzS = etaS * ( x[j][i].w - x[j-1][i].w )/dz;
577: if (i==ilim) dwdxE = etaE * ( x[j][i].u - x[j+1][i].u )/dz;
578: else dwdxE = etaE * ( x[j][i+1].w - x[j][i].w )/dx;
579: dwdxW = 2.0*etaW * ( x[j][i].w - x[j][i-1].w )/dx;
580:
581: /* Z-MOMENTUM */
582: residual = -dPdz /* constant viscosity terms */
583: +( dwdzN - dwdzS )/dz
584: +( dwdxE - dwdxW )/dx;
586: if ( param->ivisc!=VISC_CONST ) {
587: dudzE = etaE * ( x[j+1][i].u - x[j][i].u )/dz;
588: dudzW = etaW * ( x[j+1][i-1].u - x[j][i-1].u )/dz;
590: residual += ( dwdzN - dwdzS )/dz + ( dudzE - dudzW )/dx;
591: }
593: return residual;
594: }
596: /*---------------------------------------------------------------------*/
599: /* computes the residual of eqn (2) above */
600: PETSC_STATIC_INLINE PetscScalar ContinuityResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
601: /*---------------------------------------------------------------------*/
602: {
603: GridInfo *grid =user->grid;
604: PetscScalar uE,uW,wN,wS,dudx,dwdz;
606: uW = x[j][i-1].u; uE = x[j][i].u; dudx = ( uE - uW )/grid->dx;
607: wS = x[j-1][i].w; wN = x[j][i].w; dwdz = ( wN - wS )/grid->dz;
609: return dudx + dwdz;
610: }
612: /*---------------------------------------------------------------------*/
615: /* computes the residual of eqn (3) above */
616: PETSC_STATIC_INLINE PetscScalar EnergyResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
617: /*---------------------------------------------------------------------*/
618: {
619: Parameter *param=user->param;
620: GridInfo *grid =user->grid;
621: PetscScalar dx = grid->dx, dz=grid->dz;
622: PetscInt ilim=grid->ni-1, jlim=grid->nj-1, jlid=grid->jlid;
623: PetscScalar TE, TN, TS, TW, residual;
624: PetscScalar uE,uW,wN,wS;
625: PetscScalar fN,fS,fE,fW,dTdxW,dTdxE,dTdzN,dTdzS;
627: dTdzN = ( x[j+1][i].T - x[j][i].T )/dz;
628: dTdzS = ( x[j][i].T - x[j-1][i].T )/dz;
629: dTdxE = ( x[j][i+1].T - x[j][i].T )/dx;
630: dTdxW = ( x[j][i].T - x[j][i-1].T )/dx;
631:
632: residual = ( ( dTdzN - dTdzS )/dz + /* diffusion term */
633: ( dTdxE - dTdxW )/dx )*dx*dz/param->peclet;
635: if (j<=jlid && i>=j) {
636: /* don't advect in the lid */
637: return residual;
639: } else if (i<j) {
640: /* beneath the slab sfc */
641: uW = uE = param->cb;
642: wS = wN = param->sb;
644: } else {
645: /* advect in the slab and wedge */
646: uW = x[j][i-1].u; uE = x[j][i].u;
647: wS = x[j-1][i].w; wN = x[j][i].w;
648: }
650: if ( param->adv_scheme==ADVECT_FV || i==ilim-1 || j==jlim-1 || i==1 || j==1 ) {
651: /* finite volume advection */
652: TS = ( x[j][i].T + x[j-1][i].T )/2.0;
653: TN = ( x[j][i].T + x[j+1][i].T )/2.0;
654: TE = ( x[j][i].T + x[j][i+1].T )/2.0;
655: TW = ( x[j][i].T + x[j][i-1].T )/2.0;
656: fN = wN*TN*dx; fS = wS*TS*dx;
657: fE = uE*TE*dz; fW = uW*TW*dz;
658:
659: } else {
660: /* Fromm advection scheme */
661: 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
662: - 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;
663: 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
664: - 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;
665: 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
666: - 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;
667: 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
668: - 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;
669: }
670:
671: residual -= ( fE - fW + fN - fS );
673: return residual;
674: }
676: /*---------------------------------------------------------------------*/
679: /* computes the shear stress---used on the boundaries */
680: PETSC_STATIC_INLINE PetscScalar ShearStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
681: /*---------------------------------------------------------------------*/
682: {
683: Parameter *param=user->param;
684: GridInfo *grid=user->grid;
685: PetscInt ilim=grid->ni-1, jlim=grid->nj-1;
686: PetscScalar uN, uS, wE, wW;
688: if ( j<=grid->jlid || i<j || i==ilim || j==jlim ) return EPS_ZERO;
690: if (ipos==CELL_CENTER) { /* on cell center */
692: wE = WInterp(x,i,j-1);
693: if (i==j) { wW = param->sb; uN = param->cb;}
694: else { wW = WInterp(x,i-1,j-1); uN = UInterp(x,i-1,j); }
695: if (j==grid->jlid+1) uS = 0.0;
696: else uS = UInterp(x,i-1,j-1);
698: } else { /* on cell corner */
700: uN = x[j+1][i].u; uS = x[j][i].u;
701: wW = x[j][i].w; wE = x[j][i+1].w;
703: }
705: return (uN-uS)/grid->dz + (wE-wW)/grid->dx;
706: }
708: /*---------------------------------------------------------------------*/
711: /* computes the normal stress---used on the boundaries */
712: PETSC_STATIC_INLINE PetscScalar XNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
713: /*---------------------------------------------------------------------*/
714: {
715: Parameter *param=user->param;
716: GridInfo *grid =user->grid;
717: PetscScalar dx = grid->dx, dz=grid->dz;
718: PetscInt ilim=grid->ni-1, jlim=grid->nj-1, ivisc;
719: PetscScalar epsC=0.0, etaC, TC, uE, uW, pC, z_scale;
720: if (i<j || j<=grid->jlid) return EPS_ZERO;
722: ivisc=param->ivisc; z_scale = param->z_scale;
724: if (ipos==CELL_CENTER) { /* on cell center */
726: TC = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*z_scale );
727: if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CENTER,user);
728: etaC = Viscosity(TC,epsC,dz*j,param);
730: uW = x[j][i-1].u; uE = x[j][i].u;
731: pC = x[j][i].p;
733: } else { /* on cell corner */
734: if ( i==ilim || j==jlim ) return EPS_ZERO;
736: TC = param->potentialT * TInterp(x,i,j) * exp( j*dz*z_scale );
737: if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CORNER,user);
738: etaC = Viscosity(TC,epsC,dz*(j+0.5),param);
740: if (i==j) uW = param->sb;
741: else uW = UInterp(x,i-1,j);
742: uE = UInterp(x,i,j); pC = PInterp(x,i,j);
743: }
744:
745: return 2.0*etaC*(uE-uW)/dx - pC;
746: }
748: /*---------------------------------------------------------------------*/
751: /* computes the normal stress---used on the boundaries */
752: PETSC_STATIC_INLINE PetscScalar ZNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
753: /*---------------------------------------------------------------------*/
754: {
755: Parameter *param=user->param;
756: GridInfo *grid =user->grid;
757: PetscScalar dz=grid->dz;
758: PetscInt ilim=grid->ni-1, jlim=grid->nj-1, ivisc;
759: PetscScalar epsC=0.0, etaC, TC;
760: PetscScalar pC, wN, wS, z_scale;
761: if (i<j || j<=grid->jlid) return EPS_ZERO;
763: ivisc=param->ivisc; z_scale = param->z_scale;
765: if (ipos==CELL_CENTER) { /* on cell center */
767: TC = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*z_scale );
768: if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CENTER,user);
769: etaC = Viscosity(TC,epsC,dz*j,param);
770: wN = x[j][i].w; wS = x[j-1][i].w; pC = x[j][i].p;
772: } else { /* on cell corner */
773: if ( (i==ilim) || (j==jlim) ) return EPS_ZERO;
775: TC = param->potentialT * TInterp(x,i,j) * exp( j*dz*z_scale );
776: if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CORNER,user);
777: etaC = Viscosity(TC,epsC,dz*(j+0.5),param);
778: if (i==j) wN = param->sb;
779: else wN = WInterp(x,i,j);
780: wS = WInterp(x,i,j-1); pC = PInterp(x,i,j);
781: }
783: return 2.0*etaC*(wN-wS)/dz - pC;
784: }
786: /*---------------------------------------------------------------------*/
788: /*=====================================================================
789: INITIALIZATION, POST-PROCESSING AND OUTPUT FUNCTIONS
790: =====================================================================*/
792: /*---------------------------------------------------------------------*/
795: /* initializes the problem parameters and checks for
796: command line changes */
797: PetscErrorCode SetParams(Parameter *param, GridInfo *grid)
798: /*---------------------------------------------------------------------*/
799: {
800: PetscErrorCode ierr, ierr_out=0;
801: PetscReal SEC_PER_YR = 3600.00*24.00*365.2500;
802: PetscReal PI = 3.14159265358979323846;
803: PetscReal alpha_g_on_cp_units_inverse_km=4.0e-5*9.8;
804:
805: /* domain geometry */
806: param->slab_dip = 45.0;
807: param->width = 320.0; /* km */
808: param->depth = 300.0; /* km */
809: param->lid_depth = 35.0; /* km */
810: param->fault_depth = 35.0; /* km */
811: PetscOptionsGetReal(PETSC_NULL,"-slab_dip",&(param->slab_dip),PETSC_NULL);
812: PetscOptionsGetReal(PETSC_NULL,"-width",&(param->width),PETSC_NULL);
813: PetscOptionsGetReal(PETSC_NULL,"-depth",&(param->depth),PETSC_NULL);
814: PetscOptionsGetReal(PETSC_NULL,"-lid_depth",&(param->lid_depth),PETSC_NULL);
815: PetscOptionsGetReal(PETSC_NULL,"-fault_depth",&(param->fault_depth),PETSC_NULL);
816: param->slab_dip = param->slab_dip*PI/180.0; /* radians */
818: /* grid information */
819: PetscOptionsGetInt(PETSC_NULL, "-jfault",&(grid->jfault),PETSC_NULL);
820: grid->ni = 82;
821: PetscOptionsGetInt(PETSC_NULL, "-ni",&(grid->ni),PETSC_NULL);
822: grid->dx = param->width/((double)(grid->ni-2)); /* km */
823: grid->dz = grid->dx*tan(param->slab_dip); /* km */
824: grid->nj = (PetscInt)(param->depth/grid->dz + 3.0); /* gridpoints*/
825: param->depth = grid->dz*(grid->nj-2); /* km */
826: grid->inose = 0; /* gridpoints*/
827: PetscOptionsGetInt(PETSC_NULL,"-inose",&(grid->inose),PETSC_NULL);
828: grid->bx = DMDA_BOUNDARY_NONE;
829: grid->by = DMDA_BOUNDARY_NONE;
830: grid->stencil = DMDA_STENCIL_BOX;
831: grid->dof = 4;
832: grid->stencil_width = 2;
833: grid->mglevels = 1;
835: /* boundary conditions */
836: param->pv_analytic = PETSC_FALSE;
837: param->ibound = BC_NOSTRESS;
838: PetscOptionsGetInt(PETSC_NULL,"-ibound",&(param->ibound),PETSC_NULL);
840: /* physical constants */
841: param->slab_velocity = 5.0; /* cm/yr */
842: param->slab_age = 50.0; /* Ma */
843: param->lid_age = 50.0; /* Ma */
844: param->kappa = 0.7272e-6; /* m^2/sec */
845: param->potentialT = 1300.0; /* degrees C */
846: PetscOptionsGetReal(PETSC_NULL,"-slab_velocity",&(param->slab_velocity),PETSC_NULL);
847: PetscOptionsGetReal(PETSC_NULL,"-slab_age",&(param->slab_age),PETSC_NULL);
848: PetscOptionsGetReal(PETSC_NULL,"-lid_age",&(param->lid_age),PETSC_NULL);
849: PetscOptionsGetReal(PETSC_NULL,"-kappa",&(param->kappa),PETSC_NULL);
850: PetscOptionsGetReal(PETSC_NULL,"-potentialT",&(param->potentialT),PETSC_NULL);
852: /* viscosity */
853: param->ivisc = 3; /* 0=isovisc, 1=difn creep, 2=disl creep, 3=full */
854: param->eta0 = 1e24; /* Pa-s */
855: param->visc_cutoff = 0.0; /* factor of eta_0 */
856: param->continuation = 1.0;
857: /* constants for diffusion creep */
858: param->diffusion.A = 1.8e7; /* Pa-s */
859: param->diffusion.n = 1.0; /* dim'less */
860: param->diffusion.Estar = 375e3; /* J/mol */
861: param->diffusion.Vstar = 5e-6; /* m^3/mol */
862: /* constants for param->dislocationocation creep */
863: param->dislocation.A = 2.8969e4; /* Pa-s */
864: param->dislocation.n = 3.5; /* dim'less */
865: param->dislocation.Estar = 530e3; /* J/mol */
866: param->dislocation.Vstar = 14e-6; /* m^3/mol */
867: PetscOptionsGetInt(PETSC_NULL, "-ivisc",&(param->ivisc),PETSC_NULL);
868: PetscOptionsGetReal(PETSC_NULL,"-visc_cutoff",&(param->visc_cutoff),PETSC_NULL);
869: param->output_ivisc = param->ivisc;
870: PetscOptionsGetInt(PETSC_NULL,"-output_ivisc",&(param->output_ivisc),PETSC_NULL);
871: PetscOptionsGetReal(PETSC_NULL,"-vstar",&(param->dislocation.Vstar),PETSC_NULL);
873: /* output options */
874: param->quiet = PETSC_FALSE;
875: param->param_test = PETSC_FALSE;
876: PetscOptionsHasName(PETSC_NULL,"-quiet",&(param->quiet));
877: PetscOptionsHasName(PETSC_NULL,"-test",&(param->param_test));
878: PetscOptionsGetString(PETSC_NULL,"-file",param->filename,PETSC_MAX_PATH_LEN,&(param->output_to_file));
880: /* advection */
881: param->adv_scheme = ADVECT_FROMM; /* advection scheme: 0=finite vol, 1=Fromm */
882: PetscOptionsGetInt(PETSC_NULL,"-adv_scheme",&(param->adv_scheme),PETSC_NULL);
884: /* misc. flags */
885: param->stop_solve = PETSC_FALSE;
886: param->interrupted = PETSC_FALSE;
887: param->kspmon = PETSC_FALSE;
888: param->toggle_kspmon = PETSC_FALSE;
890: /* derived parameters for slab angle */
891: param->sb = sin(param->slab_dip);
892: param->cb = cos(param->slab_dip);
893: param->c = param->slab_dip*param->sb/(param->slab_dip*param->slab_dip-param->sb*param->sb);
894: param->d = (param->slab_dip*param->cb-param->sb)/(param->slab_dip*param->slab_dip-param->sb*param->sb);
896: /* length, velocity and time scale for non-dimensionalization */
897: param->L = PetscMin(param->width,param->depth); /* km */
898: param->V = param->slab_velocity/100.0/SEC_PER_YR; /* m/sec */
900: /* other unit conversions and derived parameters */
901: param->scaled_width = param->width/param->L; /* dim'less */
902: param->scaled_depth = param->depth/param->L; /* dim'less */
903: param->lid_depth = param->lid_depth/param->L; /* dim'less */
904: param->fault_depth = param->fault_depth/param->L; /* dim'less */
905: grid->dx = grid->dx/param->L; /* dim'less */
906: grid->dz = grid->dz/param->L; /* dim'less */
907: grid->jlid = (PetscInt)(param->lid_depth/grid->dz); /* gridcells */
908: grid->jfault = (PetscInt)(param->fault_depth/grid->dz); /* gridcells */
909: param->lid_depth = grid->jlid*grid->dz; /* dim'less */
910: param->fault_depth = grid->jfault*grid->dz; /* dim'less */
911: grid->corner = grid->jlid+1; /* gridcells */
912: param->peclet = param->V /* m/sec */
913: * param->L*1000.0 /* m */
914: / param->kappa; /* m^2/sec */
915: param->z_scale = param->L * alpha_g_on_cp_units_inverse_km;
916: param->skt = sqrt(param->kappa*param->slab_age*SEC_PER_YR);
917: PetscOptionsGetReal(PETSC_NULL,"-peclet",&(param->peclet),PETSC_NULL);
918:
919: return ierr_out;
920: }
922: /*---------------------------------------------------------------------*/
925: /* prints a report of the problem parameters to stdout */
926: PetscErrorCode ReportParams(Parameter *param, GridInfo *grid)
927: /*---------------------------------------------------------------------*/
928: {
929: PetscErrorCode ierr, ierr_out=0;
930: char date[30];
931: PetscReal PI = 3.14159265358979323846;
933: PetscGetDate(date,30);
935: if ( !(param->quiet) ) {
936: PetscPrintf(PETSC_COMM_WORLD,"---------------------BEGIN ex30 PARAM REPORT-------------------\n");
937: /* PetscPrintf(PETSC_COMM_WORLD," %s\n",&(date[0]));*/
939: PetscPrintf(PETSC_COMM_WORLD,"Domain: \n");
940: PetscPrintf(PETSC_COMM_WORLD," Width = %G km, Depth = %G km\n",param->width,param->depth);
941: PetscPrintf(PETSC_COMM_WORLD," Slab dip = %G degrees, Slab velocity = %G cm/yr\n",param->slab_dip*180.0/PI,param->slab_velocity);
942: 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);
944: PetscPrintf(PETSC_COMM_WORLD,"\nGrid: \n");
945: 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);
946: PetscPrintf(PETSC_COMM_WORLD," jlid = %3D jfault = %3D \n",grid->jlid,grid->jfault);
947: PetscPrintf(PETSC_COMM_WORLD," Pe = %G\n",param->peclet);
949: PetscPrintf(PETSC_COMM_WORLD,"\nRheology:");
950: if (param->ivisc==VISC_CONST) {
951: PetscPrintf(PETSC_COMM_WORLD," Isoviscous \n");
952: if (param->pv_analytic)
953: PetscPrintf(PETSC_COMM_WORLD," Pressure and Velocity prescribed! \n");
954: } else if (param->ivisc==VISC_DIFN) {
955: PetscPrintf(PETSC_COMM_WORLD," Diffusion Creep (T-Dependent Newtonian) \n");
956: PetscPrintf(PETSC_COMM_WORLD," Viscosity range: %G--%G Pa-sec \n",param->eta0,param->visc_cutoff*param->eta0);
957: } else if (param->ivisc==VISC_DISL ) {
958: PetscPrintf(PETSC_COMM_WORLD," Dislocation Creep (T-Dependent Non-Newtonian) \n");
959: PetscPrintf(PETSC_COMM_WORLD," Viscosity range: %G--%G Pa-sec \n",param->eta0,param->visc_cutoff*param->eta0);
960: } else if (param->ivisc==VISC_FULL ) {
961: PetscPrintf(PETSC_COMM_WORLD," Full Rheology \n");
962: PetscPrintf(PETSC_COMM_WORLD," Viscosity range: %G--%G Pa-sec \n",param->eta0,param->visc_cutoff*param->eta0);
963: } else {
964: PetscPrintf(PETSC_COMM_WORLD," Invalid! \n");
965: ierr_out=1;
966: }
968: PetscPrintf(PETSC_COMM_WORLD,"Boundary condition:");
969: if ( param->ibound==BC_ANALYTIC ) {
970: PetscPrintf(PETSC_COMM_WORLD," Isoviscous Analytic Dirichlet \n");
971: } else if ( param->ibound==BC_NOSTRESS ) {
972: PetscPrintf(PETSC_COMM_WORLD," Stress-Free (normal & shear stress)\n");
973: } else if ( param->ibound==BC_EXPERMNT ) {
974: PetscPrintf(PETSC_COMM_WORLD," Experimental boundary condition \n");
975: } else {
976: PetscPrintf(PETSC_COMM_WORLD," Invalid! \n");
977: ierr_out=1;
978: }
980: if (param->output_to_file)
981: #if defined(PETSC_HAVE_MATLAB_ENGINE)
982: PetscPrintf(PETSC_COMM_WORLD,"Output Destination: Mat file \"%s\"\n",param->filename);
983: #else
984: PetscPrintf(PETSC_COMM_WORLD,"Output Destination: PETSc binary file \"%s\"\n",param->filename);
985: #endif
986: if ( param->output_ivisc != param->ivisc )
987: PetscPrintf(PETSC_COMM_WORLD," Output viscosity: -ivisc %D\n",param->output_ivisc);
989: PetscPrintf(PETSC_COMM_WORLD,"---------------------END ex30 PARAM REPORT---------------------\n");
990: }
991: if ( param->param_test ) PetscEnd();
992: return ierr_out;
993: }
995: /* ------------------------------------------------------------------- */
998: /* generates an inital guess using the analytic solution for isoviscous
999: corner flow */
1000: PetscErrorCode Initialize(DM da)
1001: /* ------------------------------------------------------------------- */
1002: {
1003: AppCtx *user;
1004: Parameter *param;
1005: GridInfo *grid;
1006: PetscInt i,j,is,js,im,jm;
1008: Field **x;
1009: Vec Xguess;
1011: /* Get the fine grid */
1012: DMGetApplicationContext(da,&user);
1013: Xguess = user->Xguess;
1014: param = user->param;
1015: grid = user->grid;
1016: DMDAGetCorners(da,&is,&js,PETSC_NULL,&im,&jm,PETSC_NULL);
1017: DMDAVecGetArray(da,Xguess,(void**)&x);
1019: /* Compute initial guess */
1020: for (j=js; j<js+jm; j++) {
1021: for (i=is; i<is+im; i++) {
1022: if (i<j) {
1023: x[j][i].u = param->cb;
1024: } else if (j<=grid->jlid) {
1025: x[j][i].u = 0.0;
1026: } else {
1027: x[j][i].u = HorizVelocity(i,j,user);
1028: }
1029: if (i<=j) {
1030: x[j][i].w = param->sb;
1031: } else if (j<=grid->jlid) {
1032: x[j][i].w = 0.0;
1033: } else {
1034: x[j][i].w = VertVelocity(i,j,user);
1035: }
1036: if (i<j || j<=grid->jlid) {
1037: x[j][i].p = 0.0;
1038: } else {
1039: x[j][i].p = Pressure(i,j,user);
1040: }
1041: x[j][i].T = PetscMin(grid->dz*(j-0.5),1.0);
1042: }
1043: }
1045: /* Restore x to Xguess */
1046: DMDAVecRestoreArray(da,Xguess,(void**)&x);
1048: return 0;
1049: }
1051: /*---------------------------------------------------------------------*/
1054: /* controls output to a file */
1055: PetscErrorCode DoOutput(SNES snes, PetscInt its)
1056: /*---------------------------------------------------------------------*/
1057: {
1058: AppCtx *user;
1059: Parameter *param;
1060: GridInfo *grid;
1061: PetscInt ivt;
1063: PetscMPIInt rank;
1064: PetscViewer viewer;
1065: Vec res, pars;
1066: MPI_Comm comm;
1067: DM da;
1069: SNESGetDM(snes,&da);
1070: DMGetApplicationContext(da,&user);
1071: param = user->param;
1072: grid = user->grid;
1073: ivt = param->ivisc;
1075: param->ivisc = param->output_ivisc;
1077: /* compute final residual and final viscosity/strain rate fields */
1078: SNESGetFunction(snes, &res, PETSC_NULL, PETSC_NULL);
1079: ViscosityField(da, user->x, user->Xguess);
1081: /* get the communicator and the rank of the processor */
1082: PetscObjectGetComm((PetscObject)snes, &comm);
1083: MPI_Comm_rank(comm, &rank);
1085: if (param->output_to_file) { /* send output to binary file */
1086: VecCreate(comm, &pars);
1087: if (rank == 0) { /* on processor 0 */
1088: VecSetSizes(pars, 20, PETSC_DETERMINE);
1089: VecSetFromOptions(pars);
1090: VecSetValue(pars,0, (PetscScalar)(grid->ni),INSERT_VALUES);
1091: VecSetValue(pars,1, (PetscScalar)(grid->nj),INSERT_VALUES);
1092: VecSetValue(pars,2, (PetscScalar)(grid->dx),INSERT_VALUES);
1093: VecSetValue(pars,3, (PetscScalar)(grid->dz),INSERT_VALUES);
1094: VecSetValue(pars,4, (PetscScalar)(param->L),INSERT_VALUES);
1095: VecSetValue(pars,5, (PetscScalar)(param->V),INSERT_VALUES);
1096: /* skipped 6 intentionally */
1097: VecSetValue(pars,7, (PetscScalar)(param->slab_dip),INSERT_VALUES);
1098: VecSetValue(pars,8, (PetscScalar)(grid->jlid),INSERT_VALUES);
1099: VecSetValue(pars,9, (PetscScalar)(param->lid_depth),INSERT_VALUES);
1100: VecSetValue(pars,10,(PetscScalar)(grid->jfault),INSERT_VALUES);
1101: VecSetValue(pars,11,(PetscScalar)(param->fault_depth),INSERT_VALUES);
1102: VecSetValue(pars,12,(PetscScalar)(param->potentialT),INSERT_VALUES);
1103: VecSetValue(pars,13,(PetscScalar)(param->ivisc),INSERT_VALUES);
1104: VecSetValue(pars,14,(PetscScalar)(param->visc_cutoff),INSERT_VALUES);
1105: VecSetValue(pars,15,(PetscScalar)(param->ibound),INSERT_VALUES);
1106: VecSetValue(pars,16,(PetscScalar)(its),INSERT_VALUES);
1107: } else { /* on some other processor */
1108: VecSetSizes(pars, 0, PETSC_DETERMINE);
1109: VecSetFromOptions(pars);
1110: }
1111: VecAssemblyBegin(pars); VecAssemblyEnd(pars);
1113: /* create viewer */
1114: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1115: PetscViewerMatlabOpen(PETSC_COMM_WORLD,param->filename,FILE_MODE_WRITE,&viewer);
1116: #else
1117: PetscViewerBinaryOpen(PETSC_COMM_WORLD,param->filename,FILE_MODE_WRITE,&viewer);
1118: #endif
1120: /* send vectors to viewer */
1121: PetscObjectSetName((PetscObject)res,"res");
1122: VecView(res,viewer);
1123: PetscObjectSetName((PetscObject)user->x,"out");
1124: VecView(user->x, viewer);
1125: PetscObjectSetName((PetscObject)(user->Xguess),"aux");
1126: VecView(user->Xguess, viewer);
1127: StressField(da); /* compute stress fields */
1128: PetscObjectSetName((PetscObject)(user->Xguess),"str");
1129: VecView(user->Xguess, viewer);
1130: PetscObjectSetName((PetscObject)pars,"par");
1131: VecView(pars, viewer);
1132:
1133: /* destroy viewer and vector */
1134: PetscViewerDestroy(&viewer);
1135: VecDestroy(&pars);
1136: }
1137:
1138: param->ivisc = ivt;
1139: return 0;
1140: }
1142: /* ------------------------------------------------------------------- */
1145: /* Compute both the second invariant of the strain rate tensor and the viscosity, at both cell centers and cell corners */
1146: PetscErrorCode ViscosityField(DM da, Vec X, Vec V)
1147: /* ------------------------------------------------------------------- */
1148: {
1149: AppCtx *user;
1150: Parameter *param;
1151: GridInfo *grid;
1152: Vec localX;
1153: Field **v, **x;
1154: PassiveReal eps, /* dx,*/ dz, T, epsC, TC;
1155: PetscInt i,j,is,js,im,jm,ilim,jlim,ivt;
1159: DMGetApplicationContext(da,&user);
1160: param = user->param;
1161: grid = user->grid;
1162: ivt = param->ivisc;
1163: param->ivisc = param->output_ivisc;
1165: DMGetLocalVector(da, &localX);
1166: DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX);
1167: DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX);
1168: DMDAVecGetArray(da,localX,(void**)&x);
1169: DMDAVecGetArray(da,V,(void**)&v);
1171: /* Parameters */
1172: /* dx = grid->dx; */ dz = grid->dz;
1173: ilim = grid->ni-1; jlim = grid->nj-1;
1175: /* Compute real temperature, strain rate and viscosity */
1176: DMDAGetCorners(da,&is,&js,PETSC_NULL,&im,&jm,PETSC_NULL);
1177: for (j=js; j<js+jm; j++) {
1178: for (i=is; i<is+im; i++) {
1179: T = param->potentialT * x[j][i].T * exp( (j-0.5)*dz*param->z_scale );
1180: if (i<ilim && j<jlim) {
1181: TC = param->potentialT * TInterp(x,i,j) * exp( j*dz*param->z_scale );
1182: } else {
1183: TC = T;
1184: }
1185: eps = CalcSecInv(x,i,j,CELL_CENTER,user);
1186: epsC = CalcSecInv(x,i,j,CELL_CORNER,user);
1187: v[j][i].u = eps;
1188: v[j][i].w = epsC;
1189: v[j][i].p = Viscosity(T,eps,dz*(j-0.5),param);
1190: v[j][i].T = Viscosity(TC,epsC,dz*j,param);
1191: }
1192: }
1193: DMDAVecRestoreArray(da,V,(void**)&v);
1194: DMDAVecRestoreArray(da,localX,(void**)&x);
1195: DMRestoreLocalVector(da, &localX);
1196: param->ivisc = ivt;
1197: return(0);
1198: }
1200: /* ------------------------------------------------------------------- */
1203: /* post-processing: compute stress everywhere */
1204: PetscErrorCode StressField(DM da)
1205: /* ------------------------------------------------------------------- */
1206: {
1207: AppCtx *user;
1208: PetscInt i,j,is,js,im,jm;
1210: Vec locVec;
1211: Field **x, **y;
1213: DMGetApplicationContext(da,&user);
1215: /* Get the fine grid of Xguess and X */
1216: DMDAGetCorners(da,&is,&js,PETSC_NULL,&im,&jm,PETSC_NULL);
1217: DMDAVecGetArray(da,user->Xguess,(void**)&x);
1219: DMGetLocalVector(da, &locVec);
1220: DMGlobalToLocalBegin(da, user->x, INSERT_VALUES, locVec);
1221: DMGlobalToLocalEnd(da, user->x, INSERT_VALUES, locVec);
1222: DMDAVecGetArray(da,locVec,(void**)&y);
1224: /* Compute stress on the corner points */
1225: for (j=js; j<js+jm; j++) {
1226: for (i=is; i<is+im; i++) {
1227:
1228: x[j][i].u = ShearStress(y,i,j,CELL_CENTER,user);
1229: x[j][i].w = ShearStress(y,i,j,CELL_CORNER,user);
1230: x[j][i].p = XNormalStress(y,i,j,CELL_CENTER,user);
1231: x[j][i].T = ZNormalStress(y,i,j,CELL_CENTER,user);
1232: }
1233: }
1235: /* Restore the fine grid of Xguess and X */
1236: DMDAVecRestoreArray(da,user->Xguess,(void**)&x);
1237: DMDAVecRestoreArray(da,locVec,(void**)&y);
1238: DMRestoreLocalVector(da, &locVec);
1239: return 0;
1240: }
1242: /*=====================================================================
1243: UTILITY FUNCTIONS
1244: =====================================================================*/
1246: /*---------------------------------------------------------------------*/
1249: /* returns the velocity of the subducting slab and handles fault nodes
1250: for BC */
1251: PETSC_STATIC_INLINE PassiveScalar SlabVel(char c, PetscInt i, PetscInt j, AppCtx *user)
1252: /*---------------------------------------------------------------------*/
1253: {
1254: Parameter *param = user->param;
1255: GridInfo *grid = user->grid;
1257: if (c=='U' || c=='u') {
1258: if (i<j-1) {
1259: return param->cb;
1260: } else if (j<=grid->jfault) {
1261: return 0.0;
1262: } else
1263: return param->cb;
1265: } else {
1266: if (i<j) {
1267: return param->sb;
1268: } else if (j<=grid->jfault) {
1269: return 0.0;
1270: } else
1271: return param->sb;
1272: }
1273: }
1275: /*---------------------------------------------------------------------*/
1278: /* solution to diffusive half-space cooling model for BC */
1279: PETSC_STATIC_INLINE PassiveScalar PlateModel(PetscInt j, PetscInt plate, AppCtx *user)
1280: /*---------------------------------------------------------------------*/
1281: {
1282: Parameter *param = user->param;
1283: PassiveScalar z;
1284: if (plate==PLATE_LID)
1285: z = (j-0.5)*user->grid->dz;
1286: else /* PLATE_SLAB */
1287: z = (j-0.5)*user->grid->dz*param->cb;
1288: #if defined (PETSC_HAVE_ERF)
1289: return erf(z*param->L/2.0/param->skt);
1290: #else
1291: SETERRQ(PETSC_COMM_SELF,1,"erf() not available on this machine");
1292: #endif
1293: }
1296: /* ------------------------------------------------------------------- */
1299: /* utility function */
1300: PetscBool OptionsHasName(const char pre[],const char name[])
1301: /* ------------------------------------------------------------------- */
1302: {
1303: PetscBool retval;
1305: PetscOptionsHasName(pre,name,&retval);CHKERRABORT(PETSC_COMM_WORLD,ierr);
1306: return retval;
1307: }
1309: /*=====================================================================
1310: INTERACTIVE SIGNAL HANDLING
1311: =====================================================================*/
1313: /* ------------------------------------------------------------------- */
1316: PetscErrorCode SNESConverged_Interactive(SNES snes, PetscInt it,PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *ctx)
1317: /* ------------------------------------------------------------------- */
1318: {
1319: AppCtx *user = (AppCtx *) ctx;
1320: Parameter *param = user->param;
1321: KSP ksp;
1325: if (param->interrupted) {
1326: param->interrupted = PETSC_FALSE;
1327: PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: exiting SNES solve. \n");
1328: *reason = SNES_CONVERGED_FNORM_ABS;
1329: return(0);
1330: } else if (param->toggle_kspmon) {
1331: param->toggle_kspmon = PETSC_FALSE;
1332: SNESGetKSP(snes, &ksp);
1333: if (param->kspmon) {
1334: KSPMonitorCancel(ksp);
1335: param->kspmon = PETSC_FALSE;
1336: PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: deactivating ksp singular value monitor. \n");
1337: } else {
1338: KSPMonitorSet(ksp,KSPMonitorSingularValue,PETSC_NULL,PETSC_NULL);
1339: param->kspmon = PETSC_TRUE;
1340: PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: activating ksp singular value monitor. \n");
1341: }
1342: }
1343: PetscFunctionReturn(SNESDefaultConverged(snes,it,xnorm,snorm,fnorm,reason,ctx));
1344: }
1346: /* ------------------------------------------------------------------- */
1347: #include <signal.h>
1350: PetscErrorCode InteractiveHandler(int signum, void *ctx)
1351: /* ------------------------------------------------------------------- */
1352: {
1353: AppCtx *user = (AppCtx *) ctx;
1354: Parameter *param = user->param;
1356: if (signum == SIGILL) {
1357: param->toggle_kspmon = PETSC_TRUE;
1358: #if !defined(PETSC_HAVE_MISSING_SIGCONT)
1359: } else if (signum == SIGCONT) {
1360: param->interrupted = PETSC_TRUE;
1361: #endif
1362: #if !defined(PETSC_HAVE_MISSING_SIGURG)
1363: } else if (signum == SIGURG) {
1364: param->stop_solve = PETSC_TRUE;
1365: #endif
1366: }
1367: return 0;
1368: }
1370: /*---------------------------------------------------------------------*/
1373: /* main call-back function that computes the processor-local piece
1374: of the residual */
1375: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
1376: /*---------------------------------------------------------------------*/
1377: {
1378: AppCtx *user = (AppCtx*)ptr;
1379: Parameter *param = user->param;
1380: GridInfo *grid = user->grid;
1381: PetscScalar mag_w, mag_u;
1382: PetscInt i,j,mx,mz,ilim,jlim;
1383: PetscInt is,ie,js,je,ibound; /* ,ivisc */
1387: /* Define global and local grid parameters */
1388: mx = info->mx; mz = info->my;
1389: ilim = mx-1; jlim = mz-1;
1390: is = info->xs; ie = info->xs+info->xm;
1391: js = info->ys; je = info->ys+info->ym;
1393: /* Define geometric and numeric parameters */
1394: /* ivisc = param->ivisc; */ ibound = param->ibound;
1396: for (j=js; j<je; j++) {
1397: for (i=is; i<ie; i++) {
1399: /************* X-MOMENTUM/VELOCITY *************/
1400: if (i<j) {
1401: f[j][i].u = x[j][i].u - SlabVel('U',i,j,user);
1403: } else if (j<=grid->jlid || (j<grid->corner+grid->inose && i<grid->corner+grid->inose)) {
1404: /* in the lithospheric lid */
1405: f[j][i].u = x[j][i].u - 0.0;
1407: } else if (i==ilim) {
1408: /* on the right side boundary */
1409: if (ibound==BC_ANALYTIC) {
1410: f[j][i].u = x[j][i].u - HorizVelocity(i,j,user);
1411: } else {
1412: f[j][i].u = XNormalStress(x,i,j,CELL_CENTER,user) - EPS_ZERO;
1413: }
1415: } else if (j==jlim) {
1416: /* on the bottom boundary */
1417: if (ibound==BC_ANALYTIC) {
1418: f[j][i].u = x[j][i].u - HorizVelocity(i,j,user);
1419: } else if (ibound==BC_NOSTRESS) {
1420: f[j][i].u = XMomentumResidual(x,i,j,user);
1421: } else {
1422: /* experimental boundary condition */
1423: }
1425: } else {
1426: /* in the mantle wedge */
1427: f[j][i].u = XMomentumResidual(x,i,j,user);
1428: }
1429:
1430: /************* Z-MOMENTUM/VELOCITY *************/
1431: if (i<=j) {
1432: f[j][i].w = x[j][i].w - SlabVel('W',i,j,user);
1434: } else if (j<=grid->jlid || (j<grid->corner+grid->inose && i<grid->corner+grid->inose)) {
1435: /* in the lithospheric lid */
1436: f[j][i].w = x[j][i].w - 0.0;
1438: } else if (j==jlim) {
1439: /* on the bottom boundary */
1440: if (ibound==BC_ANALYTIC) {
1441: f[j][i].w = x[j][i].w - VertVelocity(i,j,user);
1442: } else {
1443: f[j][i].w = ZNormalStress(x,i,j,CELL_CENTER,user) - EPS_ZERO;
1444: }
1446: } else if (i==ilim) {
1447: /* on the right side boundary */
1448: if (ibound==BC_ANALYTIC) {
1449: f[j][i].w = x[j][i].w - VertVelocity(i,j,user);
1450: } else if (ibound==BC_NOSTRESS) {
1451: f[j][i].w = ZMomentumResidual(x,i,j,user);
1452: } else {
1453: /* experimental boundary condition */
1454: }
1456: } else {
1457: /* in the mantle wedge */
1458: f[j][i].w = ZMomentumResidual(x,i,j,user);
1459: }
1461: /************* CONTINUITY/PRESSURE *************/
1462: if (i<j || j<=grid->jlid || (j<grid->corner+grid->inose && i<grid->corner+grid->inose)) {
1463: /* in the lid or slab */
1464: f[j][i].p = x[j][i].p;
1465:
1466: } else if ((i==ilim || j==jlim) && ibound==BC_ANALYTIC) {
1467: /* on an analytic boundary */
1468: f[j][i].p = x[j][i].p - Pressure(i,j,user);
1470: } else {
1471: /* in the mantle wedge */
1472: f[j][i].p = ContinuityResidual(x,i,j,user);
1473: }
1475: /************* TEMPERATURE *************/
1476: if (j==0) {
1477: /* on the surface */
1478: f[j][i].T = x[j][i].T + x[j+1][i].T + PetscMax(PetscRealPart(x[j][i].T),0.0);
1480: } else if (i==0) {
1481: /* slab inflow boundary */
1482: f[j][i].T = x[j][i].T - PlateModel(j,PLATE_SLAB,user);
1484: } else if (i==ilim) {
1485: /* right side boundary */
1486: mag_u = 1.0 - pow( (1.0-PetscMax(PetscMin(PetscRealPart(x[j][i-1].u)/param->cb,1.0),0.0)), 5.0 );
1487: 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);
1489: } else if (j==jlim) {
1490: /* bottom boundary */
1491: mag_w = 1.0 - pow( (1.0-PetscMax(PetscMin(PetscRealPart(x[j-1][i].w)/param->sb,1.0),0.0)), 5.0 );
1492: f[j][i].T = x[j][i].T - mag_w*x[j-1][i-1].T - (1.0-mag_w);
1494: } else {
1495: /* in the mantle wedge */
1496: f[j][i].T = EnergyResidual(x,i,j,user);
1497: }
1498: }
1499: }
1500: return(0);
1501: }