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