Actual source code: snesut.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

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

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

 10:    Collective on SNES

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

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

 21:    Level: intermediate

 23: .keywords: SNES, nonlinear, vector, monitor, view

 25: .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
 26: @*/
 27: PetscErrorCode  SNESMonitorSolution(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
 28: {
 30:   Vec            x;
 31:   PetscViewer    viewer = vf->viewer;

 35:   SNESGetSolution(snes,&x);
 36:   PetscViewerPushFormat(viewer,vf->format);
 37:   VecView(x,viewer);
 38:   PetscViewerPopFormat(viewer);
 39:   return(0);
 40: }

 42: /*@C
 43:    SNESMonitorResidual - Monitors progress of the SNES solvers by calling
 44:    VecView() for the residual at each iteration.

 46:    Collective on SNES

 48:    Input Parameters:
 49: +  snes - the SNES context
 50: .  its - iteration number
 51: .  fgnorm - 2-norm of residual
 52: -  dummy -  a viewer

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

 57:    Level: intermediate

 59: .keywords: SNES, nonlinear, vector, monitor, view

 61: .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
 62: @*/
 63: PetscErrorCode  SNESMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
 64: {
 66:   Vec            x;
 67:   PetscViewer    viewer = vf->viewer;

 71:   SNESGetFunction(snes,&x,0,0);
 72:   PetscViewerPushFormat(viewer,vf->format);
 73:   VecView(x,viewer);
 74:   PetscViewerPopFormat(viewer);
 75:   return(0);
 76: }

 78: /*@C
 79:    SNESMonitorSolutionUpdate - Monitors progress of the SNES solvers by calling
 80:    VecView() for the UPDATE to the solution at each iteration.

 82:    Collective on SNES

 84:    Input Parameters:
 85: +  snes - the SNES context
 86: .  its - iteration number
 87: .  fgnorm - 2-norm of residual
 88: -  dummy - a viewer

 90:    Options Database Keys:
 91: .  -snes_monitor_solution_update [ascii binary draw][:filename][:viewer format] - plots update to solution at each iteration

 93:    Level: intermediate

 95: .keywords: SNES, nonlinear, vector, monitor, view

 97: .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
 98: @*/
 99: PetscErrorCode  SNESMonitorSolutionUpdate(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
100: {
102:   Vec            x;
103:   PetscViewer    viewer = vf->viewer;

107:   SNESGetSolutionUpdate(snes,&x);
108:   PetscViewerPushFormat(viewer,vf->format);
109:   VecView(x,viewer);
110:   PetscViewerPopFormat(viewer);
111:   return(0);
112: }

114: /*@C
115:    KSPMonitorSNES - Print the residual norm of the nonlinear function at each iteration of the linear iterative solver.

117:    Collective on KSP

119:    Input Parameters:
120: +  ksp   - iterative context
121: .  n     - iteration number
122: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
123: -  dummy - unused monitor context

125:    Level: intermediate

127: .keywords: KSP, default, monitor, residual

129: .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGResidualNormCreate()
130: @*/
131: PetscErrorCode  KSPMonitorSNES(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
132: {
134:   PetscViewer    viewer;
135:   SNES           snes = (SNES) dummy;
136:   Vec            snes_solution,work1,work2;
137:   PetscReal      snorm;

140:   SNESGetSolution(snes,&snes_solution);
141:   VecDuplicate(snes_solution,&work1);
142:   VecDuplicate(snes_solution,&work2);
143:   KSPBuildSolution(ksp,work1,NULL);
144:   VecAYPX(work1,-1.0,snes_solution);
145:   SNESComputeFunction(snes,work1,work2);
146:   VecNorm(work2,NORM_2,&snorm);
147:   VecDestroy(&work1);
148:   VecDestroy(&work2);

150:   PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
151:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
152:   if (n == 0 && ((PetscObject)ksp)->prefix) {
153:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
154:   }
155:   PetscViewerASCIIPrintf(viewer,"%3D SNES Residual norm %5.3e KSP Residual norm %5.3e \n",n,(double)snorm,(double)rnorm);
156:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
157:   return(0);
158: }

160:  #include <petscdraw.h>

162: /*@C
163:    KSPMonitorSNESLGResidualNormCreate - Creates a line graph context for use with
164:    KSP to monitor convergence of preconditioned residual norms.

166:    Collective on KSP

168:    Input Parameters:
169: +  comm - communicator context
170: .  host - the X display to open, or null for the local machine
171: .  label - the title to put in the title bar
172: .  x, y - the screen coordinates of the upper left coordinate of
173:           the window
174: -  m, n - the screen width and height in pixels

176:    Output Parameter:
177: .  draw - the drawing context

179:    Options Database Key:
180: .  -ksp_monitor_lg_residualnorm - Sets line graph monitor

182:    Notes:
183:    Use KSPMonitorSNESLGResidualNormDestroy() to destroy this line graph; do not use PetscDrawLGDestroy().

185:    Level: intermediate

187: .keywords: KSP, monitor, line graph, residual, create

189: .seealso: KSPMonitorSNESLGResidualNormDestroy(), KSPMonitorSet(), KSPMonitorSNESLGTrueResidualCreate()
190: @*/
191: PetscErrorCode  KSPMonitorSNESLGResidualNormCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscObject **objs)
192: {
193:   PetscDraw      draw;
195:   PetscDrawAxis  axis;
196:   PetscDrawLG    lg;
197:   const char     *names[] = {"Linear residual","Nonlinear residual"};

200:   PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
201:   PetscDrawSetFromOptions(draw);
202:   PetscDrawLGCreate(draw,2,&lg);
203:   PetscDrawLGSetLegend(lg,names);
204:   PetscDrawLGSetFromOptions(lg);
205:   PetscDrawLGGetAxis(lg,&axis);
206:   PetscDrawAxisSetLabels(axis,"Convergence of Residual Norm","Iteration","Residual Norm");
207:   PetscDrawDestroy(&draw);

209:   PetscMalloc1(2,objs);
210:   (*objs)[1] = (PetscObject)lg;
211:   return(0);
212: }

214: PetscErrorCode  KSPMonitorSNESLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscObject *objs)
215: {
216:   SNES           snes = (SNES) objs[0];
217:   PetscDrawLG    lg   = (PetscDrawLG) objs[1];
219:   PetscReal      y[2];
220:   Vec            snes_solution,work1,work2;

223:   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
224:   else y[0] = -15.0;

226:   SNESGetSolution(snes,&snes_solution);
227:   VecDuplicate(snes_solution,&work1);
228:   VecDuplicate(snes_solution,&work2);
229:   KSPBuildSolution(ksp,work1,NULL);
230:   VecAYPX(work1,-1.0,snes_solution);
231:   SNESComputeFunction(snes,work1,work2);
232:   VecNorm(work2,NORM_2,y+1);
233:   if (y[1] > 0.0) y[1] = PetscLog10Real(y[1]);
234:   else y[1] = -15.0;
235:   VecDestroy(&work1);
236:   VecDestroy(&work2);

238:   PetscDrawLGAddPoint(lg,NULL,y);
239:   if (n < 20 || !(n % 5) || snes->reason) {
240:     PetscDrawLGDraw(lg);
241:     PetscDrawLGSave(lg);
242:   }
243:   return(0);
244: }

246: /*@
247:    KSPMonitorSNESLGResidualNormDestroy - Destroys a line graph context that was created
248:    with KSPMonitorSNESLGResidualNormCreate().

250:    Collective on KSP

252:    Input Parameter:
253: .  draw - the drawing context

255:    Level: intermediate

257: .keywords: KSP, monitor, line graph, destroy

259: .seealso: KSPMonitorSNESLGResidualNormCreate(), KSPMonitorSNESLGTrueResidualDestroy(), KSPMonitorSet()
260: @*/
261: PetscErrorCode  KSPMonitorSNESLGResidualNormDestroy(PetscObject **objs)
262: {
264:   PetscDrawLG    lg = (PetscDrawLG) (*objs)[1];

267:   PetscDrawLGDestroy(&lg);
268:   PetscFree(*objs);
269:   return(0);
270: }

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

275:    Collective on SNES

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:    Notes:
284:    This routine prints the residual norm at each iteration.

286:    Level: intermediate

288: .keywords: SNES, nonlinear, default, monitor, norm

290: .seealso: SNESMonitorSet(), SNESMonitorSolution()
291: @*/
292: PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
293: {
295:   PetscViewer    viewer = vf->viewer;

299:   PetscViewerPushFormat(viewer,vf->format);
300:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
301:   PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
302:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
303:   PetscViewerPopFormat(viewer);CHKERRQ(ierr) ;
304:   return(0);
305: }

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

310:    Collective on SNES

312:    Input Parameters:
313: +  snes - the SNES context
314: .  its - iteration number
315: .  fgnorm - 2-norm of residual
316: -  vf - viewer and format structure

318:    Notes:
319:    This routine prints the largest value in each row of the Jacobian 

321:    Level: intermediate

323: .keywords: SNES, nonlinear, default, monitor, norm

325: .seealso: SNESMonitorSet(), SNESMonitorSolution()
326: @*/
327: PetscErrorCode  SNESMonitorScaling(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
328: {
330:   PetscViewer    viewer = vf->viewer;
331:   KSP            ksp;
332:   Mat            J;
333:   Vec            v;

337:   SNESGetKSP(snes,&ksp);
338:   KSPGetOperators(ksp,&J,NULL);
339:   MatCreateVecs(J,&v,NULL);
340:   MatGetRowMaxAbs(J,v,NULL);
341:   PetscViewerPushFormat(viewer,vf->format);
342:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
343:   PetscViewerASCIIPrintf(viewer,"%3D SNES Jacobian maximum row entries \n");
344:   VecView(v,viewer);
345:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
346:   PetscViewerPopFormat(viewer);CHKERRQ(ierr) ;
347:   VecDestroy(&v);
348:   return(0);
349: }

351: PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,PetscViewerAndFormat *vf)
352: {
353: #if defined(PETSC_MISSING_LAPACK_GEEV)
354:   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
355: #elif defined(PETSC_HAVE_ESSL)
356:   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - No support for ESSL Lapack Routines");
357: #else
358:   Vec            X;
359:   Mat            J,dJ,dJdense;
361:   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
362:   PetscInt       n,i;
363:   PetscBLASInt   nb,lwork;
364:   PetscReal      *eigr,*eigi;
365:   PetscScalar    *work;
366:   PetscScalar    *a;

369:   if (it == 0) return(0);
370:   /* create the difference between the current update and the current jacobian */
371:   SNESGetSolution(snes,&X);
372:   SNESGetJacobian(snes,NULL,&J,&func,NULL);
373:   MatDuplicate(J,MAT_COPY_VALUES,&dJ);
374:   SNESComputeJacobian(snes,X,dJ,dJ);
375:   MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);

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

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

413: PetscErrorCode  SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
414: {
416:   Vec            resid;
417:   PetscReal      rmax,pwork;
418:   PetscInt       i,n,N;
419:   PetscScalar    *r;

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

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

440:    Collective on SNES

442:    Input Parameters:
443: +  snes   - iterative context
444: .  it    - iteration number
445: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
446: -  dummy - unused monitor context

448:    Options Database Key:
449: .  -snes_monitor_range - Activates SNESMonitorRange()

451:    Level: intermediate

453: .keywords: SNES, default, monitor, residual

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

467:   if (!it) prev = rnorm;
468:   SNESMonitorRange_Private(snes,it,&perc);

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

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

484:    Collective on SNES

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

492:    Level: intermediate

494:    Notes:
495:     Insure that SNESMonitorRatio() is called when you set this monitor
496: .keywords: SNES, nonlinear, monitor, norm

498: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorRatio()
499: @*/
500: PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
501: {
502:   PetscErrorCode          ierr;
503:   PetscInt                len;
504:   PetscReal               *history;
505:   PetscViewer             viewer = vf->viewer;
506: 
508:   SNESGetConvergenceHistory(snes,&history,NULL,&len);
509:   PetscViewerPushFormat(viewer,vf->format);
510:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
511:   if (!its || !history || its > len) {
512:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
513:   } else {
514:     PetscReal ratio = fgnorm/history[its-1];
515:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);
516:   }
517:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
518:   PetscViewerPopFormat(viewer);
519:   return(0);
520: }

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

525:    Collective on SNES

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

531:    Level: intermediate

533: .keywords: SNES, nonlinear, monitor, norm

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

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

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

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

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

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

584:   Collective on SNES

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

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

595:   Level: intermediate

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

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

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

643:    Collective on SNES

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

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

662:    where
663: +    maxf - maximum number of function evaluations,
664:             set with SNESSetTolerances()
665: .    nfct - number of function evaluations,
666: .    abstol - absolute function norm tolerance,
667:             set with SNESSetTolerances()
668: -    rtol - relative function norm tolerance, set with SNESSetTolerances()

670:    Level: intermediate

672: .keywords: SNES, nonlinear, default, converged, convergence

674: .seealso: SNESSetConvergenceTest()
675: @*/
676: PetscErrorCode  SNESConvergedDefault(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
677: {


684:   *reason = SNES_CONVERGED_ITERATING;

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

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

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

722:    Logically Collective on SNES

724:    Input Parameters:
725: +  snes - the SNES context
726: .  it - the iteration (0 indicates before any Newton steps)
727: .  xnorm - 2-norm of current iterate
728: .  snorm - 2-norm of current step
729: .  fnorm - 2-norm of function at current iterate
730: -  dummy - unused context

732:    Output Parameter:
733: .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN

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

738:    Level: advanced

740: .keywords: SNES, nonlinear, skip, converged, convergence

742: .seealso: SNESSetConvergenceTest()
743: @*/
744: PetscErrorCode  SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
745: {


752:   *reason = SNES_CONVERGED_ITERATING;

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

763: /*@C
764:   SNESSetWorkVecs - Gets a number of work vectors.

766:   Input Parameters:
767: . snes  - the SNES context
768: . nw - number of work vectors to allocate

770:   Level: developer

772: @*/
773: PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw)
774: {
775:   DM             dm;
776:   Vec            v;

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

783:   SNESGetDM(snes, &dm);
784:   DMGetGlobalVector(dm, &v);
785:   VecDuplicateVecs(v,snes->nwork,&snes->work);
786:   DMRestoreGlobalVector(dm, &v);
787:   PetscLogObjectParents(snes,nw,snes->work);
788:   return(0);
789: }