1 static const char help[] = "p-Bratu nonlinear PDE in 2d.\n\ 2 We solve the p-Laplacian (nonlinear diffusion) combined with\n\ 3 the Bratu (solid fuel ignition) nonlinearity in a 2D rectangular\n\ 4 domain, using distributed arrays (DAs) to partition the parallel grid.\n\ 5 The command line options include:\n\ 6 -p <2>: `p' in p-Laplacian term\n\ 7 -epsilon <1e-05>: Strain-regularization in p-Laplacian\n\ 8 -lambda <6>: Bratu parameter\n\ 9 -blocks <bx,by>: number of coefficient interfaces in x and y direction\n\ 10 -kappa <1e-3>: diffusivity in odd regions\n\ 11 \n"; 12 13 /*F 14 The $p$-Bratu problem is a combination of the $p$-Laplacian (nonlinear diffusion) and the Brutu solid fuel ignition problem. 15 This problem is modeled by the partial differential equation 16 17 \begin{equation*} 18 -\nabla\cdot (\eta \nabla u) - \lambda \exp(u) = 0 19 \end{equation*} 20 21 on $\Omega = (-1,1)^2$ with closure 22 23 \begin{align*} 24 \eta(\gamma) &= (\epsilon^2 + \gamma)^{(p-2)/2} & \gamma &= \frac 1 2 |\nabla u|^2 25 \end{align*} 26 27 and boundary conditions $u = 0$ for $(x,y) \in \partial \Omega$ 28 29 A 9-point finite difference stencil is used to discretize 30 the boundary value problem to obtain a nonlinear system of equations. 31 This would be a 5-point stencil if not for the $p$-Laplacian's nonlinearity. 32 F*/ 33 34 /* 35 mpiexec -n 2 ./ex15 -snes_monitor -ksp_monitor log_summary 36 */ 37 38 /* 39 Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 40 Include "petscsnes.h" so that we can use SNES solvers. Note that this 41 file automatically includes: 42 petsc.h - base PETSc routines petscvec.h - vectors 43 petscsys.h - system routines petscmat.h - matrices 44 petscis.h - index sets petscksp.h - Krylov subspace methods 45 petscviewer.h - viewers petscpc.h - preconditioners 46 petscksp.h - linear solvers 47 */ 48 49 #include <petscdm.h> 50 #include <petscdmda.h> 51 #include <petscsnes.h> 52 53 typedef enum { 54 JAC_BRATU, 55 JAC_PICARD, 56 JAC_STAR, 57 JAC_NEWTON 58 } JacType; 59 static const char *const JacTypes[] = {"BRATU", "PICARD", "STAR", "NEWTON", "JacType", "JAC_", 0}; 60 61 /* 62 User-defined application context - contains data needed by the 63 application-provided call-back routines, FormJacobianLocal() and 64 FormFunctionLocal(). 65 */ 66 typedef struct { 67 PetscReal lambda; /* Bratu parameter */ 68 PetscReal p; /* Exponent in p-Laplacian */ 69 PetscReal epsilon; /* Regularization */ 70 PetscReal source; /* Source term */ 71 JacType jtype; /* What type of Jacobian to assemble */ 72 PetscBool picard; 73 PetscInt blocks[2]; 74 PetscReal kappa; 75 PetscInt initial; /* initial conditions type */ 76 } AppCtx; 77 78 /* 79 User-defined routines 80 */ 81 static PetscErrorCode FormRHS(AppCtx *, DM, Vec); 82 static PetscErrorCode FormInitialGuess(AppCtx *, DM, Vec); 83 static PetscErrorCode FormFunctionLocal(DMDALocalInfo *, PetscScalar **, PetscScalar **, AppCtx *); 84 static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *, PetscScalar **, PetscScalar **, AppCtx *); 85 static PetscErrorCode FormJacobianLocal(DMDALocalInfo *, PetscScalar **, Mat, Mat, AppCtx *); 86 static PetscErrorCode NonlinearGS(SNES, Vec, Vec, void *); 87 88 typedef struct _n_PreCheck *PreCheck; 89 struct _n_PreCheck { 90 MPI_Comm comm; 91 PetscReal angle; 92 Vec Ylast; 93 PetscViewer monitor; 94 }; 95 PetscErrorCode PreCheckCreate(MPI_Comm, PreCheck *); 96 PetscErrorCode PreCheckDestroy(PreCheck *); 97 PetscErrorCode PreCheckFunction(SNESLineSearch, Vec, Vec, PetscBool *, void *); 98 PetscErrorCode PreCheckSetFromOptions(PreCheck); 99 100 int main(int argc, char **argv) { 101 SNES snes; /* nonlinear solver */ 102 Vec x, r, b; /* solution, residual, rhs vectors */ 103 AppCtx user; /* user-defined work context */ 104 PetscInt its; /* iterations for convergence */ 105 SNESConvergedReason reason; /* Check convergence */ 106 PetscBool alloc_star; /* Only allocate for the STAR stencil */ 107 PetscBool write_output; 108 char filename[PETSC_MAX_PATH_LEN] = "ex15.vts"; 109 PetscReal bratu_lambda_max = 6.81, bratu_lambda_min = 0.; 110 DM da; /* distributed array data structure */ 111 PreCheck precheck = NULL; /* precheck context for version in this file */ 112 PetscInt use_precheck; /* 0=none, 1=version in this file, 2=SNES-provided version */ 113 PetscReal precheck_angle; /* When manually setting the SNES-provided precheck function */ 114 SNESLineSearch linesearch; 115 116 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117 Initialize program 118 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 119 120 PetscFunctionBeginUser; 121 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 122 123 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124 Initialize problem parameters 125 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 126 user.lambda = 0.0; 127 user.p = 2.0; 128 user.epsilon = 1e-5; 129 user.source = 0.1; 130 user.jtype = JAC_NEWTON; 131 user.initial = -1; 132 user.blocks[0] = 1; 133 user.blocks[1] = 1; 134 user.kappa = 1e-3; 135 alloc_star = PETSC_FALSE; 136 use_precheck = 0; 137 precheck_angle = 10.; 138 user.picard = PETSC_FALSE; 139 PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "p-Bratu options", __FILE__); 140 { 141 PetscInt two = 2; 142 PetscCall(PetscOptionsReal("-lambda", "Bratu parameter", "", user.lambda, &user.lambda, NULL)); 143 PetscCall(PetscOptionsReal("-p", "Exponent `p' in p-Laplacian", "", user.p, &user.p, NULL)); 144 PetscCall(PetscOptionsReal("-epsilon", "Strain-regularization in p-Laplacian", "", user.epsilon, &user.epsilon, NULL)); 145 PetscCall(PetscOptionsReal("-source", "Constant source term", "", user.source, &user.source, NULL)); 146 PetscCall(PetscOptionsEnum("-jtype", "Jacobian approximation to assemble", "", JacTypes, (PetscEnum)user.jtype, (PetscEnum *)&user.jtype, NULL)); 147 PetscCall(PetscOptionsName("-picard", "Solve with defect-correction Picard iteration", "", &user.picard)); 148 if (user.picard) { 149 user.jtype = JAC_PICARD; 150 PetscCheck(user.p == 3, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Picard iteration is only supported for p == 3"); 151 /* the Picard linearization only requires a 5 point stencil, while the Newton linearization requires a 9 point stencil */ 152 /* hence allocating the 5 point stencil gives the same convergence as the 9 point stencil since the extra stencil points are not used */ 153 PetscCall(PetscOptionsBool("-alloc_star", "Allocate for STAR stencil (5-point)", "", alloc_star, &alloc_star, NULL)); 154 } 155 PetscCall(PetscOptionsInt("-precheck", "Use a pre-check correction intended for use with Picard iteration 1=this version, 2=library", "", use_precheck, &use_precheck, NULL)); 156 PetscCall(PetscOptionsInt("-initial", "Initial conditions type (-1: default, 0: zero-valued, 1: peaked guess)", "", user.initial, &user.initial, NULL)); 157 if (use_precheck == 2) { /* Using library version, get the angle */ 158 PetscCall(PetscOptionsReal("-precheck_angle", "Angle in degrees between successive search directions necessary to activate step correction", "", precheck_angle, &precheck_angle, NULL)); 159 } 160 PetscCall(PetscOptionsIntArray("-blocks", "number of coefficient interfaces in x and y direction", "", user.blocks, &two, NULL)); 161 PetscCall(PetscOptionsReal("-kappa", "diffusivity in odd regions", "", user.kappa, &user.kappa, NULL)); 162 PetscCall(PetscOptionsString("-o", "Output solution in vts format", "", filename, filename, sizeof(filename), &write_output)); 163 } 164 PetscOptionsEnd(); 165 if (user.lambda > bratu_lambda_max || user.lambda < bratu_lambda_min) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "WARNING: lambda %g out of range for p=2\n", (double)user.lambda)); 166 167 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 168 Create nonlinear solver context 169 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 170 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 171 172 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 173 Create distributed array (DMDA) to manage parallel grid and vectors 174 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 175 PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, alloc_star ? DMDA_STENCIL_STAR : DMDA_STENCIL_BOX, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da)); 176 PetscCall(DMSetFromOptions(da)); 177 PetscCall(DMSetUp(da)); 178 179 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 180 Extract global vectors from DM; then duplicate for remaining 181 vectors that are the same types 182 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 183 PetscCall(DMCreateGlobalVector(da, &x)); 184 PetscCall(VecDuplicate(x, &r)); 185 PetscCall(VecDuplicate(x, &b)); 186 187 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 188 User can override with: 189 -snes_mf : matrix-free Newton-Krylov method with no preconditioning 190 (unless user explicitly sets preconditioner) 191 -snes_mf_operator : form preconditioning matrix as set by the user, 192 but use matrix-free approx for Jacobian-vector 193 products within Newton-Krylov method 194 195 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 196 197 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 198 Set local function evaluation routine 199 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 200 PetscCall(DMSetApplicationContext(da, &user)); 201 PetscCall(SNESSetDM(snes, da)); 202 if (user.picard) { 203 /* 204 This is not really right requiring the user to call SNESSetFunction/Jacobian but the DMDASNESSetPicardLocal() cannot access 205 the SNES to set it 206 */ 207 PetscCall(DMDASNESSetPicardLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionPicardLocal, (PetscErrorCode(*)(DMDALocalInfo *, void *, Mat, Mat, void *))FormJacobianLocal, &user)); 208 PetscCall(SNESSetFunction(snes, NULL, SNESPicardComputeFunction, &user)); 209 PetscCall(SNESSetJacobian(snes, NULL, NULL, SNESPicardComputeJacobian, &user)); 210 } else { 211 PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, &user)); 212 PetscCall(DMDASNESSetJacobianLocal(da, (PetscErrorCode(*)(DMDALocalInfo *, void *, Mat, Mat, void *))FormJacobianLocal, &user)); 213 } 214 215 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 216 Customize nonlinear solver; set runtime options 217 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 218 PetscCall(SNESSetFromOptions(snes)); 219 PetscCall(SNESSetNGS(snes, NonlinearGS, &user)); 220 PetscCall(SNESGetLineSearch(snes, &linesearch)); 221 /* Set up the precheck context if requested */ 222 if (use_precheck == 1) { /* Use the precheck routines in this file */ 223 PetscCall(PreCheckCreate(PETSC_COMM_WORLD, &precheck)); 224 PetscCall(PreCheckSetFromOptions(precheck)); 225 PetscCall(SNESLineSearchSetPreCheck(linesearch, PreCheckFunction, precheck)); 226 } else if (use_precheck == 2) { /* Use the version provided by the library */ 227 PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &precheck_angle)); 228 } 229 230 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 231 Evaluate initial guess 232 Note: The user should initialize the vector, x, with the initial guess 233 for the nonlinear solver prior to calling SNESSolve(). In particular, 234 to employ an initial guess of zero, the user should explicitly set 235 this vector to zero by calling VecSet(). 236 */ 237 238 PetscCall(FormInitialGuess(&user, da, x)); 239 PetscCall(FormRHS(&user, da, b)); 240 241 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 242 Solve nonlinear system 243 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 244 PetscCall(SNESSolve(snes, b, x)); 245 PetscCall(SNESGetIterationNumber(snes, &its)); 246 PetscCall(SNESGetConvergedReason(snes, &reason)); 247 248 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s Number of nonlinear iterations = %" PetscInt_FMT "\n", SNESConvergedReasons[reason], its)); 249 250 if (write_output) { 251 PetscViewer viewer; 252 PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD, filename, FILE_MODE_WRITE, &viewer)); 253 PetscCall(VecView(x, viewer)); 254 PetscCall(PetscViewerDestroy(&viewer)); 255 } 256 257 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 258 Free work space. All PETSc objects should be destroyed when they 259 are no longer needed. 260 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 261 262 PetscCall(VecDestroy(&x)); 263 PetscCall(VecDestroy(&r)); 264 PetscCall(VecDestroy(&b)); 265 PetscCall(SNESDestroy(&snes)); 266 PetscCall(DMDestroy(&da)); 267 PetscCall(PreCheckDestroy(&precheck)); 268 PetscCall(PetscFinalize()); 269 return 0; 270 } 271 272 /* ------------------------------------------------------------------- */ 273 /* 274 FormInitialGuess - Forms initial approximation. 275 276 Input Parameters: 277 user - user-defined application context 278 X - vector 279 280 Output Parameter: 281 X - vector 282 */ 283 static PetscErrorCode FormInitialGuess(AppCtx *user, DM da, Vec X) { 284 PetscInt i, j, Mx, My, xs, ys, xm, ym; 285 PetscReal temp1, temp, hx, hy; 286 PetscScalar **x; 287 288 PetscFunctionBeginUser; 289 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)); 290 291 hx = 1.0 / (PetscReal)(Mx - 1); 292 hy = 1.0 / (PetscReal)(My - 1); 293 temp1 = user->lambda / (user->lambda + 1.); 294 295 /* 296 Get a pointer to vector data. 297 - For default PETSc vectors, VecGetArray() returns a pointer to 298 the data array. Otherwise, the routine is implementation dependent. 299 - You MUST call VecRestoreArray() when you no longer need access to 300 the array. 301 */ 302 PetscCall(DMDAVecGetArray(da, X, &x)); 303 304 /* 305 Get local grid boundaries (for 2-dimensional DA): 306 xs, ys - starting grid indices (no ghost points) 307 xm, ym - widths of local grid (no ghost points) 308 309 */ 310 PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 311 312 /* 313 Compute initial guess over the locally owned part of the grid 314 */ 315 for (j = ys; j < ys + ym; j++) { 316 temp = (PetscReal)(PetscMin(j, My - j - 1)) * hy; 317 for (i = xs; i < xs + xm; i++) { 318 if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 319 /* boundary conditions are all zero Dirichlet */ 320 x[j][i] = 0.0; 321 } else { 322 if (user->initial == -1) { 323 if (user->lambda != 0) { 324 x[j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp)); 325 } else { 326 /* The solution above is an exact solution for lambda=0, this avoids "accidentally" starting 327 * with an exact solution. */ 328 const PetscReal xx = 2 * (PetscReal)i / (Mx - 1) - 1, yy = 2 * (PetscReal)j / (My - 1) - 1; 329 x[j][i] = (1 - xx * xx) * (1 - yy * yy) * xx * yy; 330 } 331 } else if (user->initial == 0) { 332 x[j][i] = 0.; 333 } else if (user->initial == 1) { 334 const PetscReal xx = 2 * (PetscReal)i / (Mx - 1) - 1, yy = 2 * (PetscReal)j / (My - 1) - 1; 335 x[j][i] = (1 - xx * xx) * (1 - yy * yy) * xx * yy; 336 } else { 337 if (user->lambda != 0) { 338 x[j][i] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp)); 339 } else { 340 x[j][i] = 0.5 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, Mx - i - 1)) * hx, temp)); 341 } 342 } 343 } 344 } 345 } 346 /* 347 Restore vector 348 */ 349 PetscCall(DMDAVecRestoreArray(da, X, &x)); 350 PetscFunctionReturn(0); 351 } 352 353 /* 354 FormRHS - Forms constant RHS for the problem. 355 356 Input Parameters: 357 user - user-defined application context 358 B - RHS vector 359 360 Output Parameter: 361 B - vector 362 */ 363 static PetscErrorCode FormRHS(AppCtx *user, DM da, Vec B) { 364 PetscInt i, j, Mx, My, xs, ys, xm, ym; 365 PetscReal hx, hy; 366 PetscScalar **b; 367 368 PetscFunctionBeginUser; 369 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)); 370 371 hx = 1.0 / (PetscReal)(Mx - 1); 372 hy = 1.0 / (PetscReal)(My - 1); 373 PetscCall(DMDAVecGetArray(da, B, &b)); 374 PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 375 for (j = ys; j < ys + ym; j++) { 376 for (i = xs; i < xs + xm; i++) { 377 if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 378 b[j][i] = 0.0; 379 } else { 380 b[j][i] = hx * hy * user->source; 381 } 382 } 383 } 384 PetscCall(DMDAVecRestoreArray(da, B, &b)); 385 PetscFunctionReturn(0); 386 } 387 388 static inline PetscReal kappa(const AppCtx *ctx, PetscReal x, PetscReal y) { 389 return (((PetscInt)(x * ctx->blocks[0])) + ((PetscInt)(y * ctx->blocks[1]))) % 2 ? ctx->kappa : 1.0; 390 } 391 /* p-Laplacian diffusivity */ 392 static inline PetscScalar eta(const AppCtx *ctx, PetscReal x, PetscReal y, PetscScalar ux, PetscScalar uy) { 393 return kappa(ctx, x, y) * PetscPowScalar(PetscSqr(ctx->epsilon) + 0.5 * (ux * ux + uy * uy), 0.5 * (ctx->p - 2.)); 394 } 395 static inline PetscScalar deta(const AppCtx *ctx, PetscReal x, PetscReal y, PetscScalar ux, PetscScalar uy) { 396 return (ctx->p == 2) ? 0 : kappa(ctx, x, y) * PetscPowScalar(PetscSqr(ctx->epsilon) + 0.5 * (ux * ux + uy * uy), 0.5 * (ctx->p - 4)) * 0.5 * (ctx->p - 2.); 397 } 398 399 /* ------------------------------------------------------------------- */ 400 /* 401 FormFunctionLocal - Evaluates nonlinear function, F(x). 402 */ 403 static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user) { 404 PetscReal hx, hy, dhx, dhy, sc; 405 PetscInt i, j; 406 PetscScalar eu; 407 408 PetscFunctionBeginUser; 409 hx = 1.0 / (PetscReal)(info->mx - 1); 410 hy = 1.0 / (PetscReal)(info->my - 1); 411 sc = hx * hy * user->lambda; 412 dhx = 1 / hx; 413 dhy = 1 / hy; 414 /* 415 Compute function over the locally owned part of the grid 416 */ 417 for (j = info->ys; j < info->ys + info->ym; j++) { 418 for (i = info->xs; i < info->xs + info->xm; i++) { 419 PetscReal xx = i * hx, yy = j * hy; 420 if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 421 f[j][i] = x[j][i]; 422 } else { 423 const PetscScalar u = x[j][i], ux_E = dhx * (x[j][i + 1] - x[j][i]), uy_E = 0.25 * dhy * (x[j + 1][i] + x[j + 1][i + 1] - x[j - 1][i] - x[j - 1][i + 1]), ux_W = dhx * (x[j][i] - x[j][i - 1]), uy_W = 0.25 * dhy * (x[j + 1][i - 1] + x[j + 1][i] - x[j - 1][i - 1] - x[j - 1][i]), ux_N = 0.25 * dhx * (x[j][i + 1] + x[j + 1][i + 1] - x[j][i - 1] - x[j + 1][i - 1]), uy_N = dhy * (x[j + 1][i] - x[j][i]), ux_S = 0.25 * dhx * (x[j - 1][i + 1] + x[j][i + 1] - x[j - 1][i - 1] - x[j][i - 1]), uy_S = dhy * (x[j][i] - x[j - 1][i]), e_E = eta(user, xx, yy, ux_E, uy_E), e_W = eta(user, xx, yy, ux_W, uy_W), e_N = eta(user, xx, yy, ux_N, uy_N), e_S = eta(user, xx, yy, ux_S, uy_S), uxx = -hy * (e_E * ux_E - e_W * ux_W), uyy = -hx * (e_N * uy_N - e_S * uy_S); 424 if (sc) eu = PetscExpScalar(u); 425 else eu = 0.; 426 /* For p=2, these terms decay to: 427 uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx 428 uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy 429 */ 430 f[j][i] = uxx + uyy - sc * eu; 431 } 432 } 433 } 434 PetscCall(PetscLogFlops(info->xm * info->ym * (72.0))); 435 PetscFunctionReturn(0); 436 } 437 438 /* 439 This is the opposite sign of the part of FormFunctionLocal that excludes the A(x) x part of the operation, 440 that is FormFunction applies A(x) x - b(x) while this applies b(x) because for Picard we think of it as solving A(x) x = b(x) 441 442 */ 443 static PetscErrorCode FormFunctionPicardLocal(DMDALocalInfo *info, PetscScalar **x, PetscScalar **f, AppCtx *user) { 444 PetscReal hx, hy, sc; 445 PetscInt i, j; 446 447 PetscFunctionBeginUser; 448 hx = 1.0 / (PetscReal)(info->mx - 1); 449 hy = 1.0 / (PetscReal)(info->my - 1); 450 sc = hx * hy * user->lambda; 451 /* 452 Compute function over the locally owned part of the grid 453 */ 454 for (j = info->ys; j < info->ys + info->ym; j++) { 455 for (i = info->xs; i < info->xs + info->xm; i++) { 456 if (!(i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1)) { 457 const PetscScalar u = x[j][i]; 458 f[j][i] = sc * PetscExpScalar(u); 459 } else { 460 f[j][i] = 0.0; /* this is zero because the A(x) x term forces the x to be zero on the boundary */ 461 } 462 } 463 } 464 PetscCall(PetscLogFlops(info->xm * info->ym * 2.0)); 465 PetscFunctionReturn(0); 466 } 467 468 /* 469 FormJacobianLocal - Evaluates Jacobian matrix. 470 */ 471 static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, PetscScalar **x, Mat J, Mat B, AppCtx *user) { 472 PetscInt i, j; 473 MatStencil col[9], row; 474 PetscScalar v[9]; 475 PetscReal hx, hy, hxdhy, hydhx, dhx, dhy, sc; 476 477 PetscFunctionBeginUser; 478 PetscCall(MatZeroEntries(B)); 479 hx = 1.0 / (PetscReal)(info->mx - 1); 480 hy = 1.0 / (PetscReal)(info->my - 1); 481 sc = hx * hy * user->lambda; 482 dhx = 1 / hx; 483 dhy = 1 / hy; 484 hxdhy = hx / hy; 485 hydhx = hy / hx; 486 487 /* 488 Compute entries for the locally owned part of the Jacobian. 489 - PETSc parallel matrix formats are partitioned by 490 contiguous chunks of rows across the processors. 491 - Each processor needs to insert only elements that it owns 492 locally (but any non-local elements will be sent to the 493 appropriate processor during matrix assembly). 494 - Here, we set all entries for a particular row at once. 495 */ 496 for (j = info->ys; j < info->ys + info->ym; j++) { 497 for (i = info->xs; i < info->xs + info->xm; i++) { 498 PetscReal xx = i * hx, yy = j * hy; 499 row.j = j; 500 row.i = i; 501 /* boundary points */ 502 if (i == 0 || j == 0 || i == info->mx - 1 || j == info->my - 1) { 503 v[0] = 1.0; 504 PetscCall(MatSetValuesStencil(B, 1, &row, 1, &row, v, INSERT_VALUES)); 505 } else { 506 /* interior grid points */ 507 const PetscScalar ux_E = dhx * (x[j][i + 1] - x[j][i]), uy_E = 0.25 * dhy * (x[j + 1][i] + x[j + 1][i + 1] - x[j - 1][i] - x[j - 1][i + 1]), ux_W = dhx * (x[j][i] - x[j][i - 1]), uy_W = 0.25 * dhy * (x[j + 1][i - 1] + x[j + 1][i] - x[j - 1][i - 1] - x[j - 1][i]), ux_N = 0.25 * dhx * (x[j][i + 1] + x[j + 1][i + 1] - x[j][i - 1] - x[j + 1][i - 1]), uy_N = dhy * (x[j + 1][i] - x[j][i]), ux_S = 0.25 * dhx * (x[j - 1][i + 1] + x[j][i + 1] - x[j - 1][i - 1] - x[j][i - 1]), uy_S = dhy * (x[j][i] - x[j - 1][i]), u = x[j][i], e_E = eta(user, xx, yy, ux_E, uy_E), e_W = eta(user, xx, yy, ux_W, uy_W), e_N = eta(user, xx, yy, ux_N, uy_N), e_S = eta(user, xx, yy, ux_S, uy_S), de_E = deta(user, xx, yy, ux_E, uy_E), de_W = deta(user, xx, yy, ux_W, uy_W), de_N = deta(user, xx, yy, ux_N, uy_N), de_S = deta(user, xx, yy, ux_S, uy_S), skew_E = de_E * ux_E * uy_E, skew_W = de_W * ux_W * uy_W, skew_N = de_N * ux_N * uy_N, skew_S = de_S * ux_S * uy_S, cross_EW = 0.25 * (skew_E - skew_W), cross_NS = 0.25 * (skew_N - skew_S), newt_E = e_E + de_E * PetscSqr(ux_E), newt_W = e_W + de_W * PetscSqr(ux_W), newt_N = e_N + de_N * PetscSqr(uy_N), newt_S = e_S + de_S * PetscSqr(uy_S); 508 /* interior grid points */ 509 switch (user->jtype) { 510 case JAC_BRATU: 511 /* Jacobian from p=2 */ 512 v[0] = -hxdhy; 513 col[0].j = j - 1; 514 col[0].i = i; 515 v[1] = -hydhx; 516 col[1].j = j; 517 col[1].i = i - 1; 518 v[2] = 2.0 * (hydhx + hxdhy) - sc * PetscExpScalar(u); 519 col[2].j = row.j; 520 col[2].i = row.i; 521 v[3] = -hydhx; 522 col[3].j = j; 523 col[3].i = i + 1; 524 v[4] = -hxdhy; 525 col[4].j = j + 1; 526 col[4].i = i; 527 PetscCall(MatSetValuesStencil(B, 1, &row, 5, col, v, INSERT_VALUES)); 528 break; 529 case JAC_PICARD: 530 /* Jacobian arising from Picard linearization */ 531 v[0] = -hxdhy * e_S; 532 col[0].j = j - 1; 533 col[0].i = i; 534 v[1] = -hydhx * e_W; 535 col[1].j = j; 536 col[1].i = i - 1; 537 v[2] = (e_W + e_E) * hydhx + (e_S + e_N) * hxdhy; 538 col[2].j = row.j; 539 col[2].i = row.i; 540 v[3] = -hydhx * e_E; 541 col[3].j = j; 542 col[3].i = i + 1; 543 v[4] = -hxdhy * e_N; 544 col[4].j = j + 1; 545 col[4].i = i; 546 PetscCall(MatSetValuesStencil(B, 1, &row, 5, col, v, INSERT_VALUES)); 547 break; 548 case JAC_STAR: 549 /* Full Jacobian, but only a star stencil */ 550 col[0].j = j - 1; 551 col[0].i = i; 552 col[1].j = j; 553 col[1].i = i - 1; 554 col[2].j = j; 555 col[2].i = i; 556 col[3].j = j; 557 col[3].i = i + 1; 558 col[4].j = j + 1; 559 col[4].i = i; 560 v[0] = -hxdhy * newt_S + cross_EW; 561 v[1] = -hydhx * newt_W + cross_NS; 562 v[2] = hxdhy * (newt_N + newt_S) + hydhx * (newt_E + newt_W) - sc * PetscExpScalar(u); 563 v[3] = -hydhx * newt_E - cross_NS; 564 v[4] = -hxdhy * newt_N - cross_EW; 565 PetscCall(MatSetValuesStencil(B, 1, &row, 5, col, v, INSERT_VALUES)); 566 break; 567 case JAC_NEWTON: 568 /* The Jacobian is 569 570 -div [ eta (grad u) + deta (grad u0 . grad u) grad u0 ] - (eE u0) u 571 572 */ 573 col[0].j = j - 1; 574 col[0].i = i - 1; 575 col[1].j = j - 1; 576 col[1].i = i; 577 col[2].j = j - 1; 578 col[2].i = i + 1; 579 col[3].j = j; 580 col[3].i = i - 1; 581 col[4].j = j; 582 col[4].i = i; 583 col[5].j = j; 584 col[5].i = i + 1; 585 col[6].j = j + 1; 586 col[6].i = i - 1; 587 col[7].j = j + 1; 588 col[7].i = i; 589 col[8].j = j + 1; 590 col[8].i = i + 1; 591 v[0] = -0.25 * (skew_S + skew_W); 592 v[1] = -hxdhy * newt_S + cross_EW; 593 v[2] = 0.25 * (skew_S + skew_E); 594 v[3] = -hydhx * newt_W + cross_NS; 595 v[4] = hxdhy * (newt_N + newt_S) + hydhx * (newt_E + newt_W) - sc * PetscExpScalar(u); 596 v[5] = -hydhx * newt_E - cross_NS; 597 v[6] = 0.25 * (skew_N + skew_W); 598 v[7] = -hxdhy * newt_N - cross_EW; 599 v[8] = -0.25 * (skew_N + skew_E); 600 PetscCall(MatSetValuesStencil(B, 1, &row, 9, col, v, INSERT_VALUES)); 601 break; 602 default: SETERRQ(PetscObjectComm((PetscObject)info->da), PETSC_ERR_SUP, "Jacobian type %d not implemented", user->jtype); 603 } 604 } 605 } 606 } 607 608 /* 609 Assemble matrix, using the 2-step process: 610 MatAssemblyBegin(), MatAssemblyEnd(). 611 */ 612 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 613 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 614 615 if (J != B) { 616 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 617 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 618 } 619 620 /* 621 Tell the matrix we will never add a new nonzero location to the 622 matrix. If we do, it will generate an error. 623 */ 624 if (user->jtype == JAC_NEWTON) PetscCall(PetscLogFlops(info->xm * info->ym * (131.0))); 625 PetscFunctionReturn(0); 626 } 627 628 /*********************************************************** 629 * PreCheck implementation 630 ***********************************************************/ 631 PetscErrorCode PreCheckSetFromOptions(PreCheck precheck) { 632 PetscBool flg; 633 634 PetscFunctionBeginUser; 635 PetscOptionsBegin(precheck->comm, NULL, "PreCheck Options", "none"); 636 PetscCall(PetscOptionsReal("-precheck_angle", "Angle in degrees between successive search directions necessary to activate step correction", "", precheck->angle, &precheck->angle, NULL)); 637 flg = PETSC_FALSE; 638 PetscCall(PetscOptionsBool("-precheck_monitor", "Monitor choices made by precheck routine", "", flg, &flg, NULL)); 639 if (flg) PetscCall(PetscViewerASCIIOpen(precheck->comm, "stdout", &precheck->monitor)); 640 PetscOptionsEnd(); 641 PetscFunctionReturn(0); 642 } 643 644 /* 645 Compare the direction of the current and previous step, modify the current step accordingly 646 */ 647 PetscErrorCode PreCheckFunction(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx) { 648 PreCheck precheck; 649 Vec Ylast; 650 PetscScalar dot; 651 PetscInt iter; 652 PetscReal ynorm, ylastnorm, theta, angle_radians; 653 SNES snes; 654 655 PetscFunctionBeginUser; 656 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 657 precheck = (PreCheck)ctx; 658 if (!precheck->Ylast) PetscCall(VecDuplicate(Y, &precheck->Ylast)); 659 Ylast = precheck->Ylast; 660 PetscCall(SNESGetIterationNumber(snes, &iter)); 661 if (iter < 1) { 662 PetscCall(VecCopy(Y, Ylast)); 663 *changed = PETSC_FALSE; 664 PetscFunctionReturn(0); 665 } 666 667 PetscCall(VecDot(Y, Ylast, &dot)); 668 PetscCall(VecNorm(Y, NORM_2, &ynorm)); 669 PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm)); 670 /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 671 theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0)); 672 angle_radians = precheck->angle * PETSC_PI / 180.; 673 if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 674 /* Modify the step Y */ 675 PetscReal alpha, ydiffnorm; 676 PetscCall(VecAXPY(Ylast, -1.0, Y)); 677 PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm)); 678 alpha = ylastnorm / ydiffnorm; 679 PetscCall(VecCopy(Y, Ylast)); 680 PetscCall(VecScale(Y, alpha)); 681 if (precheck->monitor) PetscCall(PetscViewerASCIIPrintf(precheck->monitor, "Angle %g degrees less than threshold %g, corrected step by alpha=%g\n", (double)(theta * 180. / PETSC_PI), (double)precheck->angle, (double)alpha)); 682 } else { 683 PetscCall(VecCopy(Y, Ylast)); 684 *changed = PETSC_FALSE; 685 if (precheck->monitor) PetscCall(PetscViewerASCIIPrintf(precheck->monitor, "Angle %g degrees exceeds threshold %g, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)precheck->angle)); 686 } 687 PetscFunctionReturn(0); 688 } 689 690 PetscErrorCode PreCheckDestroy(PreCheck *precheck) { 691 PetscFunctionBeginUser; 692 if (!*precheck) PetscFunctionReturn(0); 693 PetscCall(VecDestroy(&(*precheck)->Ylast)); 694 PetscCall(PetscViewerDestroy(&(*precheck)->monitor)); 695 PetscCall(PetscFree(*precheck)); 696 PetscFunctionReturn(0); 697 } 698 699 PetscErrorCode PreCheckCreate(MPI_Comm comm, PreCheck *precheck) { 700 PetscFunctionBeginUser; 701 PetscCall(PetscNew(precheck)); 702 703 (*precheck)->comm = comm; 704 (*precheck)->angle = 10.; /* only active if angle is less than 10 degrees */ 705 PetscFunctionReturn(0); 706 } 707 708 /* 709 Applies some sweeps on nonlinear Gauss-Seidel on each process 710 711 */ 712 PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ctx) { 713 PetscInt i, j, k, xs, ys, xm, ym, its, tot_its, sweeps, l, m; 714 PetscReal hx, hy, hxdhy, hydhx, dhx, dhy, sc; 715 PetscScalar **x, **b, bij, F, F0 = 0, J, y, u, eu; 716 PetscReal atol, rtol, stol; 717 DM da; 718 AppCtx *user = (AppCtx *)ctx; 719 Vec localX, localB; 720 DMDALocalInfo info; 721 722 PetscFunctionBeginUser; 723 PetscCall(SNESGetDM(snes, &da)); 724 PetscCall(DMDAGetLocalInfo(da, &info)); 725 726 hx = 1.0 / (PetscReal)(info.mx - 1); 727 hy = 1.0 / (PetscReal)(info.my - 1); 728 sc = hx * hy * user->lambda; 729 dhx = 1 / hx; 730 dhy = 1 / hy; 731 hxdhy = hx / hy; 732 hydhx = hy / hx; 733 734 tot_its = 0; 735 PetscCall(SNESNGSGetSweeps(snes, &sweeps)); 736 PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its)); 737 PetscCall(DMGetLocalVector(da, &localX)); 738 if (B) PetscCall(DMGetLocalVector(da, &localB)); 739 if (B) { 740 PetscCall(DMGlobalToLocalBegin(da, B, INSERT_VALUES, localB)); 741 PetscCall(DMGlobalToLocalEnd(da, B, INSERT_VALUES, localB)); 742 } 743 if (B) PetscCall(DMDAVecGetArrayRead(da, localB, &b)); 744 PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 745 PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 746 PetscCall(DMDAVecGetArray(da, localX, &x)); 747 for (l = 0; l < sweeps; l++) { 748 /* 749 Get local grid boundaries (for 2-dimensional DMDA): 750 xs, ys - starting grid indices (no ghost points) 751 xm, ym - widths of local grid (no ghost points) 752 */ 753 PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 754 for (m = 0; m < 2; m++) { 755 for (j = ys; j < ys + ym; j++) { 756 for (i = xs + (m + j) % 2; i < xs + xm; i += 2) { 757 PetscReal xx = i * hx, yy = j * hy; 758 if (B) bij = b[j][i]; 759 else bij = 0.; 760 761 if (i == 0 || j == 0 || i == info.mx - 1 || j == info.my - 1) { 762 /* boundary conditions are all zero Dirichlet */ 763 x[j][i] = 0.0 + bij; 764 } else { 765 const PetscScalar u_E = x[j][i + 1], u_W = x[j][i - 1], u_N = x[j + 1][i], u_S = x[j - 1][i]; 766 const PetscScalar uy_E = 0.25 * dhy * (x[j + 1][i] + x[j + 1][i + 1] - x[j - 1][i] - x[j - 1][i + 1]), uy_W = 0.25 * dhy * (x[j + 1][i - 1] + x[j + 1][i] - x[j - 1][i - 1] - x[j - 1][i]), ux_N = 0.25 * dhx * (x[j][i + 1] + x[j + 1][i + 1] - x[j][i - 1] - x[j + 1][i - 1]), ux_S = 0.25 * dhx * (x[j - 1][i + 1] + x[j][i + 1] - x[j - 1][i - 1] - x[j][i - 1]); 767 u = x[j][i]; 768 for (k = 0; k < its; k++) { 769 const PetscScalar ux_E = dhx * (u_E - u), ux_W = dhx * (u - u_W), uy_N = dhy * (u_N - u), uy_S = dhy * (u - u_S), e_E = eta(user, xx, yy, ux_E, uy_E), e_W = eta(user, xx, yy, ux_W, uy_W), e_N = eta(user, xx, yy, ux_N, uy_N), e_S = eta(user, xx, yy, ux_S, uy_S), de_E = deta(user, xx, yy, ux_E, uy_E), de_W = deta(user, xx, yy, ux_W, uy_W), de_N = deta(user, xx, yy, ux_N, uy_N), de_S = deta(user, xx, yy, ux_S, uy_S), newt_E = e_E + de_E * PetscSqr(ux_E), newt_W = e_W + de_W * PetscSqr(ux_W), newt_N = e_N + de_N * PetscSqr(uy_N), newt_S = e_S + de_S * PetscSqr(uy_S), uxx = -hy * (e_E * ux_E - e_W * ux_W), uyy = -hx * (e_N * uy_N - e_S * uy_S); 770 771 if (sc) eu = PetscExpScalar(u); 772 else eu = 0; 773 774 F = uxx + uyy - sc * eu - bij; 775 if (k == 0) F0 = F; 776 J = hxdhy * (newt_N + newt_S) + hydhx * (newt_E + newt_W) - sc * eu; 777 y = F / J; 778 u -= y; 779 tot_its++; 780 if (atol > PetscAbsReal(PetscRealPart(F)) || rtol * PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) || stol * PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) break; 781 } 782 x[j][i] = u; 783 } 784 } 785 } 786 } 787 /* 788 x Restore vector 789 */ 790 } 791 PetscCall(DMDAVecRestoreArray(da, localX, &x)); 792 PetscCall(DMLocalToGlobalBegin(da, localX, INSERT_VALUES, X)); 793 PetscCall(DMLocalToGlobalEnd(da, localX, INSERT_VALUES, X)); 794 PetscCall(PetscLogFlops(tot_its * (118.0))); 795 PetscCall(DMRestoreLocalVector(da, &localX)); 796 if (B) { 797 PetscCall(DMDAVecRestoreArrayRead(da, localB, &b)); 798 PetscCall(DMRestoreLocalVector(da, &localB)); 799 } 800 PetscFunctionReturn(0); 801 } 802 803 /*TEST 804 805 test: 806 nsize: 2 807 args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype NEWTON 808 requires: !single 809 810 test: 811 suffix: 2 812 nsize: 2 813 args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -precheck 1 814 requires: !single 815 816 test: 817 suffix: 3 818 nsize: 2 819 args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -picard -precheck 1 -p 3 820 requires: !single 821 822 test: 823 suffix: 3_star 824 nsize: 2 825 args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -jtype PICARD -picard -precheck 1 -p 3 -alloc_star 826 output_file: output/ex15_3.out 827 requires: !single 828 829 test: 830 suffix: 4 831 args: -snes_monitor_short -snes_type newtonls -npc_snes_type ngs -snes_npc_side left -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -ksp_monitor_short -pc_type none 832 requires: !single 833 834 test: 835 suffix: lag_jac 836 nsize: 4 837 args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 6.0 -lambda 0 -jtype NEWTON -snes_type ngmres -npc_snes_type newtonls -npc_snes_lag_jacobian 5 -npc_pc_type asm -npc_ksp_converged_reason -npc_snes_lag_jacobian_persists 838 requires: !single 839 840 test: 841 suffix: lag_pc 842 nsize: 4 843 args: -snes_monitor_short -da_grid_x 20 -da_grid_y 20 -p 6.0 -lambda 0 -jtype NEWTON -snes_type ngmres -npc_snes_type newtonls -npc_snes_lag_preconditioner 5 -npc_pc_type asm -npc_ksp_converged_reason -npc_snes_lag_preconditioner_persists 844 requires: !single 845 846 test: 847 suffix: nleqerr 848 args: -snes_monitor_short -snes_type newtonls -da_grid_x 20 -da_grid_y 20 -p 1.3 -lambda 1 -snes_linesearch_monitor -pc_type lu -snes_linesearch_type nleqerr 849 requires: !single 850 851 test: 852 suffix: mf 853 args: -snes_monitor_short -pc_type lu -da_refine 4 -p 3 -ksp_rtol 1.e-12 -snes_mf_operator 854 requires: !single 855 856 test: 857 suffix: mf_picard 858 args: -snes_monitor_short -pc_type lu -da_refine 4 -p 3 -ksp_rtol 1.e-12 -snes_mf_operator -picard 859 requires: !single 860 output_file: output/ex15_mf.out 861 862 test: 863 suffix: fd_picard 864 args: -snes_monitor_short -pc_type lu -da_refine 2 -p 3 -ksp_rtol 1.e-12 -snes_fd -picard 865 requires: !single 866 867 test: 868 suffix: fd_color_picard 869 args: -snes_monitor_short -pc_type lu -da_refine 4 -p 3 -ksp_rtol 1.e-12 -snes_fd_color -picard 870 requires: !single 871 output_file: output/ex15_mf.out 872 873 TEST*/ 874