1 static char help[] = "Newton methods to solve u'' + u^{2} = f in parallel.\n\ 2 This example employs a user-defined monitoring routine and optionally a user-defined\n\ 3 routine to check candidate iterates produced by line search routines.\n\ 4 The command line options include:\n\ 5 -pre_check_iterates : activate checking of iterates\n\ 6 -post_check_iterates : activate checking of iterates\n\ 7 -check_tol <tol>: set tolerance for iterate checking\n\ 8 -user_precond : activate a (trivial) user-defined preconditioner\n\n"; 9 10 /* 11 Include "petscdm.h" so that we can use data management objects (DMs) 12 Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 13 Include "petscsnes.h" so that we can use SNES solvers. Note that this 14 file automatically includes: 15 petscsys.h - base PETSc routines 16 petscvec.h - vectors 17 petscmat.h - matrices 18 petscis.h - index sets 19 petscksp.h - Krylov subspace methods 20 petscviewer.h - viewers 21 petscpc.h - preconditioners 22 petscksp.h - linear solvers 23 */ 24 25 #include <petscdm.h> 26 #include <petscdmda.h> 27 #include <petscsnes.h> 28 29 /* 30 User-defined routines. 31 */ 32 PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *); 33 PetscErrorCode FormFunction(SNES, Vec, Vec, void *); 34 PetscErrorCode FormInitialGuess(Vec); 35 PetscErrorCode Monitor(SNES, PetscInt, PetscReal, void *); 36 PetscErrorCode PreCheck(SNESLineSearch, Vec, Vec, PetscBool *, void *); 37 PetscErrorCode PostCheck(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *); 38 PetscErrorCode PostSetSubKSP(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *); 39 PetscErrorCode MatrixFreePreconditioner(PC, Vec, Vec); 40 41 /* 42 User-defined application context 43 */ 44 typedef struct { 45 DM da; /* distributed array */ 46 Vec F; /* right-hand-side of PDE */ 47 PetscMPIInt rank; /* rank of processor */ 48 PetscMPIInt size; /* size of communicator */ 49 PetscReal h; /* mesh spacing */ 50 PetscBool sjerr; /* if or not to test jacobian domain error */ 51 } ApplicationCtx; 52 53 /* 54 User-defined context for monitoring 55 */ 56 typedef struct { 57 PetscViewer viewer; 58 } MonitorCtx; 59 60 /* 61 User-defined context for checking candidate iterates that are 62 determined by line search methods 63 */ 64 typedef struct { 65 Vec last_step; /* previous iterate */ 66 PetscReal tolerance; /* tolerance for changes between successive iterates */ 67 ApplicationCtx *user; 68 } StepCheckCtx; 69 70 typedef struct { 71 PetscInt its0; /* num of prevous outer KSP iterations */ 72 } SetSubKSPCtx; 73 74 int main(int argc, char **argv) { 75 SNES snes; /* SNES context */ 76 SNESLineSearch linesearch; /* SNESLineSearch context */ 77 Mat J; /* Jacobian matrix */ 78 ApplicationCtx ctx; /* user-defined context */ 79 Vec x, r, U, F; /* vectors */ 80 MonitorCtx monP; /* monitoring context */ 81 StepCheckCtx checkP; /* step-checking context */ 82 SetSubKSPCtx checkP1; 83 PetscBool pre_check, post_check, post_setsubksp; /* flag indicating whether we're checking candidate iterates */ 84 PetscScalar xp, *FF, *UU, none = -1.0; 85 PetscInt its, N = 5, i, maxit, maxf, xs, xm; 86 PetscReal abstol, rtol, stol, norm; 87 PetscBool flg, viewinitial = PETSC_FALSE; 88 89 PetscFunctionBeginUser; 90 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 91 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &ctx.rank)); 92 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &ctx.size)); 93 PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &N, NULL)); 94 ctx.h = 1.0 / (N - 1); 95 ctx.sjerr = PETSC_FALSE; 96 PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_jacobian_domain_error", &ctx.sjerr, NULL)); 97 PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_initial", &viewinitial, NULL)); 98 99 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 100 Create nonlinear solver context 101 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 102 103 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 104 105 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 106 Create vector data structures; set function evaluation routine 107 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 108 109 /* 110 Create distributed array (DMDA) to manage parallel grid and vectors 111 */ 112 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, N, 1, 1, NULL, &ctx.da)); 113 PetscCall(DMSetFromOptions(ctx.da)); 114 PetscCall(DMSetUp(ctx.da)); 115 116 /* 117 Extract global and local vectors from DMDA; then duplicate for remaining 118 vectors that are the same types 119 */ 120 PetscCall(DMCreateGlobalVector(ctx.da, &x)); 121 PetscCall(VecDuplicate(x, &r)); 122 PetscCall(VecDuplicate(x, &F)); 123 ctx.F = F; 124 PetscCall(VecDuplicate(x, &U)); 125 126 /* 127 Set function evaluation routine and vector. Whenever the nonlinear 128 solver needs to compute the nonlinear function, it will call this 129 routine. 130 - Note that the final routine argument is the user-defined 131 context that provides application-specific data for the 132 function evaluation routine. 133 */ 134 PetscCall(SNESSetFunction(snes, r, FormFunction, &ctx)); 135 136 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 137 Create matrix data structure; set Jacobian evaluation routine 138 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 139 140 PetscCall(MatCreate(PETSC_COMM_WORLD, &J)); 141 PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, N, N)); 142 PetscCall(MatSetFromOptions(J)); 143 PetscCall(MatSeqAIJSetPreallocation(J, 3, NULL)); 144 PetscCall(MatMPIAIJSetPreallocation(J, 3, NULL, 3, NULL)); 145 146 /* 147 Set Jacobian matrix data structure and default Jacobian evaluation 148 routine. Whenever the nonlinear solver needs to compute the 149 Jacobian matrix, it will call this routine. 150 - Note that the final routine argument is the user-defined 151 context that provides application-specific data for the 152 Jacobian evaluation routine. 153 */ 154 PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &ctx)); 155 156 /* 157 Optionally allow user-provided preconditioner 158 */ 159 PetscCall(PetscOptionsHasName(NULL, NULL, "-user_precond", &flg)); 160 if (flg) { 161 KSP ksp; 162 PC pc; 163 PetscCall(SNESGetKSP(snes, &ksp)); 164 PetscCall(KSPGetPC(ksp, &pc)); 165 PetscCall(PCSetType(pc, PCSHELL)); 166 PetscCall(PCShellSetApply(pc, MatrixFreePreconditioner)); 167 } 168 169 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 170 Customize nonlinear solver; set runtime options 171 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 172 173 /* 174 Set an optional user-defined monitoring routine 175 */ 176 PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, 0, 0, 0, 0, 400, 400, &monP.viewer)); 177 PetscCall(SNESMonitorSet(snes, Monitor, &monP, 0)); 178 179 /* 180 Set names for some vectors to facilitate monitoring (optional) 181 */ 182 PetscCall(PetscObjectSetName((PetscObject)x, "Approximate Solution")); 183 PetscCall(PetscObjectSetName((PetscObject)U, "Exact Solution")); 184 185 /* 186 Set SNES/KSP/KSP/PC runtime options, e.g., 187 -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 188 */ 189 PetscCall(SNESSetFromOptions(snes)); 190 191 /* 192 Set an optional user-defined routine to check the validity of candidate 193 iterates that are determined by line search methods 194 */ 195 PetscCall(SNESGetLineSearch(snes, &linesearch)); 196 PetscCall(PetscOptionsHasName(NULL, NULL, "-post_check_iterates", &post_check)); 197 198 if (post_check) { 199 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activating post step checking routine\n")); 200 PetscCall(SNESLineSearchSetPostCheck(linesearch, PostCheck, &checkP)); 201 PetscCall(VecDuplicate(x, &(checkP.last_step))); 202 203 checkP.tolerance = 1.0; 204 checkP.user = &ctx; 205 206 PetscCall(PetscOptionsGetReal(NULL, NULL, "-check_tol", &checkP.tolerance, NULL)); 207 } 208 209 PetscCall(PetscOptionsHasName(NULL, NULL, "-post_setsubksp", &post_setsubksp)); 210 if (post_setsubksp) { 211 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activating post setsubksp\n")); 212 PetscCall(SNESLineSearchSetPostCheck(linesearch, PostSetSubKSP, &checkP1)); 213 } 214 215 PetscCall(PetscOptionsHasName(NULL, NULL, "-pre_check_iterates", &pre_check)); 216 if (pre_check) { 217 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activating pre step checking routine\n")); 218 PetscCall(SNESLineSearchSetPreCheck(linesearch, PreCheck, &checkP)); 219 } 220 221 /* 222 Print parameters used for convergence testing (optional) ... just 223 to demonstrate this routine; this information is also printed with 224 the option -snes_view 225 */ 226 PetscCall(SNESGetTolerances(snes, &abstol, &rtol, &stol, &maxit, &maxf)); 227 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "atol=%g, rtol=%g, stol=%g, maxit=%" PetscInt_FMT ", maxf=%" PetscInt_FMT "\n", (double)abstol, (double)rtol, (double)stol, maxit, maxf)); 228 229 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 230 Initialize application: 231 Store right-hand-side of PDE and exact solution 232 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 233 234 /* 235 Get local grid boundaries (for 1-dimensional DMDA): 236 xs, xm - starting grid index, width of local grid (no ghost points) 237 */ 238 PetscCall(DMDAGetCorners(ctx.da, &xs, NULL, NULL, &xm, NULL, NULL)); 239 240 /* 241 Get pointers to vector data 242 */ 243 PetscCall(DMDAVecGetArray(ctx.da, F, &FF)); 244 PetscCall(DMDAVecGetArray(ctx.da, U, &UU)); 245 246 /* 247 Compute local vector entries 248 */ 249 xp = ctx.h * xs; 250 for (i = xs; i < xs + xm; i++) { 251 FF[i] = 6.0 * xp + PetscPowScalar(xp + 1.e-12, 6.0); /* +1.e-12 is to prevent 0^6 */ 252 UU[i] = xp * xp * xp; 253 xp += ctx.h; 254 } 255 256 /* 257 Restore vectors 258 */ 259 PetscCall(DMDAVecRestoreArray(ctx.da, F, &FF)); 260 PetscCall(DMDAVecRestoreArray(ctx.da, U, &UU)); 261 if (viewinitial) { 262 PetscCall(VecView(U, 0)); 263 PetscCall(VecView(F, 0)); 264 } 265 266 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 267 Evaluate initial guess; then solve nonlinear system 268 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 269 270 /* 271 Note: The user should initialize the vector, x, with the initial guess 272 for the nonlinear solver prior to calling SNESSolve(). In particular, 273 to employ an initial guess of zero, the user should explicitly set 274 this vector to zero by calling VecSet(). 275 */ 276 PetscCall(FormInitialGuess(x)); 277 PetscCall(SNESSolve(snes, NULL, x)); 278 PetscCall(SNESGetIterationNumber(snes, &its)); 279 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 280 281 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 282 Check solution and clean up 283 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 284 285 /* 286 Check the error 287 */ 288 PetscCall(VecAXPY(x, none, U)); 289 PetscCall(VecNorm(x, NORM_2, &norm)); 290 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error %g Iterations %" PetscInt_FMT "\n", (double)norm, its)); 291 if (ctx.sjerr) { 292 SNESType snestype; 293 PetscCall(SNESGetType(snes, &snestype)); 294 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SNES Type %s\n", snestype)); 295 } 296 297 /* 298 Free work space. All PETSc objects should be destroyed when they 299 are no longer needed. 300 */ 301 PetscCall(PetscViewerDestroy(&monP.viewer)); 302 if (post_check) PetscCall(VecDestroy(&checkP.last_step)); 303 PetscCall(VecDestroy(&x)); 304 PetscCall(VecDestroy(&r)); 305 PetscCall(VecDestroy(&U)); 306 PetscCall(VecDestroy(&F)); 307 PetscCall(MatDestroy(&J)); 308 PetscCall(SNESDestroy(&snes)); 309 PetscCall(DMDestroy(&ctx.da)); 310 PetscCall(PetscFinalize()); 311 return 0; 312 } 313 314 /* ------------------------------------------------------------------- */ 315 /* 316 FormInitialGuess - Computes initial guess. 317 318 Input/Output Parameter: 319 . x - the solution vector 320 */ 321 PetscErrorCode FormInitialGuess(Vec x) { 322 PetscScalar pfive = .50; 323 324 PetscFunctionBeginUser; 325 PetscCall(VecSet(x, pfive)); 326 PetscFunctionReturn(0); 327 } 328 329 /* ------------------------------------------------------------------- */ 330 /* 331 FormFunction - Evaluates nonlinear function, F(x). 332 333 Input Parameters: 334 . snes - the SNES context 335 . x - input vector 336 . ctx - optional user-defined context, as set by SNESSetFunction() 337 338 Output Parameter: 339 . f - function vector 340 341 Note: 342 The user-defined context can contain any application-specific 343 data needed for the function evaluation. 344 */ 345 PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, void *ctx) { 346 ApplicationCtx *user = (ApplicationCtx *)ctx; 347 DM da = user->da; 348 PetscScalar *ff, d; 349 const PetscScalar *xx, *FF; 350 PetscInt i, M, xs, xm; 351 Vec xlocal; 352 353 PetscFunctionBeginUser; 354 PetscCall(DMGetLocalVector(da, &xlocal)); 355 /* 356 Scatter ghost points to local vector, using the 2-step process 357 DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 358 By placing code between these two statements, computations can 359 be done while messages are in transition. 360 */ 361 PetscCall(DMGlobalToLocalBegin(da, x, INSERT_VALUES, xlocal)); 362 PetscCall(DMGlobalToLocalEnd(da, x, INSERT_VALUES, xlocal)); 363 364 /* 365 Get pointers to vector data. 366 - The vector xlocal includes ghost point; the vectors x and f do 367 NOT include ghost points. 368 - Using DMDAVecGetArray() allows accessing the values using global ordering 369 */ 370 PetscCall(DMDAVecGetArrayRead(da, xlocal, (void *)&xx)); 371 PetscCall(DMDAVecGetArray(da, f, &ff)); 372 PetscCall(DMDAVecGetArrayRead(da, user->F, (void *)&FF)); 373 374 /* 375 Get local grid boundaries (for 1-dimensional DMDA): 376 xs, xm - starting grid index, width of local grid (no ghost points) 377 */ 378 PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 379 PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 380 381 /* 382 Set function values for boundary points; define local interior grid point range: 383 xsi - starting interior grid index 384 xei - ending interior grid index 385 */ 386 if (xs == 0) { /* left boundary */ 387 ff[0] = xx[0]; 388 xs++; 389 xm--; 390 } 391 if (xs + xm == M) { /* right boundary */ 392 ff[xs + xm - 1] = xx[xs + xm - 1] - 1.0; 393 xm--; 394 } 395 396 /* 397 Compute function over locally owned part of the grid (interior points only) 398 */ 399 d = 1.0 / (user->h * user->h); 400 for (i = xs; i < xs + xm; i++) ff[i] = d * (xx[i - 1] - 2.0 * xx[i] + xx[i + 1]) + xx[i] * xx[i] - FF[i]; 401 402 /* 403 Restore vectors 404 */ 405 PetscCall(DMDAVecRestoreArrayRead(da, xlocal, (void *)&xx)); 406 PetscCall(DMDAVecRestoreArray(da, f, &ff)); 407 PetscCall(DMDAVecRestoreArrayRead(da, user->F, (void *)&FF)); 408 PetscCall(DMRestoreLocalVector(da, &xlocal)); 409 PetscFunctionReturn(0); 410 } 411 412 /* ------------------------------------------------------------------- */ 413 /* 414 FormJacobian - Evaluates Jacobian matrix. 415 416 Input Parameters: 417 . snes - the SNES context 418 . x - input vector 419 . dummy - optional user-defined context (not used here) 420 421 Output Parameters: 422 . jac - Jacobian matrix 423 . B - optionally different preconditioning matrix 424 . flag - flag indicating matrix structure 425 */ 426 PetscErrorCode FormJacobian(SNES snes, Vec x, Mat jac, Mat B, void *ctx) { 427 ApplicationCtx *user = (ApplicationCtx *)ctx; 428 PetscScalar *xx, d, A[3]; 429 PetscInt i, j[3], M, xs, xm; 430 DM da = user->da; 431 432 PetscFunctionBeginUser; 433 if (user->sjerr) { 434 PetscCall(SNESSetJacobianDomainError(snes)); 435 PetscFunctionReturn(0); 436 } 437 /* 438 Get pointer to vector data 439 */ 440 PetscCall(DMDAVecGetArrayRead(da, x, &xx)); 441 PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 442 443 /* 444 Get range of locally owned matrix 445 */ 446 PetscCall(DMDAGetInfo(da, NULL, &M, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 447 448 /* 449 Determine starting and ending local indices for interior grid points. 450 Set Jacobian entries for boundary points. 451 */ 452 453 if (xs == 0) { /* left boundary */ 454 i = 0; 455 A[0] = 1.0; 456 457 PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES)); 458 xs++; 459 xm--; 460 } 461 if (xs + xm == M) { /* right boundary */ 462 i = M - 1; 463 A[0] = 1.0; 464 PetscCall(MatSetValues(jac, 1, &i, 1, &i, A, INSERT_VALUES)); 465 xm--; 466 } 467 468 /* 469 Interior grid points 470 - Note that in this case we set all elements for a particular 471 row at once. 472 */ 473 d = 1.0 / (user->h * user->h); 474 for (i = xs; i < xs + xm; i++) { 475 j[0] = i - 1; 476 j[1] = i; 477 j[2] = i + 1; 478 A[0] = A[2] = d; 479 A[1] = -2.0 * d + 2.0 * xx[i]; 480 PetscCall(MatSetValues(jac, 1, &i, 3, j, A, INSERT_VALUES)); 481 } 482 483 /* 484 Assemble matrix, using the 2-step process: 485 MatAssemblyBegin(), MatAssemblyEnd(). 486 By placing code between these two statements, computations can be 487 done while messages are in transition. 488 489 Also, restore vector. 490 */ 491 492 PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 493 PetscCall(DMDAVecRestoreArrayRead(da, x, &xx)); 494 PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 495 496 PetscFunctionReturn(0); 497 } 498 499 /* ------------------------------------------------------------------- */ 500 /* 501 Monitor - Optional user-defined monitoring routine that views the 502 current iterate with an x-window plot. Set by SNESMonitorSet(). 503 504 Input Parameters: 505 snes - the SNES context 506 its - iteration number 507 norm - 2-norm function value (may be estimated) 508 ctx - optional user-defined context for private data for the 509 monitor routine, as set by SNESMonitorSet() 510 511 Note: 512 See the manpage for PetscViewerDrawOpen() for useful runtime options, 513 such as -nox to deactivate all x-window output. 514 */ 515 PetscErrorCode Monitor(SNES snes, PetscInt its, PetscReal fnorm, void *ctx) { 516 MonitorCtx *monP = (MonitorCtx *)ctx; 517 Vec x; 518 519 PetscFunctionBeginUser; 520 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "iter = %" PetscInt_FMT ",SNES Function norm %g\n", its, (double)fnorm)); 521 PetscCall(SNESGetSolution(snes, &x)); 522 PetscCall(VecView(x, monP->viewer)); 523 PetscFunctionReturn(0); 524 } 525 526 /* ------------------------------------------------------------------- */ 527 /* 528 PreCheck - Optional user-defined routine that checks the validity of 529 candidate steps of a line search method. Set by SNESLineSearchSetPreCheck(). 530 531 Input Parameters: 532 snes - the SNES context 533 xcurrent - current solution 534 y - search direction and length 535 536 Output Parameters: 537 y - proposed step (search direction and length) (possibly changed) 538 changed_y - tells if the step has changed or not 539 */ 540 PetscErrorCode PreCheck(SNESLineSearch linesearch, Vec xcurrent, Vec y, PetscBool *changed_y, void *ctx) { 541 PetscFunctionBeginUser; 542 *changed_y = PETSC_FALSE; 543 PetscFunctionReturn(0); 544 } 545 546 /* ------------------------------------------------------------------- */ 547 /* 548 PostCheck - Optional user-defined routine that checks the validity of 549 candidate steps of a line search method. Set by SNESLineSearchSetPostCheck(). 550 551 Input Parameters: 552 snes - the SNES context 553 ctx - optional user-defined context for private data for the 554 monitor routine, as set by SNESLineSearchSetPostCheck() 555 xcurrent - current solution 556 y - search direction and length 557 x - the new candidate iterate 558 559 Output Parameters: 560 y - proposed step (search direction and length) (possibly changed) 561 x - current iterate (possibly modified) 562 563 */ 564 PetscErrorCode PostCheck(SNESLineSearch linesearch, Vec xcurrent, Vec y, Vec x, PetscBool *changed_y, PetscBool *changed_x, void *ctx) { 565 PetscInt i, iter, xs, xm; 566 StepCheckCtx *check; 567 ApplicationCtx *user; 568 PetscScalar *xa, *xa_last, tmp; 569 PetscReal rdiff; 570 DM da; 571 SNES snes; 572 573 PetscFunctionBeginUser; 574 *changed_x = PETSC_FALSE; 575 *changed_y = PETSC_FALSE; 576 577 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 578 check = (StepCheckCtx *)ctx; 579 user = check->user; 580 PetscCall(SNESGetIterationNumber(snes, &iter)); 581 582 /* iteration 1 indicates we are working on the second iteration */ 583 if (iter > 0) { 584 da = user->da; 585 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Checking candidate step at iteration %" PetscInt_FMT " with tolerance %g\n", iter, (double)check->tolerance)); 586 587 /* Access local array data */ 588 PetscCall(DMDAVecGetArray(da, check->last_step, &xa_last)); 589 PetscCall(DMDAVecGetArray(da, x, &xa)); 590 PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 591 592 /* 593 If we fail the user-defined check for validity of the candidate iterate, 594 then modify the iterate as we like. (Note that the particular modification 595 below is intended simply to demonstrate how to manipulate this data, not 596 as a meaningful or appropriate choice.) 597 */ 598 for (i = xs; i < xs + xm; i++) { 599 if (!PetscAbsScalar(xa[i])) rdiff = 2 * check->tolerance; 600 else rdiff = PetscAbsScalar((xa[i] - xa_last[i]) / xa[i]); 601 if (rdiff > check->tolerance) { 602 tmp = xa[i]; 603 xa[i] = .5 * (xa[i] + xa_last[i]); 604 *changed_x = PETSC_TRUE; 605 PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Altering entry %" PetscInt_FMT ": x=%g, x_last=%g, diff=%g, x_new=%g\n", i, (double)PetscAbsScalar(tmp), (double)PetscAbsScalar(xa_last[i]), (double)rdiff, (double)PetscAbsScalar(xa[i]))); 606 } 607 } 608 PetscCall(DMDAVecRestoreArray(da, check->last_step, &xa_last)); 609 PetscCall(DMDAVecRestoreArray(da, x, &xa)); 610 } 611 PetscCall(VecCopy(x, check->last_step)); 612 PetscFunctionReturn(0); 613 } 614 615 /* ------------------------------------------------------------------- */ 616 /* 617 PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used 618 e.g, 619 mpiexec -n 8 ./ex3 -nox -n 10000 -ksp_type fgmres -pc_type bjacobi -pc_bjacobi_blocks 4 -sub_ksp_type gmres -sub_ksp_max_it 3 -post_setsubksp -sub_ksp_rtol 1.e-16 620 Set by SNESLineSearchSetPostCheck(). 621 622 Input Parameters: 623 linesearch - the LineSearch context 624 xcurrent - current solution 625 y - search direction and length 626 x - the new candidate iterate 627 628 Output Parameters: 629 y - proposed step (search direction and length) (possibly changed) 630 x - current iterate (possibly modified) 631 632 */ 633 PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch, Vec xcurrent, Vec y, Vec x, PetscBool *changed_y, PetscBool *changed_x, void *ctx) { 634 SetSubKSPCtx *check; 635 PetscInt iter, its, sub_its, maxit; 636 KSP ksp, sub_ksp, *sub_ksps; 637 PC pc; 638 PetscReal ksp_ratio; 639 SNES snes; 640 641 PetscFunctionBeginUser; 642 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 643 check = (SetSubKSPCtx *)ctx; 644 PetscCall(SNESGetIterationNumber(snes, &iter)); 645 PetscCall(SNESGetKSP(snes, &ksp)); 646 PetscCall(KSPGetPC(ksp, &pc)); 647 PetscCall(PCBJacobiGetSubKSP(pc, NULL, NULL, &sub_ksps)); 648 sub_ksp = sub_ksps[0]; 649 PetscCall(KSPGetIterationNumber(ksp, &its)); /* outer KSP iteration number */ 650 PetscCall(KSPGetIterationNumber(sub_ksp, &sub_its)); /* inner KSP iteration number */ 651 652 if (iter) { 653 PetscCall(PetscPrintf(PETSC_COMM_WORLD, " ...PostCheck snes iteration %" PetscInt_FMT ", ksp_it %" PetscInt_FMT " %" PetscInt_FMT ", subksp_it %" PetscInt_FMT "\n", iter, check->its0, its, sub_its)); 654 ksp_ratio = ((PetscReal)(its)) / check->its0; 655 maxit = (PetscInt)(ksp_ratio * sub_its + 0.5); 656 if (maxit < 2) maxit = 2; 657 PetscCall(KSPSetTolerances(sub_ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, maxit)); 658 PetscCall(PetscPrintf(PETSC_COMM_WORLD, " ...ksp_ratio %g, new maxit %" PetscInt_FMT "\n\n", (double)ksp_ratio, maxit)); 659 } 660 check->its0 = its; /* save current outer KSP iteration number */ 661 PetscFunctionReturn(0); 662 } 663 664 /* ------------------------------------------------------------------- */ 665 /* 666 MatrixFreePreconditioner - This routine demonstrates the use of a 667 user-provided preconditioner. This code implements just the null 668 preconditioner, which of course is not recommended for general use. 669 670 Input Parameters: 671 + pc - preconditioner 672 - x - input vector 673 674 Output Parameter: 675 . y - preconditioned vector 676 */ 677 PetscErrorCode MatrixFreePreconditioner(PC pc, Vec x, Vec y) { 678 PetscCall(VecCopy(x, y)); 679 return 0; 680 } 681 682 /*TEST 683 684 test: 685 args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 686 687 test: 688 suffix: 2 689 nsize: 3 690 args: -nox -pc_type asm -mat_type mpiaij -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 691 692 test: 693 suffix: 3 694 nsize: 2 695 args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 696 697 test: 698 suffix: 4 699 args: -nox -pre_check_iterates -post_check_iterates 700 701 test: 702 suffix: 5 703 requires: double !complex !single 704 nsize: 2 705 args: -nox -snes_test_jacobian -snes_test_jacobian_view 706 707 test: 708 suffix: 6 709 requires: double !complex !single 710 nsize: 4 711 args: -test_jacobian_domain_error -snes_converged_reason -snes_check_jacobian_domain_error 1 712 713 test: 714 suffix: 7 715 requires: double !complex !single 716 nsize: 4 717 args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontr -snes_check_jacobian_domain_error 1 718 719 test: 720 suffix: 8 721 requires: double !complex !single 722 nsize: 4 723 args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonrsls -snes_check_jacobian_domain_error 1 724 725 test: 726 suffix: 9 727 requires: double !complex !single 728 nsize: 4 729 args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonssls -snes_check_jacobian_domain_error 1 730 731 test: 732 suffix: 10 733 requires: double !complex !single 734 nsize: 4 735 args: -test_jacobian_domain_error -snes_converged_reason -snes_type qn -snes_qn_scale_type jacobian -snes_check_jacobian_domain_error 1 736 737 test: 738 suffix: 11 739 requires: double !complex !single 740 nsize: 4 741 args: -test_jacobian_domain_error -snes_converged_reason -snes_type ms -snes_ms_type m62 -snes_ms_damping 0.9 -snes_check_jacobian_domain_error 1 742 743 test: 744 suffix: 12 745 args: -view_initial 746 filter: grep -v "type:" 747 748 test: 749 suffix: 13 750 requires: double !complex !single 751 nsize: 4 752 args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontrdc -snes_check_jacobian_domain_error 1 753 754 TEST*/ 755