1c4762a1bSJed Brown static char help[] = "Bratu nonlinear PDE in 2d.\n\ 2c4762a1bSJed Brown We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\ 3c4762a1bSJed Brown domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\ 4c4762a1bSJed Brown The command line options include:\n\ 5c4762a1bSJed Brown -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\ 6c4762a1bSJed Brown problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n\ 7c4762a1bSJed Brown -m_par/n_par <parameter>, where <parameter> indicates an integer\n \ 8c4762a1bSJed Brown that MMS3 will be evaluated with 2^m_par, 2^n_par"; 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* ------------------------------------------------------------------------ 11c4762a1bSJed Brown 12c4762a1bSJed Brown Solid Fuel Ignition (SFI) problem. This problem is modeled by 13c4762a1bSJed Brown the partial differential equation 14c4762a1bSJed Brown 15c4762a1bSJed Brown -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1, 16c4762a1bSJed Brown 17c4762a1bSJed Brown with boundary conditions 18c4762a1bSJed Brown 19c4762a1bSJed Brown u = 0 for x = 0, x = 1, y = 0, y = 1. 20c4762a1bSJed Brown 21c4762a1bSJed Brown A finite difference approximation with the usual 5-point stencil 22c4762a1bSJed Brown is used to discretize the boundary value problem to obtain a nonlinear 23c4762a1bSJed Brown system of equations. 24c4762a1bSJed Brown 25c4762a1bSJed Brown This example shows how geometric multigrid can be run transparently with a nonlinear solver so long 26c4762a1bSJed Brown as SNESSetDM() is provided. Example usage 27c4762a1bSJed Brown 28c4762a1bSJed Brown ./ex5 -pc_type mg -ksp_monitor -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17 29c4762a1bSJed Brown -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full 30c4762a1bSJed Brown 31c4762a1bSJed Brown or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of 32c4762a1bSJed Brown multigrid levels, it will be determined automatically based on the number of refinements done) 33c4762a1bSJed Brown 34c4762a1bSJed Brown ./ex5 -pc_type mg -ksp_monitor -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3 35c4762a1bSJed Brown -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full 36c4762a1bSJed Brown 37c4762a1bSJed Brown ------------------------------------------------------------------------- */ 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* 40c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 41c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 42c4762a1bSJed Brown */ 43c4762a1bSJed Brown #include <petscdm.h> 44c4762a1bSJed Brown #include <petscdmda.h> 45c4762a1bSJed Brown #include <petscsnes.h> 46c4762a1bSJed Brown #include <petscmatlab.h> 47c4762a1bSJed Brown #include <petsc/private/snesimpl.h> /* For SNES_Solve event */ 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown User-defined application context - contains data needed by the 51c4762a1bSJed Brown application-provided call-back routines, FormJacobianLocal() and 52c4762a1bSJed Brown FormFunctionLocal(). 53c4762a1bSJed Brown */ 54c4762a1bSJed Brown typedef struct AppCtx AppCtx; 55c4762a1bSJed Brown struct AppCtx { 56c4762a1bSJed Brown PetscReal param; /* test problem parameter */ 57c4762a1bSJed Brown PetscInt m, n; /* MMS3 parameters */ 58c4762a1bSJed Brown PetscErrorCode (*mms_solution)(AppCtx *, const DMDACoor2d *, PetscScalar *); 59c4762a1bSJed Brown PetscErrorCode (*mms_forcing)(AppCtx *, const DMDACoor2d *, PetscScalar *); 60c4762a1bSJed Brown }; 61c4762a1bSJed Brown 62227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */ 63c4762a1bSJed Brown /* 64227f9e03SJunchao Zhang FormInitialGuess - Forms initial approximation. 65227f9e03SJunchao Zhang 66227f9e03SJunchao Zhang Input Parameters: 67227f9e03SJunchao Zhang da - The DM 68227f9e03SJunchao Zhang user - user-defined application context 69227f9e03SJunchao Zhang 70227f9e03SJunchao Zhang Output Parameter: 71227f9e03SJunchao Zhang X - vector 72c4762a1bSJed Brown */ 739371c9d4SSatish Balay static PetscErrorCode FormInitialGuess(DM da, AppCtx *user, Vec X) { 74227f9e03SJunchao Zhang PetscInt i, j, Mx, My, xs, ys, xm, ym; 75227f9e03SJunchao Zhang PetscReal lambda, temp1, temp, hx, hy; 76227f9e03SJunchao Zhang PetscScalar **x; 77227f9e03SJunchao Zhang 78227f9e03SJunchao Zhang PetscFunctionBeginUser; 79227f9e03SJunchao Zhang PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 80227f9e03SJunchao Zhang 81227f9e03SJunchao Zhang lambda = user->param; 82227f9e03SJunchao Zhang hx = 1.0 / (PetscReal)(Mx - 1); 83227f9e03SJunchao Zhang hy = 1.0 / (PetscReal)(My - 1); 84227f9e03SJunchao Zhang temp1 = lambda / (lambda + 1.0); 85227f9e03SJunchao Zhang 86227f9e03SJunchao Zhang /* 87227f9e03SJunchao Zhang Get a pointer to vector data. 88227f9e03SJunchao Zhang - For default PETSc vectors, VecGetArray() returns a pointer to 89227f9e03SJunchao Zhang the data array. Otherwise, the routine is implementation dependent. 90227f9e03SJunchao Zhang - You MUST call VecRestoreArray() when you no longer need access to 91227f9e03SJunchao Zhang the array. 92227f9e03SJunchao Zhang */ 93227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(da, X, &x)); 94227f9e03SJunchao Zhang 95227f9e03SJunchao Zhang /* 96227f9e03SJunchao Zhang Get local grid boundaries (for 2-dimensional DMDA): 97227f9e03SJunchao Zhang xs, ys - starting grid indices (no ghost points) 98227f9e03SJunchao Zhang xm, ym - widths of local grid (no ghost points) 99227f9e03SJunchao Zhang 100227f9e03SJunchao Zhang */ 101227f9e03SJunchao Zhang PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 102227f9e03SJunchao Zhang 103227f9e03SJunchao Zhang /* 104227f9e03SJunchao Zhang Compute initial guess over the locally owned part of the grid 105227f9e03SJunchao Zhang */ 106227f9e03SJunchao Zhang for (j = ys; j < ys + ym; j++) { 107227f9e03SJunchao Zhang temp = (PetscReal)(PetscMin(j, My - j - 1)) * hy; 108227f9e03SJunchao Zhang for (i = xs; i < xs + xm; i++) { 109227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 110227f9e03SJunchao Zhang /* boundary conditions are all zero Dirichlet */ 111227f9e03SJunchao Zhang x[j][i] = 0.0; 112227f9e03SJunchao Zhang } else { 113227f9e03SJunchao Zhang x[j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp)); 114227f9e03SJunchao Zhang } 115227f9e03SJunchao Zhang } 116227f9e03SJunchao Zhang } 117227f9e03SJunchao Zhang 118227f9e03SJunchao Zhang /* 119227f9e03SJunchao Zhang Restore vector 120227f9e03SJunchao Zhang */ 121227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(da, X, &x)); 122227f9e03SJunchao Zhang PetscFunctionReturn(0); 123227f9e03SJunchao Zhang } 124227f9e03SJunchao Zhang 125227f9e03SJunchao Zhang /* 126227f9e03SJunchao Zhang FormExactSolution - Forms MMS solution 127227f9e03SJunchao Zhang 128227f9e03SJunchao Zhang Input Parameters: 129227f9e03SJunchao Zhang da - The DM 130227f9e03SJunchao Zhang user - user-defined application context 131227f9e03SJunchao Zhang 132227f9e03SJunchao Zhang Output Parameter: 133227f9e03SJunchao Zhang X - vector 134227f9e03SJunchao Zhang */ 1359371c9d4SSatish Balay static PetscErrorCode FormExactSolution(DM da, AppCtx *user, Vec U) { 136227f9e03SJunchao Zhang DM coordDA; 137227f9e03SJunchao Zhang Vec coordinates; 138227f9e03SJunchao Zhang DMDACoor2d **coords; 139227f9e03SJunchao Zhang PetscScalar **u; 140227f9e03SJunchao Zhang PetscInt xs, ys, xm, ym, i, j; 141227f9e03SJunchao Zhang 142227f9e03SJunchao Zhang PetscFunctionBeginUser; 143227f9e03SJunchao Zhang PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 144227f9e03SJunchao Zhang PetscCall(DMGetCoordinateDM(da, &coordDA)); 145227f9e03SJunchao Zhang PetscCall(DMGetCoordinates(da, &coordinates)); 146227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords)); 147227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(da, U, &u)); 148227f9e03SJunchao Zhang for (j = ys; j < ys + ym; ++j) { 1499371c9d4SSatish Balay for (i = xs; i < xs + xm; ++i) { user->mms_solution(user, &coords[j][i], &u[j][i]); } 150227f9e03SJunchao Zhang } 151227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(da, U, &u)); 152227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords)); 153227f9e03SJunchao Zhang PetscFunctionReturn(0); 154227f9e03SJunchao Zhang } 155227f9e03SJunchao Zhang 1569371c9d4SSatish Balay static PetscErrorCode ZeroBCSolution(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) { 157227f9e03SJunchao Zhang u[0] = 0.; 158227f9e03SJunchao Zhang return 0; 159227f9e03SJunchao Zhang } 160227f9e03SJunchao Zhang 161227f9e03SJunchao Zhang /* The functions below evaluate the MMS solution u(x,y) and associated forcing 162227f9e03SJunchao Zhang 163227f9e03SJunchao Zhang f(x,y) = -u_xx - u_yy - lambda exp(u) 164227f9e03SJunchao Zhang 165227f9e03SJunchao Zhang such that u(x,y) is an exact solution with f(x,y) as the right hand side forcing term. 166227f9e03SJunchao Zhang */ 1679371c9d4SSatish Balay static PetscErrorCode MMSSolution1(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) { 168227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 169227f9e03SJunchao Zhang u[0] = x * (1 - x) * y * (1 - y); 170227f9e03SJunchao Zhang PetscLogFlops(5); 171227f9e03SJunchao Zhang return 0; 172227f9e03SJunchao Zhang } 1739371c9d4SSatish Balay static PetscErrorCode MMSForcing1(AppCtx *user, const DMDACoor2d *c, PetscScalar *f) { 174227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 175227f9e03SJunchao Zhang f[0] = 2 * x * (1 - x) + 2 * y * (1 - y) - user->param * PetscExpReal(x * (1 - x) * y * (1 - y)); 176227f9e03SJunchao Zhang return 0; 177227f9e03SJunchao Zhang } 178227f9e03SJunchao Zhang 1799371c9d4SSatish Balay static PetscErrorCode MMSSolution2(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) { 180227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 181227f9e03SJunchao Zhang u[0] = PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y); 182227f9e03SJunchao Zhang PetscLogFlops(5); 183227f9e03SJunchao Zhang return 0; 184227f9e03SJunchao Zhang } 1859371c9d4SSatish Balay static PetscErrorCode MMSForcing2(AppCtx *user, const DMDACoor2d *c, PetscScalar *f) { 186227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 187227f9e03SJunchao Zhang f[0] = 2 * PetscSqr(PETSC_PI) * PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y) - user->param * PetscExpReal(PetscSinReal(PETSC_PI * x) * PetscSinReal(PETSC_PI * y)); 188227f9e03SJunchao Zhang return 0; 189227f9e03SJunchao Zhang } 190227f9e03SJunchao Zhang 1919371c9d4SSatish Balay static PetscErrorCode MMSSolution3(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) { 192227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 193227f9e03SJunchao Zhang u[0] = PetscSinReal(user->m * PETSC_PI * x * (1 - y)) * PetscSinReal(user->n * PETSC_PI * y * (1 - x)); 194227f9e03SJunchao Zhang PetscLogFlops(5); 195227f9e03SJunchao Zhang return 0; 196227f9e03SJunchao Zhang } 1979371c9d4SSatish Balay static PetscErrorCode MMSForcing3(AppCtx *user, const DMDACoor2d *c, PetscScalar *f) { 198227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 199227f9e03SJunchao Zhang PetscReal m = user->m, n = user->n, lambda = user->param; 2009371c9d4SSatish Balay f[0] = (-(PetscExpReal(PetscSinReal(m * PETSC_PI * x * (1 - y)) * PetscSinReal(n * PETSC_PI * (1 - x) * y)) * lambda) + PetscSqr(PETSC_PI) * (-2 * m * n * ((-1 + x) * x + (-1 + y) * y) * PetscCosReal(m * PETSC_PI * x * (-1 + y)) * PetscCosReal(n * PETSC_PI * (-1 + x) * y) + (PetscSqr(m) * (PetscSqr(x) + PetscSqr(-1 + y)) + PetscSqr(n) * (PetscSqr(-1 + x) + PetscSqr(y))) * PetscSinReal(m * PETSC_PI * x * (-1 + y)) * PetscSinReal(n * PETSC_PI * (-1 + x) * y))); 201227f9e03SJunchao Zhang return 0; 202227f9e03SJunchao Zhang } 203227f9e03SJunchao Zhang 2049371c9d4SSatish Balay static PetscErrorCode MMSSolution4(AppCtx *user, const DMDACoor2d *c, PetscScalar *u) { 205227f9e03SJunchao Zhang const PetscReal Lx = 1., Ly = 1.; 206227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 207227f9e03SJunchao Zhang u[0] = (PetscPowReal(x, 4) - PetscSqr(Lx) * PetscSqr(x)) * (PetscPowReal(y, 4) - PetscSqr(Ly) * PetscSqr(y)); 208227f9e03SJunchao Zhang PetscLogFlops(9); 209227f9e03SJunchao Zhang return 0; 210227f9e03SJunchao Zhang } 2119371c9d4SSatish Balay static PetscErrorCode MMSForcing4(AppCtx *user, const DMDACoor2d *c, PetscScalar *f) { 212227f9e03SJunchao Zhang const PetscReal Lx = 1., Ly = 1.; 213227f9e03SJunchao Zhang PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y); 2149371c9d4SSatish Balay f[0] = (2 * PetscSqr(x) * (PetscSqr(x) - PetscSqr(Lx)) * (PetscSqr(Ly) - 6 * PetscSqr(y)) + 2 * PetscSqr(y) * (PetscSqr(Lx) - 6 * PetscSqr(x)) * (PetscSqr(y) - PetscSqr(Ly)) - user->param * PetscExpReal((PetscPowReal(x, 4) - PetscSqr(Lx) * PetscSqr(x)) * (PetscPowReal(y, 4) - PetscSqr(Ly) * PetscSqr(y)))); 215227f9e03SJunchao Zhang return 0; 216227f9e03SJunchao Zhang } 217227f9e03SJunchao Zhang 218227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */ 219227f9e03SJunchao Zhang /* 220227f9e03SJunchao Zhang FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch 221227f9e03SJunchao Zhang 222227f9e03SJunchao Zhang */ 2239371c9d4SSatish Balay static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user) { 224227f9e03SJunchao Zhang PetscInt i, j; 225227f9e03SJunchao Zhang PetscReal lambda, hx, hy, hxdhy, hydhx; 226227f9e03SJunchao Zhang PetscScalar u, ue, uw, un, us, uxx, uyy, mms_solution, mms_forcing; 227227f9e03SJunchao Zhang DMDACoor2d c; 228227f9e03SJunchao Zhang 229227f9e03SJunchao Zhang PetscFunctionBeginUser; 230227f9e03SJunchao Zhang lambda = user->param; 231227f9e03SJunchao Zhang hx = 1.0 / (PetscReal)(info->mx - 1); 232227f9e03SJunchao Zhang hy = 1.0 / (PetscReal)(info->my - 1); 233227f9e03SJunchao Zhang hxdhy = hx / hy; 234227f9e03SJunchao Zhang hydhx = hy / hx; 235227f9e03SJunchao Zhang /* 236227f9e03SJunchao Zhang Compute function over the locally owned part of the grid 237227f9e03SJunchao Zhang */ 238227f9e03SJunchao Zhang for (j = info->ys; j < info->ys + info->ym; j++) { 239227f9e03SJunchao Zhang for (i = info->xs; i < info->xs + info->xm; i++) { 240227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 2419371c9d4SSatish Balay c.x = i * hx; 2429371c9d4SSatish Balay c.y = j * hy; 243227f9e03SJunchao Zhang PetscCall(user->mms_solution(user, &c, &mms_solution)); 244227f9e03SJunchao Zhang f[j][i] = 2.0 * (hydhx + hxdhy) * (x[j][i] - mms_solution); 245227f9e03SJunchao Zhang } else { 246227f9e03SJunchao Zhang u = x[j][i]; 247227f9e03SJunchao Zhang uw = x[j][i - 1]; 248227f9e03SJunchao Zhang ue = x[j][i + 1]; 249227f9e03SJunchao Zhang un = x[j - 1][i]; 250227f9e03SJunchao Zhang us = x[j + 1][i]; 251227f9e03SJunchao Zhang 252227f9e03SJunchao Zhang /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */ 2539371c9d4SSatish Balay if (i - 1 == 0) { 2549371c9d4SSatish Balay c.x = (i - 1) * hx; 2559371c9d4SSatish Balay c.y = j * hy; 2569371c9d4SSatish Balay PetscCall(user->mms_solution(user, &c, &uw)); 2579371c9d4SSatish Balay } 2589371c9d4SSatish Balay if (i + 1 == info->mx - 1) { 2599371c9d4SSatish Balay c.x = (i + 1) * hx; 2609371c9d4SSatish Balay c.y = j * hy; 2619371c9d4SSatish Balay PetscCall(user->mms_solution(user, &c, &ue)); 2629371c9d4SSatish Balay } 2639371c9d4SSatish Balay if (j - 1 == 0) { 2649371c9d4SSatish Balay c.x = i * hx; 2659371c9d4SSatish Balay c.y = (j - 1) * hy; 2669371c9d4SSatish Balay PetscCall(user->mms_solution(user, &c, &un)); 2679371c9d4SSatish Balay } 2689371c9d4SSatish Balay if (j + 1 == info->my - 1) { 2699371c9d4SSatish Balay c.x = i * hx; 2709371c9d4SSatish Balay c.y = (j + 1) * hy; 2719371c9d4SSatish Balay PetscCall(user->mms_solution(user, &c, &us)); 2729371c9d4SSatish Balay } 273227f9e03SJunchao Zhang 274227f9e03SJunchao Zhang uxx = (2.0 * u - uw - ue) * hydhx; 275227f9e03SJunchao Zhang uyy = (2.0 * u - un - us) * hxdhy; 276227f9e03SJunchao Zhang mms_forcing = 0; 2779371c9d4SSatish Balay c.x = i * hx; 2789371c9d4SSatish Balay c.y = j * hy; 279227f9e03SJunchao Zhang if (user->mms_forcing) PetscCall(user->mms_forcing(user, &c, &mms_forcing)); 280227f9e03SJunchao Zhang f[j][i] = uxx + uyy - hx * hy * (lambda * PetscExpScalar(u) + mms_forcing); 281227f9e03SJunchao Zhang } 282227f9e03SJunchao Zhang } 283227f9e03SJunchao Zhang } 284227f9e03SJunchao Zhang PetscCall(PetscLogFlops(11.0 * info->ym * info->xm)); 285227f9e03SJunchao Zhang PetscFunctionReturn(0); 286227f9e03SJunchao Zhang } 287227f9e03SJunchao Zhang 288227f9e03SJunchao Zhang /* FormObjectiveLocal - Evaluates nonlinear function, F(x) on local process patch */ 2899371c9d4SSatish Balay static PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info, PetscScalar **x, PetscReal *obj, AppCtx *user) { 290227f9e03SJunchao Zhang PetscInt i, j; 291227f9e03SJunchao Zhang PetscReal lambda, hx, hy, hxdhy, hydhx, sc, lobj = 0; 292227f9e03SJunchao Zhang PetscScalar u, ue, uw, un, us, uxux, uyuy; 293227f9e03SJunchao Zhang MPI_Comm comm; 294227f9e03SJunchao Zhang 295227f9e03SJunchao Zhang PetscFunctionBeginUser; 296227f9e03SJunchao Zhang *obj = 0; 297227f9e03SJunchao Zhang PetscCall(PetscObjectGetComm((PetscObject)info->da, &comm)); 298227f9e03SJunchao Zhang lambda = user->param; 299227f9e03SJunchao Zhang hx = 1.0 / (PetscReal)(info->mx - 1); 300227f9e03SJunchao Zhang hy = 1.0 / (PetscReal)(info->my - 1); 301227f9e03SJunchao Zhang sc = hx * hy * lambda; 302227f9e03SJunchao Zhang hxdhy = hx / hy; 303227f9e03SJunchao Zhang hydhx = hy / hx; 304227f9e03SJunchao Zhang /* 305227f9e03SJunchao Zhang Compute function over the locally owned part of the grid 306227f9e03SJunchao Zhang */ 307227f9e03SJunchao Zhang for (j = info->ys; j < info->ys + info->ym; j++) { 308227f9e03SJunchao Zhang for (i = info->xs; i < info->xs + info->xm; i++) { 309227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 310227f9e03SJunchao Zhang lobj += PetscRealPart((hydhx + hxdhy) * x[j][i] * x[j][i]); 311227f9e03SJunchao Zhang } else { 312227f9e03SJunchao Zhang u = x[j][i]; 313227f9e03SJunchao Zhang uw = x[j][i - 1]; 314227f9e03SJunchao Zhang ue = x[j][i + 1]; 315227f9e03SJunchao Zhang un = x[j - 1][i]; 316227f9e03SJunchao Zhang us = x[j + 1][i]; 317227f9e03SJunchao Zhang 318227f9e03SJunchao Zhang if (i - 1 == 0) uw = 0.; 319227f9e03SJunchao Zhang if (i + 1 == info->mx - 1) ue = 0.; 320227f9e03SJunchao Zhang if (j - 1 == 0) un = 0.; 321227f9e03SJunchao Zhang if (j + 1 == info->my - 1) us = 0.; 322227f9e03SJunchao Zhang 323227f9e03SJunchao Zhang /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */ 324227f9e03SJunchao Zhang 325227f9e03SJunchao Zhang uxux = u * (2. * u - ue - uw) * hydhx; 326227f9e03SJunchao Zhang uyuy = u * (2. * u - un - us) * hxdhy; 327227f9e03SJunchao Zhang 328227f9e03SJunchao Zhang lobj += PetscRealPart(0.5 * (uxux + uyuy) - sc * PetscExpScalar(u)); 329227f9e03SJunchao Zhang } 330227f9e03SJunchao Zhang } 331227f9e03SJunchao Zhang } 332227f9e03SJunchao Zhang PetscCall(PetscLogFlops(12.0 * info->ym * info->xm)); 333227f9e03SJunchao Zhang PetscCallMPI(MPI_Allreduce(&lobj, obj, 1, MPIU_REAL, MPIU_SUM, comm)); 334227f9e03SJunchao Zhang PetscFunctionReturn(0); 335227f9e03SJunchao Zhang } 336227f9e03SJunchao Zhang 337227f9e03SJunchao Zhang /* 338227f9e03SJunchao Zhang FormJacobianLocal - Evaluates Jacobian matrix on local process patch 339227f9e03SJunchao Zhang */ 3409371c9d4SSatish Balay static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat jac, Mat jacpre, AppCtx *user) { 341227f9e03SJunchao Zhang PetscInt i, j, k; 342227f9e03SJunchao Zhang MatStencil col[5], row; 343227f9e03SJunchao Zhang PetscScalar lambda, v[5], hx, hy, hxdhy, hydhx, sc; 344227f9e03SJunchao Zhang DM coordDA; 345227f9e03SJunchao Zhang Vec coordinates; 346227f9e03SJunchao Zhang DMDACoor2d **coords; 347227f9e03SJunchao Zhang 348227f9e03SJunchao Zhang PetscFunctionBeginUser; 349227f9e03SJunchao Zhang lambda = user->param; 350227f9e03SJunchao Zhang /* Extract coordinates */ 351227f9e03SJunchao Zhang PetscCall(DMGetCoordinateDM(info->da, &coordDA)); 352227f9e03SJunchao Zhang PetscCall(DMGetCoordinates(info->da, &coordinates)); 353227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords)); 354227f9e03SJunchao Zhang hx = info->xm > 1 ? PetscRealPart(coords[info->ys][info->xs + 1].x) - PetscRealPart(coords[info->ys][info->xs].x) : 1.0; 355227f9e03SJunchao Zhang hy = info->ym > 1 ? PetscRealPart(coords[info->ys + 1][info->xs].y) - PetscRealPart(coords[info->ys][info->xs].y) : 1.0; 356227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords)); 357227f9e03SJunchao Zhang hxdhy = hx / hy; 358227f9e03SJunchao Zhang hydhx = hy / hx; 359227f9e03SJunchao Zhang sc = hx * hy * lambda; 360227f9e03SJunchao Zhang 361227f9e03SJunchao Zhang /* 362227f9e03SJunchao Zhang Compute entries for the locally owned part of the Jacobian. 363227f9e03SJunchao Zhang - Currently, all PETSc parallel matrix formats are partitioned by 364227f9e03SJunchao Zhang contiguous chunks of rows across the processors. 365227f9e03SJunchao Zhang - Each processor needs to insert only elements that it owns 366227f9e03SJunchao Zhang locally (but any non-local elements will be sent to the 367227f9e03SJunchao Zhang appropriate processor during matrix assembly). 368227f9e03SJunchao Zhang - Here, we set all entries for a particular row at once. 369227f9e03SJunchao Zhang - We can set matrix entries either using either 370227f9e03SJunchao Zhang MatSetValuesLocal() or MatSetValues(), as discussed above. 371227f9e03SJunchao Zhang */ 372227f9e03SJunchao Zhang for (j = info->ys; j < info->ys + info->ym; j++) { 373227f9e03SJunchao Zhang for (i = info->xs; i < info->xs + info->xm; i++) { 3749371c9d4SSatish Balay row.j = j; 3759371c9d4SSatish Balay row.i = i; 376227f9e03SJunchao Zhang /* boundary points */ 377227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 378227f9e03SJunchao Zhang v[0] = 2.0 * (hydhx + hxdhy); 379227f9e03SJunchao Zhang PetscCall(MatSetValuesStencil(jacpre, 1, &row, 1, &row, v, INSERT_VALUES)); 380227f9e03SJunchao Zhang } else { 381227f9e03SJunchao Zhang k = 0; 382227f9e03SJunchao Zhang /* interior grid points */ 383227f9e03SJunchao Zhang if (j - 1 != 0) { 384227f9e03SJunchao Zhang v[k] = -hxdhy; 3859371c9d4SSatish Balay col[k].j = j - 1; 3869371c9d4SSatish Balay col[k].i = i; 387227f9e03SJunchao Zhang k++; 388227f9e03SJunchao Zhang } 389227f9e03SJunchao Zhang if (i - 1 != 0) { 390227f9e03SJunchao Zhang v[k] = -hydhx; 3919371c9d4SSatish Balay col[k].j = j; 3929371c9d4SSatish Balay col[k].i = i - 1; 393227f9e03SJunchao Zhang k++; 394227f9e03SJunchao Zhang } 395227f9e03SJunchao Zhang 3969371c9d4SSatish Balay v[k] = 2.0 * (hydhx + hxdhy) - sc * PetscExpScalar(x[j][i]); 3979371c9d4SSatish Balay col[k].j = row.j; 3989371c9d4SSatish Balay col[k].i = row.i; 3999371c9d4SSatish Balay k++; 400227f9e03SJunchao Zhang 401227f9e03SJunchao Zhang if (i + 1 != info->mx - 1) { 402227f9e03SJunchao Zhang v[k] = -hydhx; 4039371c9d4SSatish Balay col[k].j = j; 4049371c9d4SSatish Balay col[k].i = i + 1; 405227f9e03SJunchao Zhang k++; 406227f9e03SJunchao Zhang } 407227f9e03SJunchao Zhang if (j + 1 != info->mx - 1) { 408227f9e03SJunchao Zhang v[k] = -hxdhy; 4099371c9d4SSatish Balay col[k].j = j + 1; 4109371c9d4SSatish Balay col[k].i = i; 411227f9e03SJunchao Zhang k++; 412227f9e03SJunchao Zhang } 413227f9e03SJunchao Zhang PetscCall(MatSetValuesStencil(jacpre, 1, &row, k, col, v, INSERT_VALUES)); 414227f9e03SJunchao Zhang } 415227f9e03SJunchao Zhang } 416227f9e03SJunchao Zhang } 417227f9e03SJunchao Zhang 418227f9e03SJunchao Zhang /* 419227f9e03SJunchao Zhang Assemble matrix, using the 2-step process: 420227f9e03SJunchao Zhang MatAssemblyBegin(), MatAssemblyEnd(). 421227f9e03SJunchao Zhang */ 422227f9e03SJunchao Zhang PetscCall(MatAssemblyBegin(jacpre, MAT_FINAL_ASSEMBLY)); 423227f9e03SJunchao Zhang PetscCall(MatAssemblyEnd(jacpre, MAT_FINAL_ASSEMBLY)); 424227f9e03SJunchao Zhang /* 425227f9e03SJunchao Zhang Tell the matrix we will never add a new nonzero location to the 426227f9e03SJunchao Zhang matrix. If we do, it will generate an error. 427227f9e03SJunchao Zhang */ 428227f9e03SJunchao Zhang PetscCall(MatSetOption(jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 429227f9e03SJunchao Zhang PetscFunctionReturn(0); 430227f9e03SJunchao Zhang } 431227f9e03SJunchao Zhang 4329371c9d4SSatish Balay static PetscErrorCode FormFunctionMatlab(SNES snes, Vec X, Vec F, void *ptr) { 433227f9e03SJunchao Zhang #if PetscDefined(HAVE_MATLAB_ENGINE) 434227f9e03SJunchao Zhang AppCtx *user = (AppCtx *)ptr; 435227f9e03SJunchao Zhang PetscInt Mx, My; 436227f9e03SJunchao Zhang PetscReal lambda, hx, hy; 437227f9e03SJunchao Zhang Vec localX, localF; 438227f9e03SJunchao Zhang MPI_Comm comm; 439227f9e03SJunchao Zhang DM da; 440227f9e03SJunchao Zhang 441227f9e03SJunchao Zhang PetscFunctionBeginUser; 442227f9e03SJunchao Zhang PetscCall(SNESGetDM(snes, &da)); 443227f9e03SJunchao Zhang PetscCall(DMGetLocalVector(da, &localX)); 444227f9e03SJunchao Zhang PetscCall(DMGetLocalVector(da, &localF)); 445227f9e03SJunchao Zhang PetscCall(PetscObjectSetName((PetscObject)localX, "localX")); 446227f9e03SJunchao Zhang PetscCall(PetscObjectSetName((PetscObject)localF, "localF")); 447227f9e03SJunchao Zhang PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 448227f9e03SJunchao Zhang 449227f9e03SJunchao Zhang lambda = user->param; 450227f9e03SJunchao Zhang hx = 1.0 / (PetscReal)(Mx - 1); 451227f9e03SJunchao Zhang hy = 1.0 / (PetscReal)(My - 1); 452227f9e03SJunchao Zhang 453227f9e03SJunchao Zhang PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 454227f9e03SJunchao Zhang /* 455227f9e03SJunchao Zhang Scatter ghost points to local vector,using the 2-step process 456227f9e03SJunchao Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 457227f9e03SJunchao Zhang By placing code between these two statements, computations can be 458227f9e03SJunchao Zhang done while messages are in transition. 459227f9e03SJunchao Zhang */ 460227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 461227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 462227f9e03SJunchao Zhang PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm), (PetscObject)localX)); 463227f9e03SJunchao Zhang PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm), "localF=ex5m(localX,%18.16e,%18.16e,%18.16e)", (double)hx, (double)hy, (double)lambda)); 464227f9e03SJunchao Zhang PetscCall(PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm), (PetscObject)localF)); 465227f9e03SJunchao Zhang 466227f9e03SJunchao Zhang /* 467227f9e03SJunchao Zhang Insert values into global vector 468227f9e03SJunchao Zhang */ 469227f9e03SJunchao Zhang PetscCall(DMLocalToGlobalBegin(da, localF, INSERT_VALUES, F)); 470227f9e03SJunchao Zhang PetscCall(DMLocalToGlobalEnd(da, localF, INSERT_VALUES, F)); 471227f9e03SJunchao Zhang PetscCall(DMRestoreLocalVector(da, &localX)); 472227f9e03SJunchao Zhang PetscCall(DMRestoreLocalVector(da, &localF)); 473227f9e03SJunchao Zhang PetscFunctionReturn(0); 474227f9e03SJunchao Zhang #else 475227f9e03SJunchao Zhang return 0; /* Never called */ 476227f9e03SJunchao Zhang #endif 477227f9e03SJunchao Zhang } 478227f9e03SJunchao Zhang 479227f9e03SJunchao Zhang /* ------------------------------------------------------------------- */ 480227f9e03SJunchao Zhang /* 481227f9e03SJunchao Zhang Applies some sweeps on nonlinear Gauss-Seidel on each process 482227f9e03SJunchao Zhang 483227f9e03SJunchao Zhang */ 4849371c9d4SSatish Balay static PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx) { 485227f9e03SJunchao Zhang PetscInt i, j, k, Mx, My, xs, ys, xm, ym, its, tot_its, sweeps, l; 486227f9e03SJunchao Zhang PetscReal lambda, hx, hy, hxdhy, hydhx, sc; 487227f9e03SJunchao Zhang PetscScalar **x, **b, bij, F, F0 = 0, J, u, un, us, ue, eu, uw, uxx, uyy, y; 488227f9e03SJunchao Zhang PetscReal atol, rtol, stol; 489227f9e03SJunchao Zhang DM da; 490227f9e03SJunchao Zhang AppCtx *user; 491227f9e03SJunchao Zhang Vec localX, localB; 492227f9e03SJunchao Zhang 493227f9e03SJunchao Zhang PetscFunctionBeginUser; 494227f9e03SJunchao Zhang tot_its = 0; 495227f9e03SJunchao Zhang PetscCall(SNESNGSGetSweeps(snes, &sweeps)); 496227f9e03SJunchao Zhang PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its)); 497227f9e03SJunchao Zhang PetscCall(SNESGetDM(snes, &da)); 498227f9e03SJunchao Zhang PetscCall(DMGetApplicationContext(da, &user)); 499227f9e03SJunchao Zhang 500227f9e03SJunchao Zhang PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 501227f9e03SJunchao Zhang 502227f9e03SJunchao Zhang lambda = user->param; 503227f9e03SJunchao Zhang hx = 1.0 / (PetscReal)(Mx - 1); 504227f9e03SJunchao Zhang hy = 1.0 / (PetscReal)(My - 1); 505227f9e03SJunchao Zhang sc = hx * hy * lambda; 506227f9e03SJunchao Zhang hxdhy = hx / hy; 507227f9e03SJunchao Zhang hydhx = hy / hx; 508227f9e03SJunchao Zhang 509227f9e03SJunchao Zhang PetscCall(DMGetLocalVector(da, &localX)); 510*48a46eb9SPierre Jolivet if (B) PetscCall(DMGetLocalVector(da, &localB)); 511227f9e03SJunchao Zhang for (l = 0; l < sweeps; l++) { 512227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 513227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 514227f9e03SJunchao Zhang if (B) { 515227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalBegin(da, B, INSERT_VALUES, localB)); 516227f9e03SJunchao Zhang PetscCall(DMGlobalToLocalEnd(da, B, INSERT_VALUES, localB)); 517227f9e03SJunchao Zhang } 518227f9e03SJunchao Zhang /* 519227f9e03SJunchao Zhang Get a pointer to vector data. 520227f9e03SJunchao Zhang - For default PETSc vectors, VecGetArray() returns a pointer to 521227f9e03SJunchao Zhang the data array. Otherwise, the routine is implementation dependent. 522227f9e03SJunchao Zhang - You MUST call VecRestoreArray() when you no longer need access to 523227f9e03SJunchao Zhang the array. 524227f9e03SJunchao Zhang */ 525227f9e03SJunchao Zhang PetscCall(DMDAVecGetArray(da, localX, &x)); 526227f9e03SJunchao Zhang if (B) PetscCall(DMDAVecGetArray(da, localB, &b)); 527227f9e03SJunchao Zhang /* 528227f9e03SJunchao Zhang Get local grid boundaries (for 2-dimensional DMDA): 529227f9e03SJunchao Zhang xs, ys - starting grid indices (no ghost points) 530227f9e03SJunchao Zhang xm, ym - widths of local grid (no ghost points) 531227f9e03SJunchao Zhang */ 532227f9e03SJunchao Zhang PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 533227f9e03SJunchao Zhang 534227f9e03SJunchao Zhang for (j = ys; j < ys + ym; j++) { 535227f9e03SJunchao Zhang for (i = xs; i < xs + xm; i++) { 536227f9e03SJunchao Zhang if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 537227f9e03SJunchao Zhang /* boundary conditions are all zero Dirichlet */ 538227f9e03SJunchao Zhang x[j][i] = 0.0; 539227f9e03SJunchao Zhang } else { 540227f9e03SJunchao Zhang if (B) bij = b[j][i]; 541227f9e03SJunchao Zhang else bij = 0.; 542227f9e03SJunchao Zhang 543227f9e03SJunchao Zhang u = x[j][i]; 544227f9e03SJunchao Zhang un = x[j - 1][i]; 545227f9e03SJunchao Zhang us = x[j + 1][i]; 546227f9e03SJunchao Zhang ue = x[j][i - 1]; 547227f9e03SJunchao Zhang uw = x[j][i + 1]; 548227f9e03SJunchao Zhang 549227f9e03SJunchao Zhang for (k = 0; k < its; k++) { 550227f9e03SJunchao Zhang eu = PetscExpScalar(u); 551227f9e03SJunchao Zhang uxx = (2.0 * u - ue - uw) * hydhx; 552227f9e03SJunchao Zhang uyy = (2.0 * u - un - us) * hxdhy; 553227f9e03SJunchao Zhang F = uxx + uyy - sc * eu - bij; 554227f9e03SJunchao Zhang if (k == 0) F0 = F; 555227f9e03SJunchao Zhang J = 2.0 * (hydhx + hxdhy) - sc * eu; 556227f9e03SJunchao Zhang y = F / J; 557227f9e03SJunchao Zhang u -= y; 558227f9e03SJunchao Zhang tot_its++; 559227f9e03SJunchao Zhang 5609371c9d4SSatish Balay if (atol > PetscAbsReal(PetscRealPart(F)) || rtol * PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || stol * PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) { break; } 561227f9e03SJunchao Zhang } 562227f9e03SJunchao Zhang x[j][i] = u; 563227f9e03SJunchao Zhang } 564227f9e03SJunchao Zhang } 565227f9e03SJunchao Zhang } 566227f9e03SJunchao Zhang /* 567227f9e03SJunchao Zhang Restore vector 568227f9e03SJunchao Zhang */ 569227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(da, localX, &x)); 570227f9e03SJunchao Zhang PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X)); 571227f9e03SJunchao Zhang PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X)); 572227f9e03SJunchao Zhang } 573227f9e03SJunchao Zhang PetscCall(PetscLogFlops(tot_its * (21.0))); 574227f9e03SJunchao Zhang PetscCall(DMRestoreLocalVector(da, &localX)); 575227f9e03SJunchao Zhang if (B) { 576227f9e03SJunchao Zhang PetscCall(DMDAVecRestoreArray(da, localB, &b)); 577227f9e03SJunchao Zhang PetscCall(DMRestoreLocalVector(da, &localB)); 578227f9e03SJunchao Zhang } 579227f9e03SJunchao Zhang PetscFunctionReturn(0); 580227f9e03SJunchao Zhang } 581c4762a1bSJed Brown 5829371c9d4SSatish Balay int main(int argc, char **argv) { 583c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 584c4762a1bSJed Brown Vec x; /* solution vector */ 585c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 586c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 587c4762a1bSJed Brown PetscReal bratu_lambda_max = 6.81; 588c4762a1bSJed Brown PetscReal bratu_lambda_min = 0.; 58916d037c0SJunchao Zhang PetscInt MMS = 1; 59016d037c0SJunchao Zhang PetscBool flg = PETSC_FALSE, setMMS; 591c4762a1bSJed Brown DM da; 592c4762a1bSJed Brown Vec r = NULL; 593c4762a1bSJed Brown KSP ksp; 594c4762a1bSJed Brown PetscInt lits, slits; 595c4762a1bSJed Brown 596c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 597c4762a1bSJed Brown Initialize program 598c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 599c4762a1bSJed Brown 600327415f7SBarry Smith PetscFunctionBeginUser; 6019566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 602c4762a1bSJed Brown 603c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 604c4762a1bSJed Brown Initialize problem parameters 605c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 606c4762a1bSJed Brown user.param = 6.0; 6079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-par", &user.param, NULL)); 608e00437b9SBarry Smith PetscCheck(user.param <= bratu_lambda_max && user.param >= bratu_lambda_min, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Lambda, %g, is out of range, [%g, %g]", (double)user.param, (double)bratu_lambda_min, (double)bratu_lambda_max); 60916d037c0SJunchao Zhang PetscCall(PetscOptionsGetInt(NULL, NULL, "-mms", &MMS, &setMMS)); 610c4762a1bSJed Brown if (MMS == 3) { 611c4762a1bSJed Brown PetscInt mPar = 2, nPar = 1; 6129566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m_par", &mPar, NULL)); 6139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n_par", &nPar, NULL)); 614c4762a1bSJed Brown user.m = PetscPowInt(2, mPar); 615c4762a1bSJed Brown user.n = PetscPowInt(2, nPar); 616c4762a1bSJed Brown } 617c4762a1bSJed Brown 618c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 619c4762a1bSJed Brown Create nonlinear solver context 620c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6219566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 6229566063dSJacob Faibussowitsch PetscCall(SNESSetCountersReset(snes, PETSC_FALSE)); 6239566063dSJacob Faibussowitsch PetscCall(SNESSetNGS(snes, NonlinearGS, NULL)); 624c4762a1bSJed Brown 625c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 626c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 627c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6289566063dSJacob 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)); 6299566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 6309566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 6319566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 6329566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da, &user)); 6339566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, da)); 634c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 635c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 636c4762a1bSJed Brown vectors that are the same types 637c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6389566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 639c4762a1bSJed Brown 640c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 641c4762a1bSJed Brown Set local function evaluation routine 642c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 643c4762a1bSJed Brown switch (MMS) { 6449371c9d4SSatish Balay case 0: 6459371c9d4SSatish Balay user.mms_solution = ZeroBCSolution; 6469371c9d4SSatish Balay user.mms_forcing = NULL; 6479371c9d4SSatish Balay break; 6489371c9d4SSatish Balay case 1: 6499371c9d4SSatish Balay user.mms_solution = MMSSolution1; 6509371c9d4SSatish Balay user.mms_forcing = MMSForcing1; 6519371c9d4SSatish Balay break; 6529371c9d4SSatish Balay case 2: 6539371c9d4SSatish Balay user.mms_solution = MMSSolution2; 6549371c9d4SSatish Balay user.mms_forcing = MMSForcing2; 6559371c9d4SSatish Balay break; 6569371c9d4SSatish Balay case 3: 6579371c9d4SSatish Balay user.mms_solution = MMSSolution3; 6589371c9d4SSatish Balay user.mms_forcing = MMSForcing3; 6599371c9d4SSatish Balay break; 6609371c9d4SSatish Balay case 4: 6619371c9d4SSatish Balay user.mms_solution = MMSSolution4; 6629371c9d4SSatish Balay user.mms_forcing = MMSForcing4; 6639371c9d4SSatish Balay break; 66463a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Unknown MMS type %" PetscInt_FMT, MMS); 665c4762a1bSJed Brown } 6669566063dSJacob Faibussowitsch PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (DMDASNESFunction)FormFunctionLocal, &user)); 6679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-fd", &flg, NULL)); 668*48a46eb9SPierre Jolivet if (!flg) PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobian)FormJacobianLocal, &user)); 669c4762a1bSJed Brown 6709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-obj", &flg, NULL)); 671*48a46eb9SPierre Jolivet if (flg) PetscCall(DMDASNESSetObjectiveLocal(da, (DMDASNESObjective)FormObjectiveLocal, &user)); 672c4762a1bSJed Brown 6734e1ad211SJed Brown if (PetscDefined(HAVE_MATLAB_ENGINE)) { 6744e1ad211SJed Brown PetscBool matlab_function = PETSC_FALSE; 6759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-matlab_function", &matlab_function, 0)); 676c4762a1bSJed Brown if (matlab_function) { 6779566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 6789566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunctionMatlab, &user)); 679c4762a1bSJed Brown } 6804e1ad211SJed Brown } 681c4762a1bSJed Brown 682c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 683c4762a1bSJed Brown Customize nonlinear solver; set runtime options 684c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6859566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 686c4762a1bSJed Brown 6879566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(da, &user, x)); 688c4762a1bSJed Brown 689c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 690c4762a1bSJed Brown Solve nonlinear system 691c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 6929566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 6939566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 694c4762a1bSJed Brown 6959566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &slits)); 6969566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 6979566063dSJacob Faibussowitsch PetscCall(KSPGetTotalIterations(ksp, &lits)); 69863a3b9bcSJacob Faibussowitsch PetscCheck(lits == slits, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Number of total linear iterations reported by SNES %" PetscInt_FMT " does not match reported by KSP %" PetscInt_FMT, slits, lits); 699c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 700c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 701c4762a1bSJed Brown 702c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 703c4762a1bSJed Brown If using MMS, check the l_2 error 704c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 70516d037c0SJunchao Zhang if (setMMS) { 706c4762a1bSJed Brown Vec e; 707c4762a1bSJed Brown PetscReal errorl2, errorinf; 708c4762a1bSJed Brown PetscInt N; 709c4762a1bSJed Brown 7109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &e)); 7119566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)x, NULL, "-sol_view")); 7129566063dSJacob Faibussowitsch PetscCall(FormExactSolution(da, &user, e)); 7139566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)e, NULL, "-exact_view")); 7149566063dSJacob Faibussowitsch PetscCall(VecAXPY(e, -1.0, x)); 7159566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)e, NULL, "-error_view")); 7169566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_2, &errorl2)); 7179566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_INFINITY, &errorinf)); 7189566063dSJacob Faibussowitsch PetscCall(VecGetSize(e, &N)); 71963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " error L2 %g inf %g\n", N, (double)(errorl2 / PetscSqrtReal((PetscReal)N)), (double)errorinf)); 7209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 7219566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetDof(SNES_Solve, 0, N)); 7229566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetError(SNES_Solve, 0, errorl2 / PetscSqrtReal(N))); 723c4762a1bSJed Brown } 724c4762a1bSJed Brown 725c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 726c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 727c4762a1bSJed Brown are no longer needed. 728c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 7299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 7309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 7319566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 7329566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 7339566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 734b122ec5aSJacob Faibussowitsch return 0; 735c4762a1bSJed Brown } 736c4762a1bSJed Brown 737c4762a1bSJed Brown /*TEST 738c4762a1bSJed Brown 739c4762a1bSJed Brown test: 740c4762a1bSJed Brown suffix: asm_0 741c4762a1bSJed Brown requires: !single 742c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu 743c4762a1bSJed Brown 744c4762a1bSJed Brown test: 745c4762a1bSJed Brown suffix: msm_0 746c4762a1bSJed Brown requires: !single 747c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu 748c4762a1bSJed Brown 749c4762a1bSJed Brown test: 750c4762a1bSJed Brown suffix: asm_1 751c4762a1bSJed Brown requires: !single 752c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 753c4762a1bSJed Brown 754c4762a1bSJed Brown test: 755c4762a1bSJed Brown suffix: msm_1 756c4762a1bSJed Brown requires: !single 757c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 758c4762a1bSJed Brown 759c4762a1bSJed Brown test: 760c4762a1bSJed Brown suffix: asm_2 761c4762a1bSJed Brown requires: !single 762c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 763c4762a1bSJed Brown 764c4762a1bSJed Brown test: 765c4762a1bSJed Brown suffix: msm_2 766c4762a1bSJed Brown requires: !single 767c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 768c4762a1bSJed Brown 769c4762a1bSJed Brown test: 770c4762a1bSJed Brown suffix: asm_3 771c4762a1bSJed Brown requires: !single 772c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 773c4762a1bSJed Brown 774c4762a1bSJed Brown test: 775c4762a1bSJed Brown suffix: msm_3 776c4762a1bSJed Brown requires: !single 777c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 778c4762a1bSJed Brown 779c4762a1bSJed Brown test: 780c4762a1bSJed Brown suffix: asm_4 781c4762a1bSJed Brown requires: !single 782c4762a1bSJed Brown nsize: 2 783c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 784c4762a1bSJed Brown 785c4762a1bSJed Brown test: 786c4762a1bSJed Brown suffix: msm_4 787c4762a1bSJed Brown requires: !single 788c4762a1bSJed Brown nsize: 2 789c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 790c4762a1bSJed Brown 791c4762a1bSJed Brown test: 792c4762a1bSJed Brown suffix: asm_5 793c4762a1bSJed Brown requires: !single 794c4762a1bSJed Brown nsize: 2 795c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8 796c4762a1bSJed Brown 797c4762a1bSJed Brown test: 798c4762a1bSJed Brown suffix: msm_5 799c4762a1bSJed Brown requires: !single 800c4762a1bSJed Brown nsize: 2 801c4762a1bSJed Brown args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8 802c4762a1bSJed Brown 803c4762a1bSJed Brown test: 804c4762a1bSJed Brown args: -snes_rtol 1.e-5 -pc_type mg -ksp_monitor_short -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17 -mg_levels_ksp_monitor_short -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full 805c4762a1bSJed Brown requires: !single 806c4762a1bSJed Brown 807c4762a1bSJed Brown test: 808c4762a1bSJed Brown suffix: 2 809c4762a1bSJed Brown args: -pc_type mg -ksp_converged_reason -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3 -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full -ksp_atol -1. 810c4762a1bSJed Brown 811c4762a1bSJed Brown test: 812c4762a1bSJed Brown suffix: 3 813c4762a1bSJed Brown nsize: 2 814c4762a1bSJed Brown args: -snes_grid_sequence 2 -snes_mf_operator -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -snes_rtol 1.e-2 815c4762a1bSJed Brown filter: grep -v "otal number of function evaluations" 816c4762a1bSJed Brown 817c4762a1bSJed Brown test: 818c4762a1bSJed Brown suffix: 4 819c4762a1bSJed Brown nsize: 2 820c4762a1bSJed Brown args: -snes_grid_sequence 2 -snes_monitor_short -ksp_converged_reason -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -ksp_atol -1 821c4762a1bSJed Brown 822c4762a1bSJed Brown test: 823c4762a1bSJed Brown suffix: 5_anderson 824c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type anderson 825c4762a1bSJed Brown 826c4762a1bSJed Brown test: 827c4762a1bSJed Brown suffix: 5_aspin 828c4762a1bSJed Brown nsize: 4 829c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type aspin -snes_view 830c4762a1bSJed Brown 831c4762a1bSJed Brown test: 832c4762a1bSJed Brown suffix: 5_broyden 833c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_qn_type broyden -snes_qn_m 10 834c4762a1bSJed Brown 835c4762a1bSJed Brown test: 836c4762a1bSJed Brown suffix: 5_fas 837c4762a1bSJed Brown args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6 838c4762a1bSJed Brown requires: !single 839c4762a1bSJed Brown 840c4762a1bSJed Brown test: 841c4762a1bSJed Brown suffix: 5_fas_additive 842c4762a1bSJed Brown args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6 -snes_fas_type additive -snes_max_it 50 843c4762a1bSJed Brown 844c4762a1bSJed Brown test: 845c4762a1bSJed Brown suffix: 5_fas_monitor 846c4762a1bSJed Brown args: -da_refine 1 -snes_type fas -snes_fas_monitor 847c4762a1bSJed Brown requires: !single 848c4762a1bSJed Brown 849c4762a1bSJed Brown test: 850c4762a1bSJed Brown suffix: 5_ls 851c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls 852c4762a1bSJed Brown 853c4762a1bSJed Brown test: 854c4762a1bSJed Brown suffix: 5_ls_sell_sor 855c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls -dm_mat_type sell -pc_type sor 856c4762a1bSJed Brown output_file: output/ex5_5_ls.out 857c4762a1bSJed Brown 858c4762a1bSJed Brown test: 859c4762a1bSJed Brown suffix: 5_nasm 860c4762a1bSJed Brown nsize: 4 861c4762a1bSJed Brown args: -snes_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type nasm -snes_nasm_type restrict -snes_max_it 10 862c4762a1bSJed Brown 863c4762a1bSJed Brown test: 864c4762a1bSJed Brown suffix: 5_ncg 865c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ncg -snes_ncg_type fr 866c4762a1bSJed Brown 867c4762a1bSJed Brown test: 868c4762a1bSJed Brown suffix: 5_newton_asm_dmda 869c4762a1bSJed Brown nsize: 4 870c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type asm -pc_asm_dm_subdomains -malloc_dump 871c4762a1bSJed Brown requires: !single 872c4762a1bSJed Brown 873c4762a1bSJed Brown test: 874c4762a1bSJed Brown suffix: 5_newton_gasm_dmda 875c4762a1bSJed Brown nsize: 4 876c4762a1bSJed Brown args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type gasm -malloc_dump 877c4762a1bSJed Brown requires: !single 878c4762a1bSJed Brown 879c4762a1bSJed Brown test: 880c4762a1bSJed Brown suffix: 5_ngmres 881c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10 882c4762a1bSJed Brown 883c4762a1bSJed Brown test: 884c4762a1bSJed Brown suffix: 5_ngmres_fas 885c4762a1bSJed Brown args: -snes_rtol 1.e-4 -snes_type ngmres -npc_fas_coarse_snes_max_it 1 -npc_fas_coarse_snes_type newtonls -npc_fas_coarse_pc_type lu -npc_fas_coarse_ksp_type preonly -snes_ngmres_m 10 -snes_monitor_short -npc_snes_max_it 1 -npc_snes_type fas -npc_fas_coarse_ksp_type richardson -da_refine 6 886c4762a1bSJed Brown 887c4762a1bSJed Brown test: 888c4762a1bSJed Brown suffix: 5_ngmres_ngs 889c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -npc_snes_type ngs -npc_snes_max_it 1 890c4762a1bSJed Brown 891c4762a1bSJed Brown test: 892c4762a1bSJed Brown suffix: 5_ngmres_nrichardson 893c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10 -npc_snes_type nrichardson -npc_snes_max_it 3 894c4762a1bSJed Brown output_file: output/ex5_5_ngmres_richardson.out 895c4762a1bSJed Brown 896c4762a1bSJed Brown test: 897c4762a1bSJed Brown suffix: 5_nrichardson 898c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type nrichardson 899c4762a1bSJed Brown 900c4762a1bSJed Brown test: 901c4762a1bSJed Brown suffix: 5_qn 902c4762a1bSJed Brown args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_linesearch_type cp -snes_qn_m 10 903c4762a1bSJed Brown 904c4762a1bSJed Brown test: 905c4762a1bSJed Brown suffix: 6 906c4762a1bSJed Brown nsize: 4 907c4762a1bSJed Brown args: -snes_converged_reason -ksp_converged_reason -da_grid_x 129 -da_grid_y 129 -pc_type mg -pc_mg_levels 8 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.5,0,1.1 -mg_levels_ksp_max_it 2 908c4762a1bSJed Brown 909c4762a1bSJed Brown test: 910c4762a1bSJed Brown requires: complex !single 911c4762a1bSJed Brown suffix: complex 912c4762a1bSJed Brown args: -snes_mf_operator -mat_mffd_complex -snes_monitor 913c4762a1bSJed Brown 914c4762a1bSJed Brown TEST*/ 915