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