Actual source code: ex30.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  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";



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}$-----

 55: #include <petscsnes.h>
 56: #include <petscdm.h>
 57: #include <petscdmda.h>

 59: #define VISC_CONST   0
 60: #define VISC_DIFN    1
 61: #define VISC_DISL    2
 62: #define VISC_FULL    3
 63: #define CELL_CENTER  0
 64: #define CELL_CORNER  1
 65: #define BC_ANALYTIC  0
 66: #define BC_NOSTRESS  1
 67: #define BC_EXPERMNT  2
 68: #define ADVECT_FV    0
 69: #define ADVECT_FROMM 1
 70: #define PLATE_SLAB   0
 71: #define PLATE_LID    1
 72: #define EPS_ZERO     0.00000001

 74: typedef struct { /* holds the variables to be solved for */
 75:   PetscScalar u,w,p,T;
 76: } Field;

 78: typedef struct { /* parameters needed to compute viscosity */
 79:   PetscReal A,n,Estar,Vstar;
 80: } ViscParam;

 82: typedef struct { /* physical and miscelaneous parameters */
 83:   PetscReal width, depth, scaled_width, scaled_depth, peclet, potentialT;
 84:   PetscReal slab_dip, slab_age, slab_velocity, kappa, z_scale;
 85:   PetscReal c, d, sb, cb, skt, visc_cutoff, lid_age, eta0, continuation;
 86:   PetscReal L, V, lid_depth, fault_depth;
 87:   ViscParam diffusion, dislocation;
 88:   PetscInt  ivisc, adv_scheme, ibound, output_ivisc;
 89:   PetscBool quiet, param_test, output_to_file, pv_analytic;
 90:   PetscBool interrupted, stop_solve, toggle_kspmon, kspmon;
 91:   char      filename[PETSC_MAX_PATH_LEN];
 92: } Parameter;

 94: typedef struct { /* grid parameters */
 95:   DMBoundaryType   bx,by;
 96:   DMDAStencilType  stencil;
 97:   PetscInt         corner,ni,nj,jlid,jfault,inose;
 98:   PetscInt         dof,stencil_width,mglevels;
 99:   PetscReal        dx,dz;
100: } GridInfo;

102: typedef struct { /* application context */
103:   Vec       x,Xguess;
104:   Parameter *param;
105:   GridInfo  *grid;
106: } AppCtx;

108: /* Callback functions (static interface) */
109: extern PetscErrorCode FormFunctionLocal(DMDALocalInfo*,Field**,Field**,void*);

111: /* Main routines */
112: extern PetscErrorCode SetParams(Parameter*, GridInfo*);
113: extern PetscErrorCode ReportParams(Parameter*, GridInfo*);
114: extern PetscErrorCode Initialize(DM);
115: extern PetscErrorCode UpdateSolution(SNES,AppCtx*, PetscInt*);
116: extern PetscErrorCode DoOutput(SNES,PetscInt);

118: /* Post-processing & misc */
119: extern PetscErrorCode ViscosityField(DM,Vec,Vec);
120: extern PetscErrorCode StressField(DM);
121: extern PetscErrorCode SNESConverged_Interactive(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason*, void*);
122: extern PetscErrorCode InteractiveHandler(int, void*);
123: extern PetscBool  OptionsHasName(const char pre[],const char name[]);

125: /*-----------------------------------------------------------------------*/
128: int main(int argc,char **argv)
129: /*-----------------------------------------------------------------------*/
130: {
131:   SNES           snes;
132:   AppCtx         *user;               /* user-defined work context */
133:   Parameter      param;
134:   GridInfo       grid;
135:   PetscInt       nits;
137:   MPI_Comm       comm;
138:   DM             da;

140:   PetscInitialize(&argc,&argv,(char*)0,help);
141:   PetscOptionsSetValue(NULL,"-file","ex30_output");
142:   PetscOptionsSetValue(NULL,"-snes_monitor_short",NULL);
143:   PetscOptionsSetValue(NULL,"-snes_max_it","20");
144:   PetscOptionsSetValue(NULL,"-ksp_max_it","1500");
145:   PetscOptionsSetValue(NULL,"-ksp_gmres_restart","300");
146:   PetscOptionsInsert(NULL,&argc,&argv,NULL);

148:   comm = PETSC_COMM_WORLD;

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151:      Set up the problem parameters.
152:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153:   SetParams(&param,&grid);
154:   ReportParams(&param,&grid);

156:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158:   SNESCreate(comm,&snes);
159:   DMDACreate2d(comm,grid.bx,grid.by,grid.stencil,grid.ni,grid.nj,PETSC_DECIDE,PETSC_DECIDE,grid.dof,grid.stencil_width,0,0,&da);
160:   SNESSetDM(snes,da);
161:   DMDASetFieldName(da,0,"x-velocity");
162:   DMDASetFieldName(da,1,"y-velocity");
163:   DMDASetFieldName(da,2,"pressure");
164:   DMDASetFieldName(da,3,"temperature");


167:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168:      Create user context, set problem data, create vector data structures.
169:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170:   PetscNew(&user);
171:   user->param = &param;
172:   user->grid  = &grid;
173:   DMSetApplicationContext(da,user);
174:   DMCreateGlobalVector(da,&(user->Xguess));


177:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178:      Set up the SNES solver with callback functions.
179:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180:   DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,(void*)user);
181:   SNESSetFromOptions(snes);


184:   SNESSetConvergenceTest(snes,SNESConverged_Interactive,(void*)user,NULL);
185:   PetscPushSignalHandler(InteractiveHandler,(void*)user);

187:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188:      Initialize and solve the nonlinear system
189:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
190:   Initialize(da);
191:   UpdateSolution(snes,user,&nits);

193:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194:      Output variables.
195:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
196:   DoOutput(snes,nits);

198:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199:      Free work space.
200:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201:   VecDestroy(&user->Xguess);
202:   VecDestroy(&user->x);
203:   PetscFree(user);
204:   SNESDestroy(&snes);
205:   DMDestroy(&da);
206:   PetscPopSignalHandler();
207:   PetscFinalize();
208:   return 0;
209: }

211: /*=====================================================================
212:   PETSc INTERACTION FUNCTIONS (initialize & call SNESSolve)
213:   =====================================================================*/

215: /*---------------------------------------------------------------------*/
218: /*  manages solve: adaptive continuation method  */
219: PetscErrorCode UpdateSolution(SNES snes, AppCtx *user, PetscInt *nits)
220: {
221:   KSP                 ksp;
222:   PC                  pc;
223:   SNESConvergedReason reason;
224:   Parameter           *param   = user->param;
225:   PetscReal           cont_incr=0.3;
226:   PetscInt            its;
227:   PetscErrorCode      ierr;
228:   PetscBool           q = PETSC_FALSE;
229:   DM                  dm;

232:   SNESGetDM(snes,&dm);
233:   DMCreateGlobalVector(dm,&user->x);
234:   SNESGetKSP(snes,&ksp);
235:   KSPGetPC(ksp, &pc);
236:   KSPSetComputeSingularValues(ksp, PETSC_TRUE);

238:   *nits=0;

240:   /* Isoviscous solve */
241:   if (param->ivisc == VISC_CONST && !param->stop_solve) {
242:     param->ivisc = VISC_CONST;

244:     SNESSolve(snes,0,user->x);
245:     SNESGetConvergedReason(snes,&reason);
246:     SNESGetIterationNumber(snes,&its);
247:     *nits += its;
248:     VecCopy(user->x,user->Xguess);
249:     if (param->stop_solve) goto done;
250:   }

252:   /* Olivine diffusion creep */
253:   if (param->ivisc >= VISC_DIFN && !param->stop_solve) {
254:     if (!q) PetscPrintf(PETSC_COMM_WORLD,"Computing Variable Viscosity Solution\n");

256:     /* continuation method on viscosity cutoff */
257:     for (param->continuation=0.0;; param->continuation+=cont_incr) {
258:       if (!q) PetscPrintf(PETSC_COMM_WORLD," Continuation parameter = %g\n", (double)param->continuation);

260:       /* solve the non-linear system */
261:       VecCopy(user->Xguess,user->x);
262:       SNESSolve(snes,0,user->x);
263:       SNESGetConvergedReason(snes,&reason);
264:       SNESGetIterationNumber(snes,&its);
265:       *nits += its;
266:       if (!q) PetscPrintf(PETSC_COMM_WORLD," SNES iterations: %D, Cumulative: %D\n", its, *nits);
267:       if (param->stop_solve) goto done;

269:       if (reason<0) {
270:         /* NOT converged */
271:         cont_incr = -fabs(cont_incr)/2.0;
272:         if (fabs(cont_incr)<0.01) goto done;

274:       } else {
275:         /* converged */
276:         VecCopy(user->x,user->Xguess);
277:         if (param->continuation >= 1.0) goto done;
278:         if (its<=3)      cont_incr = 0.30001;
279:         else if (its<=8) cont_incr = 0.15001;
280:         else             cont_incr = 0.10001;

282:         if (param->continuation+cont_incr > 1.0) cont_incr = 1.0 - param->continuation;
283:       } /* endif reason<0 */
284:     }
285:   }
286: done:
287:   if (param->stop_solve && !q) PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: stopping solve.\n");
288:   if (reason<0 && !q) PetscPrintf(PETSC_COMM_WORLD,"FAILED TO CONVERGE: stopping solve.\n");
289:   return(0);
290: }


293: /*=====================================================================
294:   PHYSICS FUNCTIONS (compute the discrete residual)
295:   =====================================================================*/


298: /*---------------------------------------------------------------------*/
301: PETSC_STATIC_INLINE PetscScalar UInterp(Field **x, PetscInt i, PetscInt j)
302: /*---------------------------------------------------------------------*/
303: {
304:   return 0.25*(x[j][i].u+x[j+1][i].u+x[j][i+1].u+x[j+1][i+1].u);
305: }

307: /*---------------------------------------------------------------------*/
310: PETSC_STATIC_INLINE PetscScalar WInterp(Field **x, PetscInt i, PetscInt j)
311: /*---------------------------------------------------------------------*/
312: {
313:   return 0.25*(x[j][i].w+x[j+1][i].w+x[j][i+1].w+x[j+1][i+1].w);
314: }

316: /*---------------------------------------------------------------------*/
319: PETSC_STATIC_INLINE PetscScalar PInterp(Field **x, PetscInt i, PetscInt j)
320: /*---------------------------------------------------------------------*/
321: {
322:   return 0.25*(x[j][i].p+x[j+1][i].p+x[j][i+1].p+x[j+1][i+1].p);
323: }

325: /*---------------------------------------------------------------------*/
328: PETSC_STATIC_INLINE PetscScalar TInterp(Field **x, PetscInt i, PetscInt j)
329: /*---------------------------------------------------------------------*/
330: {
331:   return 0.25*(x[j][i].T+x[j+1][i].T+x[j][i+1].T+x[j+1][i+1].T);
332: }

334: /*---------------------------------------------------------------------*/
337: /*  isoviscous analytic solution for IC */
338: PETSC_STATIC_INLINE PetscScalar HorizVelocity(PetscInt i, PetscInt j, AppCtx *user)
339: /*---------------------------------------------------------------------*/
340: {
341:   Parameter   *param = user->param;
342:   GridInfo    *grid  = user->grid;
343:   PetscScalar st, ct, th, c=param->c, d=param->d;
344:   PetscReal   x, z,r;

346:   x  = (i - grid->jlid)*grid->dx;  z = (j - grid->jlid - 0.5)*grid->dz;
347:   r  = PetscSqrtReal(x*x+z*z);
348:   st = z/r;
349:   ct = x/r;
350:   th = atan(z/x);
351:   return ct*(c*th*st+d*(st+th*ct)) + st*(c*(st-th*ct)+d*th*st);
352: }

354: /*---------------------------------------------------------------------*/
357: /*  isoviscous analytic solution for IC */
358: PETSC_STATIC_INLINE PetscScalar VertVelocity(PetscInt i, PetscInt j, AppCtx *user)
359: /*---------------------------------------------------------------------*/
360: {
361:   Parameter   *param = user->param;
362:   GridInfo    *grid  = user->grid;
363:   PetscScalar st, ct, th, c=param->c, d=param->d;
364:   PetscReal   x, z, r;

366:   x = (i - grid->jlid - 0.5)*grid->dx;  z = (j - grid->jlid)*grid->dz;
367:   r = PetscSqrtReal(x*x+z*z); st = z/r;  ct = x/r;  th = atan(z/x);
368:   return st*(c*th*st+d*(st+th*ct)) - ct*(c*(st-th*ct)+d*th*st);
369: }

371: /*---------------------------------------------------------------------*/
374: /*  isoviscous analytic solution for IC */
375: PETSC_STATIC_INLINE PetscScalar Pressure(PetscInt i, PetscInt j, AppCtx *user)
376: /*---------------------------------------------------------------------*/
377: {
378:   Parameter   *param = user->param;
379:   GridInfo    *grid  = user->grid;
380:   PetscScalar x, z, r, st, ct, c=param->c, d=param->d;

382:   x = (i - grid->jlid - 0.5)*grid->dx;  z = (j - grid->jlid - 0.5)*grid->dz;
383:   r = PetscSqrtReal(x*x+z*z);  st = z/r;  ct = x/r;
384:   return (-2.0*(c*ct-d*st)/r);
385: }

389: /*  computes the second invariant of the strain rate tensor */
390: PETSC_STATIC_INLINE PetscScalar CalcSecInv(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
391: /*---------------------------------------------------------------------*/
392: {
393:   Parameter   *param = user->param;
394:   GridInfo    *grid  = user->grid;
395:   PetscInt    ilim   =grid->ni-1, jlim=grid->nj-1;
396:   PetscScalar uN,uS,uE,uW,wN,wS,wE,wW;
397:   PetscScalar eps11, eps12, eps22;

399:   if (i<j) return EPS_ZERO;
400:   if (i==ilim) i--;
401:   if (j==jlim) j--;

403:   if (ipos==CELL_CENTER) { /* on cell center */
404:     if (j<=grid->jlid) return EPS_ZERO;

406:     uE = x[j][i].u; uW = x[j][i-1].u;
407:     wN = x[j][i].w; wS = x[j-1][i].w;
408:     wE = WInterp(x,i,j-1);
409:     if (i==j) {
410:       uN = param->cb; wW = param->sb;
411:     } else {
412:       uN = UInterp(x,i-1,j); wW = WInterp(x,i-1,j-1);
413:     }

415:     if (j==grid->jlid+1) uS = 0.0;
416:     else                 uS = UInterp(x,i-1,j-1);

418:   } else {       /* on CELL_CORNER */
419:     if (j<grid->jlid) return EPS_ZERO;

421:     uN = x[j+1][i].u;  uS = x[j][i].u;
422:     wE = x[j][i+1].w;  wW = x[j][i].w;
423:     if (i==j) {
424:       wN = param->sb;
425:       uW = param->cb;
426:     } else {
427:       wN = WInterp(x,i,j);
428:       uW = UInterp(x,i-1,j);
429:     }

431:     if (j==grid->jlid) {
432:       uE = 0.0;  uW = 0.0;
433:       uS = -uN;
434:       wS = -wN;
435:     } else {
436:       uE = UInterp(x,i,j);
437:       wS = WInterp(x,i,j-1);
438:     }
439:   }

441:   eps11 = (uE-uW)/grid->dx;  eps22 = (wN-wS)/grid->dz;
442:   eps12 = 0.5*((uN-uS)/grid->dz + (wE-wW)/grid->dx);

444:   return PetscSqrtReal(0.5*(eps11*eps11 + 2.0*eps12*eps12 + eps22*eps22));
445: }

447: /*---------------------------------------------------------------------*/
450: /*  computes the shear viscosity */
451: PETSC_STATIC_INLINE PetscScalar Viscosity(PetscScalar T, PetscScalar eps, PetscScalar z, Parameter *param)
452: /*---------------------------------------------------------------------*/
453: {
454:   PetscReal   result   =0.0;
455:   ViscParam   difn     =param->diffusion, disl=param->dislocation;
456:   PetscInt    iVisc    =param->ivisc;
457:   PetscScalar eps_scale=param->V/(param->L*1000.0);
458:   PetscScalar strain_power, v1, v2, P;
459:   PetscScalar rho_g = 32340.0, R=8.3144;

461:   P = rho_g*(z*param->L*1000.0); /* Pa */

463:   if (iVisc==VISC_CONST) {
464:     /* constant viscosity */
465:     return 1.0;
466:   } else if (iVisc==VISC_DIFN) {
467:     /* diffusion creep rheology */
468:     result = PetscRealPart((difn.A*PetscExpScalar((difn.Estar + P*difn.Vstar)/R/(T+273.0))/param->eta0));
469:   } else if (iVisc==VISC_DISL) {
470:     /* dislocation creep rheology */
471:     strain_power = PetscPowScalar(eps*eps_scale, (1.0-disl.n)/disl.n);

473:     result = PetscRealPart(disl.A*PetscExpScalar((disl.Estar + P*disl.Vstar)/disl.n/R/(T+273.0))*strain_power/param->eta0);
474:   } else if (iVisc==VISC_FULL) {
475:     /* dislocation/diffusion creep rheology */
476:     strain_power = PetscPowScalar(eps*eps_scale, (1.0-disl.n)/disl.n);

478:     v1 = difn.A*PetscExpScalar((difn.Estar + P*difn.Vstar)/R/(T+273.0))/param->eta0;
479:     v2 = disl.A*PetscExpScalar((disl.Estar + P*disl.Vstar)/disl.n/R/(T+273.0))*strain_power/param->eta0;

481:     result = PetscRealPart(1.0/(1.0/v1 + 1.0/v2));
482:   }

484:   /* max viscosity is param->eta0 */
485:   result = PetscMin(result, 1.0);
486:   /* min viscosity is param->visc_cutoff */
487:   result = PetscMax(result, param->visc_cutoff);
488:   /* continuation method */
489:   result = PetscPowReal(result,param->continuation);
490:   return result;
491: }

493: /*---------------------------------------------------------------------*/
496: /*  computes the residual of the x-component of eqn (1) above */
497: PETSC_STATIC_INLINE PetscScalar XMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
498: /*---------------------------------------------------------------------*/
499: {
500:   Parameter   *param=user->param;
501:   GridInfo    *grid =user->grid;
502:   PetscScalar dx    = grid->dx, dz=grid->dz;
503:   PetscScalar etaN,etaS,etaE,etaW,epsN=0.0,epsS=0.0,epsE=0.0,epsW=0.0;
504:   PetscScalar TE=0.0,TN=0.0,TS=0.0,TW=0.0, dPdx, residual, z_scale;
505:   PetscScalar dudxW,dudxE,dudzN,dudzS,dwdxN,dwdxS;
506:   PetscInt    jlim = grid->nj-1;

508:   z_scale = param->z_scale;

510:   if (param->ivisc==VISC_DIFN || param->ivisc>=VISC_DISL) { /* viscosity is T-dependent */
511:     TS = param->potentialT * TInterp(x,i,j-1) * PetscExpScalar((j-1.0)*dz*z_scale);
512:     if (j==jlim) TN = TS;
513:     else         TN = param->potentialT * TInterp(x,i,j) * PetscExpScalar(j*dz*z_scale);
514:     TW = param->potentialT * x[j][i].T        * PetscExpScalar((j-0.5)*dz*z_scale);
515:     TE = param->potentialT * x[j][i+1].T      * PetscExpScalar((j-0.5)*dz*z_scale);
516:     if (param->ivisc>=VISC_DISL) { /* olivine dislocation creep */
517:       epsN = CalcSecInv(x,i,j,  CELL_CORNER,user);
518:       epsS = CalcSecInv(x,i,j-1,CELL_CORNER,user);
519:       epsE = CalcSecInv(x,i+1,j,CELL_CENTER,user);
520:       epsW = CalcSecInv(x,i,j,  CELL_CENTER,user);
521:     }
522:   }
523:   etaN = Viscosity(TN,epsN,dz*(j+0.5),param);
524:   etaS = Viscosity(TS,epsS,dz*(j-0.5),param);
525:   etaW = Viscosity(TW,epsW,dz*j,param);
526:   etaE = Viscosity(TE,epsE,dz*j,param);

528:   dPdx = (x[j][i+1].p - x[j][i].p)/dx;
529:   if (j==jlim) dudzN = etaN * (x[j][i].w   - x[j][i+1].w)/dx;
530:   else         dudzN = etaN * (x[j+1][i].u - x[j][i].u)  /dz;
531:   dudzS = etaS * (x[j][i].u    - x[j-1][i].u)/dz;
532:   dudxE = etaE * (x[j][i+1].u  - x[j][i].u)  /dx;
533:   dudxW = etaW * (x[j][i].u    - x[j][i-1].u)/dx;

535:   residual = -dPdx                          /* X-MOMENTUM EQUATION*/
536:              +(dudxE - dudxW)/dx
537:              +(dudzN - dudzS)/dz;

539:   if (param->ivisc!=VISC_CONST) {
540:     dwdxN = etaN * (x[j][i+1].w   - x[j][i].w)  /dx;
541:     dwdxS = etaS * (x[j-1][i+1].w - x[j-1][i].w)/dx;

543:     residual += (dudxE - dudxW)/dx + (dwdxN - dwdxS)/dz;
544:   }

546:   return residual;
547: }

549: /*---------------------------------------------------------------------*/
552: /*  computes the residual of the z-component of eqn (1) above */
553: PETSC_STATIC_INLINE PetscScalar ZMomentumResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
554: /*---------------------------------------------------------------------*/
555: {
556:   Parameter   *param=user->param;
557:   GridInfo    *grid =user->grid;
558:   PetscScalar dx    = grid->dx, dz=grid->dz;
559:   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;
560:   PetscScalar TE    =0.0,TN=0.0,TS=0.0,TW=0.0, dPdz, residual,z_scale;
561:   PetscScalar dudzE,dudzW,dwdxW,dwdxE,dwdzN,dwdzS;
562:   PetscInt    ilim = grid->ni-1;

564:   /* geometric and other parameters */
565:   z_scale = param->z_scale;

567:   /* viscosity */
568:   if (param->ivisc==VISC_DIFN || param->ivisc>=VISC_DISL) { /* viscosity is T-dependent */
569:     TN = param->potentialT * x[j+1][i].T      * PetscExpScalar((j+0.5)*dz*z_scale);
570:     TS = param->potentialT * x[j][i].T        * PetscExpScalar((j-0.5)*dz*z_scale);
571:     TW = param->potentialT * TInterp(x,i-1,j) * PetscExpScalar(j*dz*z_scale);
572:     if (i==ilim) TE = TW;
573:     else         TE = param->potentialT * TInterp(x,i,j) * PetscExpScalar(j*dz*z_scale);
574:     if (param->ivisc>=VISC_DISL) { /* olivine dislocation creep */
575:       epsN = CalcSecInv(x,i,j+1,CELL_CENTER,user);
576:       epsS = CalcSecInv(x,i,j,  CELL_CENTER,user);
577:       epsE = CalcSecInv(x,i,j,  CELL_CORNER,user);
578:       epsW = CalcSecInv(x,i-1,j,CELL_CORNER,user);
579:     }
580:   }
581:   etaN = Viscosity(TN,epsN,dz*(j+1.0),param);
582:   etaS = Viscosity(TS,epsS,dz*(j+0.0),param);
583:   etaW = Viscosity(TW,epsW,dz*(j+0.5),param);
584:   etaE = Viscosity(TE,epsE,dz*(j+0.5),param);

586:   dPdz  = (x[j+1][i].p - x[j][i].p)/dz;
587:   dwdzN = etaN * (x[j+1][i].w - x[j][i].w)/dz;
588:   dwdzS = etaS * (x[j][i].w - x[j-1][i].w)/dz;
589:   if (i==ilim) dwdxE = etaE * (x[j][i].u   - x[j+1][i].u)/dz;
590:   else         dwdxE = etaE * (x[j][i+1].w - x[j][i].w)  /dx;
591:   dwdxW = 2.0*etaW * (x[j][i].w - x[j][i-1].w)/dx;

593:   /* Z-MOMENTUM */
594:   residual = -dPdz                 /* constant viscosity terms */
595:              +(dwdzN - dwdzS)/dz
596:              +(dwdxE - dwdxW)/dx;

598:   if (param->ivisc!=VISC_CONST) {
599:     dudzE = etaE * (x[j+1][i].u - x[j][i].u)/dz;
600:     dudzW = etaW * (x[j+1][i-1].u - x[j][i-1].u)/dz;

602:     residual += (dwdzN - dwdzS)/dz + (dudzE - dudzW)/dx;
603:   }

605:   return residual;
606: }

608: /*---------------------------------------------------------------------*/
611: /*  computes the residual of eqn (2) above */
612: PETSC_STATIC_INLINE PetscScalar ContinuityResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
613: /*---------------------------------------------------------------------*/
614: {
615:   GridInfo    *grid =user->grid;
616:   PetscScalar uE,uW,wN,wS,dudx,dwdz;

618:   uW = x[j][i-1].u; uE = x[j][i].u; dudx = (uE - uW)/grid->dx;
619:   wS = x[j-1][i].w; wN = x[j][i].w; dwdz = (wN - wS)/grid->dz;

621:   return dudx + dwdz;
622: }

624: /*---------------------------------------------------------------------*/
627: /*  computes the residual of eqn (3) above */
628: PETSC_STATIC_INLINE PetscScalar EnergyResidual(Field **x, PetscInt i, PetscInt j, AppCtx *user)
629: /*---------------------------------------------------------------------*/
630: {
631:   Parameter   *param=user->param;
632:   GridInfo    *grid =user->grid;
633:   PetscScalar dx    = grid->dx, dz=grid->dz;
634:   PetscInt    ilim  =grid->ni-1, jlim=grid->nj-1, jlid=grid->jlid;
635:   PetscScalar TE, TN, TS, TW, residual;
636:   PetscScalar uE,uW,wN,wS;
637:   PetscScalar fN,fS,fE,fW,dTdxW,dTdxE,dTdzN,dTdzS;

639:   dTdzN = (x[j+1][i].T - x[j][i].T)  /dz;
640:   dTdzS = (x[j][i].T   - x[j-1][i].T)/dz;
641:   dTdxE = (x[j][i+1].T - x[j][i].T)  /dx;
642:   dTdxW = (x[j][i].T   - x[j][i-1].T)/dx;

644:   residual = ((dTdzN - dTdzS)/dz + /* diffusion term */
645:               (dTdxE - dTdxW)/dx)*dx*dz/param->peclet;

647:   if (j<=jlid && i>=j) {
648:     /* don't advect in the lid */
649:     return residual;
650:   } else if (i<j) {
651:     /* beneath the slab sfc */
652:     uW = uE = param->cb;
653:     wS = wN = param->sb;
654:   } else {
655:     /* advect in the slab and wedge */
656:     uW = x[j][i-1].u; uE = x[j][i].u;
657:     wS = x[j-1][i].w; wN = x[j][i].w;
658:   }

660:   if (param->adv_scheme==ADVECT_FV || i==ilim-1 || j==jlim-1 || i==1 || j==1) {
661:     /* finite volume advection */
662:     TS = (x[j][i].T + x[j-1][i].T)/2.0;
663:     TN = (x[j][i].T + x[j+1][i].T)/2.0;
664:     TE = (x[j][i].T + x[j][i+1].T)/2.0;
665:     TW = (x[j][i].T + x[j][i-1].T)/2.0;
666:     fN = wN*TN*dx; fS = wS*TS*dx;
667:     fE = uE*TE*dz; fW = uW*TW*dz;

669:   } else {
670:     /* Fromm advection scheme */
671:     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
672:               - PetscAbsScalar(uE)*(-x[j][i+2].T + 3.0*(x[j][i+1].T-x[j][i].T)+x[j][i-1].T)/8.0)*dz;
673:     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
674:               - PetscAbsScalar(uW)*(-x[j][i+1].T + 3.0*(x[j][i].T-x[j][i-1].T)+x[j][i-2].T)/8.0)*dz;
675:     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
676:               - PetscAbsScalar(wN)*(-x[j+2][i].T + 3.0*(x[j+1][i].T-x[j][i].T)+x[j-1][i].T)/8.0)*dx;
677:     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
678:               - PetscAbsScalar(wS)*(-x[j+1][i].T + 3.0*(x[j][i].T-x[j-1][i].T)+x[j-2][i].T)/8.0)*dx;
679:   }

681:   residual -= (fE - fW + fN - fS);

683:   return residual;
684: }

686: /*---------------------------------------------------------------------*/
689: /*  computes the shear stress---used on the boundaries */
690: PETSC_STATIC_INLINE PetscScalar ShearStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
691: /*---------------------------------------------------------------------*/
692: {
693:   Parameter   *param=user->param;
694:   GridInfo    *grid =user->grid;
695:   PetscInt    ilim  =grid->ni-1, jlim=grid->nj-1;
696:   PetscScalar uN, uS, wE, wW;

698:   if (j<=grid->jlid || i<j || i==ilim || j==jlim) return EPS_ZERO;

700:   if (ipos==CELL_CENTER) { /* on cell center */

702:     wE = WInterp(x,i,j-1);
703:     if (i==j) {
704:       wW = param->sb;
705:       uN = param->cb;
706:     } else {
707:       wW = WInterp(x,i-1,j-1);
708:       uN = UInterp(x,i-1,j);
709:     }
710:     if (j==grid->jlid+1) uS = 0.0;
711:     else                 uS = UInterp(x,i-1,j-1);

713:   } else { /* on cell corner */

715:     uN = x[j+1][i].u;         uS = x[j][i].u;
716:     wW = x[j][i].w;           wE = x[j][i+1].w;

718:   }

720:   return (uN-uS)/grid->dz + (wE-wW)/grid->dx;
721: }

723: /*---------------------------------------------------------------------*/
726: /*  computes the normal stress---used on the boundaries */
727: PETSC_STATIC_INLINE PetscScalar XNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
728: /*---------------------------------------------------------------------*/
729: {
730:   Parameter   *param=user->param;
731:   GridInfo    *grid =user->grid;
732:   PetscScalar dx    = grid->dx, dz=grid->dz;
733:   PetscInt    ilim  =grid->ni-1, jlim=grid->nj-1, ivisc;
734:   PetscScalar epsC  =0.0, etaC, TC, uE, uW, pC, z_scale;
735:   if (i<j || j<=grid->jlid) return EPS_ZERO;

737:   ivisc=param->ivisc;  z_scale = param->z_scale;

739:   if (ipos==CELL_CENTER) { /* on cell center */

741:     TC = param->potentialT * x[j][i].T * PetscExpScalar((j-0.5)*dz*z_scale);
742:     if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CENTER,user);
743:     etaC = Viscosity(TC,epsC,dz*j,param);

745:     uW = x[j][i-1].u;   uE = x[j][i].u;
746:     pC = x[j][i].p;

748:   } else { /* on cell corner */
749:     if (i==ilim || j==jlim) return EPS_ZERO;

751:     TC = param->potentialT * TInterp(x,i,j) * PetscExpScalar(j*dz*z_scale);
752:     if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CORNER,user);
753:     etaC = Viscosity(TC,epsC,dz*(j+0.5),param);

755:     if (i==j) uW = param->sb;
756:     else      uW = UInterp(x,i-1,j);
757:     uE = UInterp(x,i,j); pC = PInterp(x,i,j);
758:   }

760:   return 2.0*etaC*(uE-uW)/dx - pC;
761: }

763: /*---------------------------------------------------------------------*/
766: /*  computes the normal stress---used on the boundaries */
767: PETSC_STATIC_INLINE PetscScalar ZNormalStress(Field **x, PetscInt i, PetscInt j, PetscInt ipos, AppCtx *user)
768: /*---------------------------------------------------------------------*/
769: {
770:   Parameter   *param=user->param;
771:   GridInfo    *grid =user->grid;
772:   PetscScalar dz    =grid->dz;
773:   PetscInt    ilim  =grid->ni-1, jlim=grid->nj-1, ivisc;
774:   PetscScalar epsC  =0.0, etaC, TC;
775:   PetscScalar pC, wN, wS, z_scale;
776:   if (i<j || j<=grid->jlid) return EPS_ZERO;

778:   ivisc=param->ivisc;  z_scale = param->z_scale;

780:   if (ipos==CELL_CENTER) { /* on cell center */

782:     TC = param->potentialT * x[j][i].T * PetscExpScalar((j-0.5)*dz*z_scale);
783:     if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CENTER,user);
784:     etaC = Viscosity(TC,epsC,dz*j,param);
785:     wN   = x[j][i].w; wS = x[j-1][i].w; pC = x[j][i].p;

787:   } else { /* on cell corner */
788:     if ((i==ilim) || (j==jlim)) return EPS_ZERO;

790:     TC = param->potentialT * TInterp(x,i,j) * PetscExpScalar(j*dz*z_scale);
791:     if (ivisc>=VISC_DISL) epsC = CalcSecInv(x,i,j,CELL_CORNER,user);
792:     etaC = Viscosity(TC,epsC,dz*(j+0.5),param);
793:     if (i==j) wN = param->sb;
794:     else      wN = WInterp(x,i,j);
795:     wS = WInterp(x,i,j-1); pC = PInterp(x,i,j);
796:   }

798:   return 2.0*etaC*(wN-wS)/dz - pC;
799: }

801: /*---------------------------------------------------------------------*/

803: /*=====================================================================
804:   INITIALIZATION, POST-PROCESSING AND OUTPUT FUNCTIONS
805:   =====================================================================*/

807: /*---------------------------------------------------------------------*/
810: /* initializes the problem parameters and checks for
811:    command line changes */
812: PetscErrorCode SetParams(Parameter *param, GridInfo *grid)
813: /*---------------------------------------------------------------------*/
814: {
815:   PetscErrorCode ierr, ierr_out=0;
816:   PetscReal      SEC_PER_YR                    = 3600.00*24.00*365.2500;
817:   PetscReal      alpha_g_on_cp_units_inverse_km=4.0e-5*9.8;

819:   /* domain geometry */
820:   param->slab_dip    = 45.0;
821:   param->width       = 320.0;                                              /* km */
822:   param->depth       = 300.0;                                              /* km */
823:   param->lid_depth   = 35.0;                                               /* km */
824:   param->fault_depth = 35.0;                                               /* km */

826:   PetscOptionsGetReal(NULL,NULL,"-slab_dip",&(param->slab_dip),NULL);
827:   PetscOptionsGetReal(NULL,NULL,"-width",&(param->width),NULL);
828:   PetscOptionsGetReal(NULL,NULL,"-depth",&(param->depth),NULL);
829:   PetscOptionsGetReal(NULL,NULL,"-lid_depth",&(param->lid_depth),NULL);
830:   PetscOptionsGetReal(NULL,NULL,"-fault_depth",&(param->fault_depth),NULL);

832:   param->slab_dip = param->slab_dip*PETSC_PI/180.0;                    /* radians */

834:   /* grid information */
835:   PetscOptionsGetInt(NULL,NULL, "-jfault",&(grid->jfault),NULL);
836:   grid->ni = 82;
837:   PetscOptionsGetInt(NULL,NULL, "-ni",&(grid->ni),NULL);

839:   grid->dx            = param->width/((PetscReal)(grid->ni-2));               /* km */
840:   grid->dz            = grid->dx*tan(param->slab_dip);                     /* km */
841:   grid->nj            = (PetscInt)(param->depth/grid->dz + 3.0);         /* gridpoints*/
842:   param->depth        = grid->dz*(grid->nj-2);                             /* km */
843:   grid->inose         = 0;                                          /* gridpoints*/
844:   PetscOptionsGetInt(NULL,NULL,"-inose",&(grid->inose),NULL);
845:   grid->bx            = DM_BOUNDARY_NONE;
846:   grid->by            = DM_BOUNDARY_NONE;
847:   grid->stencil       = DMDA_STENCIL_BOX;
848:   grid->dof           = 4;
849:   grid->stencil_width = 2;
850:   grid->mglevels      = 1;

852:   /* boundary conditions */
853:   param->pv_analytic = PETSC_FALSE;
854:   param->ibound      = BC_NOSTRESS;
855:   PetscOptionsGetInt(NULL,NULL,"-ibound",&(param->ibound),NULL);

857:   /* physical constants */
858:   param->slab_velocity = 5.0;               /* cm/yr */
859:   param->slab_age      = 50.0;              /* Ma */
860:   param->lid_age       = 50.0;              /* Ma */
861:   param->kappa         = 0.7272e-6;         /* m^2/sec */
862:   param->potentialT    = 1300.0;            /* degrees C */

864:   PetscOptionsGetReal(NULL,NULL,"-slab_velocity",&(param->slab_velocity),NULL);
865:   PetscOptionsGetReal(NULL,NULL,"-slab_age",&(param->slab_age),NULL);
866:   PetscOptionsGetReal(NULL,NULL,"-lid_age",&(param->lid_age),NULL);
867:   PetscOptionsGetReal(NULL,NULL,"-kappa",&(param->kappa),NULL);
868:   PetscOptionsGetReal(NULL,NULL,"-potentialT",&(param->potentialT),NULL);

870:   /* viscosity */
871:   param->ivisc        = 3;                  /* 0=isovisc, 1=difn creep, 2=disl creep, 3=full */
872:   param->eta0         = 1e24;               /* Pa-s */
873:   param->visc_cutoff  = 0.0;                /* factor of eta_0 */
874:   param->continuation = 1.0;

876:   /* constants for diffusion creep */
877:   param->diffusion.A     = 1.8e7;             /* Pa-s */
878:   param->diffusion.n     = 1.0;               /* dim'less */
879:   param->diffusion.Estar = 375e3;             /* J/mol */
880:   param->diffusion.Vstar = 5e-6;              /* m^3/mol */

882:   /* constants for param->dislocationocation creep */
883:   param->dislocation.A     = 2.8969e4;        /* Pa-s */
884:   param->dislocation.n     = 3.5;             /* dim'less */
885:   param->dislocation.Estar = 530e3;           /* J/mol */
886:   param->dislocation.Vstar = 14e-6;           /* m^3/mol */

888:   PetscOptionsGetInt(NULL,NULL, "-ivisc",&(param->ivisc),NULL);
889:   PetscOptionsGetReal(NULL,NULL,"-visc_cutoff",&(param->visc_cutoff),NULL);

891:   param->output_ivisc = param->ivisc;

893:   PetscOptionsGetInt(NULL,NULL,"-output_ivisc",&(param->output_ivisc),NULL);
894:   PetscOptionsGetReal(NULL,NULL,"-vstar",&(param->dislocation.Vstar),NULL);

896:   /* output options */
897:   param->quiet      = PETSC_FALSE;
898:   param->param_test = PETSC_FALSE;

900:   PetscOptionsHasName(NULL,NULL,"-quiet",&(param->quiet));
901:   PetscOptionsHasName(NULL,NULL,"-test",&(param->param_test));
902:   PetscOptionsGetString(NULL,NULL,"-file",param->filename,PETSC_MAX_PATH_LEN,&(param->output_to_file));

904:   /* advection */
905:   param->adv_scheme = ADVECT_FROMM;       /* advection scheme: 0=finite vol, 1=Fromm */

907:   PetscOptionsGetInt(NULL,NULL,"-adv_scheme",&(param->adv_scheme),NULL);

909:   /* misc. flags */
910:   param->stop_solve    = PETSC_FALSE;
911:   param->interrupted   = PETSC_FALSE;
912:   param->kspmon        = PETSC_FALSE;
913:   param->toggle_kspmon = PETSC_FALSE;

915:   /* derived parameters for slab angle */
916:   param->sb = PetscSinReal(param->slab_dip);
917:   param->cb = PetscCosReal(param->slab_dip);
918:   param->c  =  param->slab_dip*param->sb/(param->slab_dip*param->slab_dip-param->sb*param->sb);
919:   param->d  = (param->slab_dip*param->cb-param->sb)/(param->slab_dip*param->slab_dip-param->sb*param->sb);

921:   /* length, velocity and time scale for non-dimensionalization */
922:   param->L = PetscMin(param->width,param->depth);               /* km */
923:   param->V = param->slab_velocity/100.0/SEC_PER_YR;             /* m/sec */

925:   /* other unit conversions and derived parameters */
926:   param->scaled_width = param->width/param->L;                  /* dim'less */
927:   param->scaled_depth = param->depth/param->L;                  /* dim'less */
928:   param->lid_depth    = param->lid_depth/param->L;              /* dim'less */
929:   param->fault_depth  = param->fault_depth/param->L;            /* dim'less */
930:   grid->dx            = grid->dx/param->L;                      /* dim'less */
931:   grid->dz            = grid->dz/param->L;                      /* dim'less */
932:   grid->jlid          = (PetscInt)(param->lid_depth/grid->dz);       /* gridcells */
933:   grid->jfault        = (PetscInt)(param->fault_depth/grid->dz);     /* gridcells */
934:   param->lid_depth    = grid->jlid*grid->dz;                    /* dim'less */
935:   param->fault_depth  = grid->jfault*grid->dz;                  /* dim'less */
936:   grid->corner        = grid->jlid+1;                           /* gridcells */
937:   param->peclet       = param->V                                /* m/sec */
938:                         * param->L*1000.0                       /* m */
939:                         / param->kappa;                         /* m^2/sec */
940:   param->z_scale = param->L * alpha_g_on_cp_units_inverse_km;
941:   param->skt     = PetscSqrtReal(param->kappa*param->slab_age*SEC_PER_YR);
942:   PetscOptionsGetReal(NULL,NULL,"-peclet",&(param->peclet),NULL);

944:   return ierr_out;
945: }

947: /*---------------------------------------------------------------------*/
950: /*  prints a report of the problem parameters to stdout */
951: PetscErrorCode ReportParams(Parameter *param, GridInfo *grid)
952: /*---------------------------------------------------------------------*/
953: {
954:   PetscErrorCode ierr, ierr_out=0;
955:   char           date[30];

957:   PetscGetDate(date,30);

959:   if (!(param->quiet)) {
960:     PetscPrintf(PETSC_COMM_WORLD,"---------------------BEGIN ex30 PARAM REPORT-------------------\n");
961:     /* PetscPrintf(PETSC_COMM_WORLD,"                   %s\n",&(date[0]));*/

963:     PetscPrintf(PETSC_COMM_WORLD,"Domain: \n");
964:     PetscPrintf(PETSC_COMM_WORLD,"  Width = %g km,         Depth = %g km\n",(double)param->width,(double)param->depth);
965:     PetscPrintf(PETSC_COMM_WORLD,"  Slab dip = %g degrees,  Slab velocity = %g cm/yr\n",(double)(param->slab_dip*180.0/PETSC_PI),(double)param->slab_velocity);
966:     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);

968:     PetscPrintf(PETSC_COMM_WORLD,"\nGrid: \n");
969:     PetscPrintf(PETSC_COMM_WORLD,"  [ni,nj] = %D, %D       [dx,dz] = %g, %g km\n",grid->ni,grid->nj,(double)grid->dx*param->L,(double)(grid->dz*param->L));
970:     PetscPrintf(PETSC_COMM_WORLD,"  jlid = %3D              jfault = %3D \n",grid->jlid,grid->jfault);
971:     PetscPrintf(PETSC_COMM_WORLD,"  Pe = %g\n",(double)param->peclet);

973:     PetscPrintf(PETSC_COMM_WORLD,"\nRheology:");
974:     if (param->ivisc==VISC_CONST) {
975:       PetscPrintf(PETSC_COMM_WORLD,"                 Isoviscous \n");
976:       if (param->pv_analytic) {
977:         PetscPrintf(PETSC_COMM_WORLD,"                          Pressure and Velocity prescribed! \n");
978:       }
979:     } else if (param->ivisc==VISC_DIFN) {
980:       PetscPrintf(PETSC_COMM_WORLD,"                 Diffusion Creep (T-Dependent Newtonian) \n");
981:       PetscPrintf(PETSC_COMM_WORLD,"                          Viscosity range: %g--%g Pa-sec \n",(double)param->eta0,(double)(param->visc_cutoff*param->eta0));
982:     } else if (param->ivisc==VISC_DISL) {
983:       PetscPrintf(PETSC_COMM_WORLD,"                 Dislocation Creep (T-Dependent Non-Newtonian) \n");
984:       PetscPrintf(PETSC_COMM_WORLD,"                          Viscosity range: %g--%g Pa-sec \n",(double)param->eta0,(double)(param->visc_cutoff*param->eta0));
985:     } else if (param->ivisc==VISC_FULL) {
986:       PetscPrintf(PETSC_COMM_WORLD,"                 Full Rheology \n");
987:       PetscPrintf(PETSC_COMM_WORLD,"                          Viscosity range: %g--%g Pa-sec \n",(double)param->eta0,(double)(param->visc_cutoff*param->eta0));
988:     } else {
989:       PetscPrintf(PETSC_COMM_WORLD,"                 Invalid! \n");
990:       ierr_out = 1;
991:     }

993:     PetscPrintf(PETSC_COMM_WORLD,"Boundary condition:");
994:     if (param->ibound==BC_ANALYTIC) {
995:       PetscPrintf(PETSC_COMM_WORLD,"       Isoviscous Analytic Dirichlet \n");
996:     } else if (param->ibound==BC_NOSTRESS) {
997:       PetscPrintf(PETSC_COMM_WORLD,"       Stress-Free (normal & shear stress)\n");
998:     } else if (param->ibound==BC_EXPERMNT) {
999:       PetscPrintf(PETSC_COMM_WORLD,"       Experimental boundary condition \n");
1000:     } else {
1001:       PetscPrintf(PETSC_COMM_WORLD,"       Invalid! \n");
1002:       ierr_out = 1;
1003:     }

1005:     if (param->output_to_file)
1006: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1007:       PetscPrintf(PETSC_COMM_WORLD,"Output Destination:       Mat file \"%s\"\n",param->filename);
1008: #else
1009:       PetscPrintf(PETSC_COMM_WORLD,"Output Destination:       PETSc binary file \"%s\"\n",param->filename);
1010: #endif
1011:     if (param->output_ivisc != param->ivisc) {
1012:       PetscPrintf(PETSC_COMM_WORLD,"                          Output viscosity: -ivisc %D\n",param->output_ivisc);
1013:     }

1015:     PetscPrintf(PETSC_COMM_WORLD,"---------------------END ex30 PARAM REPORT---------------------\n");
1016:   }
1017:   if (param->param_test) PetscEnd();
1018:   return ierr_out;
1019: }

1021: /* ------------------------------------------------------------------- */
1024: /*  generates an inital guess using the analytic solution for isoviscous
1025:     corner flow */
1026: PetscErrorCode Initialize(DM da)
1027: /* ------------------------------------------------------------------- */
1028: {
1029:   AppCtx         *user;
1030:   Parameter      *param;
1031:   GridInfo       *grid;
1032:   PetscInt       i,j,is,js,im,jm;
1034:   Field          **x;
1035:   Vec            Xguess;

1037:   /* Get the fine grid */
1038:   DMGetApplicationContext(da,&user);
1039:   Xguess = user->Xguess;
1040:   param  = user->param;
1041:   grid   = user->grid;
1042:   DMDAGetCorners(da,&is,&js,NULL,&im,&jm,NULL);
1043:   DMDAVecGetArray(da,Xguess,(void**)&x);

1045:   /* Compute initial guess */
1046:   for (j=js; j<js+jm; j++) {
1047:     for (i=is; i<is+im; i++) {
1048:       if (i<j)                x[j][i].u = param->cb;
1049:       else if (j<=grid->jlid) x[j][i].u = 0.0;
1050:       else                    x[j][i].u = HorizVelocity(i,j,user);

1052:       if (i<=j)               x[j][i].w = param->sb;
1053:       else if (j<=grid->jlid) x[j][i].w = 0.0;
1054:       else                    x[j][i].w = VertVelocity(i,j,user);

1056:       if (i<j || j<=grid->jlid) x[j][i].p = 0.0;
1057:       else                      x[j][i].p = Pressure(i,j,user);

1059:       x[j][i].T = PetscMin(grid->dz*(j-0.5),1.0);
1060:     }
1061:   }

1063:   /* Restore x to Xguess */
1064:   DMDAVecRestoreArray(da,Xguess,(void**)&x);

1066:   return 0;
1067: }

1069: /*---------------------------------------------------------------------*/
1072: /*  controls output to a file */
1073: PetscErrorCode DoOutput(SNES snes, PetscInt its)
1074: /*---------------------------------------------------------------------*/
1075: {
1076:   AppCtx         *user;
1077:   Parameter      *param;
1078:   GridInfo       *grid;
1079:   PetscInt       ivt;
1081:   PetscMPIInt    rank;
1082:   PetscViewer    viewer;
1083:   Vec            res, pars;
1084:   MPI_Comm       comm;
1085:   DM             da;

1087:   SNESGetDM(snes,&da);
1088:   DMGetApplicationContext(da,&user);
1089:   param = user->param;
1090:   grid  = user->grid;
1091:   ivt   = param->ivisc;

1093:   param->ivisc = param->output_ivisc;

1095:   /* compute final residual and final viscosity/strain rate fields */
1096:   SNESGetFunction(snes, &res, NULL, NULL);
1097:   ViscosityField(da, user->x, user->Xguess);

1099:   /* get the communicator and the rank of the processor */
1100:   PetscObjectGetComm((PetscObject)snes, &comm);
1101:   MPI_Comm_rank(comm, &rank);

1103:   if (param->output_to_file) { /* send output to binary file */
1104:     VecCreate(comm, &pars);
1105:     if (!rank) { /* on processor 0 */
1106:       VecSetSizes(pars, 20, PETSC_DETERMINE);
1107:       VecSetFromOptions(pars);
1108:       VecSetValue(pars,0, (PetscScalar)(grid->ni),INSERT_VALUES);
1109:       VecSetValue(pars,1, (PetscScalar)(grid->nj),INSERT_VALUES);
1110:       VecSetValue(pars,2, (PetscScalar)(grid->dx),INSERT_VALUES);
1111:       VecSetValue(pars,3, (PetscScalar)(grid->dz),INSERT_VALUES);
1112:       VecSetValue(pars,4, (PetscScalar)(param->L),INSERT_VALUES);
1113:       VecSetValue(pars,5, (PetscScalar)(param->V),INSERT_VALUES);
1114:       /* skipped 6 intentionally */
1115:       VecSetValue(pars,7, (PetscScalar)(param->slab_dip),INSERT_VALUES);
1116:       VecSetValue(pars,8, (PetscScalar)(grid->jlid),INSERT_VALUES);
1117:       VecSetValue(pars,9, (PetscScalar)(param->lid_depth),INSERT_VALUES);
1118:       VecSetValue(pars,10,(PetscScalar)(grid->jfault),INSERT_VALUES);
1119:       VecSetValue(pars,11,(PetscScalar)(param->fault_depth),INSERT_VALUES);
1120:       VecSetValue(pars,12,(PetscScalar)(param->potentialT),INSERT_VALUES);
1121:       VecSetValue(pars,13,(PetscScalar)(param->ivisc),INSERT_VALUES);
1122:       VecSetValue(pars,14,(PetscScalar)(param->visc_cutoff),INSERT_VALUES);
1123:       VecSetValue(pars,15,(PetscScalar)(param->ibound),INSERT_VALUES);
1124:       VecSetValue(pars,16,(PetscScalar)(its),INSERT_VALUES);
1125:     } else { /* on some other processor */
1126:       VecSetSizes(pars, 0, PETSC_DETERMINE);
1127:       VecSetFromOptions(pars);
1128:     }
1129:     VecAssemblyBegin(pars); VecAssemblyEnd(pars);

1131:     /* create viewer */
1132: #if defined(PETSC_HAVE_MATLAB_ENGINE)
1133:     PetscViewerMatlabOpen(PETSC_COMM_WORLD,param->filename,FILE_MODE_WRITE,&viewer);
1134: #else
1135:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,param->filename,FILE_MODE_WRITE,&viewer);
1136: #endif

1138:     /* send vectors to viewer */
1139:     PetscObjectSetName((PetscObject)res,"res");
1140:     VecView(res,viewer);
1141:     PetscObjectSetName((PetscObject)user->x,"out");
1142:     VecView(user->x, viewer);
1143:     PetscObjectSetName((PetscObject)(user->Xguess),"aux");
1144:     VecView(user->Xguess, viewer);
1145:     StressField(da); /* compute stress fields */
1146:     PetscObjectSetName((PetscObject)(user->Xguess),"str");
1147:     VecView(user->Xguess, viewer);
1148:     PetscObjectSetName((PetscObject)pars,"par");
1149:     VecView(pars, viewer);

1151:     /* destroy viewer and vector */
1152:     PetscViewerDestroy(&viewer);
1153:     VecDestroy(&pars);
1154:   }

1156:   param->ivisc = ivt;
1157:   return 0;
1158: }

1160: /* ------------------------------------------------------------------- */
1163: /* Compute both the second invariant of the strain rate tensor and the viscosity, at both cell centers and cell corners */
1164: PetscErrorCode ViscosityField(DM da, Vec X, Vec V)
1165: /* ------------------------------------------------------------------- */
1166: {
1167:   AppCtx         *user;
1168:   Parameter      *param;
1169:   GridInfo       *grid;
1170:   Vec            localX;
1171:   Field          **v, **x;
1172:   PetscReal      eps, /* dx,*/ dz, T, epsC, TC;
1173:   PetscInt       i,j,is,js,im,jm,ilim,jlim,ivt;

1177:   DMGetApplicationContext(da,&user);
1178:   param        = user->param;
1179:   grid         = user->grid;
1180:   ivt          = param->ivisc;
1181:   param->ivisc = param->output_ivisc;

1183:   DMGetLocalVector(da, &localX);
1184:   DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX);
1185:   DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX);
1186:   DMDAVecGetArray(da,localX,(void**)&x);
1187:   DMDAVecGetArray(da,V,(void**)&v);

1189:   /* Parameters */
1190:   /* dx = grid->dx; */ dz = grid->dz;

1192:   ilim = grid->ni-1; jlim = grid->nj-1;

1194:   /* Compute real temperature, strain rate and viscosity */
1195:   DMDAGetCorners(da,&is,&js,NULL,&im,&jm,NULL);
1196:   for (j=js; j<js+jm; j++) {
1197:     for (i=is; i<is+im; i++) {
1198:       T = PetscRealPart(param->potentialT * x[j][i].T * PetscExpScalar((j-0.5)*dz*param->z_scale));
1199:       if (i<ilim && j<jlim) {
1200:         TC = PetscRealPart(param->potentialT * TInterp(x,i,j) * PetscExpScalar(j*dz*param->z_scale));
1201:       } else {
1202:         TC = T;
1203:       }
1204:       eps  = PetscRealPart((CalcSecInv(x,i,j,CELL_CENTER,user)));
1205:       epsC = PetscRealPart(CalcSecInv(x,i,j,CELL_CORNER,user));

1207:       v[j][i].u = eps;
1208:       v[j][i].w = epsC;
1209:       v[j][i].p = Viscosity(T,eps,dz*(j-0.5),param);
1210:       v[j][i].T = Viscosity(TC,epsC,dz*j,param);
1211:     }
1212:   }
1213:   DMDAVecRestoreArray(da,V,(void**)&v);
1214:   DMDAVecRestoreArray(da,localX,(void**)&x);
1215:   DMRestoreLocalVector(da, &localX);

1217:   param->ivisc = ivt;
1218:   return(0);
1219: }

1221: /* ------------------------------------------------------------------- */
1224: /* post-processing: compute stress everywhere */
1225: PetscErrorCode StressField(DM da)
1226: /* ------------------------------------------------------------------- */
1227: {
1228:   AppCtx         *user;
1229:   PetscInt       i,j,is,js,im,jm;
1231:   Vec            locVec;
1232:   Field          **x, **y;

1234:   DMGetApplicationContext(da,&user);

1236:   /* Get the fine grid of Xguess and X */
1237:   DMDAGetCorners(da,&is,&js,NULL,&im,&jm,NULL);
1238:   DMDAVecGetArray(da,user->Xguess,(void**)&x);

1240:   DMGetLocalVector(da, &locVec);
1241:   DMGlobalToLocalBegin(da, user->x, INSERT_VALUES, locVec);
1242:   DMGlobalToLocalEnd(da, user->x, INSERT_VALUES, locVec);
1243:   DMDAVecGetArray(da,locVec,(void**)&y);

1245:   /* Compute stress on the corner points */
1246:   for (j=js; j<js+jm; j++) {
1247:     for (i=is; i<is+im; i++) {
1248:       x[j][i].u = ShearStress(y,i,j,CELL_CENTER,user);
1249:       x[j][i].w = ShearStress(y,i,j,CELL_CORNER,user);
1250:       x[j][i].p = XNormalStress(y,i,j,CELL_CENTER,user);
1251:       x[j][i].T = ZNormalStress(y,i,j,CELL_CENTER,user);
1252:     }
1253:   }

1255:   /* Restore the fine grid of Xguess and X */
1256:   DMDAVecRestoreArray(da,user->Xguess,(void**)&x);
1257:   DMDAVecRestoreArray(da,locVec,(void**)&y);
1258:   DMRestoreLocalVector(da, &locVec);
1259:   return 0;
1260: }

1262: /*=====================================================================
1263:   UTILITY FUNCTIONS
1264:   =====================================================================*/

1266: /*---------------------------------------------------------------------*/
1269: /* returns the velocity of the subducting slab and handles fault nodes
1270:    for BC */
1271: PETSC_STATIC_INLINE PetscScalar SlabVel(char c, PetscInt i, PetscInt j, AppCtx *user)
1272: /*---------------------------------------------------------------------*/
1273: {
1274:   Parameter *param = user->param;
1275:   GridInfo  *grid  = user->grid;

1277:   if (c=='U' || c=='u') {
1278:     if (i<j-1) return param->cb;
1279:     else if (j<=grid->jfault) return 0.0;
1280:     else return param->cb;

1282:   } else {
1283:     if (i<j) return param->sb;
1284:     else if (j<=grid->jfault) return 0.0;
1285:     else return param->sb;
1286:   }
1287: }

1289: /*---------------------------------------------------------------------*/
1292: /*  solution to diffusive half-space cooling model for BC */
1293: PETSC_STATIC_INLINE PetscScalar PlateModel(PetscInt j, PetscInt plate, AppCtx *user)
1294: /*---------------------------------------------------------------------*/
1295: {
1296:   Parameter     *param = user->param;
1297:   PetscScalar   z;
1298:   if (plate==PLATE_LID) z = (j-0.5)*user->grid->dz;
1299:   else z = (j-0.5)*user->grid->dz*param->cb;  /* PLATE_SLAB */
1300: #if defined(PETSC_HAVE_ERF)
1301:   return (erf(PetscRealPart(z*param->L/2.0/param->skt)));
1302: #else
1303:   SETERRQ(PETSC_COMM_SELF,1,"erf() not available on this machine");
1304: #endif
1305: }


1308: /* ------------------------------------------------------------------- */
1311: /*  utility function */
1312: PetscBool  OptionsHasName(const char pre[],const char name[])
1313: /* ------------------------------------------------------------------- */
1314: {
1315:   PetscBool      retval;
1317:   PetscOptionsHasName(NULL,pre,name,&retval);CHKERRABORT(PETSC_COMM_WORLD,ierr);
1318:   return retval;
1319: }

1321: /*=====================================================================
1322:   INTERACTIVE SIGNAL HANDLING
1323:   =====================================================================*/

1325: /* ------------------------------------------------------------------- */
1328: PetscErrorCode SNESConverged_Interactive(SNES snes, PetscInt it,PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *ctx)
1329: /* ------------------------------------------------------------------- */
1330: {
1331:   AppCtx         *user  = (AppCtx*) ctx;
1332:   Parameter      *param = user->param;
1333:   KSP            ksp;

1337:   if (param->interrupted) {
1338:     param->interrupted = PETSC_FALSE;
1339:     PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: exiting SNES solve. \n");
1340:     *reason = SNES_CONVERGED_FNORM_ABS;
1341:     return(0);
1342:   } else if (param->toggle_kspmon) {
1343:     param->toggle_kspmon = PETSC_FALSE;

1345:     SNESGetKSP(snes, &ksp);

1347:     if (param->kspmon) {
1348:       KSPMonitorCancel(ksp);

1350:       param->kspmon = PETSC_FALSE;
1351:       PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: deactivating ksp singular value monitor. \n");
1352:     } else {
1353:       PetscViewerAndFormat *vf;
1354:       PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
1355:       KSPMonitorSet(ksp,(PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*))KSPMonitorSingularValue,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);

1357:       param->kspmon = PETSC_TRUE;
1358:       PetscPrintf(PETSC_COMM_WORLD,"USER SIGNAL: activating ksp singular value monitor. \n");
1359:     }
1360:   }
1361:   PetscFunctionReturn(SNESConvergedDefault(snes,it,xnorm,snorm,fnorm,reason,ctx));
1362: }

1364: /* ------------------------------------------------------------------- */
1365: #include <signal.h>
1368: PetscErrorCode InteractiveHandler(int signum, void *ctx)
1369: /* ------------------------------------------------------------------- */
1370: {
1371:   AppCtx    *user  = (AppCtx*) ctx;
1372:   Parameter *param = user->param;

1374:   if (signum == SIGILL) {
1375:     param->toggle_kspmon = PETSC_TRUE;
1376: #if !defined(PETSC_MISSING_SIGCONT)
1377:   } else if (signum == SIGCONT) {
1378:     param->interrupted = PETSC_TRUE;
1379: #endif
1380: #if !defined(PETSC_MISSING_SIGURG)
1381:   } else if (signum == SIGURG) {
1382:     param->stop_solve = PETSC_TRUE;
1383: #endif
1384:   }
1385:   return 0;
1386: }

1388: /*---------------------------------------------------------------------*/
1391: /*  main call-back function that computes the processor-local piece
1392:     of the residual */
1393: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
1394: /*---------------------------------------------------------------------*/
1395: {
1396:   AppCtx      *user  = (AppCtx*)ptr;
1397:   Parameter   *param = user->param;
1398:   GridInfo    *grid  = user->grid;
1399:   PetscScalar mag_w, mag_u;
1400:   PetscInt    i,j,mx,mz,ilim,jlim;
1401:   PetscInt    is,ie,js,je,ibound;    /* ,ivisc */

1404:   /* Define global and local grid parameters */
1405:   mx   = info->mx;     mz   = info->my;
1406:   ilim = mx-1;         jlim = mz-1;
1407:   is   = info->xs;     ie   = info->xs+info->xm;
1408:   js   = info->ys;     je   = info->ys+info->ym;

1410:   /* Define geometric and numeric parameters */
1411:   /* ivisc = param->ivisc; */ ibound = param->ibound;

1413:   for (j=js; j<je; j++) {
1414:     for (i=is; i<ie; i++) {

1416:       /************* X-MOMENTUM/VELOCITY *************/
1417:       if (i<j) f[j][i].u = x[j][i].u - SlabVel('U',i,j,user);
1418:       else if (j<=grid->jlid || (j<grid->corner+grid->inose && i<grid->corner+grid->inose)) {
1419:         /* in the lithospheric lid */
1420:         f[j][i].u = x[j][i].u - 0.0;
1421:       } else if (i==ilim) {
1422:         /* on the right side boundary */
1423:         if (ibound==BC_ANALYTIC) {
1424:           f[j][i].u = x[j][i].u - HorizVelocity(i,j,user);
1425:         } else {
1426:           f[j][i].u = XNormalStress(x,i,j,CELL_CENTER,user) - EPS_ZERO;
1427:         }

1429:       } else if (j==jlim) {
1430:         /* on the bottom boundary */
1431:         if (ibound==BC_ANALYTIC) {
1432:           f[j][i].u = x[j][i].u - HorizVelocity(i,j,user);
1433:         } else if (ibound==BC_NOSTRESS) {
1434:           f[j][i].u = XMomentumResidual(x,i,j,user);
1435:         } else {
1436:           /* experimental boundary condition */
1437:         }

1439:       } else {
1440:         /* in the mantle wedge */
1441:         f[j][i].u = XMomentumResidual(x,i,j,user);
1442:       }

1444:       /************* Z-MOMENTUM/VELOCITY *************/
1445:       if (i<=j) {
1446:         f[j][i].w = x[j][i].w - SlabVel('W',i,j,user);

1448:       } else if (j<=grid->jlid || (j<grid->corner+grid->inose && i<grid->corner+grid->inose)) {
1449:         /* in the lithospheric lid */
1450:         f[j][i].w = x[j][i].w - 0.0;

1452:       } else if (j==jlim) {
1453:         /* on the bottom boundary */
1454:         if (ibound==BC_ANALYTIC) {
1455:           f[j][i].w = x[j][i].w - VertVelocity(i,j,user);
1456:         } else {
1457:           f[j][i].w = ZNormalStress(x,i,j,CELL_CENTER,user) - EPS_ZERO;
1458:         }

1460:       } else if (i==ilim) {
1461:         /* on the right side boundary */
1462:         if (ibound==BC_ANALYTIC) {
1463:           f[j][i].w = x[j][i].w - VertVelocity(i,j,user);
1464:         } else if (ibound==BC_NOSTRESS) {
1465:           f[j][i].w = ZMomentumResidual(x,i,j,user);
1466:         } else {
1467:           /* experimental boundary condition */
1468:         }

1470:       } else {
1471:         /* in the mantle wedge */
1472:         f[j][i].w =  ZMomentumResidual(x,i,j,user);
1473:       }

1475:       /************* CONTINUITY/PRESSURE *************/
1476:       if (i<j || j<=grid->jlid || (j<grid->corner+grid->inose && i<grid->corner+grid->inose)) {
1477:         /* in the lid or slab */
1478:         f[j][i].p = x[j][i].p;

1480:       } else if ((i==ilim || j==jlim) && ibound==BC_ANALYTIC) {
1481:         /* on an analytic boundary */
1482:         f[j][i].p = x[j][i].p - Pressure(i,j,user);

1484:       } else {
1485:         /* in the mantle wedge */
1486:         f[j][i].p = ContinuityResidual(x,i,j,user);
1487:       }

1489:       /************* TEMPERATURE *************/
1490:       if (j==0) {
1491:         /* on the surface */
1492:         f[j][i].T = x[j][i].T + x[j+1][i].T + PetscMax(PetscRealPart(x[j][i].T),0.0);

1494:       } else if (i==0) {
1495:         /* slab inflow boundary */
1496:         f[j][i].T = x[j][i].T - PlateModel(j,PLATE_SLAB,user);

1498:       } else if (i==ilim) {
1499:         /* right side boundary */
1500:         mag_u = 1.0 - PetscPowRealInt((1.0-PetscMax(PetscMin(PetscRealPart(x[j][i-1].u)/param->cb,1.0),0.0)), 5);
1501:         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);

1503:       } else if (j==jlim) {
1504:         /* bottom boundary */
1505:         mag_w = 1.0 - PetscPowRealInt((1.0-PetscMax(PetscMin(PetscRealPart(x[j-1][i].w)/param->sb,1.0),0.0)), 5);
1506:         f[j][i].T = x[j][i].T - mag_w*x[j-1][i-1].T - (1.0-mag_w);

1508:       } else {
1509:         /* in the mantle wedge */
1510:         f[j][i].T = EnergyResidual(x,i,j,user);
1511:       }
1512:     }
1513:   }
1514:   return(0);
1515: }