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 { 76 SNES snes; /* SNES context */ 77 SNESLineSearch linesearch; /* SNESLineSearch context */ 78 Mat J; /* Jacobian matrix */ 79 ApplicationCtx ctx; /* user-defined context */ 80 Vec x,r,U,F; /* vectors */ 81 MonitorCtx monP; /* monitoring context */ 82 StepCheckCtx checkP; /* step-checking context */ 83 SetSubKSPCtx checkP1; 84 PetscBool pre_check,post_check,post_setsubksp; /* flag indicating whether we're checking candidate iterates */ 85 PetscScalar xp,*FF,*UU,none = -1.0; 86 PetscInt its,N = 5,i,maxit,maxf,xs,xm; 87 PetscReal abstol,rtol,stol,norm; 88 PetscBool flg,viewinitial = PETSC_FALSE; 89 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)); ctx.F = F; 123 PetscCall(VecDuplicate(x,&U)); 124 125 /* 126 Set function evaluation routine and vector. Whenever the nonlinear 127 solver needs to compute the nonlinear function, it will call this 128 routine. 129 - Note that the final routine argument is the user-defined 130 context that provides application-specific data for the 131 function evaluation routine. 132 */ 133 PetscCall(SNESSetFunction(snes,r,FormFunction,&ctx)); 134 135 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136 Create matrix data structure; set Jacobian evaluation routine 137 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138 139 PetscCall(MatCreate(PETSC_COMM_WORLD,&J)); 140 PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,N,N)); 141 PetscCall(MatSetFromOptions(J)); 142 PetscCall(MatSeqAIJSetPreallocation(J,3,NULL)); 143 PetscCall(MatMPIAIJSetPreallocation(J,3,NULL,3,NULL)); 144 145 /* 146 Set Jacobian matrix data structure and default Jacobian evaluation 147 routine. Whenever the nonlinear solver needs to compute the 148 Jacobian matrix, it will call this routine. 149 - Note that the final routine argument is the user-defined 150 context that provides application-specific data for the 151 Jacobian evaluation routine. 152 */ 153 PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,&ctx)); 154 155 /* 156 Optionally allow user-provided preconditioner 157 */ 158 PetscCall(PetscOptionsHasName(NULL,NULL,"-user_precond",&flg)); 159 if (flg) { 160 KSP ksp; 161 PC pc; 162 PetscCall(SNESGetKSP(snes,&ksp)); 163 PetscCall(KSPGetPC(ksp,&pc)); 164 PetscCall(PCSetType(pc,PCSHELL)); 165 PetscCall(PCShellSetApply(pc,MatrixFreePreconditioner)); 166 } 167 168 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 169 Customize nonlinear solver; set runtime options 170 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 171 172 /* 173 Set an optional user-defined monitoring routine 174 */ 175 PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,0,0,0,400,400,&monP.viewer)); 176 PetscCall(SNESMonitorSet(snes,Monitor,&monP,0)); 177 178 /* 179 Set names for some vectors to facilitate monitoring (optional) 180 */ 181 PetscCall(PetscObjectSetName((PetscObject)x,"Approximate Solution")); 182 PetscCall(PetscObjectSetName((PetscObject)U,"Exact Solution")); 183 184 /* 185 Set SNES/KSP/KSP/PC runtime options, e.g., 186 -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 187 */ 188 PetscCall(SNESSetFromOptions(snes)); 189 190 /* 191 Set an optional user-defined routine to check the validity of candidate 192 iterates that are determined by line search methods 193 */ 194 PetscCall(SNESGetLineSearch(snes, &linesearch)); 195 PetscCall(PetscOptionsHasName(NULL,NULL,"-post_check_iterates",&post_check)); 196 197 if (post_check) { 198 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Activating post step checking routine\n")); 199 PetscCall(SNESLineSearchSetPostCheck(linesearch,PostCheck,&checkP)); 200 PetscCall(VecDuplicate(x,&(checkP.last_step))); 201 202 checkP.tolerance = 1.0; 203 checkP.user = &ctx; 204 205 PetscCall(PetscOptionsGetReal(NULL,NULL,"-check_tol",&checkP.tolerance,NULL)); 206 } 207 208 PetscCall(PetscOptionsHasName(NULL,NULL,"-post_setsubksp",&post_setsubksp)); 209 if (post_setsubksp) { 210 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Activating post setsubksp\n")); 211 PetscCall(SNESLineSearchSetPostCheck(linesearch,PostSetSubKSP,&checkP1)); 212 } 213 214 PetscCall(PetscOptionsHasName(NULL,NULL,"-pre_check_iterates",&pre_check)); 215 if (pre_check) { 216 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Activating pre step checking routine\n")); 217 PetscCall(SNESLineSearchSetPreCheck(linesearch,PreCheck,&checkP)); 218 } 219 220 /* 221 Print parameters used for convergence testing (optional) ... just 222 to demonstrate this routine; this information is also printed with 223 the option -snes_view 224 */ 225 PetscCall(SNESGetTolerances(snes,&abstol,&rtol,&stol,&maxit,&maxf)); 226 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)); 227 228 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 229 Initialize application: 230 Store right-hand-side of PDE and exact solution 231 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 232 233 /* 234 Get local grid boundaries (for 1-dimensional DMDA): 235 xs, xm - starting grid index, width of local grid (no ghost points) 236 */ 237 PetscCall(DMDAGetCorners(ctx.da,&xs,NULL,NULL,&xm,NULL,NULL)); 238 239 /* 240 Get pointers to vector data 241 */ 242 PetscCall(DMDAVecGetArray(ctx.da,F,&FF)); 243 PetscCall(DMDAVecGetArray(ctx.da,U,&UU)); 244 245 /* 246 Compute local vector entries 247 */ 248 xp = ctx.h*xs; 249 for (i=xs; i<xs+xm; i++) { 250 FF[i] = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */ 251 UU[i] = xp*xp*xp; 252 xp += ctx.h; 253 } 254 255 /* 256 Restore vectors 257 */ 258 PetscCall(DMDAVecRestoreArray(ctx.da,F,&FF)); 259 PetscCall(DMDAVecRestoreArray(ctx.da,U,&UU)); 260 if (viewinitial) { 261 PetscCall(VecView(U,0)); 262 PetscCall(VecView(F,0)); 263 } 264 265 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 266 Evaluate initial guess; then solve nonlinear system 267 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 268 269 /* 270 Note: The user should initialize the vector, x, with the initial guess 271 for the nonlinear solver prior to calling SNESSolve(). In particular, 272 to employ an initial guess of zero, the user should explicitly set 273 this vector to zero by calling VecSet(). 274 */ 275 PetscCall(FormInitialGuess(x)); 276 PetscCall(SNESSolve(snes,NULL,x)); 277 PetscCall(SNESGetIterationNumber(snes,&its)); 278 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %" PetscInt_FMT "\n",its)); 279 280 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 281 Check solution and clean up 282 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 283 284 /* 285 Check the error 286 */ 287 PetscCall(VecAXPY(x,none,U)); 288 PetscCall(VecNorm(x,NORM_2,&norm)); 289 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g Iterations %" PetscInt_FMT "\n",(double)norm,its)); 290 if (ctx.sjerr) { 291 SNESType snestype; 292 PetscCall(SNESGetType(snes,&snestype)); 293 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"SNES Type %s\n",snestype)); 294 } 295 296 /* 297 Free work space. All PETSc objects should be destroyed when they 298 are no longer needed. 299 */ 300 PetscCall(PetscViewerDestroy(&monP.viewer)); 301 if (post_check) PetscCall(VecDestroy(&checkP.last_step)); 302 PetscCall(VecDestroy(&x)); 303 PetscCall(VecDestroy(&r)); 304 PetscCall(VecDestroy(&U)); 305 PetscCall(VecDestroy(&F)); 306 PetscCall(MatDestroy(&J)); 307 PetscCall(SNESDestroy(&snes)); 308 PetscCall(DMDestroy(&ctx.da)); 309 PetscCall(PetscFinalize()); 310 return 0; 311 } 312 313 /* ------------------------------------------------------------------- */ 314 /* 315 FormInitialGuess - Computes initial guess. 316 317 Input/Output Parameter: 318 . x - the solution vector 319 */ 320 PetscErrorCode FormInitialGuess(Vec x) 321 { 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 { 347 ApplicationCtx *user = (ApplicationCtx*) ctx; 348 DM da = user->da; 349 PetscScalar *ff,d; 350 const PetscScalar *xx,*FF; 351 PetscInt i,M,xs,xm; 352 Vec xlocal; 353 354 PetscFunctionBeginUser; 355 PetscCall(DMGetLocalVector(da,&xlocal)); 356 /* 357 Scatter ghost points to local vector, using the 2-step process 358 DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). 359 By placing code between these two statements, computations can 360 be done while messages are in transition. 361 */ 362 PetscCall(DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal)); 363 PetscCall(DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal)); 364 365 /* 366 Get pointers to vector data. 367 - The vector xlocal includes ghost point; the vectors x and f do 368 NOT include ghost points. 369 - Using DMDAVecGetArray() allows accessing the values using global ordering 370 */ 371 PetscCall(DMDAVecGetArrayRead(da,xlocal,(void*)&xx)); 372 PetscCall(DMDAVecGetArray(da,f,&ff)); 373 PetscCall(DMDAVecGetArrayRead(da,user->F,(void*)&FF)); 374 375 /* 376 Get local grid boundaries (for 1-dimensional DMDA): 377 xs, xm - starting grid index, width of local grid (no ghost points) 378 */ 379 PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 380 PetscCall(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 381 382 /* 383 Set function values for boundary points; define local interior grid point range: 384 xsi - starting interior grid index 385 xei - ending interior grid index 386 */ 387 if (xs == 0) { /* left boundary */ 388 ff[0] = xx[0]; 389 xs++;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 { 428 ApplicationCtx *user = (ApplicationCtx*) ctx; 429 PetscScalar *xx,d,A[3]; 430 PetscInt i,j[3],M,xs,xm; 431 DM da = user->da; 432 433 PetscFunctionBeginUser; 434 if (user->sjerr) { 435 PetscCall(SNESSetJacobianDomainError(snes)); 436 PetscFunctionReturn(0); 437 } 438 /* 439 Get pointer to vector data 440 */ 441 PetscCall(DMDAVecGetArrayRead(da,x,&xx)); 442 PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 443 444 /* 445 Get range of locally owned matrix 446 */ 447 PetscCall(DMDAGetInfo(da,NULL,&M,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 448 449 /* 450 Determine starting and ending local indices for interior grid points. 451 Set Jacobian entries for boundary points. 452 */ 453 454 if (xs == 0) { /* left boundary */ 455 i = 0; A[0] = 1.0; 456 457 PetscCall(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES)); 458 xs++;xm--; 459 } 460 if (xs+xm == M) { /* right boundary */ 461 i = M-1; 462 A[0] = 1.0; 463 PetscCall(MatSetValues(jac,1,&i,1,&i,A,INSERT_VALUES)); 464 xm--; 465 } 466 467 /* 468 Interior grid points 469 - Note that in this case we set all elements for a particular 470 row at once. 471 */ 472 d = 1.0/(user->h*user->h); 473 for (i=xs; i<xs+xm; i++) { 474 j[0] = i - 1; j[1] = i; j[2] = i + 1; 475 A[0] = A[2] = d; A[1] = -2.0*d + 2.0*xx[i]; 476 PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES)); 477 } 478 479 /* 480 Assemble matrix, using the 2-step process: 481 MatAssemblyBegin(), MatAssemblyEnd(). 482 By placing code between these two statements, computations can be 483 done while messages are in transition. 484 485 Also, restore vector. 486 */ 487 488 PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 489 PetscCall(DMDAVecRestoreArrayRead(da,x,&xx)); 490 PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 491 492 PetscFunctionReturn(0); 493 } 494 495 /* ------------------------------------------------------------------- */ 496 /* 497 Monitor - Optional user-defined monitoring routine that views the 498 current iterate with an x-window plot. Set by SNESMonitorSet(). 499 500 Input Parameters: 501 snes - the SNES context 502 its - iteration number 503 norm - 2-norm function value (may be estimated) 504 ctx - optional user-defined context for private data for the 505 monitor routine, as set by SNESMonitorSet() 506 507 Note: 508 See the manpage for PetscViewerDrawOpen() for useful runtime options, 509 such as -nox to deactivate all x-window output. 510 */ 511 PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *ctx) 512 { 513 MonitorCtx *monP = (MonitorCtx*) ctx; 514 Vec x; 515 516 PetscFunctionBeginUser; 517 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"iter = %" PetscInt_FMT ",SNES Function norm %g\n",its,(double)fnorm)); 518 PetscCall(SNESGetSolution(snes,&x)); 519 PetscCall(VecView(x,monP->viewer)); 520 PetscFunctionReturn(0); 521 } 522 523 /* ------------------------------------------------------------------- */ 524 /* 525 PreCheck - Optional user-defined routine that checks the validity of 526 candidate steps of a line search method. Set by SNESLineSearchSetPreCheck(). 527 528 Input Parameters: 529 snes - the SNES context 530 xcurrent - current solution 531 y - search direction and length 532 533 Output Parameters: 534 y - proposed step (search direction and length) (possibly changed) 535 changed_y - tells if the step has changed or not 536 */ 537 PetscErrorCode PreCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y, PetscBool *changed_y, void * ctx) 538 { 539 PetscFunctionBeginUser; 540 *changed_y = PETSC_FALSE; 541 PetscFunctionReturn(0); 542 } 543 544 /* ------------------------------------------------------------------- */ 545 /* 546 PostCheck - Optional user-defined routine that checks the validity of 547 candidate steps of a line search method. Set by SNESLineSearchSetPostCheck(). 548 549 Input Parameters: 550 snes - the SNES context 551 ctx - optional user-defined context for private data for the 552 monitor routine, as set by SNESLineSearchSetPostCheck() 553 xcurrent - current solution 554 y - search direction and length 555 x - the new candidate iterate 556 557 Output Parameters: 558 y - proposed step (search direction and length) (possibly changed) 559 x - current iterate (possibly modified) 560 561 */ 562 PetscErrorCode PostCheck(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool *changed_y,PetscBool *changed_x, void * ctx) 563 { 564 PetscInt i,iter,xs,xm; 565 StepCheckCtx *check; 566 ApplicationCtx *user; 567 PetscScalar *xa,*xa_last,tmp; 568 PetscReal rdiff; 569 DM da; 570 SNES snes; 571 572 PetscFunctionBeginUser; 573 *changed_x = PETSC_FALSE; 574 *changed_y = PETSC_FALSE; 575 576 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 577 check = (StepCheckCtx*)ctx; 578 user = check->user; 579 PetscCall(SNESGetIterationNumber(snes,&iter)); 580 581 /* iteration 1 indicates we are working on the second iteration */ 582 if (iter > 0) { 583 da = user->da; 584 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Checking candidate step at iteration %" PetscInt_FMT " with tolerance %g\n",iter,(double)check->tolerance)); 585 586 /* Access local array data */ 587 PetscCall(DMDAVecGetArray(da,check->last_step,&xa_last)); 588 PetscCall(DMDAVecGetArray(da,x,&xa)); 589 PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 590 591 /* 592 If we fail the user-defined check for validity of the candidate iterate, 593 then modify the iterate as we like. (Note that the particular modification 594 below is intended simply to demonstrate how to manipulate this data, not 595 as a meaningful or appropriate choice.) 596 */ 597 for (i=xs; i<xs+xm; i++) { 598 if (!PetscAbsScalar(xa[i])) rdiff = 2*check->tolerance; 599 else rdiff = PetscAbsScalar((xa[i] - xa_last[i])/xa[i]); 600 if (rdiff > check->tolerance) { 601 tmp = xa[i]; 602 xa[i] = .5*(xa[i] + xa_last[i]); 603 *changed_x = PETSC_TRUE; 604 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]))); 605 } 606 } 607 PetscCall(DMDAVecRestoreArray(da,check->last_step,&xa_last)); 608 PetscCall(DMDAVecRestoreArray(da,x,&xa)); 609 } 610 PetscCall(VecCopy(x,check->last_step)); 611 PetscFunctionReturn(0); 612 } 613 614 /* ------------------------------------------------------------------- */ 615 /* 616 PostSetSubKSP - Optional user-defined routine that reset SubKSP options when hierarchical bjacobi PC is used 617 e.g, 618 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 619 Set by SNESLineSearchSetPostCheck(). 620 621 Input Parameters: 622 linesearch - the LineSearch context 623 xcurrent - current solution 624 y - search direction and length 625 x - the new candidate iterate 626 627 Output Parameters: 628 y - proposed step (search direction and length) (possibly changed) 629 x - current iterate (possibly modified) 630 631 */ 632 PetscErrorCode PostSetSubKSP(SNESLineSearch linesearch,Vec xcurrent,Vec y,Vec x,PetscBool *changed_y,PetscBool *changed_x, void * ctx) 633 { 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 { 679 PetscCall(VecCopy(x,y)); 680 return 0; 681 } 682 683 /*TEST 684 685 test: 686 args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 687 688 test: 689 suffix: 2 690 nsize: 3 691 args: -nox -pc_type asm -mat_type mpiaij -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 692 693 test: 694 suffix: 3 695 nsize: 2 696 args: -nox -snes_monitor_cancel -snes_monitor_short -ksp_gmres_cgs_refinement_type refine_always 697 698 test: 699 suffix: 4 700 args: -nox -pre_check_iterates -post_check_iterates 701 702 test: 703 suffix: 5 704 requires: double !complex !single 705 nsize: 2 706 args: -nox -snes_test_jacobian -snes_test_jacobian_view 707 708 test: 709 suffix: 6 710 requires: double !complex !single 711 nsize: 4 712 args: -test_jacobian_domain_error -snes_converged_reason -snes_check_jacobian_domain_error 1 713 714 test: 715 suffix: 7 716 requires: double !complex !single 717 nsize: 4 718 args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontr -snes_check_jacobian_domain_error 1 719 720 test: 721 suffix: 8 722 requires: double !complex !single 723 nsize: 4 724 args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonrsls -snes_check_jacobian_domain_error 1 725 726 test: 727 suffix: 9 728 requires: double !complex !single 729 nsize: 4 730 args: -test_jacobian_domain_error -snes_converged_reason -snes_type vinewtonssls -snes_check_jacobian_domain_error 1 731 732 test: 733 suffix: 10 734 requires: double !complex !single 735 nsize: 4 736 args: -test_jacobian_domain_error -snes_converged_reason -snes_type qn -snes_qn_scale_type jacobian -snes_check_jacobian_domain_error 1 737 738 test: 739 suffix: 11 740 requires: double !complex !single 741 nsize: 4 742 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 743 744 test: 745 suffix: 12 746 args: -view_initial 747 filter: grep -v "type:" 748 749 test: 750 suffix: 13 751 requires: double !complex !single 752 nsize: 4 753 args: -test_jacobian_domain_error -snes_converged_reason -snes_type newtontrdc -snes_check_jacobian_domain_error 1 754 755 TEST*/ 756