Actual source code: snesut.c
2: #include <petsc/private/snesimpl.h>
3: #include <petscdm.h>
4: #include <petscsection.h>
5: #include <petscblaslapack.h>
7: /*@C
8: SNESMonitorSolution - Monitors progress of the SNES solvers by calling
9: VecView() for the approximate solution at each iteration.
11: Collective on SNES
13: Input Parameters:
14: + snes - the SNES context
15: . its - iteration number
16: . fgnorm - 2-norm of residual
17: - dummy - a viewer
19: Options Database Keys:
20: . -snes_monitor_solution [ascii binary draw][:filename][:viewer format] - plots solution at each iteration
22: Level: intermediate
24: .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
25: @*/
26: PetscErrorCode SNESMonitorSolution(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
27: {
29: Vec x;
30: PetscViewer viewer = vf->viewer;
34: SNESGetSolution(snes,&x);
35: PetscViewerPushFormat(viewer,vf->format);
36: VecView(x,viewer);
37: PetscViewerPopFormat(viewer);
38: return(0);
39: }
41: /*@C
42: SNESMonitorResidual - Monitors progress of the SNES solvers by calling
43: VecView() for the residual at each iteration.
45: Collective on SNES
47: Input Parameters:
48: + snes - the SNES context
49: . its - iteration number
50: . fgnorm - 2-norm of residual
51: - dummy - a viewer
53: Options Database Keys:
54: . -snes_monitor_residual [ascii binary draw][:filename][:viewer format] - plots residual (not its norm) at each iteration
56: Level: intermediate
58: .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
59: @*/
60: PetscErrorCode SNESMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
61: {
63: Vec x;
64: PetscViewer viewer = vf->viewer;
68: SNESGetFunction(snes,&x,NULL,NULL);
69: PetscViewerPushFormat(viewer,vf->format);
70: VecView(x,viewer);
71: PetscViewerPopFormat(viewer);
72: return(0);
73: }
75: /*@C
76: SNESMonitorSolutionUpdate - Monitors progress of the SNES solvers by calling
77: VecView() for the UPDATE to the solution at each iteration.
79: Collective on SNES
81: Input Parameters:
82: + snes - the SNES context
83: . its - iteration number
84: . fgnorm - 2-norm of residual
85: - dummy - a viewer
87: Options Database Keys:
88: . -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration
90: Level: intermediate
92: .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
93: @*/
94: PetscErrorCode SNESMonitorSolutionUpdate(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
95: {
97: Vec x;
98: PetscViewer viewer = vf->viewer;
102: SNESGetSolutionUpdate(snes,&x);
103: PetscViewerPushFormat(viewer,vf->format);
104: VecView(x,viewer);
105: PetscViewerPopFormat(viewer);
106: return(0);
107: }
109: #include <petscdraw.h>
111: /*@C
112: KSPMonitorSNESResidual - Prints the SNES residual norm, as well as the linear residual norm, at each iteration of an iterative solver.
114: Collective on ksp
116: Input Parameters:
117: + ksp - iterative context
118: . n - iteration number
119: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
120: - vf - The viewer context
122: Options Database Key:
123: . -snes_monitor_ksp - Activates KSPMonitorSNESResidual()
125: Level: intermediate
127: .seealso: KSPMonitorSet(), KSPMonitorResidual(),KSPMonitorTrueResidualMaxNorm()
128: @*/
129: PetscErrorCode KSPMonitorSNESResidual(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
130: {
131: PetscViewer viewer = vf->viewer;
132: PetscViewerFormat format = vf->format;
133: SNES snes = (SNES) vf->data;
134: Vec snes_solution, work1, work2;
135: PetscReal snorm;
136: PetscInt tablevel;
137: const char *prefix;
138: PetscErrorCode ierr;
142: SNESGetSolution(snes, &snes_solution);
143: VecDuplicate(snes_solution, &work1);
144: VecDuplicate(snes_solution, &work2);
145: KSPBuildSolution(ksp, work1, NULL);
146: VecAYPX(work1, -1.0, snes_solution);
147: SNESComputeFunction(snes, work1, work2);
148: VecNorm(work2, NORM_2, &snorm);
149: VecDestroy(&work1);
150: VecDestroy(&work2);
152: PetscObjectGetTabLevel((PetscObject) ksp, &tablevel);
153: PetscObjectGetOptionsPrefix((PetscObject) ksp, &prefix);
154: PetscViewerPushFormat(viewer, format);
155: PetscViewerASCIIAddTab(viewer, tablevel);
156: if (n == 0 && prefix) {PetscViewerASCIIPrintf(viewer, " Residual norms for %s solve.\n", prefix);}
157: PetscViewerASCIIPrintf(viewer, "%3D SNES Residual norm %5.3e KSP Residual norm %5.3e \n", n, (double) snorm, (double) rnorm);
158: PetscViewerASCIISubtractTab(viewer, tablevel);
159: PetscViewerPopFormat(viewer);
160: return(0);
161: }
163: /*@C
164: KSPMonitorSNESResidualDrawLG - Plots the linear SNES residual norm at each iteration of an iterative solver.
166: Collective on ksp
168: Input Parameters:
169: + ksp - iterative context
170: . n - iteration number
171: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
172: - vf - The viewer context
174: Options Database Key:
175: . -snes_monitor_ksp draw::draw_lg - Activates KSPMonitorSNESResidualDrawLG()
177: Level: intermediate
179: .seealso: KSPMonitorSet(), KSPMonitorTrueResidual()
180: @*/
181: PetscErrorCode KSPMonitorSNESResidualDrawLG(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf)
182: {
183: PetscViewer viewer = vf->viewer;
184: PetscViewerFormat format = vf->format;
185: PetscDrawLG lg = vf->lg;
186: SNES snes = (SNES) vf->data;
187: Vec snes_solution, work1, work2;
188: PetscReal snorm;
189: KSPConvergedReason reason;
190: PetscReal x[2], y[2];
191: PetscErrorCode ierr;
196: SNESGetSolution(snes, &snes_solution);
197: VecDuplicate(snes_solution, &work1);
198: VecDuplicate(snes_solution, &work2);
199: KSPBuildSolution(ksp, work1, NULL);
200: VecAYPX(work1, -1.0, snes_solution);
201: SNESComputeFunction(snes, work1, work2);
202: VecNorm(work2, NORM_2, &snorm);
203: VecDestroy(&work1);
204: VecDestroy(&work2);
206: PetscViewerPushFormat(viewer, format);
207: if (!n) {PetscDrawLGReset(lg);}
208: x[0] = (PetscReal) n;
209: if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
210: else y[0] = -15.0;
211: x[1] = (PetscReal) n;
212: if (snorm > 0.0) y[1] = PetscLog10Real(snorm);
213: else y[1] = -15.0;
214: PetscDrawLGAddPoint(lg, x, y);
215: KSPGetConvergedReason(ksp, &reason);
216: if (n <= 20 || !(n % 5) || reason) {
217: PetscDrawLGDraw(lg);
218: PetscDrawLGSave(lg);
219: }
220: PetscViewerPopFormat(viewer);
221: return(0);
222: }
224: /*@C
225: KSPMonitorSNESResidualDrawLGCreate - Creates the plotter for the linear SNES residual.
227: Collective on ksp
229: Input Parameters:
230: + viewer - The PetscViewer
231: . format - The viewer format
232: - ctx - An optional user context
234: Output Parameter:
235: . vf - The viewer context
237: Level: intermediate
239: .seealso: KSPMonitorSet(), KSPMonitorTrueResidual()
240: @*/
241: PetscErrorCode KSPMonitorSNESResidualDrawLGCreate(PetscViewer viewer, PetscViewerFormat format, void *ctx, PetscViewerAndFormat **vf)
242: {
243: const char *names[] = {"linear", "nonlinear"};
247: PetscViewerAndFormatCreate(viewer, format, vf);
248: (*vf)->data = ctx;
249: KSPMonitorLGCreate(PetscObjectComm((PetscObject) viewer), NULL, NULL, "Log Residual Norm", 2, names, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &(*vf)->lg);
250: return(0);
251: }
253: /*@C
254: SNESMonitorDefault - Monitors progress of the SNES solvers (default).
256: Collective on SNES
258: Input Parameters:
259: + snes - the SNES context
260: . its - iteration number
261: . fgnorm - 2-norm of residual
262: - vf - viewer and format structure
264: Notes:
265: This routine prints the residual norm at each iteration.
267: Level: intermediate
269: .seealso: SNESMonitorSet(), SNESMonitorSolution()
270: @*/
271: PetscErrorCode SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
272: {
273: PetscViewer viewer = vf->viewer;
274: PetscViewerFormat format = vf->format;
275: PetscBool isascii, isdraw;
276: PetscErrorCode ierr;
280: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
281: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);
282: PetscViewerPushFormat(viewer,format);
283: if (isascii) {
284: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
285: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
286: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
287: } else if (isdraw) {
288: if (format == PETSC_VIEWER_DRAW_LG) {
289: PetscDrawLG lg = (PetscDrawLG) vf->lg;
290: PetscReal x, y;
293: if (!its) {PetscDrawLGReset(lg);}
294: x = (PetscReal) its;
295: if (fgnorm > 0.0) y = PetscLog10Real(fgnorm);
296: else y = -15.0;
297: PetscDrawLGAddPoint(lg,&x,&y);
298: if (its <= 20 || !(its % 5) || snes->reason) {
299: PetscDrawLGDraw(lg);
300: PetscDrawLGSave(lg);
301: }
302: }
303: }
304: PetscViewerPopFormat(viewer);
305: return(0);
306: }
308: /*@C
309: SNESMonitorScaling - Monitors the largest value in each row of the Jacobian.
311: Collective on SNES
313: Input Parameters:
314: + snes - the SNES context
315: . its - iteration number
316: . fgnorm - 2-norm of residual
317: - vf - viewer and format structure
319: Notes:
320: This routine prints the largest value in each row of the Jacobian
322: Level: intermediate
324: .seealso: SNESMonitorSet(), SNESMonitorSolution()
325: @*/
326: PetscErrorCode SNESMonitorScaling(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
327: {
329: PetscViewer viewer = vf->viewer;
330: KSP ksp;
331: Mat J;
332: Vec v;
336: SNESGetKSP(snes,&ksp);
337: KSPGetOperators(ksp,&J,NULL);
338: MatCreateVecs(J,&v,NULL);
339: MatGetRowMaxAbs(J,v,NULL);
340: PetscViewerPushFormat(viewer,vf->format);
341: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
342: PetscViewerASCIIPrintf(viewer,"%3D SNES Jacobian maximum row entries \n");
343: VecView(v,viewer);
344: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
345: PetscViewerPopFormat(viewer);
346: VecDestroy(&v);
347: return(0);
348: }
350: PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,PetscViewerAndFormat *vf)
351: {
352: Vec X;
353: Mat J,dJ,dJdense;
355: PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
356: PetscInt n;
357: PetscBLASInt nb = 0,lwork;
358: PetscReal *eigr,*eigi;
359: PetscScalar *work;
360: PetscScalar *a;
363: if (it == 0) return(0);
364: /* create the difference between the current update and the current jacobian */
365: SNESGetSolution(snes,&X);
366: SNESGetJacobian(snes,NULL,&J,&func,NULL);
367: MatDuplicate(J,MAT_COPY_VALUES,&dJ);
368: SNESComputeJacobian(snes,X,dJ,dJ);
369: MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);
371: /* compute the spectrum directly */
372: MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);
373: MatGetSize(dJ,&n,NULL);
374: PetscBLASIntCast(n,&nb);
375: lwork = 3*nb;
376: PetscMalloc1(n,&eigr);
377: PetscMalloc1(n,&eigi);
378: PetscMalloc1(lwork,&work);
379: MatDenseGetArray(dJdense,&a);
380: #if !defined(PETSC_USE_COMPLEX)
381: {
382: PetscBLASInt lierr;
383: PetscInt i;
384: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
385: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,NULL,&nb,NULL,&nb,work,&lwork,&lierr));
386: if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr);
387: PetscFPTrapPop();
388: PetscPrintf(PetscObjectComm((PetscObject)snes),"Eigenvalues of J_%d - J_%d:\n",it,it-1);
389: for (i=0;i<n;i++) {
390: PetscPrintf(PetscObjectComm((PetscObject)snes),"%5d: %20.5g + %20.5gi\n",i,(double)eigr[i],(double)eigi[i]);
391: }
392: }
393: MatDenseRestoreArray(dJdense,&a);
394: MatDestroy(&dJ);
395: MatDestroy(&dJdense);
396: PetscFree(eigr);
397: PetscFree(eigi);
398: PetscFree(work);
399: return(0);
400: #else
401: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
402: #endif
403: }
405: PETSC_INTERN PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
407: PetscErrorCode SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
408: {
410: Vec resid;
411: PetscReal rmax,pwork;
412: PetscInt i,n,N;
413: PetscScalar *r;
416: SNESGetFunction(snes,&resid,NULL,NULL);
417: VecNorm(resid,NORM_INFINITY,&rmax);
418: VecGetLocalSize(resid,&n);
419: VecGetSize(resid,&N);
420: VecGetArray(resid,&r);
421: pwork = 0.0;
422: for (i=0; i<n; i++) {
423: pwork += (PetscAbsScalar(r[i]) > .20*rmax);
424: }
425: MPIU_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));
426: VecRestoreArray(resid,&r);
427: *per = *per/N;
428: return(0);
429: }
431: /*@C
432: SNESMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value.
434: Collective on SNES
436: Input Parameters:
437: + snes - iterative context
438: . it - iteration number
439: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
440: - dummy - unused monitor context
442: Options Database Key:
443: . -snes_monitor_range - Activates SNESMonitorRange()
445: Level: intermediate
447: .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
448: @*/
449: PetscErrorCode SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,PetscViewerAndFormat *vf)
450: {
452: PetscReal perc,rel;
453: PetscViewer viewer = vf->viewer;
454: /* should be in a MonitorRangeContext */
455: static PetscReal prev;
459: if (!it) prev = rnorm;
460: SNESMonitorRange_Private(snes,it,&perc);
462: rel = (prev - rnorm)/prev;
463: prev = rnorm;
464: PetscViewerPushFormat(viewer,vf->format);
465: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
466: PetscViewerASCIIPrintf(viewer,"%3D SNES preconditioned resid norm %14.12e Percent values above 20 percent of maximum %5.2f relative decrease %5.2e ratio %5.2e \n",it,(double)rnorm,(double)(100.0*perc),(double)rel,(double)(rel/perc));
467: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
468: PetscViewerPopFormat(viewer);
469: return(0);
470: }
472: /*@C
473: SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio
474: of residual norm at each iteration to the previous.
476: Collective on SNES
478: Input Parameters:
479: + snes - the SNES context
480: . its - iteration number
481: . fgnorm - 2-norm of residual (or gradient)
482: - dummy - context of monitor
484: Level: intermediate
486: Notes:
487: Insure that SNESMonitorRatio() is called when you set this monitor
488: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorRatio()
489: @*/
490: PetscErrorCode SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
491: {
492: PetscErrorCode ierr;
493: PetscInt len;
494: PetscReal *history;
495: PetscViewer viewer = vf->viewer;
498: SNESGetConvergenceHistory(snes,&history,NULL,&len);
499: PetscViewerPushFormat(viewer,vf->format);
500: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
501: if (!its || !history || its > len) {
502: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
503: } else {
504: PetscReal ratio = fgnorm/history[its-1];
505: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);
506: }
507: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
508: PetscViewerPopFormat(viewer);
509: return(0);
510: }
512: /*@C
513: SNESMonitorRatioSetUp - Insures the SNES object is saving its history since this monitor needs access to it
515: Collective on SNES
517: Input Parameters:
518: + snes - the SNES context
519: - viewer - the PetscViewer object (ignored)
521: Level: intermediate
523: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorRatio()
524: @*/
525: PetscErrorCode SNESMonitorRatioSetUp(SNES snes,PetscViewerAndFormat *vf)
526: {
527: PetscErrorCode ierr;
528: PetscReal *history;
531: SNESGetConvergenceHistory(snes,&history,NULL,NULL);
532: if (!history) {
533: SNESSetConvergenceHistory(snes,NULL,NULL,100,PETSC_TRUE);
534: }
535: return(0);
536: }
538: /* ---------------------------------------------------------------- */
539: /*
540: Default (short) SNES Monitor, same as SNESMonitorDefault() except
541: it prints fewer digits of the residual as the residual gets smaller.
542: This is because the later digits are meaningless and are often
543: different on different machines; by using this routine different
544: machines will usually generate the same output.
546: Deprecated: Intentionally has no manual page
547: */
548: PetscErrorCode SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
549: {
551: PetscViewer viewer = vf->viewer;
555: PetscViewerPushFormat(viewer,vf->format);
556: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
557: if (fgnorm > 1.e-9) {
558: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %g \n",its,(double)fgnorm);
559: } else if (fgnorm > 1.e-11) {
560: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,(double)fgnorm);
561: } else {
562: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);
563: }
564: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
565: PetscViewerPopFormat(viewer);
566: return(0);
567: }
569: /*@C
570: SNESMonitorDefaultField - Monitors progress of the SNES solvers, separated into fields.
572: Collective on SNES
574: Input Parameters:
575: + snes - the SNES context
576: . its - iteration number
577: . fgnorm - 2-norm of residual
578: - ctx - the PetscViewer
580: Notes:
581: This routine uses the DM attached to the residual vector
583: Level: intermediate
585: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault()
586: @*/
587: PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
588: {
589: PetscViewer viewer = vf->viewer;
590: Vec r;
591: DM dm;
592: PetscReal res[256];
593: PetscInt tablevel;
598: SNESGetFunction(snes, &r, NULL, NULL);
599: VecGetDM(r, &dm);
600: if (!dm) {SNESMonitorDefault(snes, its, fgnorm, vf);}
601: else {
602: PetscSection s, gs;
603: PetscInt Nf, f;
605: DMGetLocalSection(dm, &s);
606: DMGetGlobalSection(dm, &gs);
607: if (!s || !gs) {SNESMonitorDefault(snes, its, fgnorm, vf);}
608: PetscSectionGetNumFields(s, &Nf);
609: if (Nf > 256) SETERRQ1(PetscObjectComm((PetscObject) snes), PETSC_ERR_SUP, "Do not support %d fields > 256", Nf);
610: PetscSectionVecNorm(s, gs, r, NORM_2, res);
611: PetscObjectGetTabLevel((PetscObject) snes, &tablevel);
612: PetscViewerPushFormat(viewer,vf->format);
613: PetscViewerASCIIAddTab(viewer, tablevel);
614: PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);
615: for (f = 0; f < Nf; ++f) {
616: if (f) {PetscViewerASCIIPrintf(viewer, ", ");}
617: PetscViewerASCIIPrintf(viewer, "%14.12e", res[f]);
618: }
619: PetscViewerASCIIPrintf(viewer, "] \n");
620: PetscViewerASCIISubtractTab(viewer, tablevel);
621: PetscViewerPopFormat(viewer);
622: }
623: return(0);
624: }
625: /* ---------------------------------------------------------------- */
626: /*@C
627: SNESConvergedDefault - Default onvergence test of the solvers for
628: systems of nonlinear equations.
630: Collective on SNES
632: Input Parameters:
633: + snes - the SNES context
634: . it - the iteration (0 indicates before any Newton steps)
635: . xnorm - 2-norm of current iterate
636: . snorm - 2-norm of current step
637: . fnorm - 2-norm of function at current iterate
638: - dummy - unused context
640: Output Parameter:
641: . reason - one of
642: $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol),
643: $ SNES_CONVERGED_SNORM_RELATIVE - (snorm < stol*xnorm),
644: $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0),
645: $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf),
646: $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN),
647: $ SNES_CONVERGED_ITERATING - (otherwise),
648: $ SNES_DIVERGED_DTOL - (fnorm > divtol*snes->fnorm0)
650: where
651: + maxf - maximum number of function evaluations, set with SNESSetTolerances()
652: . nfct - number of function evaluations,
653: . abstol - absolute function norm tolerance, set with SNESSetTolerances()
654: . rtol - relative function norm tolerance, set with SNESSetTolerances()
655: . divtol - divergence tolerance, set with SNESSetDivergenceTolerance()
656: - fnorm0 - 2-norm of the function at the initial solution (initial guess; zeroth iteration)
658: Options Database Keys:
659: + -snes_stol - convergence tolerance in terms of the norm of the change in the solution between steps
660: . -snes_atol <abstol> - absolute tolerance of residual norm
661: . -snes_rtol <rtol> - relative decrease in tolerance norm from the initial 2-norm of the solution
662: . -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
663: . -snes_max_funcs <max_funcs> - maximum number of function evaluations
664: . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
665: - -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
667: Level: intermediate
669: .seealso: SNESSetConvergenceTest(), SNESConvergedSkip(), SNESSetTolerances(), SNESSetDivergenceTolerance()
670: @*/
671: PetscErrorCode SNESConvergedDefault(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
672: {
679: *reason = SNES_CONVERGED_ITERATING;
681: if (!it) {
682: /* set parameter for default relative tolerance convergence test */
683: snes->ttol = fnorm*snes->rtol;
684: snes->rnorm0 = fnorm;
685: }
686: if (PetscIsInfOrNanReal(fnorm)) {
687: PetscInfo(snes,"Failed to converged, function norm is NaN\n");
688: *reason = SNES_DIVERGED_FNORM_NAN;
689: } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
690: PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);
691: *reason = SNES_CONVERGED_FNORM_ABS;
692: } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
693: PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);
694: *reason = SNES_DIVERGED_FUNCTION_COUNT;
695: }
697: if (it && !*reason) {
698: if (fnorm <= snes->ttol) {
699: PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);
700: *reason = SNES_CONVERGED_FNORM_RELATIVE;
701: } else if (snorm < snes->stol*xnorm) {
702: PetscInfo3(snes,"Converged due to small update length: %14.12e < %14.12e * %14.12e\n",(double)snorm,(double)snes->stol,(double)xnorm);
703: *reason = SNES_CONVERGED_SNORM_RELATIVE;
704: } else if (snes->divtol > 0 && (fnorm > snes->divtol*snes->rnorm0)) {
705: PetscInfo3(snes,"Diverged due to increase in function norm: %14.12e > %14.12e * %14.12e\n",(double)fnorm,(double)snes->divtol,(double)snes->rnorm0);
706: *reason = SNES_DIVERGED_DTOL;
707: }
709: }
710: return(0);
711: }
713: /*@C
714: SNESConvergedSkip - Convergence test for SNES that NEVER returns as
715: converged, UNLESS the maximum number of iteration have been reached.
717: Logically Collective on SNES
719: Input Parameters:
720: + snes - the SNES context
721: . it - the iteration (0 indicates before any Newton steps)
722: . xnorm - 2-norm of current iterate
723: . snorm - 2-norm of current step
724: . fnorm - 2-norm of function at current iterate
725: - dummy - unused context
727: Output Parameter:
728: . reason - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
730: Notes:
731: Convergence is then declared after a fixed number of iterations have been used.
733: Level: advanced
735: .seealso: SNESConvergedDefault(), SNESSetConvergenceTest()
736: @*/
737: PetscErrorCode SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
738: {
745: *reason = SNES_CONVERGED_ITERATING;
747: if (fnorm != fnorm) {
748: PetscInfo(snes,"Failed to converged, function norm is NaN\n");
749: *reason = SNES_DIVERGED_FNORM_NAN;
750: } else if (it == snes->max_its) {
751: *reason = SNES_CONVERGED_ITS;
752: }
753: return(0);
754: }
756: /*@C
757: SNESSetWorkVecs - Gets a number of work vectors.
759: Input Parameters:
760: + snes - the SNES context
761: - nw - number of work vectors to allocate
763: Level: developer
765: @*/
766: PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw)
767: {
768: DM dm;
769: Vec v;
773: if (snes->work) {VecDestroyVecs(snes->nwork,&snes->work);}
774: snes->nwork = nw;
776: SNESGetDM(snes, &dm);
777: DMGetGlobalVector(dm, &v);
778: VecDuplicateVecs(v,snes->nwork,&snes->work);
779: DMRestoreGlobalVector(dm, &v);
780: PetscLogObjectParents(snes,nw,snes->work);
781: return(0);
782: }