1c4762a1bSJed Brown static char help[] = "Surface processes in geophysics.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /*T 4c4762a1bSJed Brown Concepts: SNES^parallel Surface process example 5c4762a1bSJed Brown Concepts: DMDA^using distributed arrays; 6c4762a1bSJed Brown Concepts: IS coloirng types; 7c4762a1bSJed Brown Processors: n 8c4762a1bSJed Brown T*/ 9c4762a1bSJed Brown 10c4762a1bSJed Brown 11c4762a1bSJed Brown 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscsnes.h> 14c4762a1bSJed Brown #include <petscdm.h> 15c4762a1bSJed Brown #include <petscdmda.h> 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* 18c4762a1bSJed Brown User-defined application context - contains data needed by the 19c4762a1bSJed Brown application-provided call-back routines, FormJacobianLocal() and 20c4762a1bSJed Brown FormFunctionLocal(). 21c4762a1bSJed Brown */ 22c4762a1bSJed Brown typedef struct { 23c4762a1bSJed Brown PetscReal D; /* The diffusion coefficient */ 24c4762a1bSJed Brown PetscReal K; /* The advection coefficient */ 25c4762a1bSJed Brown PetscInt m; /* Exponent for A */ 26c4762a1bSJed Brown } AppCtx; 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* 29c4762a1bSJed Brown User-defined routines 30c4762a1bSJed Brown */ 31c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*); 32c4762a1bSJed Brown extern PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscScalar**,Mat,AppCtx*); 33c4762a1bSJed Brown 34c4762a1bSJed Brown int main(int argc,char **argv) 35c4762a1bSJed Brown { 36c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 37c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 38c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 39c4762a1bSJed Brown PetscErrorCode ierr; 40c4762a1bSJed Brown DM da; 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 43c4762a1bSJed Brown Initialize program 44c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 45c4762a1bSJed Brown 46c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 49c4762a1bSJed Brown Initialize problem parameters 50c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 51c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Surface Process Problem Options", "SNES");CHKERRQ(ierr); 52c4762a1bSJed Brown user.D = 1.0; 53c4762a1bSJed Brown ierr = PetscOptionsReal("-D", "The diffusion coefficient D", __FILE__, user.D, &user.D, NULL);CHKERRQ(ierr); 54c4762a1bSJed Brown user.K = 1.0; 55c4762a1bSJed Brown ierr = PetscOptionsReal("-K", "The advection coefficient K", __FILE__, user.K, &user.K, NULL);CHKERRQ(ierr); 56c4762a1bSJed Brown user.m = 1; 57c4762a1bSJed Brown ierr = PetscOptionsInt("-m", "The exponent for A", __FILE__, user.m, &user.m, NULL);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 61c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 62c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 63c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 66c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = DMSetApplicationContext(da,&user);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = SNESSetDM(snes, da);CHKERRQ(ierr); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 72c4762a1bSJed Brown Set local function evaluation routine 73c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 74c4762a1bSJed Brown ierr = DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user);CHKERRQ(ierr); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 77c4762a1bSJed Brown Customize solver; set runtime options 78c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 79c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 80c4762a1bSJed Brown 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 83c4762a1bSJed Brown Solve nonlinear system 84c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 85c4762a1bSJed Brown ierr = SNESSolve(snes,0,0);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 89c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 90c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 93c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 94c4762a1bSJed Brown are no longer needed. 95c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 96c4762a1bSJed Brown 97c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 99c4762a1bSJed Brown 100c4762a1bSJed Brown ierr = PetscFinalize(); 101c4762a1bSJed Brown return ierr; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscScalar funcU(DMDACoor2d *coords) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown return coords->x + coords->y; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown PetscScalar funcA(PetscScalar z, AppCtx *user) 110c4762a1bSJed Brown { 111c4762a1bSJed Brown PetscScalar v = 1.0; 112c4762a1bSJed Brown PetscInt i; 113c4762a1bSJed Brown 114c4762a1bSJed Brown for (i = 0; i < user->m; ++i) v *= z; 115c4762a1bSJed Brown return v; 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown PetscScalar funcADer(PetscScalar z, AppCtx *user) 119c4762a1bSJed Brown { 120c4762a1bSJed Brown PetscScalar v = 1.0; 121c4762a1bSJed Brown PetscInt i; 122c4762a1bSJed Brown 123c4762a1bSJed Brown for (i = 0; i < user->m-1; ++i) v *= z; 124c4762a1bSJed Brown return (PetscScalar)user->m*v; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* 128c4762a1bSJed Brown FormFunctionLocal - Evaluates nonlinear function, F(x). 129c4762a1bSJed Brown */ 130c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user) 131c4762a1bSJed Brown { 132c4762a1bSJed Brown DM coordDA; 133c4762a1bSJed Brown Vec coordinates; 134c4762a1bSJed Brown DMDACoor2d **coords; 135c4762a1bSJed Brown PetscScalar u, ux, uy, uxx, uyy; 136c4762a1bSJed Brown PetscReal D, K, hx, hy, hxdhy, hydhx; 137c4762a1bSJed Brown PetscInt i,j; 138c4762a1bSJed Brown PetscErrorCode ierr; 139c4762a1bSJed Brown 140c4762a1bSJed Brown PetscFunctionBeginUser; 141c4762a1bSJed Brown D = user->D; 142c4762a1bSJed Brown K = user->K; 143c4762a1bSJed Brown hx = 1.0/(PetscReal)(info->mx-1); 144c4762a1bSJed Brown hy = 1.0/(PetscReal)(info->my-1); 145c4762a1bSJed Brown hxdhy = hx/hy; 146c4762a1bSJed Brown hydhx = hy/hx; 147c4762a1bSJed Brown /* 148c4762a1bSJed Brown Compute function over the locally owned part of the grid 149c4762a1bSJed Brown */ 150c4762a1bSJed Brown ierr = DMGetCoordinateDM(info->da, &coordDA);CHKERRQ(ierr); 151c4762a1bSJed Brown ierr = DMGetCoordinates(info->da, &coordinates);CHKERRQ(ierr); 152c4762a1bSJed Brown ierr = DMDAVecGetArray(coordDA, coordinates, &coords);CHKERRQ(ierr); 153c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 154c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 155c4762a1bSJed Brown if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) f[j][i] = x[j][i]; 156c4762a1bSJed Brown else { 157c4762a1bSJed Brown u = x[j][i]; 158c4762a1bSJed Brown ux = (x[j][i+1] - x[j][i])/hx; 159c4762a1bSJed Brown uy = (x[j+1][i] - x[j][i])/hy; 160c4762a1bSJed Brown uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx; 161c4762a1bSJed Brown uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy; 162c4762a1bSJed Brown f[j][i] = D*(uxx + uyy) - (K*funcA(x[j][i], user)*PetscSqrtScalar(ux*ux + uy*uy) + funcU(&coords[j][i]))*hx*hy; 163c4762a1bSJed Brown if (PetscIsInfOrNanScalar(f[j][i])) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FP, "Invalid residual: %g", (double)PetscRealPart(f[j][i])); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 167c4762a1bSJed Brown ierr = DMDAVecRestoreArray(coordDA, coordinates, &coords);CHKERRQ(ierr); 168*ca0c957dSBarry Smith ierr = PetscLogFlops(11.0*info->ym*info->xm);CHKERRQ(ierr); 169c4762a1bSJed Brown PetscFunctionReturn(0); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* 173c4762a1bSJed Brown FormJacobianLocal - Evaluates Jacobian matrix. 174c4762a1bSJed Brown */ 175c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat jac,AppCtx *user) 176c4762a1bSJed Brown { 177c4762a1bSJed Brown MatStencil col[5], row; 178c4762a1bSJed Brown PetscScalar D, K, A, v[5], hx, hy, hxdhy, hydhx, ux, uy; 179c4762a1bSJed Brown PetscReal normGradZ; 180c4762a1bSJed Brown PetscInt i, j,k; 181c4762a1bSJed Brown PetscErrorCode ierr; 182c4762a1bSJed Brown 183c4762a1bSJed Brown PetscFunctionBeginUser; 184c4762a1bSJed Brown D = user->D; 185c4762a1bSJed Brown K = user->K; 186c4762a1bSJed Brown hx = 1.0/(PetscReal)(info->mx-1); 187c4762a1bSJed Brown hy = 1.0/(PetscReal)(info->my-1); 188c4762a1bSJed Brown hxdhy = hx/hy; 189c4762a1bSJed Brown hydhx = hy/hx; 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* 192c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 193c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 194c4762a1bSJed Brown contiguous chunks of rows across the processors. 195c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 196c4762a1bSJed Brown locally (but any non-local elements will be sent to the 197c4762a1bSJed Brown appropriate processor during matrix assembly). 198c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 199c4762a1bSJed Brown - We can set matrix entries either using either 200c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(), as discussed above. 201c4762a1bSJed Brown */ 202c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 203c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 204c4762a1bSJed Brown row.j = j; row.i = i; 205c4762a1bSJed Brown if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) { 206c4762a1bSJed Brown /* boundary points */ 207c4762a1bSJed Brown v[0] = 1.0; 208c4762a1bSJed Brown ierr = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr); 209c4762a1bSJed Brown } else { 210c4762a1bSJed Brown /* interior grid points */ 211c4762a1bSJed Brown ux = (x[j][i+1] - x[j][i])/hx; 212c4762a1bSJed Brown uy = (x[j+1][i] - x[j][i])/hy; 213c4762a1bSJed Brown normGradZ = PetscRealPart(PetscSqrtScalar(ux*ux + uy*uy)); 214c4762a1bSJed Brown if (normGradZ < 1.0e-8) normGradZ = 1.0e-8; 215c4762a1bSJed Brown A = funcA(x[j][i], user); 216c4762a1bSJed Brown 217c4762a1bSJed Brown v[0] = -D*hxdhy; col[0].j = j - 1; col[0].i = i; 218c4762a1bSJed Brown v[1] = -D*hydhx; col[1].j = j; col[1].i = i-1; 219c4762a1bSJed Brown v[2] = D*2.0*(hydhx + hxdhy) + K*(funcADer(x[j][i], user)*normGradZ - A/normGradZ)*hx*hy; col[2].j = row.j; col[2].i = row.i; 220c4762a1bSJed Brown v[3] = -D*hydhx + K*A*hx*hy/(2.0*normGradZ); col[3].j = j; col[3].i = i+1; 221c4762a1bSJed Brown v[4] = -D*hxdhy + K*A*hx*hy/(2.0*normGradZ); col[4].j = j + 1; col[4].i = i; 222c4762a1bSJed Brown for (k = 0; k < 5; ++k) { 223c4762a1bSJed Brown if (PetscIsInfOrNanScalar(v[k])) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FP, "Invalid residual: %g", (double)PetscRealPart(v[k])); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown ierr = MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown } 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* 231c4762a1bSJed Brown Assemble matrix, using the 2-step process: 232c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 233c4762a1bSJed Brown */ 234c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 236c4762a1bSJed Brown /* 237c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 238c4762a1bSJed Brown matrix. If we do, it will generate an error. 239c4762a1bSJed Brown */ 240c4762a1bSJed Brown ierr = MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 241c4762a1bSJed Brown PetscFunctionReturn(0); 242c4762a1bSJed Brown } 243c4762a1bSJed Brown 244c4762a1bSJed Brown 245c4762a1bSJed Brown /*TEST 246c4762a1bSJed Brown 247c4762a1bSJed Brown test: 248c4762a1bSJed Brown args: -snes_view -snes_monitor_short -da_refine 1 -pc_type mg -ksp_type fgmres -pc_mg_type full -mg_levels_ksp_chebyshev_esteig 0.5,1.1 249c4762a1bSJed Brown 250c4762a1bSJed Brown test: 251c4762a1bSJed Brown suffix: ew_1 252c4762a1bSJed Brown args: -snes_monitor_short -ksp_converged_reason -da_grid_x 20 -da_grid_y 20 -snes_ksp_ew -snes_ksp_ew_version 1 253c4762a1bSJed Brown requires: !single 254c4762a1bSJed Brown 255c4762a1bSJed Brown test: 256c4762a1bSJed Brown suffix: ew_2 257c4762a1bSJed Brown args: -snes_monitor_short -ksp_converged_reason -da_grid_x 20 -da_grid_y 20 -snes_ksp_ew -snes_ksp_ew_version 2 258c4762a1bSJed Brown 259c4762a1bSJed Brown test: 260c4762a1bSJed Brown suffix: ew_3 261c4762a1bSJed Brown args: -snes_monitor_short -ksp_converged_reason -da_grid_x 20 -da_grid_y 20 -snes_ksp_ew -snes_ksp_ew_version 3 262c4762a1bSJed Brown requires: !single 263c4762a1bSJed Brown 264c4762a1bSJed Brown test: 265c4762a1bSJed Brown suffix: fm_rise_2 266c4762a1bSJed Brown args: -K 3 -m 1 -D 0.2 -snes_monitor_short -snes_type ngmres -snes_npc_side right -npc_snes_type newtonls -npc_snes_linesearch_type basic -snes_ngmres_restart_it 1 -snes_ngmres_restart_fm_rise 267c4762a1bSJed Brown requires: !single 268c4762a1bSJed Brown 269c4762a1bSJed Brown test: 270c4762a1bSJed Brown suffix: fm_rise_4 271c4762a1bSJed Brown args: -K 3 -m 1 -D 0.2 -snes_monitor_short -snes_type ngmres -snes_npc_side right -npc_snes_type newtonls -npc_snes_linesearch_type basic -snes_ngmres_restart_it 2 -snes_ngmres_restart_fm_rise -snes_rtol 1.e-2 -snes_max_it 5 272c4762a1bSJed Brown 273c4762a1bSJed Brown TEST*/ 274