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