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