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 /*T 11 Concepts: SNES^basic parallel example 12 Concepts: SNES^setting a user-defined monitoring routine 13 Processors: n 14 T*/ 15 16 /* 17 Include "petscdm.h" so that we can use data management objects (DMs) 18 Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 19 Include "petscsnes.h" so that we can use SNES solvers. Note that this 20 file automatically includes: 21 petscsys.h - base PETSc routines 22 petscvec.h - vectors 23 petscmat.h - matrices 24 petscis.h - index sets 25 petscksp.h - Krylov subspace methods 26 petscviewer.h - viewers 27 petscpc.h - preconditioners 28 petscksp.h - linear solvers 29 */ 30 31 #include <petscdm.h> 32 #include <petscdmda.h> 33 #include <petscsnes.h> 34 35 /* 36 User-defined routines. 37 */ 38 PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 39 PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 40 PetscErrorCode FormInitialGuess(Vec); 41 PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*); 42 PetscErrorCode PreCheck(SNESLineSearch,Vec,Vec,PetscBool*,void*); 43 PetscErrorCode PostCheck(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 44 PetscErrorCode PostSetSubKSP(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 45 PetscErrorCode MatrixFreePreconditioner(PC,Vec,Vec); 46 47 /* 48 User-defined application context 49 */ 50 typedef struct { 51 DM da; /* distributed array */ 52 Vec F; /* right-hand-side of PDE */ 53 PetscMPIInt rank; /* rank of processor */ 54 PetscMPIInt size; /* size of communicator */ 55 PetscReal h; /* mesh spacing */ 56 PetscBool sjerr; /* if or not to test jacobian domain error */ 57 } ApplicationCtx; 58 59 /* 60 User-defined context for monitoring 61 */ 62 typedef struct { 63 PetscViewer viewer; 64 } MonitorCtx; 65 66 /* 67 User-defined context for checking candidate iterates that are 68 determined by line search methods 69 */ 70 typedef struct { 71 Vec last_step; /* previous iterate */ 72 PetscReal tolerance; /* tolerance for changes between successive iterates */ 73 ApplicationCtx *user; 74 } StepCheckCtx; 75 76 typedef struct { 77 PetscInt its0; /* num of prevous outer KSP iterations */ 78 } SetSubKSPCtx; 79 80 int main(int argc,char **argv) 81 { 82 SNES snes; /* SNES context */ 83 SNESLineSearch linesearch; /* SNESLineSearch context */ 84 Mat J; /* Jacobian matrix */ 85 ApplicationCtx ctx; /* user-defined context */ 86 Vec x,r,U,F; /* vectors */ 87 MonitorCtx monP; /* monitoring context */ 88 StepCheckCtx checkP; /* step-checking context */ 89 SetSubKSPCtx checkP1; 90 PetscBool pre_check,post_check,post_setsubksp; /* flag indicating whether we're checking candidate iterates */ 91 PetscScalar xp,*FF,*UU,none = -1.0; 92 PetscInt its,N = 5,i,maxit,maxf,xs,xm; 93 PetscReal abstol,rtol,stol,norm; 94 PetscBool flg,viewinitial = PETSC_FALSE; 95 96 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 97 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&ctx.rank)); 98 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&ctx.size)); 99 PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL)); 100 ctx.h = 1.0/(N-1); 101 ctx.sjerr = PETSC_FALSE; 102 PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_jacobian_domain_error",&ctx.sjerr,NULL)); 103 PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_initial",&viewinitial,NULL)); 104 105 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 106 Create nonlinear solver context 107 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 108 109 PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 110 111 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 112 Create vector data structures; set function evaluation routine 113 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 114 115 /* 116 Create distributed array (DMDA) to manage parallel grid and vectors 117 */ 118 PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,1,1,NULL,&ctx.da)); 119 PetscCall(DMSetFromOptions(ctx.da)); 120 PetscCall(DMSetUp(ctx.da)); 121 122 /* 123 Extract global and local vectors from DMDA; then duplicate for remaining 124 vectors that are the same types 125 */ 126 PetscCall(DMCreateGlobalVector(ctx.da,&x)); 127 PetscCall(VecDuplicate(x,&r)); 128 PetscCall(VecDuplicate(x,&F)); 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=%D, maxf=%D\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 = %D\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 %D\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 *xx,*ff,*FF,d; 356 PetscInt i,M,xs,xm; 357 Vec xlocal; 358 359 PetscFunctionBeginUser; 360 PetscCall(DMGetLocalVector(da,&xlocal)); 361 /* 362 Scatter ghost points to local vector, using the 2-step process 363 DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 364 By placing code between these two statements, computations can 365 be done while messages are in transition. 366 */ 367 PetscCall(DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal)); 368 PetscCall(DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal)); 369 370 /* 371 Get pointers to vector data. 372 - The vector xlocal includes ghost point; the vectors x and f do 373 NOT include ghost points. 374 - Using DMDAVecGetArray() allows accessing the values using global ordering 375 */ 376 PetscCall(DMDAVecGetArray(da,xlocal,&xx)); 377 PetscCall(DMDAVecGetArray(da,f,&ff)); 378 PetscCall(DMDAVecGetArray(da,user->F,&FF)); 379 380 /* 381 Get local grid boundaries (for 1-dimensional DMDA): 382 xs, xm - starting grid index, width of local grid (no ghost points) 383 */ 384 PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 385 PetscCall(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 386 387 /* 388 Set function values for boundary points; define local interior grid point range: 389 xsi - starting interior grid index 390 xei - ending interior grid index 391 */ 392 if (xs == 0) { /* left boundary */ 393 ff[0] = xx[0]; 394 xs++;xm--; 395 } 396 if (xs+xm == M) { /* right boundary */ 397 ff[xs+xm-1] = xx[xs+xm-1] - 1.0; 398 xm--; 399 } 400 401 /* 402 Compute function over locally owned part of the grid (interior points only) 403 */ 404 d = 1.0/(user->h*user->h); 405 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]; 406 407 /* 408 Restore vectors 409 */ 410 PetscCall(DMDAVecRestoreArray(da,xlocal,&xx)); 411 PetscCall(DMDAVecRestoreArray(da,f,&ff)); 412 PetscCall(DMDAVecRestoreArray(da,user->F,&FF)); 413 PetscCall(DMRestoreLocalVector(da,&xlocal)); 414 PetscFunctionReturn(0); 415 } 416 417 /* ------------------------------------------------------------------- */ 418 /* 419 FormJacobian - Evaluates Jacobian matrix. 420 421 Input Parameters: 422 . snes - the SNES context 423 . x - input vector 424 . dummy - optional user-defined context (not used here) 425 426 Output Parameters: 427 . jac - Jacobian matrix 428 . B - optionally different preconditioning matrix 429 . flag - flag indicating matrix structure 430 */ 431 PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *ctx) 432 { 433 ApplicationCtx *user = (ApplicationCtx*) ctx; 434 PetscScalar *xx,d,A[3]; 435 PetscInt i,j[3],M,xs,xm; 436 DM da = user->da; 437 438 PetscFunctionBeginUser; 439 if (user->sjerr) { 440 PetscCall(SNESSetJacobianDomainError(snes)); 441 PetscFunctionReturn(0); 442 } 443 /* 444 Get pointer to vector data 445 */ 446 PetscCall(DMDAVecGetArrayRead(da,x,&xx)); 447 PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 448 449 /* 450 Get range of locally owned matrix 451 */ 452 PetscCall(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 453 454 /* 455 Determine starting and ending local indices for interior grid points. 456 Set Jacobian entries for boundary points. 457 */ 458 459 if (xs == 0) { /* left boundary */ 460 i = 0; A[0] = 1.0; 461 462 PetscCall(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES)); 463 xs++;xm--; 464 } 465 if (xs+xm == M) { /* right boundary */ 466 i = M-1; 467 A[0] = 1.0; 468 PetscCall(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES)); 469 xm--; 470 } 471 472 /* 473 Interior grid points 474 - Note that in this case we set all elements for a particular 475 row at once. 476 */ 477 d = 1.0/(user->h*user->h); 478 for (i=xs; i<xs+xm; i++) { 479 j[0] = i - 1; j[1] = i; j[2] = i + 1; 480 A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i]; 481 PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES)); 482 } 483 484 /* 485 Assemble matrix, using the 2-step process: 486 MatAssemblyBegin(), MatAssemblyEnd(). 487 By placing code between these two statements, computations can be 488 done while messages are in transition. 489 490 Also, restore vector. 491 */ 492 493 PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 494 PetscCall(DMDAVecRestoreArrayRead(da,x,&xx)); 495 PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 496 497 PetscFunctionReturn(0); 498 } 499 500 /* ------------------------------------------------------------------- */ 501 /* 502 Monitor - Optional user-defined monitoring routine that views the 503 current iterate with an x-window plot. Set by SNESMonitorSet(). 504 505 Input Parameters: 506 snes - the SNES context 507 its - iteration number 508 norm - 2-norm function value (may be estimated) 509 ctx - optional user-defined context for private data for the 510 monitor routine, as set by SNESMonitorSet() 511 512 Note: 513 See the manpage for PetscViewerDrawOpen() for useful runtime options, 514 such as -nox to deactivate all x-window output. 515 */ 516 PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx) 517 { 518 MonitorCtx *monP = (MonitorCtx*) ctx; 519 Vec x; 520 521 PetscFunctionBeginUser; 522 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"iter = %D,SNES Function norm %g\n",its,(double)fnorm)); 523 PetscCall(SNESGetSolution(snes,&x)); 524 PetscCall(VecView(x,monP->viewer)); 525 PetscFunctionReturn(0); 526 } 527 528 /* ------------------------------------------------------------------- */ 529 /* 530 PreCheck - Optional user-defined routine that checks the validity of 531 candidate steps of a line search method. Set by SNESLineSearchSetPreCheck(). 532 533 Input Parameters: 534 snes - the SNES context 535 xcurrent - current solution 536 y - search direction and length 537 538 Output Parameters: 539 y - proposed step (search direction and length) (possibly changed) 540 changed_y - tells if the step has changed or not 541 */ 542 PetscErrorCode PreCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y, PetscBool *changed_y, void * ctx) 543 { 544 PetscFunctionBeginUser; 545 *changed_y = PETSC_FALSE; 546 PetscFunctionReturn(0); 547 } 548 549 /* ------------------------------------------------------------------- */ 550 /* 551 PostCheck - Optional user-defined routine that checks the validity of 552 candidate steps of a line search method. Set by SNESLineSearchSetPostCheck(). 553 554 Input Parameters: 555 snes - the SNES context 556 ctx - optional user-defined context for private data for the 557 monitor routine, as set by SNESLineSearchSetPostCheck() 558 xcurrent - current solution 559 y - search direction and length 560 x - the new candidate iterate 561 562 Output Parameters: 563 y - proposed step (search direction and length) (possibly changed) 564 x - current iterate (possibly modified) 565 566 */ 567 PetscErrorCode PostCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool *changed_y,PetscBool *changed_x, void * ctx) 568 { 569 PetscInt i,iter,xs,xm; 570 StepCheckCtx *check; 571 ApplicationCtx *user; 572 PetscScalar *xa,*xa_last,tmp; 573 PetscReal rdiff; 574 DM da; 575 SNES snes; 576 577 PetscFunctionBeginUser; 578 *changed_x = PETSC_FALSE; 579 *changed_y = PETSC_FALSE; 580 581 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 582 check = (StepCheckCtx*)ctx; 583 user = check->user; 584 PetscCall(SNESGetIterationNumber(snes,&iter)); 585 586 /* iteration 1 indicates we are working on the second iteration */ 587 if (iter > 0) { 588 da = user->da; 589 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Checking candidate step at iteration %D with tolerance %g\n",iter,(double)check->tolerance)); 590 591 /* Access local array data */ 592 PetscCall(DMDAVecGetArray(da,check->last_step,&xa_last)); 593 PetscCall(DMDAVecGetArray(da,x,&xa)); 594 PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 595 596 /* 597 If we fail the user-defined check for validity of the candidate iterate, 598 then modify the iterate as we like. (Note that the particular modification 599 below is intended simply to demonstrate how to manipulate this data, not 600 as a meaningful or appropriate choice.) 601 */ 602 for (i=xs; i<xs+xm; i++) { 603 if (!PetscAbsScalar(xa[i])) rdiff = 2*check->tolerance; 604 else rdiff = PetscAbsScalar((xa[i] - xa_last[i])/xa[i]); 605 if (rdiff > check->tolerance) { 606 tmp = xa[i]; 607 xa[i] = .5*(xa[i] + xa_last[i]); 608 *changed_x = PETSC_TRUE; 609 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Altering entry %D: 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]))); 610 } 611 } 612 PetscCall(DMDAVecRestoreArray(da,check->last_step,&xa_last)); 613 PetscCall(DMDAVecRestoreArray(da,x,&xa)); 614 } 615 PetscCall(VecCopy(x,check->last_step)); 616 PetscFunctionReturn(0); 617 } 618 619 /* ------------------------------------------------------------------- */ 620 /* 621 PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used 622 e.g, 623 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 624 Set by SNESLineSearchSetPostCheck(). 625 626 Input Parameters: 627 linesearch - the LineSearch context 628 xcurrent - current solution 629 y - search direction and length 630 x - the new candidate iterate 631 632 Output Parameters: 633 y - proposed step (search direction and length) (possibly changed) 634 x - current iterate (possibly modified) 635 636 */ 637 PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool *changed_y,PetscBool *changed_x, void * ctx) 638 { 639 SetSubKSPCtx *check; 640 PetscInt iter,its,sub_its,maxit; 641 KSP ksp,sub_ksp,*sub_ksps; 642 PC pc; 643 PetscReal ksp_ratio; 644 SNES snes; 645 646 PetscFunctionBeginUser; 647 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 648 check = (SetSubKSPCtx*)ctx; 649 PetscCall(SNESGetIterationNumber(snes,&iter)); 650 PetscCall(SNESGetKSP(snes,&ksp)); 651 PetscCall(KSPGetPC(ksp,&pc)); 652 PetscCall(PCBJacobiGetSubKSP(pc,NULL,NULL,&sub_ksps)); 653 sub_ksp = sub_ksps[0]; 654 PetscCall(KSPGetIterationNumber(ksp,&its)); /* outer KSP iteration number */ 655 PetscCall(KSPGetIterationNumber(sub_ksp,&sub_its)); /* inner KSP iteration number */ 656 657 if (iter) { 658 PetscCall(PetscPrintf(PETSC_COMM_WORLD," ...PostCheck snes iteration %D, ksp_it %D %D, subksp_it %D\n",iter,check->its0,its,sub_its)); 659 ksp_ratio = ((PetscReal)(its))/check->its0; 660 maxit = (PetscInt)(ksp_ratio*sub_its + 0.5); 661 if (maxit < 2) maxit = 2; 662 PetscCall(KSPSetTolerances(sub_ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,maxit)); 663 PetscCall(PetscPrintf(PETSC_COMM_WORLD," ...ksp_ratio %g, new maxit %D\n\n",(double)ksp_ratio,maxit)); 664 } 665 check->its0 = its; /* save current outer KSP iteration number */ 666 PetscFunctionReturn(0); 667 } 668 669 /* ------------------------------------------------------------------- */ 670 /* 671 MatrixFreePreconditioner - This routine demonstrates the use of a 672 user-provided preconditioner. This code implements just the null 673 preconditioner, which of course is not recommended for general use. 674 675 Input Parameters: 676 + pc - preconditioner 677 - x - input vector 678 679 Output Parameter: 680 . y - preconditioned vector 681 */ 682 PetscErrorCode MatrixFreePreconditioner(PC pc,Vec x,Vec y) 683 { 684 PetscCall(VecCopy(x,y)); 685 return 0; 686 } 687 688 /*TEST 689 690 test: 691 args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 692 693 test: 694 suffix: 2 695 nsize: 3 696 args: -nox -pc_type asm -mat_type mpiaij -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 697 698 test: 699 suffix: 3 700 nsize: 2 701 args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 702 703 test: 704 suffix: 4 705 args: -nox -pre_check_iterates -post_check_iterates 706 707 test: 708 suffix: 5 709 requires: double !complex !single 710 nsize: 2 711 args: -nox -snes_test_jacobian -snes_test_jacobian_view 712 713 test: 714 suffix: 6 715 requires: double !complex !single 716 nsize: 4 717 args: -test_jacobian_domain_error -snes_converged_reason -snes_check_jacobian_domain_error 1 718 719 test: 720 suffix: 7 721 requires: double !complex !single 722 nsize: 4 723 args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontr -snes_check_jacobian_domain_error 1 724 725 test: 726 suffix: 8 727 requires: double !complex !single 728 nsize: 4 729 args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonrsls -snes_check_jacobian_domain_error 1 730 731 test: 732 suffix: 9 733 requires: double !complex !single 734 nsize: 4 735 args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonssls -snes_check_jacobian_domain_error 1 736 737 test: 738 suffix: 10 739 requires: double !complex !single 740 nsize: 4 741 args: -test_jacobian_domain_error -snes_converged_reason -snes_type qn -snes_qn_scale_type jacobian -snes_check_jacobian_domain_error 1 742 743 test: 744 suffix: 11 745 requires: double !complex !single 746 nsize: 4 747 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 748 749 test: 750 suffix: 12 751 args: -view_initial 752 filter: grep -v "type:" 753 754 test: 755 suffix: 13 756 requires: double !complex !single 757 nsize: 4 758 args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontrdc -snes_check_jacobian_domain_error 1 759 760 TEST*/ 761