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