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