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