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