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