1 #define PETSCSNES_DLL 2 3 #include "include/private/snesimpl.h" /*I "petscsnes.h" I*/ 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "SNESMonitorSolution" 7 /*@C 8 SNESMonitorSolution - Monitors progress of the SNES solvers by calling 9 VecView() for the approximate solution at each iteration. 10 11 Collective on SNES 12 13 Input Parameters: 14 + snes - the SNES context 15 . its - iteration number 16 . fgnorm - 2-norm of residual 17 - dummy - either a viewer or PETSC_NULL 18 19 Level: intermediate 20 21 .keywords: SNES, nonlinear, vector, monitor, view 22 23 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView() 24 @*/ 25 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSolution(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 26 { 27 PetscErrorCode ierr; 28 Vec x; 29 PetscViewer viewer = (PetscViewer) dummy; 30 31 PetscFunctionBegin; 32 ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr); 33 if (!viewer) { 34 MPI_Comm comm; 35 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 36 viewer = PETSC_VIEWER_DRAW_(comm); 37 } 38 ierr = VecView(x,viewer);CHKERRQ(ierr); 39 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "SNESMonitorResidual" 45 /*@C 46 SNESMonitorResidual - Monitors progress of the SNES solvers by calling 47 VecView() for the residual at each iteration. 48 49 Collective on SNES 50 51 Input Parameters: 52 + snes - the SNES context 53 . its - iteration number 54 . fgnorm - 2-norm of residual 55 - dummy - either a viewer or PETSC_NULL 56 57 Level: intermediate 58 59 .keywords: SNES, nonlinear, vector, monitor, view 60 61 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView() 62 @*/ 63 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 64 { 65 PetscErrorCode ierr; 66 Vec x; 67 PetscViewer viewer = (PetscViewer) dummy; 68 69 PetscFunctionBegin; 70 ierr = SNESGetFunction(snes,&x,0,0);CHKERRQ(ierr); 71 if (!viewer) { 72 MPI_Comm comm; 73 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 74 viewer = PETSC_VIEWER_DRAW_(comm); 75 } 76 ierr = VecView(x,viewer);CHKERRQ(ierr); 77 78 PetscFunctionReturn(0); 79 } 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "SNESMonitorSolutionUpdate" 83 /*@C 84 SNESMonitorSolutionUpdate - Monitors progress of the SNES solvers by calling 85 VecView() for the UPDATE to the solution at each iteration. 86 87 Collective on SNES 88 89 Input Parameters: 90 + snes - the SNES context 91 . its - iteration number 92 . fgnorm - 2-norm of residual 93 - dummy - either a viewer or PETSC_NULL 94 95 Level: intermediate 96 97 .keywords: SNES, nonlinear, vector, monitor, view 98 99 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView() 100 @*/ 101 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSolutionUpdate(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 102 { 103 PetscErrorCode ierr; 104 Vec x; 105 PetscViewer viewer = (PetscViewer) dummy; 106 107 PetscFunctionBegin; 108 ierr = SNESGetSolutionUpdate(snes,&x);CHKERRQ(ierr); 109 if (!viewer) { 110 MPI_Comm comm; 111 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 112 viewer = PETSC_VIEWER_DRAW_(comm); 113 } 114 ierr = VecView(x,viewer);CHKERRQ(ierr); 115 116 PetscFunctionReturn(0); 117 } 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "SNESMonitorDefault" 121 /*@C 122 SNESMonitorDefault - Monitors progress of the SNES solvers (default). 123 124 Collective on SNES 125 126 Input Parameters: 127 + snes - the SNES context 128 . its - iteration number 129 . fgnorm - 2-norm of residual 130 - dummy - unused context 131 132 Notes: 133 This routine prints the residual norm at each iteration. 134 135 Level: intermediate 136 137 .keywords: SNES, nonlinear, default, monitor, norm 138 139 .seealso: SNESMonitorSet(), SNESMonitorSolution() 140 @*/ 141 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 142 { 143 PetscErrorCode ierr; 144 PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy; 145 146 PetscFunctionBegin; 147 if (!dummy) { 148 ierr = PetscViewerASCIIMonitorCreate(snes->comm,"stdout",0,&viewer);CHKERRQ(ierr); 149 } 150 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,fgnorm);CHKERRQ(ierr); 151 if (!dummy) { 152 ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 153 } 154 PetscFunctionReturn(0); 155 } 156 157 typedef struct { 158 PetscViewerASCIIMonitor viewer; 159 PetscReal *history; 160 } SNESMonitorRatioContext; 161 162 #undef __FUNCT__ 163 #define __FUNCT__ "SNESMonitorRatio" 164 /*@C 165 SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio 166 of residual norm at each iteration to the previous. 167 168 Collective on SNES 169 170 Input Parameters: 171 + snes - the SNES context 172 . its - iteration number 173 . fgnorm - 2-norm of residual (or gradient) 174 - dummy - context of monitor 175 176 Level: intermediate 177 178 .keywords: SNES, nonlinear, monitor, norm 179 180 .seealso: SNESMonitorSet(), SNESMonitorSolution() 181 @*/ 182 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 183 { 184 PetscErrorCode ierr; 185 PetscInt len; 186 PetscReal *history; 187 SNESMonitorRatioContext *ctx = (SNESMonitorRatioContext*)dummy; 188 189 PetscFunctionBegin; 190 ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,&len);CHKERRQ(ierr); 191 if (!its || !history || its > len) { 192 ierr = PetscViewerASCIIMonitorPrintf(ctx->viewer,"%3D SNES Function norm %14.12e \n",its,fgnorm);CHKERRQ(ierr); 193 } else { 194 PetscReal ratio = fgnorm/history[its-1]; 195 ierr = PetscViewerASCIIMonitorPrintf(ctx->viewer,"%3D SNES Function norm %14.12e %G \n",its,fgnorm,ratio);CHKERRQ(ierr); 196 } 197 PetscFunctionReturn(0); 198 } 199 200 /* 201 If the we set the history monitor space then we must destroy it 202 */ 203 #undef __FUNCT__ 204 #define __FUNCT__ "SNESMonitorRatioDestroy" 205 PetscErrorCode SNESMonitorRatioDestroy(void *ct) 206 { 207 PetscErrorCode ierr; 208 SNESMonitorRatioContext *ctx = (SNESMonitorRatioContext*)ct; 209 210 PetscFunctionBegin; 211 ierr = PetscFree(ctx->history);CHKERRQ(ierr); 212 ierr = PetscViewerASCIIMonitorDestroy(ctx->viewer);CHKERRQ(ierr); 213 ierr = PetscFree(ctx);CHKERRQ(ierr); 214 PetscFunctionReturn(0); 215 } 216 217 #undef __FUNCT__ 218 #define __FUNCT__ "SNESMonitorSetRatio" 219 /*@C 220 SNESMonitorSetRatio - Sets SNES to use a monitor that prints the 221 ratio of the function norm at each iteration. 222 223 Collective on SNES 224 225 Input Parameters: 226 + snes - the SNES context 227 - viewer - ASCII viewer to print output 228 229 Level: intermediate 230 231 .keywords: SNES, nonlinear, monitor, norm 232 233 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault() 234 @*/ 235 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSetRatio(SNES snes,PetscViewerASCIIMonitor viewer) 236 { 237 PetscErrorCode ierr; 238 SNESMonitorRatioContext *ctx; 239 PetscReal *history; 240 241 PetscFunctionBegin; 242 if (!viewer) { 243 ierr = PetscViewerASCIIMonitorCreate(snes->comm,"stdout",0,&viewer);CHKERRQ(ierr); 244 ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr); 245 } 246 ierr = PetscNew(SNESMonitorRatioContext,&ctx); 247 ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 248 if (!history) { 249 ierr = PetscMalloc(100*sizeof(PetscReal),&ctx->history);CHKERRQ(ierr); 250 ierr = SNESSetConvergenceHistory(snes,ctx->history,0,100,PETSC_TRUE);CHKERRQ(ierr); 251 } 252 ctx->viewer = viewer; 253 ierr = SNESMonitorSet(snes,SNESMonitorRatio,ctx,SNESMonitorRatioDestroy);CHKERRQ(ierr); 254 PetscFunctionReturn(0); 255 } 256 257 /* ---------------------------------------------------------------- */ 258 #undef __FUNCT__ 259 #define __FUNCT__ "SNESMonitorDefaultShort" 260 /* 261 Default (short) SNES Monitor, same as SNESMonitorDefault() except 262 it prints fewer digits of the residual as the residual gets smaller. 263 This is because the later digits are meaningless and are often 264 different on different machines; by using this routine different 265 machines will usually generate the same output. 266 */ 267 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 268 { 269 PetscErrorCode ierr; 270 PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy; 271 272 PetscFunctionBegin; 273 if (!dummy) { 274 ierr = PetscViewerASCIIMonitorCreate(snes->comm,"stdout",0,&viewer);CHKERRQ(ierr); 275 } 276 if (fgnorm > 1.e-9) { 277 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES Function norm %G \n",its,fgnorm);CHKERRQ(ierr); 278 } else if (fgnorm > 1.e-11){ 279 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,fgnorm);CHKERRQ(ierr); 280 } else { 281 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);CHKERRQ(ierr); 282 } 283 if (!dummy) { 284 ierr = PetscViewerASCIIMonitorDestroy(viewer);CHKERRQ(ierr); 285 } 286 PetscFunctionReturn(0); 287 } 288 /* ---------------------------------------------------------------- */ 289 #undef __FUNCT__ 290 #define __FUNCT__ "SNESConverged_LS" 291 /*@C 292 SNESConverged_LS - Monitors the convergence of the solvers for 293 systems of nonlinear equations (default). 294 295 Collective on SNES 296 297 Input Parameters: 298 + snes - the SNES context 299 . it - the iteration (0 indicates before any Newton steps) 300 . xnorm - 2-norm of current iterate 301 . pnorm - 2-norm of current step 302 . fnorm - 2-norm of function at current iterate 303 - dummy - unused context 304 305 Output Parameter: 306 . reason - one of 307 $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol), 308 $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 309 $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 310 $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 311 $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 312 $ SNES_CONVERGED_ITERATING - (otherwise), 313 314 where 315 + maxf - maximum number of function evaluations, 316 set with SNESSetTolerances() 317 . nfct - number of function evaluations, 318 . abstol - absolute function norm tolerance, 319 set with SNESSetTolerances() 320 - rtol - relative function norm tolerance, set with SNESSetTolerances() 321 322 Level: intermediate 323 324 .keywords: SNES, nonlinear, default, converged, convergence 325 326 .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged() 327 @*/ 328 PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_LS(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 329 { 330 PetscErrorCode ierr; 331 332 PetscFunctionBegin; 333 *reason = SNES_CONVERGED_ITERATING; 334 335 if (!it) { 336 /* set parameter for default relative tolerance convergence test */ 337 snes->ttol = fnorm*snes->rtol; 338 } 339 if (fnorm != fnorm) { 340 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 341 *reason = SNES_DIVERGED_FNORM_NAN; 342 } else if (fnorm < snes->abstol) { 343 ierr = PetscInfo2(snes,"Converged due to function norm %G < %G\n",fnorm,snes->abstol);CHKERRQ(ierr); 344 *reason = SNES_CONVERGED_FNORM_ABS; 345 } else if (snes->nfuncs >= snes->max_funcs) { 346 ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 347 *reason = SNES_DIVERGED_FUNCTION_COUNT; 348 } 349 350 if (it && !*reason) { 351 if (fnorm <= snes->ttol) { 352 ierr = PetscInfo2(snes,"Converged due to function norm %G < %G (relative tolerance)\n",fnorm,snes->ttol);CHKERRQ(ierr); 353 *reason = SNES_CONVERGED_FNORM_RELATIVE; 354 } else if (pnorm < snes->xtol*xnorm) { 355 ierr = PetscInfo3(snes,"Converged due to small update length: %G < %G * %G\n",pnorm,snes->xtol,xnorm);CHKERRQ(ierr); 356 *reason = SNES_CONVERGED_PNORM_RELATIVE; 357 } 358 } 359 PetscFunctionReturn(0); 360 } 361 /* ------------------------------------------------------------ */ 362 #undef __FUNCT__ 363 #define __FUNCT__ "SNES_KSP_SetConvergenceTestEW" 364 /*@ 365 SNES_KSP_SetConvergenceTestEW - Sets alternative convergence test 366 for the linear solvers within an inexact Newton method. 367 368 Collective on SNES 369 370 Input Parameter: 371 . snes - SNES context 372 373 Notes: 374 Currently, the default is to use a constant relative tolerance for 375 the inner linear solvers. Alternatively, one can use the 376 Eisenstat-Walker method, where the relative convergence tolerance 377 is reset at each Newton iteration according progress of the nonlinear 378 solver. 379 380 Level: advanced 381 382 Reference: 383 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 384 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 385 386 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 387 @*/ 388 PetscErrorCode PETSCSNES_DLLEXPORT SNES_KSP_SetConvergenceTestEW(SNES snes) 389 { 390 PetscFunctionBegin; 391 snes->ksp_ewconv = PETSC_TRUE; 392 PetscFunctionReturn(0); 393 } 394 395 #undef __FUNCT__ 396 #define __FUNCT__ "SNES_KSP_SetParametersEW" 397 /*@ 398 SNES_KSP_SetParametersEW - Sets parameters for Eisenstat-Walker 399 convergence criteria for the linear solvers within an inexact 400 Newton method. 401 402 Collective on SNES 403 404 Input Parameters: 405 + snes - SNES context 406 . version - version 1, 2 (default is 2) or 3 407 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 408 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 409 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 410 . alpha2 - power for safeguard 411 . gamma2 - multiplicative factor for version 2 rtol computation 412 (0 <= gamma2 <= 1) 413 - threshold - threshold for imposing safeguard (0 < threshold < 1) 414 415 Note: 416 Version 3 was contributed by Luis Chacon, June 2006. 417 418 Use PETSC_DEFAULT to retain the default for any of the parameters. 419 420 Level: advanced 421 422 Reference: 423 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 424 inexact Newton method", Utah State University Math. Stat. Dept. Res. 425 Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 426 427 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters 428 429 .seealso: SNES_KSP_SetConvergenceTestEW() 430 @*/ 431 PetscErrorCode PETSCSNES_DLLEXPORT SNES_KSP_SetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,PetscReal gamma2,PetscReal alpha, 432 PetscReal alpha2,PetscReal threshold) 433 { 434 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx; 435 436 PetscFunctionBegin; 437 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 438 if (version != PETSC_DEFAULT) kctx->version = version; 439 if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0; 440 if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max; 441 if (gamma2 != PETSC_DEFAULT) kctx->gamma = gamma2; 442 if (alpha != PETSC_DEFAULT) kctx->alpha = alpha; 443 if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2; 444 if (threshold != PETSC_DEFAULT) kctx->threshold = threshold; 445 if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) { 446 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0); 447 } 448 if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) { 449 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max); 450 } 451 if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) { 452 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold); 453 } 454 if (kctx->gamma < 0.0 || kctx->gamma > 1.0) { 455 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma); 456 } 457 if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) { 458 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha); 459 } 460 if (kctx->version < 1 || kctx->version > 3) { 461 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version); 462 } 463 PetscFunctionReturn(0); 464 } 465 466 #undef __FUNCT__ 467 #define __FUNCT__ "SNES_KSP_EW_ComputeRelativeTolerance_Private" 468 PetscErrorCode SNES_KSP_EW_ComputeRelativeTolerance_Private(SNES snes,KSP ksp) 469 { 470 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx; 471 PetscReal rtol = 0.0,stol; 472 PetscErrorCode ierr; 473 474 PetscFunctionBegin; 475 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 476 if (!snes->iter) { /* first time in, so use the original user rtol */ 477 rtol = kctx->rtol_0; 478 } else { 479 if (kctx->version == 1) { 480 rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last; 481 if (rtol < 0.0) rtol = -rtol; 482 stol = pow(kctx->rtol_last,kctx->alpha2); 483 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 484 } else if (kctx->version == 2) { 485 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 486 stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha); 487 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 488 } else if (kctx->version == 3) { 489 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 490 /* safeguard: avoid sharp decrease of rtol */ 491 rtol = PetscMin(kctx->rtol_0,PetscMax(rtol,kctx->gamma*pow(kctx->rtol_last,kctx->alpha))); 492 /* safeguard: avoid oversolving */ 493 rtol = PetscMin(kctx->rtol_0,PetscMax(rtol,kctx->gamma*(snes->ttol)/snes->norm)); 494 495 } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version); 496 } 497 rtol = PetscMin(rtol,kctx->rtol_max); 498 kctx->rtol_last = rtol; 499 ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol = %G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr); 500 ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 501 kctx->norm_last = snes->norm; 502 PetscFunctionReturn(0); 503 } 504 505 #undef __FUNCT__ 506 #define __FUNCT__ "SNES_KSP_EW_Converged_Private" 507 PetscErrorCode SNES_KSP_EW_Converged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx) 508 { 509 SNES snes = (SNES)ctx; 510 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx; 511 PetscErrorCode ierr; 512 513 PetscFunctionBegin; 514 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context set"); 515 if (!n) {ierr = SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);CHKERRQ(ierr);} 516 ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr); 517 kctx->lresid_last = rnorm; 518 if (*reason) { 519 ierr = PetscInfo2(snes,"KSP iterations=%D, rnorm=%G\n",n,rnorm);CHKERRQ(ierr); 520 } 521 PetscFunctionReturn(0); 522 } 523 524 525