Actual source code: snesut.c
1: #include <petsc/private/snesimpl.h>
2: #include <petscdm.h>
3: #include <petscsection.h>
4: #include <petscblaslapack.h>
6: /*@C
7: SNESMonitorSolution - Monitors progress of a `SNES` `SNESSolve()` by calling
8: `VecView()` for the approximate solution at each iteration.
10: Collective
12: Input Parameters:
13: + snes - the `SNES` context
14: . its - iteration number
15: . fgnorm - 2-norm of residual
16: - vf - a viewer
18: Options Database Key:
19: . -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration
21: Level: intermediate
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()`
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;
34: PetscFunctionBegin;
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: }
43: /*@C
44: SNESMonitorResidual - Monitors progress of a `SNESSolve()` by calling
45: `VecView()` for the residual at each iteration.
47: Collective
49: Input Parameters:
50: + snes - the `SNES` context
51: . its - iteration number
52: . fgnorm - 2-norm of residual
53: - vf - a viewer
55: Options Database Key:
56: . -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration
58: Level: intermediate
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.
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;
70: PetscFunctionBegin;
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: }
79: /*@C
80: SNESMonitorSolutionUpdate - Monitors progress of a `SNESSolve()` by calling
81: `VecView()` for the UPDATE to the solution at each iteration.
83: Collective
85: Input Parameters:
86: + snes - the `SNES` context
87: . its - iteration number
88: . fgnorm - 2-norm of residual
89: - vf - a viewer
91: Options Database Key:
92: . -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration
94: Level: intermediate
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.
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;
107: PetscFunctionBegin;
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: }
116: #include <petscdraw.h>
118: /*@C
119: KSPMonitorSNESResidual - Prints the `SNES` residual norm, as well as the `KSP` residual norm, at each iteration of a `KSPSolve()` called within a `SNESSolve()`.
121: Collective
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
129: Options Database Key:
130: . -snes_monitor_ksp - Activates `KSPMonitorSNESResidual()`
132: Level: intermediate
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.
138: .seealso: [](ch_snes), `SNES`, `KSPMonitorSet()`, `KSPMonitorResidual()`, `KSPMonitorTrueResidualMaxNorm()`, `KSPMonitor()`, `SNESMonitor()`, `PetscViewerAndFormat()`
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;
150: PetscFunctionBegin;
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));
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: }
173: /*@C
174: KSPMonitorSNESResidualDrawLG - Plots the linear `KSP` residual norm and the `SNES` residual norm of a `KSPSolve()` called within a `SNESSolve()`.
176: Collective
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()`
184: Options Database Key:
185: . -snes_monitor_ksp draw::draw_lg - Activates `KSPMonitorSNESResidualDrawLG()`
187: Level: intermediate
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()`
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];
206: PetscFunctionBegin;
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));
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: }
237: /*@C
238: KSPMonitorSNESResidualDrawLGCreate - Creates the `PetscViewer` used by `KSPMonitorSNESResidualDrawLG()`
240: Collective
242: Input Parameters:
243: + viewer - The `PetscViewer`
244: . format - The viewer format
245: - ctx - An optional user context
247: Output Parameter:
248: . vf - The viewer context
250: Level: intermediate
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"};
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: }
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: }
272: /*@C
273: SNESMonitorDefault - Monitors progress of a `SNESSolve()` (default).
275: Collective
277: Input Parameters:
278: + snes - the `SNES` context
279: . its - iteration number
280: . fgnorm - 2-norm of residual
281: - vf - viewer and format structure
283: Options Database Key:
284: . -snes_monitor - use this function to monitor the convergence of the nonlinear solver
286: Level: intermediate
288: Notes:
289: Prints the residual norm at each iteration.
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.
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;
304: PetscFunctionBegin;
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: SNESObjectiveFn *objective;
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;
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;
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: }
354: /*@C
355: SNESMonitorScaling - Monitors the largest value in each row of the Jacobian of a `SNESSolve()`
357: Collective
359: Input Parameters:
360: + snes - the `SNES` context
361: . its - iteration number
362: . fgnorm - 2-norm of residual
363: - vf - viewer and format structure
365: Level: intermediate
367: Notes:
368: This routine prints the largest value in each row of the Jacobian
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.
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;
383: PetscFunctionBegin;
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: }
399: /*@C
400: SNESMonitorJacUpdateSpectrum - Monitors the spectrun of the change in the Jacobian from the last Jacobian evaluation of a `SNESSolve()`
402: Collective
404: Input Parameters:
405: + snes - the `SNES` context
406: . it - iteration number
407: . fnorm - 2-norm of residual
408: - vf - viewer and format structure
410: Options Database Key:
411: . -snes_monitor_jacupdate_spectrum - activates this monitor
413: Level: intermediate
415: Notes:
416: This routine prints the eigenvalues of the difference in the Jacobians
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.
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;
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));
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: }
475: PETSC_INTERN PetscErrorCode SNESMonitorRange_Private(SNES, PetscInt, PetscReal *);
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;
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: PetscCallMPI(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: }
498: /*@C
499: 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()`
501: Collective
503: Input Parameters:
504: + snes - `SNES` iterative context
505: . it - iteration number
506: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
507: - vf - unused monitor context
509: Options Database Key:
510: . -snes_monitor_range - Activates `SNESMonitorRange()`
512: Level: intermediate
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.
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;
527: PetscFunctionBegin;
529: if (!it) prev = rnorm;
530: PetscCall(SNESMonitorRange_Private(snes, it, &perc));
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: }
542: /*@C
543: SNESMonitorRatio - Monitors progress of a `SNESSolve()` by printing the ratio of residual norm at each iteration to the previous.
545: Collective
547: Input Parameters:
548: + snes - the `SNES` context
549: . its - iteration number
550: . fgnorm - 2-norm of residual (or gradient)
551: - vf - context of monitor
553: Options Database Key:
554: . -snes_monitor_ratio - activate this monitor
556: Level: intermediate
558: Notes:
559: This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
560: to be used during the `SNES` solve.
562: Be sure to call `SNESMonitorRationSetUp()` before using this monitor.
564: .seealso: [](ch_snes), `SNESMonitorRationSetUp()`, `SNESMonitorSet()`, `SNESMonitorSolution()`, `SNESMonitorDefault()`, `PetscViewerFormat`, `PetscViewerAndFormat`
565: @*/
566: PetscErrorCode SNESMonitorRatio(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
567: {
568: PetscInt len;
569: PetscReal *history;
570: PetscViewer viewer = vf->viewer;
572: PetscFunctionBegin;
573: PetscCall(SNESGetConvergenceHistory(snes, &history, NULL, &len));
574: PetscCall(PetscViewerPushFormat(viewer, vf->format));
575: PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
576: if (!its || !history || its > len) {
577: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e\n", its, (double)fgnorm));
578: } else {
579: PetscReal ratio = fgnorm / history[its - 1];
580: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e %14.12e\n", its, (double)fgnorm, (double)ratio));
581: }
582: PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
583: PetscCall(PetscViewerPopFormat(viewer));
584: PetscFunctionReturn(PETSC_SUCCESS);
585: }
587: /*@C
588: SNESMonitorRatioSetUp - Insures the `SNES` object is saving its history since this monitor needs access to it
590: Collective
592: Input Parameters:
593: + snes - the `SNES` context
594: - vf - `PetscViewerAndFormat` (ignored)
596: Level: intermediate
598: .seealso: [](ch_snes), `SNESMonitorSet()`, `SNESMonitorSolution()`, `SNESMonitorDefault()`, `SNESMonitorRatio()`, `PetscViewerFormat`, `PetscViewerAndFormat`
599: @*/
600: PetscErrorCode SNESMonitorRatioSetUp(SNES snes, PetscViewerAndFormat *vf)
601: {
602: PetscReal *history;
604: PetscFunctionBegin;
605: PetscCall(SNESGetConvergenceHistory(snes, &history, NULL, NULL));
606: if (!history) PetscCall(SNESSetConvergenceHistory(snes, NULL, NULL, 100, PETSC_TRUE));
607: PetscFunctionReturn(PETSC_SUCCESS);
608: }
610: /*
611: Default (short) SNES Monitor, same as SNESMonitorDefault() except
612: it prints fewer digits of the residual as the residual gets smaller.
613: This is because the later digits are meaningless and are often
614: different on different machines; by using this routine different
615: machines will usually generate the same output.
617: Deprecated: Intentionally has no manual page
618: */
619: PetscErrorCode SNESMonitorDefaultShort(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
620: {
621: PetscViewer viewer = vf->viewer;
623: PetscFunctionBegin;
625: PetscCall(PetscViewerPushFormat(viewer, vf->format));
626: PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel));
627: if (fgnorm > 1.e-9) {
628: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %g\n", its, (double)fgnorm));
629: } else if (fgnorm > 1.e-11) {
630: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %5.3e\n", its, (double)fgnorm));
631: } else {
632: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm < 1.e-11\n", its));
633: }
634: PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel));
635: PetscCall(PetscViewerPopFormat(viewer));
636: PetscFunctionReturn(PETSC_SUCCESS);
637: }
639: /*@C
640: SNESMonitorDefaultField - Monitors progress of a `SNESSolve()`, separated into fields.
642: Collective
644: Input Parameters:
645: + snes - the `SNES` context
646: . its - iteration number
647: . fgnorm - 2-norm of residual
648: - vf - the PetscViewer
650: Options Database Key:
651: . -snes_monitor_field - activate this monitor
653: Level: intermediate
655: Notes:
656: This routine uses the `DM` attached to the residual vector to define the fields.
658: This is not called directly by users, rather one calls `SNESMonitorSet()`, with this function as an argument, to cause the monitor
659: to be used during the `SNES` solve.
661: .seealso: [](ch_snes), `SNESMonitorSet()`, `SNESMonitorSolution()`, `SNESMonitorDefault()`, `PetscViewerFormat`, `PetscViewerAndFormat`
662: @*/
663: PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
664: {
665: PetscViewer viewer = vf->viewer;
666: Vec r;
667: DM dm;
668: PetscReal res[256];
669: PetscInt tablevel;
671: PetscFunctionBegin;
673: PetscCall(SNESGetFunction(snes, &r, NULL, NULL));
674: PetscCall(VecGetDM(r, &dm));
675: if (!dm) PetscCall(SNESMonitorDefault(snes, its, fgnorm, vf));
676: else {
677: PetscSection s, gs;
678: PetscInt Nf, f;
680: PetscCall(DMGetLocalSection(dm, &s));
681: PetscCall(DMGetGlobalSection(dm, &gs));
682: if (!s || !gs) PetscCall(SNESMonitorDefault(snes, its, fgnorm, vf));
683: PetscCall(PetscSectionGetNumFields(s, &Nf));
684: PetscCheck(Nf <= 256, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Do not support %" PetscInt_FMT " fields > 256", Nf);
685: PetscCall(PetscSectionVecNorm(s, gs, r, NORM_2, res));
686: PetscCall(PetscObjectGetTabLevel((PetscObject)snes, &tablevel));
687: PetscCall(PetscViewerPushFormat(viewer, vf->format));
688: PetscCall(PetscViewerASCIIAddTab(viewer, tablevel));
689: PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e [", its, (double)fgnorm));
690: for (f = 0; f < Nf; ++f) {
691: if (f) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
692: PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)res[f]));
693: }
694: PetscCall(PetscViewerASCIIPrintf(viewer, "] \n"));
695: PetscCall(PetscViewerASCIISubtractTab(viewer, tablevel));
696: PetscCall(PetscViewerPopFormat(viewer));
697: }
698: PetscFunctionReturn(PETSC_SUCCESS);
699: }
701: /*@C
702: SNESConvergedDefault - Default convergence test for `SNESSolve()`.
704: Collective
706: Input Parameters:
707: + snes - the `SNES` context
708: . it - the iteration (0 indicates before any Newton steps)
709: . xnorm - 2-norm of current iterate
710: . snorm - 2-norm of current step
711: . fnorm - 2-norm of function at current iterate
712: - dummy - unused context
714: Output Parameter:
715: . reason - converged reason, see `SNESConvergedReason`
717: Options Database Keys:
718: + -snes_convergence_test default - see `SNESSetFromOptions()`
719: . -snes_stol - convergence tolerance in terms of the norm of the change in the solution between steps
720: . -snes_atol <abstol> - absolute tolerance of residual norm
721: . -snes_rtol <rtol> - relative decrease in tolerance norm from the initial 2-norm of the solution
722: . -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
723: . -snes_max_funcs <max_funcs> - maximum number of function evaluations, use `unlimited` for no maximum
724: . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
725: - -snes_max_linear_solve_fail - number of linear solver failures before `SNESSolve()` stops
727: Level: developer
729: Notes:
730: This routine is not generally called directly. It is set with `SNESSetConvergenceTest()` automatically before the `SNESSolve()`.
732: It can be called within a custom convergence test that should also apply the standard convergence tests
734: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESSetConvergenceTest()`, `SNESConvergedSkip()`, `SNESSetTolerances()`, `SNESSetDivergenceTolerance()`,
735: `SNESConvergedReason`
736: @*/
737: PetscErrorCode SNESConvergedDefault(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
738: {
739: PetscFunctionBegin;
741: PetscAssertPointer(reason, 6);
743: *reason = SNES_CONVERGED_ITERATING;
744: if (!it) {
745: /* set parameter for default relative tolerance convergence test */
746: snes->ttol = fnorm * snes->rtol;
747: snes->rnorm0 = fnorm;
748: }
749: if (PetscIsInfOrNanReal(fnorm)) {
750: PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
751: *reason = SNES_DIVERGED_FNORM_NAN;
752: } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
753: PetscCall(PetscInfo(snes, "Converged due to function norm %14.12e < %14.12e\n", (double)fnorm, (double)snes->abstol));
754: *reason = SNES_CONVERGED_FNORM_ABS;
755: } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
756: PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", snes->nfuncs, snes->max_funcs));
757: *reason = SNES_DIVERGED_FUNCTION_COUNT;
758: }
760: if (it && !*reason) {
761: if (fnorm <= snes->ttol) {
762: PetscCall(PetscInfo(snes, "Converged due to function norm %14.12e < %14.12e (relative tolerance)\n", (double)fnorm, (double)snes->ttol));
763: *reason = SNES_CONVERGED_FNORM_RELATIVE;
764: } else if (snorm < snes->stol * xnorm) {
765: PetscCall(PetscInfo(snes, "Converged due to small update length: %14.12e < %14.12e * %14.12e\n", (double)snorm, (double)snes->stol, (double)xnorm));
766: *reason = SNES_CONVERGED_SNORM_RELATIVE;
767: } else if (snes->divtol != PETSC_UNLIMITED && (fnorm > snes->divtol * snes->rnorm0)) {
768: 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));
769: *reason = SNES_DIVERGED_DTOL;
770: }
771: }
772: PetscFunctionReturn(PETSC_SUCCESS);
773: }
775: /*@C
776: SNESConvergedSkip - Convergence test for `SNES` that NEVER returns as
777: converged, UNLESS the maximum number of iteration have been reached.
779: Logically Collective
781: Input Parameters:
782: + snes - the `SNES` context
783: . it - the iteration (0 indicates before any Newton steps)
784: . xnorm - 2-norm of current iterate
785: . snorm - 2-norm of current step
786: . fnorm - 2-norm of function at current iterate
787: - dummy - unused context
789: Output Parameter:
790: . reason - `SNES_CONVERGED_ITERATING`, `SNES_CONVERGED_ITS`, or `SNES_DIVERGED_FNORM_NAN`
792: Options Database Key:
793: . -snes_convergence_test skip - see `SNESSetFromOptions()`
795: Level: advanced
797: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESConvergedDefault()`, `SNESSetConvergenceTest()`, `SNESConvergedReason`
798: @*/
799: PetscErrorCode SNESConvergedSkip(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
800: {
801: PetscFunctionBegin;
803: PetscAssertPointer(reason, 6);
805: *reason = SNES_CONVERGED_ITERATING;
807: if (fnorm != fnorm) {
808: PetscCall(PetscInfo(snes, "Failed to converged, function norm is NaN\n"));
809: *reason = SNES_DIVERGED_FNORM_NAN;
810: } else if (it == snes->max_its) {
811: *reason = SNES_CONVERGED_ITS;
812: }
813: PetscFunctionReturn(PETSC_SUCCESS);
814: }
816: /*@
817: SNESSetWorkVecs - Allocates a number of work vectors to be used internally by `SNES` solvers
819: Input Parameters:
820: + snes - the `SNES` context
821: - nw - number of work vectors to allocate
823: Level: developer
825: Note:
826: Each `SNESType` calls this with the number of work vectors that particular type needs.
828: .seealso: [](ch_snes), `SNES`
829: @*/
830: PetscErrorCode SNESSetWorkVecs(SNES snes, PetscInt nw)
831: {
832: DM dm;
833: Vec v;
835: PetscFunctionBegin;
836: if (snes->work) PetscCall(VecDestroyVecs(snes->nwork, &snes->work));
837: snes->nwork = nw;
839: PetscCall(SNESGetDM(snes, &dm));
840: PetscCall(DMGetGlobalVector(dm, &v));
841: PetscCall(VecDuplicateVecs(v, snes->nwork, &snes->work));
842: PetscCall(DMRestoreGlobalVector(dm, &v));
843: PetscFunctionReturn(PETSC_SUCCESS);
844: }