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