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