1 2 #include <petsc/private/snesimpl.h> /*I "petsc/private/snesimpl.h" I*/ 3 #include <petscdm.h> 4 #include <petscblaslapack.h> 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "SNESMonitorSolution" 8 /*@C 9 SNESMonitorSolution - Monitors progress of the SNES solvers by calling 10 VecView() for the approximate solution at each iteration. 11 12 Collective on SNES 13 14 Input Parameters: 15 + snes - the SNES context 16 . its - iteration number 17 . fgnorm - 2-norm of residual 18 - dummy - a viewer 19 20 Level: intermediate 21 22 .keywords: SNES, nonlinear, vector, monitor, view 23 24 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView() 25 @*/ 26 PetscErrorCode SNESMonitorSolution(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 27 { 28 PetscErrorCode ierr; 29 Vec x; 30 PetscViewer viewer = (PetscViewer) dummy; 31 32 PetscFunctionBegin; 33 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 34 ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr); 35 ierr = VecView(x,viewer);CHKERRQ(ierr); 36 PetscFunctionReturn(0); 37 } 38 39 #undef __FUNCT__ 40 #define __FUNCT__ "SNESMonitorResidual" 41 /*@C 42 SNESMonitorResidual - Monitors progress of the SNES solvers by calling 43 VecView() for the residual at each iteration. 44 45 Collective on SNES 46 47 Input Parameters: 48 + snes - the SNES context 49 . its - iteration number 50 . fgnorm - 2-norm of residual 51 - dummy - a viewer 52 53 Level: intermediate 54 55 .keywords: SNES, nonlinear, vector, monitor, view 56 57 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView() 58 @*/ 59 PetscErrorCode SNESMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 60 { 61 PetscErrorCode ierr; 62 Vec x; 63 PetscViewer viewer = (PetscViewer) dummy; 64 65 PetscFunctionBegin; 66 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 67 ierr = SNESGetFunction(snes,&x,0,0);CHKERRQ(ierr); 68 ierr = VecView(x,viewer);CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 #undef __FUNCT__ 73 #define __FUNCT__ "SNESMonitorSolutionUpdate" 74 /*@C 75 SNESMonitorSolutionUpdate - Monitors progress of the SNES solvers by calling 76 VecView() for the UPDATE to the solution at each iteration. 77 78 Collective on SNES 79 80 Input Parameters: 81 + snes - the SNES context 82 . its - iteration number 83 . fgnorm - 2-norm of residual 84 - dummy - a viewer 85 86 Level: intermediate 87 88 .keywords: SNES, nonlinear, vector, monitor, view 89 90 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView() 91 @*/ 92 PetscErrorCode SNESMonitorSolutionUpdate(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 93 { 94 PetscErrorCode ierr; 95 Vec x; 96 PetscViewer viewer = (PetscViewer) dummy; 97 98 PetscFunctionBegin; 99 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 100 ierr = SNESGetSolutionUpdate(snes,&x);CHKERRQ(ierr); 101 ierr = VecView(x,viewer);CHKERRQ(ierr); 102 PetscFunctionReturn(0); 103 } 104 105 #undef __FUNCT__ 106 #define __FUNCT__ "KSPMonitorSNES" 107 /*@C 108 KSPMonitorSNES - Print the residual norm of the nonlinear function at each iteration of the linear iterative solver. 109 110 Collective on KSP 111 112 Input Parameters: 113 + ksp - iterative context 114 . n - iteration number 115 . rnorm - 2-norm (preconditioned) residual value (may be estimated). 116 - dummy - unused monitor context 117 118 Level: intermediate 119 120 .keywords: KSP, default, monitor, residual 121 122 .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGResidualNormCreate() 123 @*/ 124 PetscErrorCode KSPMonitorSNES(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy) 125 { 126 PetscErrorCode ierr; 127 PetscViewer viewer; 128 SNES snes = (SNES) dummy; 129 Vec snes_solution,work1,work2; 130 PetscReal snorm; 131 132 PetscFunctionBegin; 133 ierr = SNESGetSolution(snes,&snes_solution);CHKERRQ(ierr); 134 ierr = VecDuplicate(snes_solution,&work1);CHKERRQ(ierr); 135 ierr = VecDuplicate(snes_solution,&work2);CHKERRQ(ierr); 136 ierr = KSPBuildSolution(ksp,work1,NULL);CHKERRQ(ierr); 137 ierr = VecAYPX(work1,-1.0,snes_solution);CHKERRQ(ierr); 138 ierr = SNESComputeFunction(snes,work1,work2);CHKERRQ(ierr); 139 ierr = VecNorm(work2,NORM_2,&snorm);CHKERRQ(ierr); 140 ierr = VecDestroy(&work1);CHKERRQ(ierr); 141 ierr = VecDestroy(&work2);CHKERRQ(ierr); 142 143 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);CHKERRQ(ierr); 144 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);CHKERRQ(ierr); 145 if (n == 0 && ((PetscObject)ksp)->prefix) { 146 ierr = PetscViewerASCIIPrintf(viewer," Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);CHKERRQ(ierr); 147 } 148 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Residual norm %5.3e KSP Residual norm %5.3e \n",n,(double)snorm,(double)rnorm);CHKERRQ(ierr); 149 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);CHKERRQ(ierr); 150 PetscFunctionReturn(0); 151 } 152 153 #include <petscdraw.h> 154 155 #undef __FUNCT__ 156 #define __FUNCT__ "KSPMonitorSNESLGResidualNormCreate" 157 /*@C 158 KSPMonitorSNESLGResidualNormCreate - Creates a line graph context for use with 159 KSP to monitor convergence of preconditioned residual norms. 160 161 Collective on KSP 162 163 Input Parameters: 164 + host - the X display to open, or null for the local machine 165 . label - the title to put in the title bar 166 . x, y - the screen coordinates of the upper left coordinate of 167 the window 168 - m, n - the screen width and height in pixels 169 170 Output Parameter: 171 . draw - the drawing context 172 173 Options Database Key: 174 . -ksp_monitor_lg_residualnorm - Sets line graph monitor 175 176 Notes: 177 Use KSPMonitorSNESLGResidualNormDestroy() to destroy this line graph; do not use PetscDrawLGDestroy(). 178 179 Level: intermediate 180 181 .keywords: KSP, monitor, line graph, residual, create 182 183 .seealso: KSPMonitorSNESLGResidualNormDestroy(), KSPMonitorSet(), KSPMonitorSNESLGTrueResidualCreate() 184 @*/ 185 PetscErrorCode KSPMonitorSNESLGResidualNormCreate(const char host[],const char label[],int x,int y,int m,int n,PetscObject **objs) 186 { 187 PetscDraw draw; 188 PetscErrorCode ierr; 189 PetscDrawAxis axis; 190 PetscDrawLG drawlg; 191 const char *names[] = {"Linear residual","Nonlinear residual"}; 192 193 PetscFunctionBegin; 194 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&draw);CHKERRQ(ierr); 195 ierr = PetscDrawSetFromOptions(draw);CHKERRQ(ierr); 196 ierr = PetscDrawLGCreate(draw,2,&drawlg);CHKERRQ(ierr); 197 ierr = PetscDrawLGSetFromOptions(drawlg);CHKERRQ(ierr); 198 ierr = PetscDrawLGGetAxis(drawlg,&axis);CHKERRQ(ierr); 199 ierr = PetscDrawAxisSetLabels(axis,"Convergence of Residual Norm","Iteration","Residual Norm");CHKERRQ(ierr); 200 ierr = PetscDrawLGSetLegend(drawlg,names);CHKERRQ(ierr); 201 ierr = PetscLogObjectParent((PetscObject)drawlg,(PetscObject)draw);CHKERRQ(ierr); 202 203 ierr = PetscMalloc1(3,objs);CHKERRQ(ierr); 204 (*objs)[1] = (PetscObject)drawlg; 205 (*objs)[2] = (PetscObject)draw; 206 PetscFunctionReturn(0); 207 } 208 209 #undef __FUNCT__ 210 #define __FUNCT__ "KSPMonitorSNESLGResidualNorm" 211 PetscErrorCode KSPMonitorSNESLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscObject *objs) 212 { 213 PetscDrawLG lg = (PetscDrawLG) objs[1]; 214 PetscErrorCode ierr; 215 PetscReal y[2]; 216 SNES snes = (SNES) objs[0]; 217 Vec snes_solution,work1,work2; 218 219 PetscFunctionBegin; 220 if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm); 221 else y[0] = -15.0; 222 223 ierr = SNESGetSolution(snes,&snes_solution);CHKERRQ(ierr); 224 ierr = VecDuplicate(snes_solution,&work1);CHKERRQ(ierr); 225 ierr = VecDuplicate(snes_solution,&work2);CHKERRQ(ierr); 226 ierr = KSPBuildSolution(ksp,work1,NULL);CHKERRQ(ierr); 227 ierr = VecAYPX(work1,-1.0,snes_solution);CHKERRQ(ierr); 228 ierr = SNESComputeFunction(snes,work1,work2);CHKERRQ(ierr); 229 ierr = VecNorm(work2,NORM_2,y+1);CHKERRQ(ierr); 230 if (y[1] > 0.0) y[1] = PetscLog10Real(y[1]); 231 else y[1] = -15.0; 232 ierr = VecDestroy(&work1);CHKERRQ(ierr); 233 ierr = VecDestroy(&work2);CHKERRQ(ierr); 234 235 ierr = PetscDrawLGAddPoint(lg,NULL,y);CHKERRQ(ierr); 236 if (n < 20 || !(n % 5)) { 237 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 238 } 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "KSPMonitorSNESLGResidualNormDestroy" 244 /*@ 245 KSPMonitorSNESLGResidualNormDestroy - Destroys a line graph context that was created 246 with KSPMonitorSNESLGResidualNormCreate(). 247 248 Collective on KSP 249 250 Input Parameter: 251 . draw - the drawing context 252 253 Level: intermediate 254 255 .keywords: KSP, monitor, line graph, destroy 256 257 .seealso: KSPMonitorSNESLGResidualNormCreate(), KSPMonitorSNESLGTrueResidualDestroy(), KSPMonitorSet() 258 @*/ 259 PetscErrorCode KSPMonitorSNESLGResidualNormDestroy(PetscObject **objs) 260 { 261 PetscErrorCode ierr; 262 PetscDrawLG drawlg = (PetscDrawLG) (*objs)[1]; 263 PetscDraw draw = (PetscDraw) (*objs)[2]; 264 265 PetscFunctionBegin; 266 ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 267 ierr = PetscDrawLGDestroy(&drawlg);CHKERRQ(ierr); 268 ierr = PetscFree(*objs);CHKERRQ(ierr); 269 PetscFunctionReturn(0); 270 } 271 272 #undef __FUNCT__ 273 #define __FUNCT__ "SNESMonitorDefault" 274 /*@C 275 SNESMonitorDefault - Monitors progress of the SNES solvers (default). 276 277 Collective on SNES 278 279 Input Parameters: 280 + snes - the SNES context 281 . its - iteration number 282 . fgnorm - 2-norm of residual 283 - dummy - unused context 284 285 Notes: 286 This routine prints the residual norm at each iteration. 287 288 Level: intermediate 289 290 .keywords: SNES, nonlinear, default, monitor, norm 291 292 .seealso: SNESMonitorSet(), SNESMonitorSolution() 293 @*/ 294 PetscErrorCode SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 295 { 296 PetscErrorCode ierr; 297 PetscViewer viewer = (PetscViewer) dummy; 298 299 PetscFunctionBegin; 300 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 301 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 302 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr); 303 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 304 PetscFunctionReturn(0); 305 } 306 307 #undef __FUNCT__ 308 #define __FUNCT__ "SNESMonitorJacUpdateSpectrum" 309 PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,void *ctx) 310 { 311 #if defined(PETSC_MISSING_LAPACK_GEEV) 312 SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values."); 313 #elif defined(PETSC_HAVE_ESSL) 314 SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - No support for ESSL Lapack Routines"); 315 #else 316 Vec X; 317 Mat J,dJ,dJdense; 318 PetscErrorCode ierr; 319 PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*); 320 PetscInt n,i; 321 PetscBLASInt nb,lwork; 322 PetscReal *eigr,*eigi; 323 PetscScalar *work; 324 PetscScalar *a; 325 326 PetscFunctionBegin; 327 if (it == 0) PetscFunctionReturn(0); 328 /* create the difference between the current update and the current jacobian */ 329 ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr); 330 ierr = SNESGetJacobian(snes,NULL,&J,&func,NULL);CHKERRQ(ierr); 331 ierr = MatDuplicate(J,MAT_COPY_VALUES,&dJ);CHKERRQ(ierr); 332 ierr = SNESComputeJacobian(snes,X,dJ,dJ);CHKERRQ(ierr); 333 ierr = MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 334 335 /* compute the spectrum directly */ 336 ierr = MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);CHKERRQ(ierr); 337 ierr = MatGetSize(dJ,&n,NULL);CHKERRQ(ierr); 338 ierr = PetscBLASIntCast(n,&nb);CHKERRQ(ierr); 339 lwork = 3*nb; 340 ierr = PetscMalloc1(n,&eigr);CHKERRQ(ierr); 341 ierr = PetscMalloc1(n,&eigi);CHKERRQ(ierr); 342 ierr = PetscMalloc1(lwork,&work);CHKERRQ(ierr); 343 ierr = MatDenseGetArray(dJdense,&a);CHKERRQ(ierr); 344 #if !defined(PETSC_USE_COMPLEX) 345 { 346 PetscBLASInt lierr; 347 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 348 PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,NULL,&nb,NULL,&nb,work,&lwork,&lierr)); 349 if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr); 350 ierr = PetscFPTrapPop();CHKERRQ(ierr); 351 } 352 #else 353 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex"); 354 #endif 355 PetscPrintf(PetscObjectComm((PetscObject)snes),"Eigenvalues of J_%d - J_%d:\n",it,it-1);CHKERRQ(ierr); 356 for (i=0;i<n;i++) { 357 PetscPrintf(PetscObjectComm((PetscObject)snes),"%5d: %20.5g + %20.5gi\n",i,(double)eigr[i],(double)eigi[i]);CHKERRQ(ierr); 358 } 359 ierr = MatDenseRestoreArray(dJdense,&a);CHKERRQ(ierr); 360 ierr = MatDestroy(&dJ);CHKERRQ(ierr); 361 ierr = MatDestroy(&dJdense);CHKERRQ(ierr); 362 ierr = PetscFree(eigr);CHKERRQ(ierr); 363 ierr = PetscFree(eigi);CHKERRQ(ierr); 364 ierr = PetscFree(work);CHKERRQ(ierr); 365 PetscFunctionReturn(0); 366 #endif 367 } 368 369 #undef __FUNCT__ 370 #define __FUNCT__ "SNESMonitorRange_Private" 371 PetscErrorCode SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per) 372 { 373 PetscErrorCode ierr; 374 Vec resid; 375 PetscReal rmax,pwork; 376 PetscInt i,n,N; 377 PetscScalar *r; 378 379 PetscFunctionBegin; 380 ierr = SNESGetFunction(snes,&resid,0,0);CHKERRQ(ierr); 381 ierr = VecNorm(resid,NORM_INFINITY,&rmax);CHKERRQ(ierr); 382 ierr = VecGetLocalSize(resid,&n);CHKERRQ(ierr); 383 ierr = VecGetSize(resid,&N);CHKERRQ(ierr); 384 ierr = VecGetArray(resid,&r);CHKERRQ(ierr); 385 pwork = 0.0; 386 for (i=0; i<n; i++) { 387 pwork += (PetscAbsScalar(r[i]) > .20*rmax); 388 } 389 ierr = MPI_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 390 ierr = VecRestoreArray(resid,&r);CHKERRQ(ierr); 391 *per = *per/N; 392 PetscFunctionReturn(0); 393 } 394 395 #undef __FUNCT__ 396 #define __FUNCT__ "SNESMonitorRange" 397 /*@C 398 SNESMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value. 399 400 Collective on SNES 401 402 Input Parameters: 403 + snes - iterative context 404 . it - iteration number 405 . rnorm - 2-norm (preconditioned) residual value (may be estimated). 406 - dummy - unused monitor context 407 408 Options Database Key: 409 . -snes_monitor_range - Activates SNESMonitorRange() 410 411 Level: intermediate 412 413 .keywords: SNES, default, monitor, residual 414 415 .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate() 416 @*/ 417 PetscErrorCode SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,void *dummy) 418 { 419 PetscErrorCode ierr; 420 PetscReal perc,rel; 421 PetscViewer viewer = (PetscViewer) dummy; 422 /* should be in a MonitorRangeContext */ 423 static PetscReal prev; 424 425 PetscFunctionBegin; 426 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 427 if (!it) prev = rnorm; 428 ierr = SNESMonitorRange_Private(snes,it,&perc);CHKERRQ(ierr); 429 430 rel = (prev - rnorm)/prev; 431 prev = rnorm; 432 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 433 ierr = PetscViewerASCIIPrintf(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); 434 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 435 PetscFunctionReturn(0); 436 } 437 438 typedef struct { 439 PetscViewer viewer; 440 PetscReal *history; 441 } SNESMonitorRatioContext; 442 443 #undef __FUNCT__ 444 #define __FUNCT__ "SNESMonitorRatio" 445 /*@C 446 SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio 447 of residual norm at each iteration to the previous. 448 449 Collective on SNES 450 451 Input Parameters: 452 + snes - the SNES context 453 . its - iteration number 454 . fgnorm - 2-norm of residual (or gradient) 455 - dummy - context of monitor 456 457 Level: intermediate 458 459 .keywords: SNES, nonlinear, monitor, norm 460 461 .seealso: SNESMonitorSet(), SNESMonitorSolution() 462 @*/ 463 PetscErrorCode SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 464 { 465 PetscErrorCode ierr; 466 PetscInt len; 467 PetscReal *history; 468 SNESMonitorRatioContext *ctx = (SNESMonitorRatioContext*)dummy; 469 470 PetscFunctionBegin; 471 ierr = SNESGetConvergenceHistory(snes,&history,NULL,&len);CHKERRQ(ierr); 472 ierr = PetscViewerASCIIAddTab(ctx->viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 473 if (!its || !history || its > len) { 474 ierr = PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr); 475 } else { 476 PetscReal ratio = fgnorm/history[its-1]; 477 ierr = PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);CHKERRQ(ierr); 478 } 479 ierr = PetscViewerASCIISubtractTab(ctx->viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 480 PetscFunctionReturn(0); 481 } 482 483 /* 484 If the we set the history monitor space then we must destroy it 485 */ 486 #undef __FUNCT__ 487 #define __FUNCT__ "SNESMonitorRatioDestroy" 488 PetscErrorCode SNESMonitorRatioDestroy(void **ct) 489 { 490 PetscErrorCode ierr; 491 SNESMonitorRatioContext *ctx = *(SNESMonitorRatioContext**)ct; 492 493 PetscFunctionBegin; 494 ierr = PetscFree(ctx->history);CHKERRQ(ierr); 495 ierr = PetscViewerDestroy(&ctx->viewer);CHKERRQ(ierr); 496 ierr = PetscFree(ctx);CHKERRQ(ierr); 497 PetscFunctionReturn(0); 498 } 499 500 #undef __FUNCT__ 501 #define __FUNCT__ "SNESMonitorSetRatio" 502 /*@C 503 SNESMonitorSetRatio - Sets SNES to use a monitor that prints the 504 ratio of the function norm at each iteration. 505 506 Collective on SNES 507 508 Input Parameters: 509 + snes - the SNES context 510 - viewer - ASCII viewer to print output 511 512 Level: intermediate 513 514 .keywords: SNES, nonlinear, monitor, norm 515 516 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault() 517 @*/ 518 PetscErrorCode SNESMonitorSetRatio(SNES snes,PetscViewer viewer) 519 { 520 PetscErrorCode ierr; 521 SNESMonitorRatioContext *ctx; 522 PetscReal *history; 523 524 PetscFunctionBegin; 525 ierr = PetscNewLog(snes,&ctx);CHKERRQ(ierr); 526 ierr = SNESGetConvergenceHistory(snes,&history,NULL,NULL);CHKERRQ(ierr); 527 if (!history) { 528 ierr = PetscMalloc1(100,&ctx->history);CHKERRQ(ierr); 529 ierr = SNESSetConvergenceHistory(snes,ctx->history,0,100,PETSC_TRUE);CHKERRQ(ierr); 530 } 531 ctx->viewer = viewer; 532 533 ierr = SNESMonitorSet(snes,SNESMonitorRatio,ctx,SNESMonitorRatioDestroy);CHKERRQ(ierr); 534 PetscFunctionReturn(0); 535 } 536 537 /* ---------------------------------------------------------------- */ 538 #undef __FUNCT__ 539 #define __FUNCT__ "SNESMonitorDefaultShort" 540 /* 541 Default (short) SNES Monitor, same as SNESMonitorDefault() except 542 it prints fewer digits of the residual as the residual gets smaller. 543 This is because the later digits are meaningless and are often 544 different on different machines; by using this routine different 545 machines will usually generate the same output. 546 */ 547 PetscErrorCode SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy) 548 { 549 PetscErrorCode ierr; 550 PetscViewer viewer = (PetscViewer) dummy; 551 552 PetscFunctionBegin; 553 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 554 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 555 if (fgnorm > 1.e-9) { 556 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %g \n",its,(double)fgnorm);CHKERRQ(ierr); 557 } else if (fgnorm > 1.e-11) { 558 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,(double)fgnorm);CHKERRQ(ierr); 559 } else { 560 ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);CHKERRQ(ierr); 561 } 562 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); 563 PetscFunctionReturn(0); 564 } 565 566 #undef __FUNCT__ 567 #define __FUNCT__ "SNESMonitorDefaultField" 568 /*@C 569 SNESMonitorDefaultField - Monitors progress of the SNES solvers, separated into fields. 570 571 Collective on SNES 572 573 Input Parameters: 574 + snes - the SNES context 575 . its - iteration number 576 . fgnorm - 2-norm of residual 577 - ctx - the PetscViewer 578 579 Notes: 580 This routine uses the DM attached to the residual vector 581 582 Level: intermediate 583 584 .keywords: SNES, nonlinear, field, monitor, norm 585 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorDefaultShort() 586 @*/ 587 PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, void *ctx) 588 { 589 PetscViewer viewer = (PetscViewer) ctx; 590 Vec r; 591 DM dm; 592 PetscReal res[256]; 593 PetscInt tablevel; 594 PetscErrorCode ierr; 595 596 PetscFunctionBegin; 597 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 598 ierr = SNESGetFunction(snes, &r, NULL, NULL);CHKERRQ(ierr); 599 ierr = VecGetDM(r, &dm);CHKERRQ(ierr); 600 if (!dm) {ierr = SNESMonitorDefault(snes, its, fgnorm, ctx);CHKERRQ(ierr);} 601 else { 602 PetscSection s, gs; 603 PetscInt Nf, f; 604 605 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 606 ierr = DMGetDefaultGlobalSection(dm, &gs);CHKERRQ(ierr); 607 if (!s || !gs) {ierr = SNESMonitorDefault(snes, its, fgnorm, ctx);CHKERRQ(ierr);} 608 ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr); 609 if (Nf > 256) SETERRQ1(PetscObjectComm((PetscObject) snes), PETSC_ERR_SUP, "Do not support %d fields > 256", Nf); 610 ierr = PetscSectionVecNorm(s, gs, r, NORM_2, res);CHKERRQ(ierr); 611 ierr = PetscObjectGetTabLevel((PetscObject) snes, &tablevel);CHKERRQ(ierr); 612 ierr = PetscViewerASCIIAddTab(viewer, tablevel);CHKERRQ(ierr); 613 ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr); 614 for (f = 0; f < Nf; ++f) { 615 if (f) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 616 ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", res[f]);CHKERRQ(ierr); 617 } 618 ierr = PetscViewerASCIIPrintf(viewer, "] \n");CHKERRQ(ierr); 619 ierr = PetscViewerASCIISubtractTab(viewer, tablevel);CHKERRQ(ierr); 620 } 621 PetscFunctionReturn(0); 622 } 623 /* ---------------------------------------------------------------- */ 624 #undef __FUNCT__ 625 #define __FUNCT__ "SNESConvergedDefault" 626 /*@C 627 SNESConvergedDefault - Convergence test of the solvers for 628 systems of nonlinear equations (default). 629 630 Collective on SNES 631 632 Input Parameters: 633 + snes - the SNES context 634 . it - the iteration (0 indicates before any Newton steps) 635 . xnorm - 2-norm of current iterate 636 . snorm - 2-norm of current step 637 . fnorm - 2-norm of function at current iterate 638 - dummy - unused context 639 640 Output Parameter: 641 . reason - one of 642 $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol), 643 $ SNES_CONVERGED_SNORM_RELATIVE - (snorm < stol*xnorm), 644 $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 645 $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 646 $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 647 $ SNES_CONVERGED_ITERATING - (otherwise), 648 649 where 650 + maxf - maximum number of function evaluations, 651 set with SNESSetTolerances() 652 . nfct - number of function evaluations, 653 . abstol - absolute function norm tolerance, 654 set with SNESSetTolerances() 655 - rtol - relative function norm tolerance, set with SNESSetTolerances() 656 657 Level: intermediate 658 659 .keywords: SNES, nonlinear, default, converged, convergence 660 661 .seealso: SNESSetConvergenceTest() 662 @*/ 663 PetscErrorCode SNESConvergedDefault(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 664 { 665 PetscErrorCode ierr; 666 667 PetscFunctionBegin; 668 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 669 PetscValidPointer(reason,6); 670 671 *reason = SNES_CONVERGED_ITERATING; 672 673 if (!it) { 674 /* set parameter for default relative tolerance convergence test */ 675 snes->ttol = fnorm*snes->rtol; 676 } 677 if (PetscIsInfOrNanReal(fnorm)) { 678 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 679 *reason = SNES_DIVERGED_FNORM_NAN; 680 } else if (fnorm < snes->abstol) { 681 ierr = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr); 682 *reason = SNES_CONVERGED_FNORM_ABS; 683 } else if (snes->nfuncs >= snes->max_funcs) { 684 ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr); 685 *reason = SNES_DIVERGED_FUNCTION_COUNT; 686 } 687 688 if (it && !*reason) { 689 if (fnorm <= snes->ttol) { 690 ierr = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr); 691 *reason = SNES_CONVERGED_FNORM_RELATIVE; 692 } else if (snorm < snes->stol*xnorm) { 693 ierr = PetscInfo3(snes,"Converged due to small update length: %14.12e < %14.12e * %14.12e\n",(double)snorm,(double)snes->stol,(double)xnorm);CHKERRQ(ierr); 694 *reason = SNES_CONVERGED_SNORM_RELATIVE; 695 } 696 } 697 PetscFunctionReturn(0); 698 } 699 700 #undef __FUNCT__ 701 #define __FUNCT__ "SNESConvergedSkip" 702 /*@C 703 SNESConvergedSkip - Convergence test for SNES that NEVER returns as 704 converged, UNLESS the maximum number of iteration have been reached. 705 706 Logically Collective on SNES 707 708 Input Parameters: 709 + snes - the SNES context 710 . it - the iteration (0 indicates before any Newton steps) 711 . xnorm - 2-norm of current iterate 712 . snorm - 2-norm of current step 713 . fnorm - 2-norm of function at current iterate 714 - dummy - unused context 715 716 Output Parameter: 717 . reason - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN 718 719 Notes: 720 Convergence is then declared after a fixed number of iterations have been used. 721 722 Level: advanced 723 724 .keywords: SNES, nonlinear, skip, converged, convergence 725 726 .seealso: SNESSetConvergenceTest() 727 @*/ 728 PetscErrorCode SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 729 { 730 PetscErrorCode ierr; 731 732 PetscFunctionBegin; 733 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 734 PetscValidPointer(reason,6); 735 736 *reason = SNES_CONVERGED_ITERATING; 737 738 if (fnorm != fnorm) { 739 ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr); 740 *reason = SNES_DIVERGED_FNORM_NAN; 741 } else if (it == snes->max_its) { 742 *reason = SNES_CONVERGED_ITS; 743 } 744 PetscFunctionReturn(0); 745 } 746 747 #undef __FUNCT__ 748 #define __FUNCT__ "SNESSetWorkVecs" 749 /*@C 750 SNESSetWorkVecs - Gets a number of work vectors. 751 752 Input Parameters: 753 . snes - the SNES context 754 . nw - number of work vectors to allocate 755 756 Level: developer 757 758 Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations 759 760 @*/ 761 PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw) 762 { 763 DM dm; 764 Vec v; 765 PetscErrorCode ierr; 766 767 PetscFunctionBegin; 768 if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 769 snes->nwork = nw; 770 771 ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 772 ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr); 773 ierr = VecDuplicateVecs(v,snes->nwork,&snes->work);CHKERRQ(ierr); 774 ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr); 775 ierr = PetscLogObjectParents(snes,nw,snes->work);CHKERRQ(ierr); 776 PetscFunctionReturn(0); 777 } 778