Actual source code: snesut.c

  1: #include <petsc/private/snesimpl.h>
  2: #include <petscdm.h>
  3: #include <petscdmshell.h>
  4: #include <petscsection.h>
  5: #include <petscblaslapack.h>

  7: /*@C
  8:   SNESMonitorSolution - Monitors progress of a `SNES` `SNESSolve()` by calling
  9:   `VecView()` for the approximate solution at each iteration.

 11:   Collective

 13:   Input Parameters:
 14: + snes   - the `SNES` context
 15: . its    - iteration number
 16: . fgnorm - 2-norm of residual
 17: - vf     - a viewer

 19:   Options Database Key:
 20: . -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration

 22:   Level: intermediate

 24:   Note:
 25:   This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
 26:   to be used during the `SNESSolve()`

 28: .seealso: [](ch_snes), `SNES`, `SNESMonitorSet()`, `SNESMonitorDefault()`, `VecView()`
 29: @*/
 30: PetscErrorCode SNESMonitorSolution(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
 31: {
 32:   Vec         x;
 33:   PetscViewer viewer = vf->viewer;

 35:   PetscFunctionBegin;
 37:   PetscCall(SNESGetSolution(snes, &x));
 38:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
 39:   PetscCall(VecView(x, viewer));
 40:   PetscCall(PetscViewerPopFormat(viewer));
 41:   PetscFunctionReturn(PETSC_SUCCESS);
 42: }

 44: /*@C
 45:   SNESMonitorResidual - Monitors progress of a `SNESSolve()` by calling
 46:   `VecView()` for the residual at each iteration.

 48:   Collective

 50:   Input Parameters:
 51: + snes   - the `SNES` context
 52: . its    - iteration number
 53: . fgnorm - 2-norm of residual
 54: - vf     - a viewer

 56:   Options Database Key:
 57: . -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration

 59:   Level: intermediate

 61:   Note:
 62:   This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
 63:   to be used during the `SNES` solve.

 65: .seealso: [](ch_snes), `SNES`, `SNESMonitorSet()`, `SNESMonitorDefault()`, `VecView()`, `SNESMonitor()`
 66: @*/
 67: PetscErrorCode SNESMonitorResidual(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
 68: {
 69:   Vec x;

 71:   PetscFunctionBegin;
 73:   PetscCall(SNESGetFunction(snes, &x, NULL, NULL));
 74:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
 75:   PetscCall(VecView(x, vf->viewer));
 76:   PetscCall(PetscViewerPopFormat(vf->viewer));
 77:   PetscFunctionReturn(PETSC_SUCCESS);
 78: }

 80: /*@C
 81:   SNESMonitorSolutionUpdate - Monitors progress of a `SNESSolve()` by calling
 82:   `VecView()` for the UPDATE to the solution at each iteration.

 84:   Collective

 86:   Input Parameters:
 87: + snes   - the `SNES` context
 88: . its    - iteration number
 89: . fgnorm - 2-norm of residual
 90: - vf     - a viewer

 92:   Options Database Key:
 93: . -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration

 95:   Level: intermediate

 97:   Note:
 98:   This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
 99:   to be used during the `SNES` solve.

101: .seealso: [](ch_snes), `SNESMonitorSet()`, `SNESMonitorDefault()`, `VecView()`, `SNESMonitor()`
102: @*/
103: PetscErrorCode SNESMonitorSolutionUpdate(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
104: {
105:   Vec         x;
106:   PetscViewer viewer = vf->viewer;

108:   PetscFunctionBegin;
110:   PetscCall(SNESGetSolutionUpdate(snes, &x));
111:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
112:   PetscCall(VecView(x, viewer));
113:   PetscCall(PetscViewerPopFormat(viewer));
114:   PetscFunctionReturn(PETSC_SUCCESS);
115: }

117: #include <petscdraw.h>

119: /*@C
120:   KSPMonitorSNESResidual - Prints the `SNES` residual norm, as well as the `KSP` residual norm, at each iteration of a `KSPSolve()` called within a `SNESSolve()`.

122:   Collective

124:   Input Parameters:
125: + ksp   - iterative context
126: . n     - iteration number
127: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
128: - vf    - The viewer context

130:   Options Database Key:
131: . -snes_monitor_ksp - Activates `KSPMonitorSNESResidual()`

133:   Level: intermediate

135:   Note:
136:   This is not called directly by users, rather one calls `KSPMonitorSet()`, with this function as an argument, to cause the monitor
137:   to be used during the `KSP` solve.

139: .seealso: [](ch_snes), `SNES`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`, `KSPMonitor()`, `SNESMonitor()`, `PetscViewerAndFormat()`
140: @*/
141: PetscErrorCode KSPMonitorSNESResidual(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
142: {
143:   PetscViewer       viewer = vf->viewer;
144:   PetscViewerFormat format = vf->format;
145:   SNES              snes   = (SNES)vf->data;
146:   Vec               snes_solution, work1, work2;
147:   PetscReal         snorm;
148:   PetscInt          tablevel;
149:   const char       *prefix;

151:   PetscFunctionBegin;
153:   PetscCall(SNESGetSolution(snes, &snes_solution));
154:   PetscCall(VecDuplicate(snes_solution, &work1));
155:   PetscCall(VecDuplicate(snes_solution, &work2));
156:   PetscCall(KSPBuildSolution(ksp, work1, NULL));
157:   PetscCall(VecAYPX(work1, -1.0, snes_solution));
158:   PetscCall(SNESComputeFunction(snes, work1, work2));
159:   PetscCall(VecNorm(work2, NORM_2, &snorm));
160:   PetscCall(VecDestroy(&work1));
161:   PetscCall(VecDestroy(&work2));

163:   PetscCall(PetscObjectGetTabLevel((PetscObject)ksp, &tablevel));
164:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ksp, &prefix));
165:   PetscCall(PetscViewerPushFormat(viewer, format));
166:   PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
167:   if (n == 0 && prefix) PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual norms for %s solve.\n", prefix));
168:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Residual norm %5.3e KSP Residual norm %5.3e\n", n, (double)snorm, (double)rnorm));
169:   PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
170:   PetscCall(PetscViewerPopFormat(viewer));
171:   PetscFunctionReturn(PETSC_SUCCESS);
172: }

174: /*@C
175:   KSPMonitorSNESResidualDrawLG - Plots the linear `KSP` residual norm and the `SNES` residual norm of a `KSPSolve()` called within a `SNESSolve()`.

177:   Collective

179:   Input Parameters:
180: + ksp   - iterative context
181: . n     - iteration number
182: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
183: - vf    - The viewer context, created with `KSPMonitorSNESResidualDrawLGCreate()`

185:   Options Database Key:
186: . -snes_monitor_ksp draw::draw_lg - Activates `KSPMonitorSNESResidualDrawLG()`

188:   Level: intermediate

190:   Note:
191:   This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
192:   to be used during the `SNESSolve()`

194: .seealso: [](ch_snes), `KSPMonitorSet()`, `KSPMonitorTrueResidual()`, `SNESMonitor()`, `KSPMonitor()`, `KSPMonitorSNESResidualDrawLGCreate()`
195: @*/
196: PetscErrorCode KSPMonitorSNESResidualDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
197: {
198:   PetscViewer        viewer = vf->viewer;
199:   PetscViewerFormat  format = vf->format;
200:   PetscDrawLG        lg;
201:   SNES               snes = (SNES)vf->data;
202:   Vec                snes_solution, work1, work2;
203:   PetscReal          snorm;
204:   KSPConvergedReason reason;
205:   PetscReal          x[2], y[2];

207:   PetscFunctionBegin;
209:   PetscCall(PetscViewerDrawGetDrawLG(viewer, 0, &lg));
210:   PetscCall(SNESGetSolution(snes, &snes_solution));
211:   PetscCall(VecDuplicate(snes_solution, &work1));
212:   PetscCall(VecDuplicate(snes_solution, &work2));
213:   PetscCall(KSPBuildSolution(ksp, work1, NULL));
214:   PetscCall(VecAYPX(work1, -1.0, snes_solution));
215:   PetscCall(SNESComputeFunction(snes, work1, work2));
216:   PetscCall(VecNorm(work2, NORM_2, &snorm));
217:   PetscCall(VecDestroy(&work1));
218:   PetscCall(VecDestroy(&work2));

220:   PetscCall(PetscViewerPushFormat(viewer, format));
221:   if (!n) PetscCall(PetscDrawLGReset(lg));
222:   x[0] = (PetscReal)n;
223:   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
224:   else y[0] = -15.0;
225:   x[1] = (PetscReal)n;
226:   if (snorm > 0.0) y[1] = PetscLog10Real(snorm);
227:   else y[1] = -15.0;
228:   PetscCall(PetscDrawLGAddPoint(lg, x, y));
229:   PetscCall(KSPGetConvergedReason(ksp, &reason));
230:   if (n <= 20 || !(n % 5) || reason) {
231:     PetscCall(PetscDrawLGDraw(lg));
232:     PetscCall(PetscDrawLGSave(lg));
233:   }
234:   PetscCall(PetscViewerPopFormat(viewer));
235:   PetscFunctionReturn(PETSC_SUCCESS);
236: }

238: /*@C
239:   KSPMonitorSNESResidualDrawLGCreate - Creates the `PetscViewer` used by `KSPMonitorSNESResidualDrawLG()`

241:   Collective

243:   Input Parameters:
244: + viewer - The `PetscViewer`
245: . format - The viewer format
246: - ctx    - An optional user context

248:   Output Parameter:
249: . vf - The viewer context

251:   Level: intermediate

253: .seealso: [](ch_snes), `KSP`, `SNES`, `PetscViewerFormat`, `PetscViewerAndFormat`, `KSPMonitorSet()`, `KSPMonitorTrueResidual()`
254: @*/
255: PetscErrorCode KSPMonitorSNESResidualDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, PetscCtx ctx, PetscViewerAndFormat **vf)
256: {
257:   const char *names[] = {"linear", "nonlinear"};

259:   PetscFunctionBegin;
260:   PetscCall(PetscViewerAndFormatCreate(viewer, format, vf));
261:   (*vf)->data = ctx;
262:   PetscCall(PetscViewerMonitorLGSetUp(viewer, NULL, NULL, "Log Residual Norm", 2, names, PETSC_DECIDE, PETSC_DECIDE, 400, 300));
263:   PetscFunctionReturn(PETSC_SUCCESS);
264: }

266: PetscErrorCode SNESMonitorDefaultSetUp(SNES snes, PetscViewerAndFormat *vf)
267: {
268:   PetscFunctionBegin;
269:   if (vf->format == PETSC_VIEWER_DRAW_LG) PetscCall(PetscViewerMonitorLGSetUp(vf->viewer, NULL, NULL, "Log Residual Norm", 1, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300));
270:   PetscFunctionReturn(PETSC_SUCCESS);
271: }

273: /*@C
274:   SNESMonitorDefault - Monitors progress of a `SNESSolve()` (default).

276:   Collective

278:   Input Parameters:
279: + snes   - the `SNES` context
280: . its    - iteration number
281: . fgnorm - 2-norm of residual
282: - vf     - viewer and format structure

284:   Options Database Key:
285: . -snes_monitor - use this function to monitor the convergence of the nonlinear solver

287:   Level: intermediate

289:   Notes:
290:   Prints the residual norm at each iteration.

292:   This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
293:   to be used during the `SNES` solve.

295: .seealso: [](ch_snes), `SNESMonitorSet()`, `SNESMonitorSolution()`, `SNESMonitorFunction()`, `SNESMonitorResidual()`,
296:           `SNESMonitorSolutionUpdate()`, `SNESMonitorScaling()`, `SNESMonitorRange()`, `SNESMonitorRatio()`,
297:           `SNESMonitorDefaultField()`, `PetscViewerFormat`, `PetscViewerAndFormat`
298: @*/
299: PetscErrorCode SNESMonitorDefault(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
300: {
301:   PetscViewer       viewer = vf->viewer;
302:   PetscViewerFormat format = vf->format;
303:   PetscBool         isascii, isdraw;

305:   PetscFunctionBegin;
307:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
308:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
309:   PetscCall(PetscViewerPushFormat(viewer, format));
310:   if (isascii) {
311:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
312:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
313:       Vec              dx;
314:       PetscReal        upnorm;
315:       SNESObjectiveFn *objective;

317:       PetscCall(SNESGetSolutionUpdate(snes, &dx));
318:       PetscCall(VecNorm(dx, NORM_2, &upnorm));
319:       PetscCall(SNESGetObjective(snes, &objective, NULL));
320:       if (objective) {
321:         Vec       x;
322:         PetscReal obj;

324:         PetscCall(SNESGetSolution(snes, &x));
325:         PetscCall(SNESComputeObjective(snes, x, &obj));
326:         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));
327:       } else {
328:         PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e, Update norm %14.12e\n", its, (double)fgnorm, (double)upnorm));
329:       }
330:     } else {
331:       PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e\n", its, (double)fgnorm));
332:     }
333:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
334:   } else if (isdraw) {
335:     if (format == PETSC_VIEWER_DRAW_LG) {
336:       PetscDrawLG lg;
337:       PetscReal   x, y;

339:       PetscCall(PetscViewerDrawGetDrawLG(viewer, 0, &lg));
340:       if (!its) PetscCall(PetscDrawLGReset(lg));
341:       x = (PetscReal)its;
342:       if (fgnorm > 0.0) y = PetscLog10Real(fgnorm);
343:       else y = -15.0;
344:       PetscCall(PetscDrawLGAddPoint(lg, &x, &y));
345:       if (its <= 20 || !(its % 5) || snes->reason) {
346:         PetscCall(PetscDrawLGDraw(lg));
347:         PetscCall(PetscDrawLGSave(lg));
348:       }
349:     }
350:   }
351:   PetscCall(PetscViewerPopFormat(viewer));
352:   PetscFunctionReturn(PETSC_SUCCESS);
353: }

355: /*@C
356:   SNESMonitorScaling - Monitors the largest value in each row of the Jacobian of a `SNESSolve()`

358:   Collective

360:   Input Parameters:
361: + snes   - the `SNES` context
362: . its    - iteration number
363: . fgnorm - 2-norm of residual
364: - vf     - viewer and format structure

366:   Level: intermediate

368:   Notes:
369:   This routine prints the largest value in each row of the Jacobian

371:   This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
372:   to be used during the `SNES` solve.

374: .seealso: [](ch_snes), `SNESMonitorSet()`, `SNESMonitorSolution()`, `SNESMonitorRange()`, `SNESMonitorJacUpdateSpectrum()`,
375:           `PetscViewerFormat`, `PetscViewerAndFormat`
376: @*/
377: PetscErrorCode SNESMonitorScaling(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
378: {
379:   PetscViewer viewer = vf->viewer;
380:   KSP         ksp;
381:   Mat         J;
382:   Vec         v;

384:   PetscFunctionBegin;
386:   PetscCall(SNESGetKSP(snes, &ksp));
387:   PetscCall(KSPGetOperators(ksp, &J, NULL));
388:   PetscCall(MatCreateVecs(J, &v, NULL));
389:   PetscCall(MatGetRowMaxAbs(J, v, NULL));
390:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
391:   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
392:   PetscCall(PetscViewerASCIIPrintf(viewer, "SNES Jacobian maximum row entries\n"));
393:   PetscCall(VecView(v, viewer));
394:   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
395:   PetscCall(PetscViewerPopFormat(viewer));
396:   PetscCall(VecDestroy(&v));
397:   PetscFunctionReturn(PETSC_SUCCESS);
398: }

400: /*@C
401:   SNESMonitorJacUpdateSpectrum - Monitors the spectrun of the change in the Jacobian from the last Jacobian evaluation of a `SNESSolve()`

403:   Collective

405:   Input Parameters:
406: + snes  - the `SNES` context
407: . it    - iteration number
408: . fnorm - 2-norm of residual
409: - vf    - viewer and format structure

411:   Options Database Key:
412: . -snes_monitor_jacupdate_spectrum - activates this monitor

414:   Level: intermediate

416:   Notes:
417:   This routine prints the eigenvalues of the difference in the Jacobians

419:   This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
420:   to be used during the `SNES` solve.

422: .seealso: [](ch_snes), `SNESMonitorSet()`, `SNESMonitorSolution()`, `SNESMonitorRange()`, `PetscViewerFormat`, `PetscViewerAndFormat`
423: @*/
424: PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes, PetscInt it, PetscReal fnorm, PetscViewerAndFormat *vf)
425: {
426:   Vec X;
427:   Mat J, dJ, dJdense;
428:   PetscErrorCode (*func)(SNES, Vec, Mat, Mat, void *);
429:   PetscInt     n;
430:   PetscBLASInt nb = 0, lwork;
431:   PetscReal   *eigr, *eigi;
432:   PetscScalar *work;
433:   PetscScalar *a;

435:   PetscFunctionBegin;
436:   if (it == 0) PetscFunctionReturn(PETSC_SUCCESS);
437:   /* create the difference between the current update and the current Jacobian */
438:   PetscCall(SNESGetSolution(snes, &X));
439:   PetscCall(SNESGetJacobian(snes, NULL, &J, &func, NULL));
440:   PetscCall(MatDuplicate(J, MAT_COPY_VALUES, &dJ));
441:   PetscCall(SNESComputeJacobian(snes, X, dJ, dJ));
442:   PetscCall(MatAXPY(dJ, -1.0, J, SAME_NONZERO_PATTERN));

444:   /* compute the spectrum directly */
445:   PetscCall(MatConvert(dJ, MATSEQDENSE, MAT_INITIAL_MATRIX, &dJdense));
446:   PetscCall(MatGetSize(dJ, &n, NULL));
447:   PetscCall(PetscBLASIntCast(n, &nb));
448:   lwork = 3 * nb;
449:   PetscCall(PetscMalloc1(n, &eigr));
450:   PetscCall(PetscMalloc1(n, &eigi));
451:   PetscCall(PetscMalloc1(lwork, &work));
452:   PetscCall(MatDenseGetArray(dJdense, &a));
453: #if !defined(PETSC_USE_COMPLEX)
454:   {
455:     PetscBLASInt lierr;
456:     PetscInt     i;
457:     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
458:     PetscCallBLAS("LAPACKgeev", LAPACKgeev_("N", "N", &nb, a, &nb, eigr, eigi, NULL, &nb, NULL, &nb, work, &lwork, &lierr));
459:     PetscCheck(!lierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "geev() error %" PetscBLASInt_FMT, lierr);
460:     PetscCall(PetscFPTrapPop());
461:     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)snes), "Eigenvalues of J_%" PetscInt_FMT " - J_%" PetscInt_FMT ":\n", it, it - 1));
462:     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]));
463:   }
464:   PetscCall(MatDenseRestoreArray(dJdense, &a));
465:   PetscCall(MatDestroy(&dJ));
466:   PetscCall(MatDestroy(&dJdense));
467:   PetscCall(PetscFree(eigr));
468:   PetscCall(PetscFree(eigi));
469:   PetscCall(PetscFree(work));
470:   PetscFunctionReturn(PETSC_SUCCESS);
471: #else
472:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not coded for complex");
473: #endif
474: }

476: PETSC_INTERN PetscErrorCode SNESMonitorRange_Private(SNES, PetscInt, PetscReal *);

478: PetscErrorCode SNESMonitorRange_Private(SNES snes, PetscInt it, PetscReal *per)
479: {
480:   Vec          resid;
481:   PetscReal    rmax, pwork;
482:   PetscInt     i, n, N;
483:   PetscScalar *r;

485:   PetscFunctionBegin;
486:   PetscCall(SNESGetFunction(snes, &resid, NULL, NULL));
487:   PetscCall(VecNorm(resid, NORM_INFINITY, &rmax));
488:   PetscCall(VecGetLocalSize(resid, &n));
489:   PetscCall(VecGetSize(resid, &N));
490:   PetscCall(VecGetArray(resid, &r));
491:   pwork = 0.0;
492:   for (i = 0; i < n; i++) pwork += (PetscAbsScalar(r[i]) > .20 * rmax);
493:   PetscCallMPI(MPIU_Allreduce(&pwork, per, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)snes)));
494:   PetscCall(VecRestoreArray(resid, &r));
495:   *per = *per / N;
496:   PetscFunctionReturn(PETSC_SUCCESS);
497: }

499: /*@C
500:   SNESMonitorRange - Prints the percentage of residual elements that are more than 10 percent of the maximum entry in the residual in each iteration of a `SNESSolve()`

502:   Collective

504:   Input Parameters:
505: + snes  - `SNES` iterative context
506: . it    - iteration number
507: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
508: - vf    - unused monitor context

510:   Options Database Key:
511: . -snes_monitor_range - Activates `SNESMonitorRange()`

513:   Level: intermediate

515:   Note:
516:   This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
517:   to be used during the `SNES` solve.

519: .seealso: [](ch_snes), `SNESMonitorSet()`, `SNESMonitorDefault()`, `SNESMonitorLGCreate()`, `SNESMonitorScaling()`, `PetscViewerFormat`, `PetscViewerAndFormat`
520: @*/
521: PetscErrorCode SNESMonitorRange(SNES snes, PetscInt it, PetscReal rnorm, PetscViewerAndFormat *vf)
522: {
523:   PetscReal   perc, rel;
524:   PetscViewer viewer = vf->viewer;
525:   /* should be in a MonitorRangeContext */
526:   static PetscReal prev;

528:   PetscFunctionBegin;
530:   if (!it) prev = rnorm;
531:   PetscCall(SNESMonitorRange_Private(snes, it, &perc));

533:   rel  = (prev - rnorm) / prev;
534:   prev = rnorm;
535:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
536:   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
537:   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 * perc), (double)rel, (double)(rel / perc)));
538:   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
539:   PetscCall(PetscViewerPopFormat(viewer));
540:   PetscFunctionReturn(PETSC_SUCCESS);
541: }

543: /*@C
544:   SNESMonitorRatio - Monitors progress of a `SNESSolve()` by printing the ratio of residual norm at each iteration to the previous.

546:   Collective

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

554:   Options Database Key:
555: . -snes_monitor_ratio - activate this monitor

557:   Level: intermediate

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.

563:   Be sure to call `SNESMonitorRationSetUp()` before using this monitor.

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;

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: }

588: /*@C
589:   SNESMonitorRatioSetUp - Insures the `SNES` object is saving its history since this monitor needs access to it

591:   Collective

593:   Input Parameters:
594: + snes - the `SNES` context
595: - vf   - `PetscViewerAndFormat` (ignored)

597:   Level: intermediate

599: .seealso: [](ch_snes), `SNESMonitorSet()`, `SNESMonitorSolution()`, `SNESMonitorDefault()`, `SNESMonitorRatio()`, `PetscViewerFormat`, `PetscViewerAndFormat`
600: @*/
601: PetscErrorCode SNESMonitorRatioSetUp(SNES snes, PetscViewerAndFormat *vf)
602: {
603:   PetscReal *history;

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: }

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.

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;

624:   PetscFunctionBegin;
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: }

640: /*@C
641:   SNESMonitorDefaultField - Monitors progress of a `SNESSolve()`, separated into fields.

643:   Collective

645:   Input Parameters:
646: + snes   - the `SNES` context
647: . its    - iteration number
648: . fgnorm - 2-norm of residual
649: - vf     - the PetscViewer

651:   Options Database Key:
652: . -snes_monitor_field - activate this monitor

654:   Level: intermediate

656:   Notes:
657:   This routine uses the `DM` attached to the residual vector to define the fields.

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.

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;

672:   PetscFunctionBegin;
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;

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 (PetscInt 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: }

702: /*@C
703:   SNESConvergedDefault - Default convergence test for `SNESSolve()`.

705:   Collective

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: - ctx   - unused context

715:   Output Parameter:
716: . reason - converged reason, see `SNESConvergedReason`

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, use `unlimited` for no maximum
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

728:   Level: developer

730:   Notes:
731:   This routine is not generally called directly. It is set with `SNESSetConvergenceTest()` automatically before the `SNESSolve()`.

733:   It can be called within a custom convergence test that should also apply the standard convergence tests

735: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESSetConvergenceTest()`, `SNESConvergedSkip()`, `SNESSetTolerances()`, `SNESSetDivergenceTolerance()`,
736:           `SNESConvergedReason`
737: @*/
738: PetscErrorCode SNESConvergedDefault(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, PetscCtx ctx)
739: {
740:   PetscFunctionBegin;
742:   PetscAssertPointer(reason, 6);

744:   *reason = SNES_CONVERGED_ITERATING;
745:   if (!it) {
746:     /* set parameter for default relative tolerance convergence test */
747:     snes->ttol   = fnorm * snes->rtol;
748:     snes->rnorm0 = fnorm;
749:   }
750:   if (PetscIsInfOrNanReal(fnorm)) {
751:     PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
752:     *reason = SNES_DIVERGED_FUNCTION_NANORINF;
753:   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
754:     PetscCall(PetscInfo(snes, "Converged due to function norm %14.12e < %14.12e\n", (double)fnorm, (double)snes->abstol));
755:     *reason = SNES_CONVERGED_FNORM_ABS;
756:   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
757:     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs));
758:     *reason = SNES_DIVERGED_FUNCTION_COUNT;
759:   }

761:   if (it && !*reason) {
762:     if (fnorm <= snes->ttol) {
763:       PetscCall(PetscInfo(snes, "Converged due to function norm %14.12e < %14.12e (relative tolerance)\n", (double)fnorm, (double)snes->ttol));
764:       *reason = SNES_CONVERGED_FNORM_RELATIVE;
765:     } else if (snorm < snes->stol * xnorm) {
766:       PetscCall(PetscInfo(snes, "Converged due to small update length: %14.12e < %14.12e * %14.12e\n", (double)snorm, (double)snes->stol, (double)xnorm));
767:       *reason = SNES_CONVERGED_SNORM_RELATIVE;
768:     } else if (snes->divtol != PETSC_UNLIMITED && (fnorm > snes->divtol * snes->rnorm0)) {
769:       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));
770:       *reason = SNES_DIVERGED_DTOL;
771:     }
772:   }
773:   PetscFunctionReturn(PETSC_SUCCESS);
774: }

776: /*@C
777:   SNESConvergedSkip - Convergence test for `SNES` that NEVER returns as
778:   converged, UNLESS the maximum number of iteration have been reached.

780:   Logically Collective

782:   Input Parameters:
783: + snes  - the `SNES` context
784: . it    - the iteration (0 indicates before any Newton steps)
785: . xnorm - 2-norm of current iterate
786: . snorm - 2-norm of current step
787: . fnorm - 2-norm of function at current iterate
788: - ctx   - unused context

790:   Output Parameter:
791: . reason - `SNES_CONVERGED_ITERATING`, `SNES_CONVERGED_ITS`, or `SNES_DIVERGED_FUNCTION_NANORINF`

793:   Options Database Key:
794: . -snes_convergence_test skip - see `SNESSetFromOptions()`

796:   Level: advanced

798:   Note:
799:   This is often used if `snes` is being used as a nonlinear smoother in `SNESFAS` or possibly other `SNESType`

801: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESConvergedDefault()`, `SNESSetConvergenceTest()`, `SNESConvergedReason`
802: @*/
803: PetscErrorCode SNESConvergedSkip(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, PetscCtx ctx)
804: {
805:   PetscFunctionBegin;
807:   PetscAssertPointer(reason, 6);

809:   *reason = SNES_CONVERGED_ITERATING;

811:   if (fnorm != fnorm) {
812:     PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
813:     *reason = SNES_DIVERGED_FUNCTION_NANORINF;
814:   } else if (it == snes->max_its) {
815:     *reason = SNES_CONVERGED_ITS;
816:   }
817:   PetscFunctionReturn(PETSC_SUCCESS);
818: }

820: /*@
821:   SNESSetWorkVecs - Allocates a number of work vectors to be used internally by the `SNES` solver

823:   Input Parameters:
824: + snes - the `SNES` context
825: - nw   - number of work vectors to allocate

827:   Level: developer

829:   Note:
830:   Each `SNESType` calls this with the number of work vectors that particular type needs.

832: .seealso: [](ch_snes), `SNES`
833: @*/
834: PetscErrorCode SNESSetWorkVecs(SNES snes, PetscInt nw)
835: {
836:   DM        dm;
837:   Vec       v;
838:   PetscBool restore = PETSC_FALSE;

840:   PetscFunctionBegin;
843:   if (snes->work) PetscCall(VecDestroyVecs(snes->nwork, &snes->work));
844:   snes->nwork = nw;
845:   if (!nw) PetscFunctionReturn(PETSC_SUCCESS);

847:   PetscCall(SNESGetDM(snes, &dm));
848:   if (snes->dmAuto) {
849:     PetscCall(DMShellGetGlobalVector(dm, &v));
850:     if (!v) v = snes->vec_sol;
851:     if (!v) v = snes->vec_func;
852:     if (!v) v = snes->vec_rhs;
853:   } else {
854:     restore = PETSC_TRUE;
855:     PetscCall(DMGetGlobalVector(dm, &v));
856:   }
857:   PetscCheck(v, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Vector to be duplicated not found");
858:   PetscCall(VecDuplicateVecs(v, snes->nwork, &snes->work));
859:   if (restore) PetscCall(DMRestoreGlobalVector(dm, &v));
860:   PetscFunctionReturn(PETSC_SUCCESS);
861: }