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: PetscErrorCode SNESMonitorDefaultSetUp(SNES snes, PetscViewerAndFormat *vf)
254: {

258:   if (vf->format == PETSC_VIEWER_DRAW_LG) {
259:     KSPMonitorLGCreate(PetscObjectComm((PetscObject) vf->viewer), NULL, NULL, "Log Residual Norm", 1, NULL, PETSC_DECIDE, PETSC_DECIDE, 400, 300, &vf->lg);
260:   }
261:   return(0);
262: }

264: /*@C
265:    SNESMonitorDefault - Monitors progress of the SNES solvers (default).

267:    Collective on SNES

269:    Input Parameters:
270: +  snes - the SNES context
271: .  its - iteration number
272: .  fgnorm - 2-norm of residual
273: -  vf - viewer and format structure

275:    Notes:
276:    This routine prints the residual norm at each iteration.

278:    Level: intermediate

280: .seealso: SNESMonitorSet(), SNESMonitorSolution()
281: @*/
282: PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
283: {
284:   PetscViewer       viewer = vf->viewer;
285:   PetscViewerFormat format = vf->format;
286:   PetscBool         isascii, isdraw;
287:   PetscErrorCode    ierr;

291:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
292:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW,  &isdraw);
293:   PetscViewerPushFormat(viewer,format);
294:   if (isascii) {
295:     PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
296:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
297:     PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
298:   } else if (isdraw) {
299:     if (format == PETSC_VIEWER_DRAW_LG) {
300:       PetscDrawLG lg = (PetscDrawLG) vf->lg;
301:       PetscReal   x, y;

304:       if (!its) {PetscDrawLGReset(lg);}
305:       x = (PetscReal) its;
306:       if (fgnorm > 0.0) y = PetscLog10Real(fgnorm);
307:       else y = -15.0;
308:       PetscDrawLGAddPoint(lg,&x,&y);
309:       if (its <= 20 || !(its % 5) || snes->reason) {
310:         PetscDrawLGDraw(lg);
311:         PetscDrawLGSave(lg);
312:       }
313:     }
314:   }
315:   PetscViewerPopFormat(viewer);
316:   return(0);
317: }

319: /*@C
320:    SNESMonitorScaling - Monitors the largest value in each row of the Jacobian.

322:    Collective on SNES

324:    Input Parameters:
325: +  snes - the SNES context
326: .  its - iteration number
327: .  fgnorm - 2-norm of residual
328: -  vf - viewer and format structure

330:    Notes:
331:    This routine prints the largest value in each row of the Jacobian

333:    Level: intermediate

335: .seealso: SNESMonitorSet(), SNESMonitorSolution()
336: @*/
337: PetscErrorCode  SNESMonitorScaling(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
338: {
340:   PetscViewer    viewer = vf->viewer;
341:   KSP            ksp;
342:   Mat            J;
343:   Vec            v;

347:   SNESGetKSP(snes,&ksp);
348:   KSPGetOperators(ksp,&J,NULL);
349:   MatCreateVecs(J,&v,NULL);
350:   MatGetRowMaxAbs(J,v,NULL);
351:   PetscViewerPushFormat(viewer,vf->format);
352:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
353:   PetscViewerASCIIPrintf(viewer,"%3D SNES Jacobian maximum row entries \n");
354:   VecView(v,viewer);
355:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
356:   PetscViewerPopFormat(viewer);
357:   VecDestroy(&v);
358:   return(0);
359: }

361: PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,PetscViewerAndFormat *vf)
362: {
363:   Vec            X;
364:   Mat            J,dJ,dJdense;
366:   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
367:   PetscInt       n;
368:   PetscBLASInt   nb = 0,lwork;
369:   PetscReal      *eigr,*eigi;
370:   PetscScalar    *work;
371:   PetscScalar    *a;

374:   if (it == 0) return(0);
375:   /* create the difference between the current update and the current jacobian */
376:   SNESGetSolution(snes,&X);
377:   SNESGetJacobian(snes,NULL,&J,&func,NULL);
378:   MatDuplicate(J,MAT_COPY_VALUES,&dJ);
379:   SNESComputeJacobian(snes,X,dJ,dJ);
380:   MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);

382:   /* compute the spectrum directly */
383:   MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);
384:   MatGetSize(dJ,&n,NULL);
385:   PetscBLASIntCast(n,&nb);
386:   lwork = 3*nb;
387:   PetscMalloc1(n,&eigr);
388:   PetscMalloc1(n,&eigi);
389:   PetscMalloc1(lwork,&work);
390:   MatDenseGetArray(dJdense,&a);
391: #if !defined(PETSC_USE_COMPLEX)
392:   {
393:     PetscBLASInt lierr;
394:     PetscInt     i;
395:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
396:     PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,NULL,&nb,NULL,&nb,work,&lwork,&lierr));
397:     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr);
398:     PetscFPTrapPop();
399:     PetscPrintf(PetscObjectComm((PetscObject)snes),"Eigenvalues of J_%d - J_%d:\n",it,it-1);
400:     for (i=0;i<n;i++) {
401:       PetscPrintf(PetscObjectComm((PetscObject)snes),"%5d: %20.5g + %20.5gi\n",i,(double)eigr[i],(double)eigi[i]);
402:     }
403:   }
404:   MatDenseRestoreArray(dJdense,&a);
405:   MatDestroy(&dJ);
406:   MatDestroy(&dJdense);
407:   PetscFree(eigr);
408:   PetscFree(eigi);
409:   PetscFree(work);
410:   return(0);
411: #else
412:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
413: #endif
414: }

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

418: PetscErrorCode  SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
419: {
421:   Vec            resid;
422:   PetscReal      rmax,pwork;
423:   PetscInt       i,n,N;
424:   PetscScalar    *r;

427:   SNESGetFunction(snes,&resid,NULL,NULL);
428:   VecNorm(resid,NORM_INFINITY,&rmax);
429:   VecGetLocalSize(resid,&n);
430:   VecGetSize(resid,&N);
431:   VecGetArray(resid,&r);
432:   pwork = 0.0;
433:   for (i=0; i<n; i++) {
434:     pwork += (PetscAbsScalar(r[i]) > .20*rmax);
435:   }
436:   MPIU_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));
437:   VecRestoreArray(resid,&r);
438:   *per = *per/N;
439:   return(0);
440: }

442: /*@C
443:    SNESMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value.

445:    Collective on SNES

447:    Input Parameters:
448: +  snes   - iterative context
449: .  it    - iteration number
450: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
451: -  dummy - unused monitor context

453:    Options Database Key:
454: .  -snes_monitor_range - Activates SNESMonitorRange()

456:    Level: intermediate

458: .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
459: @*/
460: PetscErrorCode  SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,PetscViewerAndFormat *vf)
461: {
463:   PetscReal      perc,rel;
464:   PetscViewer    viewer = vf->viewer;
465:   /* should be in a MonitorRangeContext */
466:   static PetscReal prev;

470:   if (!it) prev = rnorm;
471:   SNESMonitorRange_Private(snes,it,&perc);

473:   rel  = (prev - rnorm)/prev;
474:   prev = rnorm;
475:   PetscViewerPushFormat(viewer,vf->format);
476:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
477:   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));
478:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
479:   PetscViewerPopFormat(viewer);
480:   return(0);
481: }

483: /*@C
484:    SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio
485:    of residual norm at each iteration to the previous.

487:    Collective on SNES

489:    Input Parameters:
490: +  snes - the SNES context
491: .  its - iteration number
492: .  fgnorm - 2-norm of residual (or gradient)
493: -  dummy -  context of monitor

495:    Level: intermediate

497:    Notes:
498:     Insure that SNESMonitorRatio() is called when you set this monitor
499: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorRatio()
500: @*/
501: PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
502: {
503:   PetscErrorCode          ierr;
504:   PetscInt                len;
505:   PetscReal               *history;
506:   PetscViewer             viewer = vf->viewer;

509:   SNESGetConvergenceHistory(snes,&history,NULL,&len);
510:   PetscViewerPushFormat(viewer,vf->format);
511:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
512:   if (!its || !history || its > len) {
513:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
514:   } else {
515:     PetscReal ratio = fgnorm/history[its-1];
516:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);
517:   }
518:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
519:   PetscViewerPopFormat(viewer);
520:   return(0);
521: }

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

526:    Collective on SNES

528:    Input Parameters:
529: +   snes - the SNES context
530: -   viewer - the PetscViewer object (ignored)

532:    Level: intermediate

534: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorRatio()
535: @*/
536: PetscErrorCode  SNESMonitorRatioSetUp(SNES snes,PetscViewerAndFormat *vf)
537: {
538:   PetscErrorCode          ierr;
539:   PetscReal               *history;

542:   SNESGetConvergenceHistory(snes,&history,NULL,NULL);
543:   if (!history) {
544:     SNESSetConvergenceHistory(snes,NULL,NULL,100,PETSC_TRUE);
545:   }
546:   return(0);
547: }

549: /* ---------------------------------------------------------------- */
550: /*
551:      Default (short) SNES Monitor, same as SNESMonitorDefault() except
552:   it prints fewer digits of the residual as the residual gets smaller.
553:   This is because the later digits are meaningless and are often
554:   different on different machines; by using this routine different
555:   machines will usually generate the same output.

557:   Deprecated: Intentionally has no manual page
558: */
559: PetscErrorCode  SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
560: {
562:   PetscViewer    viewer = vf->viewer;

566:   PetscViewerPushFormat(viewer,vf->format);
567:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
568:   if (fgnorm > 1.e-9) {
569:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %g \n",its,(double)fgnorm);
570:   } else if (fgnorm > 1.e-11) {
571:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,(double)fgnorm);
572:   } else {
573:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);
574:   }
575:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
576:   PetscViewerPopFormat(viewer);
577:   return(0);
578: }

580: /*@C
581:   SNESMonitorDefaultField - Monitors progress of the SNES solvers, separated into fields.

583:   Collective on SNES

585:   Input Parameters:
586: + snes   - the SNES context
587: . its    - iteration number
588: . fgnorm - 2-norm of residual
589: - ctx    - the PetscViewer

591:   Notes:
592:   This routine uses the DM attached to the residual vector

594:   Level: intermediate

596: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault()
597: @*/
598: PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
599: {
600:   PetscViewer    viewer = vf->viewer;
601:   Vec            r;
602:   DM             dm;
603:   PetscReal      res[256];
604:   PetscInt       tablevel;

609:   SNESGetFunction(snes, &r, NULL, NULL);
610:   VecGetDM(r, &dm);
611:   if (!dm) {SNESMonitorDefault(snes, its, fgnorm, vf);}
612:   else {
613:     PetscSection s, gs;
614:     PetscInt     Nf, f;

616:     DMGetLocalSection(dm, &s);
617:     DMGetGlobalSection(dm, &gs);
618:     if (!s || !gs) {SNESMonitorDefault(snes, its, fgnorm, vf);}
619:     PetscSectionGetNumFields(s, &Nf);
620:     if (Nf > 256) SETERRQ1(PetscObjectComm((PetscObject) snes), PETSC_ERR_SUP, "Do not support %d fields > 256", Nf);
621:     PetscSectionVecNorm(s, gs, r, NORM_2, res);
622:     PetscObjectGetTabLevel((PetscObject) snes, &tablevel);
623:     PetscViewerPushFormat(viewer,vf->format);
624:     PetscViewerASCIIAddTab(viewer, tablevel);
625:     PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);
626:     for (f = 0; f < Nf; ++f) {
627:       if (f) {PetscViewerASCIIPrintf(viewer, ", ");}
628:       PetscViewerASCIIPrintf(viewer, "%14.12e", res[f]);
629:     }
630:     PetscViewerASCIIPrintf(viewer, "] \n");
631:     PetscViewerASCIISubtractTab(viewer, tablevel);
632:     PetscViewerPopFormat(viewer);
633:   }
634:   return(0);
635: }
636: /* ---------------------------------------------------------------- */
637: /*@C
638:    SNESConvergedDefault - Default onvergence test of the solvers for
639:    systems of nonlinear equations.

641:    Collective on SNES

643:    Input Parameters:
644: +  snes - the SNES context
645: .  it - the iteration (0 indicates before any Newton steps)
646: .  xnorm - 2-norm of current iterate
647: .  snorm - 2-norm of current step
648: .  fnorm - 2-norm of function at current iterate
649: -  dummy - unused context

651:    Output Parameter:
652: .   reason  - one of
653: $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
654: $  SNES_CONVERGED_SNORM_RELATIVE  - (snorm < stol*xnorm),
655: $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
656: $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
657: $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
658: $  SNES_CONVERGED_ITERATING       - (otherwise),
659: $  SNES_DIVERGED_DTOL             - (fnorm > divtol*snes->fnorm0)

661:    where
662: +    maxf - maximum number of function evaluations,  set with SNESSetTolerances()
663: .    nfct - number of function evaluations,
664: .    abstol - absolute function norm tolerance, set with SNESSetTolerances()
665: .    rtol - relative function norm tolerance, set with SNESSetTolerances()
666: .    divtol - divergence tolerance, set with SNESSetDivergenceTolerance()
667: -    fnorm0 - 2-norm of the function at the initial solution (initial guess; zeroth iteration)

669:   Options Database Keys:
670: +  -snes_stol - convergence tolerance in terms of the norm  of the change in the solution between steps
671: .  -snes_atol <abstol> - absolute tolerance of residual norm
672: .  -snes_rtol <rtol> - relative decrease in tolerance norm from the initial 2-norm of the solution
673: .  -snes_divergence_tolerance <divtol> - if the residual goes above divtol*rnorm0, exit with divergence
674: .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
675: .  -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
676: -  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops

678:    Level: intermediate

680: .seealso: SNESSetConvergenceTest(), SNESConvergedSkip(), SNESSetTolerances(), SNESSetDivergenceTolerance()
681: @*/
682: PetscErrorCode  SNESConvergedDefault(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
683: {


690:   *reason = SNES_CONVERGED_ITERATING;

692:   if (!it) {
693:     /* set parameter for default relative tolerance convergence test */
694:     snes->ttol   = fnorm*snes->rtol;
695:     snes->rnorm0 = fnorm;
696:   }
697:   if (PetscIsInfOrNanReal(fnorm)) {
698:     PetscInfo(snes,"Failed to converged, function norm is NaN\n");
699:     *reason = SNES_DIVERGED_FNORM_NAN;
700:   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
701:     PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);
702:     *reason = SNES_CONVERGED_FNORM_ABS;
703:   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
704:     PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);
705:     *reason = SNES_DIVERGED_FUNCTION_COUNT;
706:   }

708:   if (it && !*reason) {
709:     if (fnorm <= snes->ttol) {
710:       PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);
711:       *reason = SNES_CONVERGED_FNORM_RELATIVE;
712:     } else if (snorm < snes->stol*xnorm) {
713:       PetscInfo3(snes,"Converged due to small update length: %14.12e < %14.12e * %14.12e\n",(double)snorm,(double)snes->stol,(double)xnorm);
714:       *reason = SNES_CONVERGED_SNORM_RELATIVE;
715:     } else if (snes->divtol > 0 && (fnorm > snes->divtol*snes->rnorm0)) {
716:       PetscInfo3(snes,"Diverged due to increase in function norm: %14.12e > %14.12e * %14.12e\n",(double)fnorm,(double)snes->divtol,(double)snes->rnorm0);
717:       *reason = SNES_DIVERGED_DTOL;
718:     }

720:   }
721:   return(0);
722: }

724: /*@C
725:    SNESConvergedSkip - Convergence test for SNES that NEVER returns as
726:    converged, UNLESS the maximum number of iteration have been reached.

728:    Logically Collective on SNES

730:    Input Parameters:
731: +  snes - the SNES context
732: .  it - the iteration (0 indicates before any Newton steps)
733: .  xnorm - 2-norm of current iterate
734: .  snorm - 2-norm of current step
735: .  fnorm - 2-norm of function at current iterate
736: -  dummy - unused context

738:    Output Parameter:
739: .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN

741:    Notes:
742:    Convergence is then declared after a fixed number of iterations have been used.

744:    Level: advanced

746: .seealso: SNESConvergedDefault(), SNESSetConvergenceTest()
747: @*/
748: PetscErrorCode  SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
749: {


756:   *reason = SNES_CONVERGED_ITERATING;

758:   if (fnorm != fnorm) {
759:     PetscInfo(snes,"Failed to converged, function norm is NaN\n");
760:     *reason = SNES_DIVERGED_FNORM_NAN;
761:   } else if (it == snes->max_its) {
762:     *reason = SNES_CONVERGED_ITS;
763:   }
764:   return(0);
765: }

767: /*@C
768:   SNESSetWorkVecs - Gets a number of work vectors.

770:   Input Parameters:
771: + snes  - the SNES context
772: - nw - number of work vectors to allocate

774:   Level: developer

776: @*/
777: PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw)
778: {
779:   DM             dm;
780:   Vec            v;

784:   if (snes->work) {VecDestroyVecs(snes->nwork,&snes->work);}
785:   snes->nwork = nw;

787:   SNESGetDM(snes, &dm);
788:   DMGetGlobalVector(dm, &v);
789:   VecDuplicateVecs(v,snes->nwork,&snes->work);
790:   DMRestoreGlobalVector(dm, &v);
791:   PetscLogObjectParents(snes,nw,snes->work);
792:   return(0);
793: }