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