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