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*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 25*d71ae5a4SJacob Faibussowitsch { 26c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 27c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 28c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 29c4762a1bSJed Brown DM da; 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 32c4762a1bSJed Brown Initialize program 33c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 34c4762a1bSJed Brown 35327415f7SBarry Smith PetscFunctionBeginUser; 369566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 39c4762a1bSJed Brown Initialize problem parameters 40c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 41d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "Surface Process Problem Options", "SNES"); 42c4762a1bSJed Brown user.D = 1.0; 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-D", "The diffusion coefficient D", __FILE__, user.D, &user.D, NULL)); 44c4762a1bSJed Brown user.K = 1.0; 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-K", "The advection coefficient K", __FILE__, user.K, &user.K, NULL)); 46c4762a1bSJed Brown user.m = 1; 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-m", "The exponent for A", __FILE__, user.m, &user.m, NULL)); 48d0609cedSBarry Smith PetscOptionsEnd(); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 51c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 52c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 539566063dSJacob 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)); 549566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 559566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 569566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 579566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da, &user)); 589566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 599566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, da)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 62c4762a1bSJed Brown Set local function evaluation routine 63c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 649566063dSJacob Faibussowitsch PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, &user)); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 67c4762a1bSJed Brown Customize solver; set runtime options 68c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 699566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 72c4762a1bSJed Brown Solve nonlinear system 73c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 749566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, 0, 0)); 759566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 78c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 7963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 82c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 83c4762a1bSJed Brown are no longer needed. 84c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 85c4762a1bSJed Brown 869566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 879566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 90b122ec5aSJacob Faibussowitsch return 0; 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93*d71ae5a4SJacob Faibussowitsch PetscScalar funcU(DMDACoor2d *coords) 94*d71ae5a4SJacob Faibussowitsch { 95c4762a1bSJed Brown return coords->x + coords->y; 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98*d71ae5a4SJacob Faibussowitsch PetscScalar funcA(PetscScalar z, AppCtx *user) 99*d71ae5a4SJacob Faibussowitsch { 100c4762a1bSJed Brown PetscScalar v = 1.0; 101c4762a1bSJed Brown PetscInt i; 102c4762a1bSJed Brown 103c4762a1bSJed Brown for (i = 0; i < user->m; ++i) v *= z; 104c4762a1bSJed Brown return v; 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107*d71ae5a4SJacob Faibussowitsch PetscScalar funcADer(PetscScalar z, AppCtx *user) 108*d71ae5a4SJacob Faibussowitsch { 109c4762a1bSJed Brown PetscScalar v = 1.0; 110c4762a1bSJed Brown PetscInt i; 111c4762a1bSJed Brown 112c4762a1bSJed Brown for (i = 0; i < user->m - 1; ++i) v *= z; 113c4762a1bSJed Brown return (PetscScalar)user->m * v; 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* 117c4762a1bSJed Brown FormFunctionLocal - Evaluates nonlinear function, F(x). 118c4762a1bSJed Brown */ 119*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user) 120*d71ae5a4SJacob Faibussowitsch { 121c4762a1bSJed Brown DM coordDA; 122c4762a1bSJed Brown Vec coordinates; 123c4762a1bSJed Brown DMDACoor2d **coords; 124c4762a1bSJed Brown PetscScalar u, ux, uy, uxx, uyy; 125c4762a1bSJed Brown PetscReal D, K, hx, hy, hxdhy, hydhx; 126c4762a1bSJed Brown PetscInt i, j; 127c4762a1bSJed Brown 128c4762a1bSJed Brown PetscFunctionBeginUser; 129c4762a1bSJed Brown D = user->D; 130c4762a1bSJed Brown K = user->K; 131c4762a1bSJed Brown hx = 1.0 / (PetscReal)(info->mx - 1); 132c4762a1bSJed Brown hy = 1.0 / (PetscReal)(info->my - 1); 133c4762a1bSJed Brown hxdhy = hx / hy; 134c4762a1bSJed Brown hydhx = hy / hx; 135c4762a1bSJed Brown /* 136c4762a1bSJed Brown Compute function over the locally owned part of the grid 137c4762a1bSJed Brown */ 1389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(info->da, &coordDA)); 1399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(info->da, &coordinates)); 1409566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords)); 141c4762a1bSJed Brown for (j = info->ys; j < info->ys + info->ym; j++) { 142c4762a1bSJed Brown for (i = info->xs; i < info->xs + info->xm; i++) { 143c4762a1bSJed Brown if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) f[j][i] = x[j][i]; 144c4762a1bSJed Brown else { 145c4762a1bSJed Brown u = x[j][i]; 146c4762a1bSJed Brown ux = (x[j][i + 1] - x[j][i]) / hx; 147c4762a1bSJed Brown uy = (x[j + 1][i] - x[j][i]) / hy; 148c4762a1bSJed Brown uxx = (2.0 * u - x[j][i - 1] - x[j][i + 1]) * hydhx; 149c4762a1bSJed Brown uyy = (2.0 * u - x[j - 1][i] - x[j + 1][i]) * hxdhy; 150c4762a1bSJed Brown f[j][i] = D * (uxx + uyy) - (K * funcA(x[j][i], user) * PetscSqrtScalar(ux * ux + uy * uy) + funcU(&coords[j][i])) * hx * hy; 151e00437b9SBarry Smith PetscCheck(!PetscIsInfOrNanScalar(f[j][i]), PETSC_COMM_SELF, PETSC_ERR_FP, "Invalid residual: %g", (double)PetscRealPart(f[j][i])); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown } 154c4762a1bSJed Brown } 1559566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords)); 1569566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(11.0 * info->ym * info->xm)); 157c4762a1bSJed Brown PetscFunctionReturn(0); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* 161c4762a1bSJed Brown FormJacobianLocal - Evaluates Jacobian matrix. 162c4762a1bSJed Brown */ 163*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat jac, AppCtx *user) 164*d71ae5a4SJacob Faibussowitsch { 165c4762a1bSJed Brown MatStencil col[5], row; 166c4762a1bSJed Brown PetscScalar D, K, A, v[5], hx, hy, hxdhy, hydhx, ux, uy; 167c4762a1bSJed Brown PetscReal normGradZ; 168c4762a1bSJed Brown PetscInt i, j, k; 169c4762a1bSJed Brown 170c4762a1bSJed Brown PetscFunctionBeginUser; 171c4762a1bSJed Brown D = user->D; 172c4762a1bSJed Brown K = user->K; 173c4762a1bSJed Brown hx = 1.0 / (PetscReal)(info->mx - 1); 174c4762a1bSJed Brown hy = 1.0 / (PetscReal)(info->my - 1); 175c4762a1bSJed Brown hxdhy = hx / hy; 176c4762a1bSJed Brown hydhx = hy / hx; 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown Compute entries for the locally owned part of the Jacobian. 180c4762a1bSJed Brown - Currently, all PETSc parallel matrix formats are partitioned by 181c4762a1bSJed Brown contiguous chunks of rows across the processors. 182c4762a1bSJed Brown - Each processor needs to insert only elements that it owns 183c4762a1bSJed Brown locally (but any non-local elements will be sent to the 184c4762a1bSJed Brown appropriate processor during matrix assembly). 185c4762a1bSJed Brown - Here, we set all entries for a particular row at once. 186c4762a1bSJed Brown - We can set matrix entries either using either 187c4762a1bSJed Brown MatSetValuesLocal() or MatSetValues(), as discussed above. 188c4762a1bSJed Brown */ 189c4762a1bSJed Brown for (j = info->ys; j < info->ys + info->ym; j++) { 190c4762a1bSJed Brown for (i = info->xs; i < info->xs + info->xm; i++) { 1919371c9d4SSatish Balay row.j = j; 1929371c9d4SSatish Balay row.i = i; 193c4762a1bSJed Brown if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 194c4762a1bSJed Brown /* boundary points */ 195c4762a1bSJed Brown v[0] = 1.0; 1969566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 1, &row, v, INSERT_VALUES)); 197c4762a1bSJed Brown } else { 198c4762a1bSJed Brown /* interior grid points */ 199c4762a1bSJed Brown ux = (x[j][i + 1] - x[j][i]) / hx; 200c4762a1bSJed Brown uy = (x[j + 1][i] - x[j][i]) / hy; 201c4762a1bSJed Brown normGradZ = PetscRealPart(PetscSqrtScalar(ux * ux + uy * uy)); 202c4762a1bSJed Brown if (normGradZ < 1.0e-8) normGradZ = 1.0e-8; 203c4762a1bSJed Brown A = funcA(x[j][i], user); 204c4762a1bSJed Brown 2059371c9d4SSatish Balay v[0] = -D * hxdhy; 2069371c9d4SSatish Balay col[0].j = j - 1; 2079371c9d4SSatish Balay col[0].i = i; 2089371c9d4SSatish Balay v[1] = -D * hydhx; 2099371c9d4SSatish Balay col[1].j = j; 2109371c9d4SSatish Balay col[1].i = i - 1; 2119371c9d4SSatish Balay v[2] = D * 2.0 * (hydhx + hxdhy) + K * (funcADer(x[j][i], user) * normGradZ - A / normGradZ) * hx * hy; 2129371c9d4SSatish Balay col[2].j = row.j; 2139371c9d4SSatish Balay col[2].i = row.i; 2149371c9d4SSatish Balay v[3] = -D * hydhx + K * A * hx * hy / (2.0 * normGradZ); 2159371c9d4SSatish Balay col[3].j = j; 2169371c9d4SSatish Balay col[3].i = i + 1; 2179371c9d4SSatish Balay v[4] = -D * hxdhy + K * A * hx * hy / (2.0 * normGradZ); 2189371c9d4SSatish Balay col[4].j = j + 1; 2199371c9d4SSatish Balay col[4].i = i; 220ad540459SPierre Jolivet for (k = 0; k < 5; ++k) PetscCheck(!PetscIsInfOrNanScalar(v[k]), PETSC_COMM_SELF, PETSC_ERR_FP, "Invalid residual: %g", (double)PetscRealPart(v[k])); 2219566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac, 1, &row, 5, col, v, INSERT_VALUES)); 222c4762a1bSJed Brown } 223c4762a1bSJed Brown } 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* 227c4762a1bSJed Brown Assemble matrix, using the 2-step process: 228c4762a1bSJed Brown MatAssemblyBegin(), MatAssemblyEnd(). 229c4762a1bSJed Brown */ 2309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 232c4762a1bSJed Brown /* 233c4762a1bSJed Brown Tell the matrix we will never add a new nonzero location to the 234c4762a1bSJed Brown matrix. If we do, it will generate an error. 235c4762a1bSJed Brown */ 2369566063dSJacob Faibussowitsch PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 237c4762a1bSJed Brown PetscFunctionReturn(0); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown /*TEST 241c4762a1bSJed Brown 242c4762a1bSJed Brown test: 243c4762a1bSJed 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 244c4762a1bSJed Brown 245c4762a1bSJed Brown test: 246c4762a1bSJed Brown suffix: ew_1 247c4762a1bSJed Brown args: -snes_monitor_short -ksp_converged_reason -da_grid_x 20 -da_grid_y 20 -snes_ksp_ew -snes_ksp_ew_version 1 248c4762a1bSJed Brown requires: !single 249c4762a1bSJed Brown 250c4762a1bSJed Brown test: 251c4762a1bSJed Brown suffix: ew_2 252c4762a1bSJed Brown args: -snes_monitor_short -ksp_converged_reason -da_grid_x 20 -da_grid_y 20 -snes_ksp_ew -snes_ksp_ew_version 2 253c4762a1bSJed Brown 254c4762a1bSJed Brown test: 255c4762a1bSJed Brown suffix: ew_3 256c4762a1bSJed Brown args: -snes_monitor_short -ksp_converged_reason -da_grid_x 20 -da_grid_y 20 -snes_ksp_ew -snes_ksp_ew_version 3 257c4762a1bSJed Brown requires: !single 258c4762a1bSJed Brown 259c4762a1bSJed Brown test: 260c4762a1bSJed Brown suffix: fm_rise_2 261c4762a1bSJed 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 262c4762a1bSJed Brown requires: !single 263c4762a1bSJed Brown 264c4762a1bSJed Brown test: 265c4762a1bSJed Brown suffix: fm_rise_4 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 2 -snes_ngmres_restart_fm_rise -snes_rtol 1.e-2 -snes_max_it 5 267c4762a1bSJed Brown 268c4762a1bSJed Brown TEST*/ 269