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: #if defined(PETSC_HAVE_ESSL)
353: SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - No support for ESSL Lapack Routines");
354: #else
355: Vec X;
356: Mat J,dJ,dJdense;
358: PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
359: PetscInt n;
360: PetscBLASInt nb = 0,lwork;
361: PetscReal *eigr,*eigi;
362: PetscScalar *work;
363: PetscScalar *a;
366: if (it == 0) return(0);
367: /* create the difference between the current update and the current jacobian */
368: SNESGetSolution(snes,&X);
369: SNESGetJacobian(snes,NULL,&J,&func,NULL);
370: MatDuplicate(J,MAT_COPY_VALUES,&dJ);
371: SNESComputeJacobian(snes,X,dJ,dJ);
372: MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);
374: /* compute the spectrum directly */
375: MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);
376: MatGetSize(dJ,&n,NULL);
377: PetscBLASIntCast(n,&nb);
378: lwork = 3*nb;
379: PetscMalloc1(n,&eigr);
380: PetscMalloc1(n,&eigi);
381: PetscMalloc1(lwork,&work);
382: MatDenseGetArray(dJdense,&a);
383: #if !defined(PETSC_USE_COMPLEX)
384: {
385: PetscBLASInt lierr;
386: PetscInt i;
387: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
388: PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,NULL,&nb,NULL,&nb,work,&lwork,&lierr));
389: if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr);
390: PetscFPTrapPop();
391: PetscPrintf(PetscObjectComm((PetscObject)snes),"Eigenvalues of J_%d - J_%d:\n",it,it-1);
392: for (i=0;i<n;i++) {
393: PetscPrintf(PetscObjectComm((PetscObject)snes),"%5d: %20.5g + %20.5gi\n",i,(double)eigr[i],(double)eigi[i]);
394: }
395: }
396: MatDenseRestoreArray(dJdense,&a);
397: MatDestroy(&dJ);
398: MatDestroy(&dJdense);
399: PetscFree(eigr);
400: PetscFree(eigi);
401: PetscFree(work);
402: return(0);
403: #else
404: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
405: #endif
406: #endif
407: }
409: PETSC_INTERN PetscErrorCode SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
411: PetscErrorCode SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
412: {
414: Vec resid;
415: PetscReal rmax,pwork;
416: PetscInt i,n,N;
417: PetscScalar *r;
420: SNESGetFunction(snes,&resid,NULL,NULL);
421: VecNorm(resid,NORM_INFINITY,&rmax);
422: VecGetLocalSize(resid,&n);
423: VecGetSize(resid,&N);
424: VecGetArray(resid,&r);
425: pwork = 0.0;
426: for (i=0; i<n; i++) {
427: pwork += (PetscAbsScalar(r[i]) > .20*rmax);
428: }
429: MPIU_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));
430: VecRestoreArray(resid,&r);
431: *per = *per/N;
432: return(0);
433: }
435: /*@C
436: SNESMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value.
438: Collective on SNES
440: Input Parameters:
441: + snes - iterative context
442: . it - iteration number
443: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
444: - dummy - unused monitor context
446: Options Database Key:
447: . -snes_monitor_range - Activates SNESMonitorRange()
449: Level: intermediate
451: .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
452: @*/
453: PetscErrorCode SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,PetscViewerAndFormat *vf)
454: {
456: PetscReal perc,rel;
457: PetscViewer viewer = vf->viewer;
458: /* should be in a MonitorRangeContext */
459: static PetscReal prev;
463: if (!it) prev = rnorm;
464: SNESMonitorRange_Private(snes,it,&perc);
466: rel = (prev - rnorm)/prev;
467: prev = rnorm;
468: PetscViewerPushFormat(viewer,vf->format);
469: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
470: 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));
471: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
472: PetscViewerPopFormat(viewer);
473: return(0);
474: }
476: /*@C
477: SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio
478: of residual norm at each iteration to the previous.
480: Collective on SNES
482: Input Parameters:
483: + snes - the SNES context
484: . its - iteration number
485: . fgnorm - 2-norm of residual (or gradient)
486: - dummy - context of monitor
488: Level: intermediate
490: Notes:
491: Insure that SNESMonitorRatio() is called when you set this monitor
492: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorRatio()
493: @*/
494: PetscErrorCode SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
495: {
496: PetscErrorCode ierr;
497: PetscInt len;
498: PetscReal *history;
499: PetscViewer viewer = vf->viewer;
502: SNESGetConvergenceHistory(snes,&history,NULL,&len);
503: PetscViewerPushFormat(viewer,vf->format);
504: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
505: if (!its || !history || its > len) {
506: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
507: } else {
508: PetscReal ratio = fgnorm/history[its-1];
509: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);
510: }
511: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
512: PetscViewerPopFormat(viewer);
513: return(0);
514: }
516: /*@C
517: SNESMonitorRatioSetUp - Insures the SNES object is saving its history since this monitor needs access to it
519: Collective on SNES
521: Input Parameters:
522: + snes - the SNES context
523: - viewer - the PetscViewer object (ignored)
525: Level: intermediate
527: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorRatio()
528: @*/
529: PetscErrorCode SNESMonitorRatioSetUp(SNES snes,PetscViewerAndFormat *vf)
530: {
531: PetscErrorCode ierr;
532: PetscReal *history;
535: SNESGetConvergenceHistory(snes,&history,NULL,NULL);
536: if (!history) {
537: SNESSetConvergenceHistory(snes,NULL,NULL,100,PETSC_TRUE);
538: }
539: return(0);
540: }
542: /* ---------------------------------------------------------------- */
543: /*
544: Default (short) SNES Monitor, same as SNESMonitorDefault() except
545: it prints fewer digits of the residual as the residual gets smaller.
546: This is because the later digits are meaningless and are often
547: different on different machines; by using this routine different
548: machines will usually generate the same output.
550: Deprecated: Intentionally has no manual page
551: */
552: PetscErrorCode SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
553: {
555: PetscViewer viewer = vf->viewer;
559: PetscViewerPushFormat(viewer,vf->format);
560: PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
561: if (fgnorm > 1.e-9) {
562: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %g \n",its,(double)fgnorm);
563: } else if (fgnorm > 1.e-11) {
564: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,(double)fgnorm);
565: } else {
566: PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);
567: }
568: PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
569: PetscViewerPopFormat(viewer);
570: return(0);
571: }
573: /*@C
574: SNESMonitorDefaultField - Monitors progress of the SNES solvers, separated into fields.
576: Collective on SNES
578: Input Parameters:
579: + snes - the SNES context
580: . its - iteration number
581: . fgnorm - 2-norm of residual
582: - ctx - the PetscViewer
584: Notes:
585: This routine uses the DM attached to the residual vector
587: Level: intermediate
589: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault()
590: @*/
591: PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
592: {
593: PetscViewer viewer = vf->viewer;
594: Vec r;
595: DM dm;
596: PetscReal res[256];
597: PetscInt tablevel;
602: SNESGetFunction(snes, &r, NULL, NULL);
603: VecGetDM(r, &dm);
604: if (!dm) {SNESMonitorDefault(snes, its, fgnorm, vf);}
605: else {
606: PetscSection s, gs;
607: PetscInt Nf, f;
609: DMGetLocalSection(dm, &s);
610: DMGetGlobalSection(dm, &gs);
611: if (!s || !gs) {SNESMonitorDefault(snes, its, fgnorm, vf);}
612: PetscSectionGetNumFields(s, &Nf);
613: if (Nf > 256) SETERRQ1(PetscObjectComm((PetscObject) snes), PETSC_ERR_SUP, "Do not support %d fields > 256", Nf);
614: PetscSectionVecNorm(s, gs, r, NORM_2, res);
615: PetscObjectGetTabLevel((PetscObject) snes, &tablevel);
616: PetscViewerPushFormat(viewer,vf->format);
617: PetscViewerASCIIAddTab(viewer, tablevel);
618: PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);
619: for (f = 0; f < Nf; ++f) {
620: if (f) {PetscViewerASCIIPrintf(viewer, ", ");}
621: PetscViewerASCIIPrintf(viewer, "%14.12e", res[f]);
622: }
623: PetscViewerASCIIPrintf(viewer, "] \n");
624: PetscViewerASCIISubtractTab(viewer, tablevel);
625: PetscViewerPopFormat(viewer);
626: }
627: return(0);
628: }
629: /* ---------------------------------------------------------------- */
630: /*@C
631: SNESConvergedDefault - Default onvergence test of the solvers for
632: systems of nonlinear equations.
634: Collective on SNES
636: Input Parameters:
637: + snes - the SNES context
638: . it - the iteration (0 indicates before any Newton steps)
639: . xnorm - 2-norm of current iterate
640: . snorm - 2-norm of current step
641: . fnorm - 2-norm of function at current iterate
642: - dummy - unused context
644: Output Parameter:
645: . reason - one of
646: $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol),
647: $ SNES_CONVERGED_SNORM_RELATIVE - (snorm < stol*xnorm),
648: $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0),
649: $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf),
650: $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN),
651: $ SNES_CONVERGED_ITERATING - (otherwise),
652: $ SNES_DIVERGED_DTOL - (fnorm > divtol*snes->fnorm0)
654: where
655: + maxf - maximum number of function evaluations, set with SNESSetTolerances()
656: . nfct - number of function evaluations,
657: . abstol - absolute function norm tolerance, set with SNESSetTolerances()
658: . rtol - relative function norm tolerance, set with SNESSetTolerances()
659: . divtol - divergence tolerance, set with SNESSetDivergenceTolerance()
660: - fnorm0 - 2-norm of the function at the initial solution (initial guess; zeroth iteration)
662: Options Database Keys:
663: + -snes_stol - convergence tolerance in terms of the norm of the change in the solution between steps
664: . -snes_atol <abstol> - absolute tolerance of residual norm
665: . -snes_rtol <rtol> - relative decrease in tolerance norm from the initial 2-norm of the solution
666: . -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
667: . -snes_max_funcs <max_funcs> - maximum number of function evaluations
668: . -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
669: - -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
671: Level: intermediate
673: .seealso: SNESSetConvergenceTest(), SNESConvergedSkip(), SNESSetTolerances(), SNESSetDivergenceTolerance()
674: @*/
675: PetscErrorCode SNESConvergedDefault(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
676: {
683: *reason = SNES_CONVERGED_ITERATING;
685: if (!it) {
686: /* set parameter for default relative tolerance convergence test */
687: snes->ttol = fnorm*snes->rtol;
688: snes->rnorm0 = fnorm;
689: }
690: if (PetscIsInfOrNanReal(fnorm)) {
691: PetscInfo(snes,"Failed to converged, function norm is NaN\n");
692: *reason = SNES_DIVERGED_FNORM_NAN;
693: } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
694: PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);
695: *reason = SNES_CONVERGED_FNORM_ABS;
696: } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
697: PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);
698: *reason = SNES_DIVERGED_FUNCTION_COUNT;
699: }
701: if (it && !*reason) {
702: if (fnorm <= snes->ttol) {
703: PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);
704: *reason = SNES_CONVERGED_FNORM_RELATIVE;
705: } else if (snorm < snes->stol*xnorm) {
706: PetscInfo3(snes,"Converged due to small update length: %14.12e < %14.12e * %14.12e\n",(double)snorm,(double)snes->stol,(double)xnorm);
707: *reason = SNES_CONVERGED_SNORM_RELATIVE;
708: } else if (snes->divtol > 0 && (fnorm > snes->divtol*snes->rnorm0)) {
709: PetscInfo3(snes,"Diverged due to increase in function norm: %14.12e > %14.12e * %14.12e\n",(double)fnorm,(double)snes->divtol,(double)snes->rnorm0);
710: *reason = SNES_DIVERGED_DTOL;
711: }
713: }
714: return(0);
715: }
717: /*@C
718: SNESConvergedSkip - Convergence test for SNES that NEVER returns as
719: converged, UNLESS the maximum number of iteration have been reached.
721: Logically Collective on SNES
723: Input Parameters:
724: + snes - the SNES context
725: . it - the iteration (0 indicates before any Newton steps)
726: . xnorm - 2-norm of current iterate
727: . snorm - 2-norm of current step
728: . fnorm - 2-norm of function at current iterate
729: - dummy - unused context
731: Output Parameter:
732: . reason - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
734: Notes:
735: Convergence is then declared after a fixed number of iterations have been used.
737: Level: advanced
739: .seealso: SNESConvergedDefault(), SNESSetConvergenceTest()
740: @*/
741: PetscErrorCode SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
742: {
749: *reason = SNES_CONVERGED_ITERATING;
751: if (fnorm != fnorm) {
752: PetscInfo(snes,"Failed to converged, function norm is NaN\n");
753: *reason = SNES_DIVERGED_FNORM_NAN;
754: } else if (it == snes->max_its) {
755: *reason = SNES_CONVERGED_ITS;
756: }
757: return(0);
758: }
760: /*@C
761: SNESSetWorkVecs - Gets a number of work vectors.
763: Input Parameters:
764: + snes - the SNES context
765: - nw - number of work vectors to allocate
767: Level: developer
769: @*/
770: PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw)
771: {
772: DM dm;
773: Vec v;
777: if (snes->work) {VecDestroyVecs(snes->nwork,&snes->work);}
778: snes->nwork = nw;
780: SNESGetDM(snes, &dm);
781: DMGetGlobalVector(dm, &v);
782: VecDuplicateVecs(v,snes->nwork,&snes->work);
783: DMRestoreGlobalVector(dm, &v);
784: PetscLogObjectParents(snes,nw,snes->work);
785: return(0);
786: }