1 2 #include <private/snesimpl.h> /*I "petscsnes.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "SNESMonitorSolution" 6 /*@C 7 SNESMonitorSolution - 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: SNESMonitorSet(), SNESMonitorDefault(), VecView() 23 @*/ 24 PetscErrorCode SNESMonitorSolution(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__ "SNESMonitorResidual" 44 /*@C 45 SNESMonitorResidual - 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: SNESMonitorSet(), SNESMonitorDefault(), VecView() 61 @*/ 62 PetscErrorCode SNESMonitorResidual(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__ "SNESMonitorSolutionUpdate" 82 /*@C 83 SNESMonitorSolutionUpdate - 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: SNESMonitorSet(), SNESMonitorDefault(), VecView() 99 @*/ 100 PetscErrorCode SNESMonitorSolutionUpdate(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__ "SNESMonitorDefault" 120 /*@C 121 SNESMonitorDefault - 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: SNESMonitorSet(), SNESMonitorSolution() 139 @*/ 140 PetscErrorCode SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 141 { 142 PetscErrorCode ierr; 143 PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy; 144 145 PetscFunctionBegin; 146 if (!dummy) { 147 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 148 } 149 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr); 150 if (!dummy) { 151 ierr = PetscViewerASCIIMonitorDestroy(&viewer);CHKERRQ(ierr); 152 } 153 PetscFunctionReturn(0); 154 } 155 156 #undef __FUNCT__ 157 #define __FUNCT__ "SNESMonitorRange_Private" 158 PetscErrorCode SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per) 159 { 160 PetscErrorCode ierr; 161 Vec resid; 162 PetscReal rmax,pwork; 163 PetscInt i,n,N; 164 PetscScalar *r; 165 166 PetscFunctionBegin; 167 ierr = SNESGetFunction(snes,&resid,0,0);CHKERRQ(ierr); 168 ierr = VecNorm(resid,NORM_INFINITY,&rmax);CHKERRQ(ierr); 169 ierr = VecGetLocalSize(resid,&n);CHKERRQ(ierr); 170 ierr = VecGetSize(resid,&N);CHKERRQ(ierr); 171 ierr = VecGetArray(resid,&r);CHKERRQ(ierr); 172 pwork = 0.0; 173 for (i=0; i<n; i++) { 174 pwork += (PetscAbsScalar(r[i]) > .20*rmax); 175 } 176 ierr = MPI_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr); 177 ierr = VecRestoreArray(resid,&r);CHKERRQ(ierr); 178 *per = *per/N; 179 PetscFunctionReturn(0); 180 } 181 182 #undef __FUNCT__ 183 #define __FUNCT__ "SNESMonitorRange" 184 /*@C 185 SNESMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value. 186 187 Collective on SNES 188 189 Input Parameters: 190 + snes - iterative context 191 . it - iteration number 192 . rnorm - 2-norm (preconditioned) residual value (may be estimated). 193 - dummy - unused monitor context 194 195 Options Database Key: 196 . -snes_monitor_range - Activates SNESMonitorRange() 197 198 Level: intermediate 199 200 .keywords: SNES, default, monitor, residual 201 202 .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate() 203 @*/ 204 PetscErrorCode SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,void *dummy) 205 { 206 PetscErrorCode ierr; 207 PetscReal perc,rel; 208 PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy; 209 /* should be in a MonitorRangeContext */ 210 static PetscReal prev; 211 212 PetscFunctionBegin; 213 if (!it) prev = rnorm; 214 if (!dummy) {ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",0,&viewer);CHKERRQ(ierr);} 215 ierr = SNESMonitorRange_Private(snes,it,&perc);CHKERRQ(ierr); 216 217 rel = (prev - rnorm)/prev; 218 prev = rnorm; 219 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES preconditioned resid norm %14.12e Percent values above 20 percent of maximum %5.2f relative decrease %5.2e ratio %5.2e \n",it,(double)rnorm,(double)100.0*perc,(double)rel,(double)rel/perc);CHKERRQ(ierr); 220 if (!dummy) {ierr = PetscViewerASCIIMonitorDestroy(&viewer);CHKERRQ(ierr);} 221 PetscFunctionReturn(0); 222 } 223 224 typedef struct { 225 PetscViewerASCIIMonitor viewer; 226 PetscReal *history; 227 } SNESMonitorRatioContext; 228 229 #undef __FUNCT__ 230 #define __FUNCT__ "SNESMonitorRatio" 231 /*@C 232 SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio 233 of residual norm at each iteration to the previous. 234 235 Collective on SNES 236 237 Input Parameters: 238 + snes - the SNES context 239 . its - iteration number 240 . fgnorm - 2-norm of residual (or gradient) 241 - dummy - context of monitor 242 243 Level: intermediate 244 245 .keywords: SNES, nonlinear, monitor, norm 246 247 .seealso: SNESMonitorSet(), SNESMonitorSolution() 248 @*/ 249 PetscErrorCode SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 250 { 251 PetscErrorCode ierr; 252 PetscInt len; 253 PetscReal *history; 254 SNESMonitorRatioContext *ctx = (SNESMonitorRatioContext*)dummy; 255 256 PetscFunctionBegin; 257 ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,&len);CHKERRQ(ierr); 258 if (!its || !history || its > len) { 259 ierr = PetscViewerASCIIMonitorPrintf(ctx->viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr); 260 } else { 261 PetscReal ratio = fgnorm/history[its-1]; 262 ierr = PetscViewerASCIIMonitorPrintf(ctx->viewer,"%3D SNES Function norm %14.12e %G \n",its,(double)fgnorm,(double)ratio);CHKERRQ(ierr); 263 } 264 PetscFunctionReturn(0); 265 } 266 267 /* 268 If the we set the history monitor space then we must destroy it 269 */ 270 #undef __FUNCT__ 271 #define __FUNCT__ "SNESMonitorRatioDestroy" 272 PetscErrorCode SNESMonitorRatioDestroy(void **ct) 273 { 274 PetscErrorCode ierr; 275 SNESMonitorRatioContext *ctx = *(SNESMonitorRatioContext**)ct; 276 277 PetscFunctionBegin; 278 ierr = PetscFree(ctx->history);CHKERRQ(ierr); 279 ierr = PetscViewerASCIIMonitorDestroy(&ctx->viewer);CHKERRQ(ierr); 280 ierr = PetscFree(ctx);CHKERRQ(ierr); 281 PetscFunctionReturn(0); 282 } 283 284 #undef __FUNCT__ 285 #define __FUNCT__ "SNESMonitorSetRatio" 286 /*@C 287 SNESMonitorSetRatio - Sets SNES to use a monitor that prints the 288 ratio of the function norm at each iteration. 289 290 Collective on SNES 291 292 Input Parameters: 293 + snes - the SNES context 294 - viewer - ASCII viewer to print output 295 296 Level: intermediate 297 298 .keywords: SNES, nonlinear, monitor, norm 299 300 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault() 301 @*/ 302 PetscErrorCode SNESMonitorSetRatio(SNES snes,PetscViewerASCIIMonitor viewer) 303 { 304 PetscErrorCode ierr; 305 SNESMonitorRatioContext *ctx; 306 PetscReal *history; 307 308 PetscFunctionBegin; 309 if (!viewer) { 310 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 311 ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr); 312 } 313 ierr = PetscNewLog(snes,SNESMonitorRatioContext,&ctx); 314 ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 315 if (!history) { 316 ierr = PetscMalloc(100*sizeof(PetscReal),&ctx->history);CHKERRQ(ierr); 317 ierr = SNESSetConvergenceHistory(snes,ctx->history,0,100,PETSC_TRUE);CHKERRQ(ierr); 318 } 319 ctx->viewer = viewer; 320 ierr = SNESMonitorSet(snes,SNESMonitorRatio,ctx,SNESMonitorRatioDestroy);CHKERRQ(ierr); 321 PetscFunctionReturn(0); 322 } 323 324 /* ---------------------------------------------------------------- */ 325 #undef __FUNCT__ 326 #define __FUNCT__ "SNESMonitorDefaultShort" 327 /* 328 Default (short) SNES Monitor, same as SNESMonitorDefault() except 329 it prints fewer digits of the residual as the residual gets smaller. 330 This is because the later digits are meaningless and are often 331 different on different machines; by using this routine different 332 machines will usually generate the same output. 333 */ 334 PetscErrorCode SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 335 { 336 PetscErrorCode ierr; 337 PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy; 338 339 PetscFunctionBegin; 340 if (!dummy) { 341 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",0,&viewer);CHKERRQ(ierr); 342 } 343 if (fgnorm > 1.e-9) { 344 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES Function norm %G \n",its,fgnorm);CHKERRQ(ierr); 345 } else if (fgnorm > 1.e-11){ 346 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,fgnorm);CHKERRQ(ierr); 347 } else { 348 ierr = PetscViewerASCIIMonitorPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);CHKERRQ(ierr); 349 } 350 if (!dummy) { 351 ierr = PetscViewerASCIIMonitorDestroy(&viewer);CHKERRQ(ierr); 352 } 353 PetscFunctionReturn(0); 354 } 355 /* ---------------------------------------------------------------- */ 356 #undef __FUNCT__ 357 #define __FUNCT__ "SNESDefaultConverged" 358 /*@C 359 SNESDefaultConverged - Convergence test of the solvers for 360 systems of nonlinear equations (default). 361 362 Collective on SNES 363 364 Input Parameters: 365 + snes - the SNES context 366 . it - the iteration (0 indicates before any Newton steps) 367 . xnorm - 2-norm of current iterate 368 . pnorm - 2-norm of current step 369 . fnorm - 2-norm of function at current iterate 370 - dummy - unused context 371 372 Output Parameter: 373 . reason - one of 374 $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol), 375 $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 376 $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 377 $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 378 $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 379 $ SNES_CONVERGED_ITERATING - (otherwise), 380 381 where 382 + maxf - maximum number of function evaluations, 383 set with SNESSetTolerances() 384 . nfct - number of function evaluations, 385 . abstol - absolute function norm tolerance, 386 set with SNESSetTolerances() 387 - rtol - relative function norm tolerance, set with SNESSetTolerances() 388 389 Level: intermediate 390 391 .keywords: SNES, nonlinear, default, converged, convergence 392 393 .seealso: SNESSetConvergenceTest() 394 @*/ 395 PetscErrorCode SNESDefaultConverged(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 396 { 397 PetscErrorCode ierr; 398 399 PetscFunctionBegin; 400 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 401 PetscValidPointer(reason,6); 402 403 *reason = SNES_CONVERGED_ITERATING; 404 405 if (!it) { 406 /* set parameter for default relative tolerance convergence test */ 407 snes->ttol = fnorm*snes->rtol; 408 } 409 if (fnorm != fnorm) { 410 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 411 *reason = SNES_DIVERGED_FNORM_NAN; 412 } else if (fnorm < snes->abstol) { 413 ierr = PetscInfo2(snes,"Converged due to function norm %G < %G\n",fnorm,snes->abstol);CHKERRQ(ierr); 414 *reason = SNES_CONVERGED_FNORM_ABS; 415 } else if (snes->nfuncs >= snes->max_funcs) { 416 ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 417 *reason = SNES_DIVERGED_FUNCTION_COUNT; 418 } 419 420 if (it && !*reason) { 421 if (fnorm <= snes->ttol) { 422 ierr = PetscInfo2(snes,"Converged due to function norm %G < %G (relative tolerance)\n",fnorm,snes->ttol);CHKERRQ(ierr); 423 *reason = SNES_CONVERGED_FNORM_RELATIVE; 424 } else if (pnorm < snes->xtol*xnorm) { 425 ierr = PetscInfo3(snes,"Converged due to small update length: %G < %G * %G\n",pnorm,snes->xtol,xnorm);CHKERRQ(ierr); 426 *reason = SNES_CONVERGED_PNORM_RELATIVE; 427 } 428 } 429 PetscFunctionReturn(0); 430 } 431 432 #undef __FUNCT__ 433 #define __FUNCT__ "SNESSkipConverged" 434 /*@C 435 SNESSkipConverged - Convergence test for SNES that NEVER returns as 436 converged, UNLESS the maximum number of iteration have been reached. 437 438 Logically Collective on SNES 439 440 Input Parameters: 441 + snes - the SNES context 442 . it - the iteration (0 indicates before any Newton steps) 443 . xnorm - 2-norm of current iterate 444 . pnorm - 2-norm of current step 445 . fnorm - 2-norm of function at current iterate 446 - dummy - unused context 447 448 Output Parameter: 449 . reason - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN 450 451 Notes: 452 Convergence is then declared after a fixed number of iterations have been used. 453 454 Level: advanced 455 456 .keywords: SNES, nonlinear, skip, converged, convergence 457 458 .seealso: SNESSetConvergenceTest() 459 @*/ 460 PetscErrorCode SNESSkipConverged(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 461 { 462 PetscErrorCode ierr; 463 464 PetscFunctionBegin; 465 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 466 PetscValidPointer(reason,6); 467 468 *reason = SNES_CONVERGED_ITERATING; 469 470 if (fnorm != fnorm) { 471 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 472 *reason = SNES_DIVERGED_FNORM_NAN; 473 } else if(it == snes->max_its) { 474 *reason = SNES_CONVERGED_ITS; 475 } 476 PetscFunctionReturn(0); 477 } 478 479 #undef __FUNCT__ 480 #define __FUNCT__ "SNESDefaultGetWork" 481 /* 482 SNESDefaultGetWork - Gets a number of work vectors. 483 484 Input Parameters: 485 . snes - the SNES context 486 . nw - number of work vectors to allocate 487 488 Level: developer 489 490 Notes: 491 Call this only if no work vectors have been allocated 492 */ 493 PetscErrorCode SNESDefaultGetWork(SNES snes,PetscInt nw) 494 { 495 PetscErrorCode ierr; 496 497 PetscFunctionBegin; 498 if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 499 snes->nwork = nw; 500 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 501 ierr = PetscLogObjectParents(snes,nw,snes->work);CHKERRQ(ierr); 502 PetscFunctionReturn(0); 503 } 504