Actual source code: snesut.c

petsc-3.12.5 2020-03-29
Report Typos and Errors

  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,0,0);
 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: /*@C
110:    KSPMonitorSNES - Print the residual norm of the nonlinear function at each iteration of the linear iterative solver.

112:    Collective on KSP

114:    Input Parameters:
115: +  ksp   - iterative context
116: .  n     - iteration number
117: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
118: -  dummy - unused monitor context

120:    Level: intermediate

122: .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGResidualNormCreate()
123: @*/
124: PetscErrorCode  KSPMonitorSNES(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
125: {
127:   PetscViewer    viewer;
128:   SNES           snes = (SNES) dummy;
129:   Vec            snes_solution,work1,work2;
130:   PetscReal      snorm;

133:   SNESGetSolution(snes,&snes_solution);
134:   VecDuplicate(snes_solution,&work1);
135:   VecDuplicate(snes_solution,&work2);
136:   KSPBuildSolution(ksp,work1,NULL);
137:   VecAYPX(work1,-1.0,snes_solution);
138:   SNESComputeFunction(snes,work1,work2);
139:   VecNorm(work2,NORM_2,&snorm);
140:   VecDestroy(&work1);
141:   VecDestroy(&work2);

143:   PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);
144:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
145:   if (n == 0 && ((PetscObject)ksp)->prefix) {
146:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
147:   }
148:   PetscViewerASCIIPrintf(viewer,"%3D SNES Residual norm %5.3e KSP Residual norm %5.3e \n",n,(double)snorm,(double)rnorm);
149:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
150:   return(0);
151: }

153:  #include <petscdraw.h>

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

159:    Collective on KSP

161:    Input Parameters:
162: +  comm - communicator context
163: .  host - the X display to open, or null for the local machine
164: .  label - the title to put in the title bar
165: .  x, y - the screen coordinates of the upper left coordinate of
166:           the window
167: -  m, n - the screen width and height in pixels

169:    Output Parameter:
170: .  draw - the drawing context

172:    Options Database Key:
173: .  -ksp_monitor_lg_residualnorm - Sets line graph monitor

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

178:    Level: intermediate

180: .seealso: KSPMonitorSNESLGResidualNormDestroy(), KSPMonitorSet(), KSPMonitorSNESLGTrueResidualCreate()
181: @*/
182: PetscErrorCode  KSPMonitorSNESLGResidualNormCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscObject **objs)
183: {
184:   PetscDraw      draw;
186:   PetscDrawAxis  axis;
187:   PetscDrawLG    lg;
188:   const char     *names[] = {"Linear residual","Nonlinear residual"};

191:   PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
192:   PetscDrawSetFromOptions(draw);
193:   PetscDrawLGCreate(draw,2,&lg);
194:   PetscDrawLGSetLegend(lg,names);
195:   PetscDrawLGSetFromOptions(lg);
196:   PetscDrawLGGetAxis(lg,&axis);
197:   PetscDrawAxisSetLabels(axis,"Convergence of Residual Norm","Iteration","Residual Norm");
198:   PetscDrawDestroy(&draw);

200:   PetscMalloc1(2,objs);
201:   (*objs)[1] = (PetscObject)lg;
202:   return(0);
203: }

205: PetscErrorCode  KSPMonitorSNESLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscObject *objs)
206: {
207:   SNES           snes = (SNES) objs[0];
208:   PetscDrawLG    lg   = (PetscDrawLG) objs[1];
210:   PetscReal      y[2];
211:   Vec            snes_solution,work1,work2;

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

217:   SNESGetSolution(snes,&snes_solution);
218:   VecDuplicate(snes_solution,&work1);
219:   VecDuplicate(snes_solution,&work2);
220:   KSPBuildSolution(ksp,work1,NULL);
221:   VecAYPX(work1,-1.0,snes_solution);
222:   SNESComputeFunction(snes,work1,work2);
223:   VecNorm(work2,NORM_2,y+1);
224:   if (y[1] > 0.0) y[1] = PetscLog10Real(y[1]);
225:   else y[1] = -15.0;
226:   VecDestroy(&work1);
227:   VecDestroy(&work2);

229:   PetscDrawLGAddPoint(lg,NULL,y);
230:   if (n < 20 || !(n % 5) || snes->reason) {
231:     PetscDrawLGDraw(lg);
232:     PetscDrawLGSave(lg);
233:   }
234:   return(0);
235: }

237: /*@
238:    KSPMonitorSNESLGResidualNormDestroy - Destroys a line graph context that was created
239:    with KSPMonitorSNESLGResidualNormCreate().

241:    Collective on KSP

243:    Input Parameter:
244: .  draw - the drawing context

246:    Level: intermediate

248: .seealso: KSPMonitorSNESLGResidualNormCreate(), KSPMonitorSNESLGTrueResidualDestroy(), KSPMonitorSet()
249: @*/
250: PetscErrorCode  KSPMonitorSNESLGResidualNormDestroy(PetscObject **objs)
251: {
253:   PetscDrawLG    lg = (PetscDrawLG) (*objs)[1];

256:   PetscDrawLGDestroy(&lg);
257:   PetscFree(*objs);
258:   return(0);
259: }

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

264:    Collective on SNES

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

272:    Notes:
273:    This routine prints the residual norm at each iteration.

275:    Level: intermediate

277: .seealso: SNESMonitorSet(), SNESMonitorSolution()
278: @*/
279: PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
280: {
282:   PetscViewer    viewer = vf->viewer;

286:   PetscViewerPushFormat(viewer,vf->format);
287:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
288:   PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
289:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
290:   PetscViewerPopFormat(viewer);
291:   return(0);
292: }

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

297:    Collective on SNES

299:    Input Parameters:
300: +  snes - the SNES context
301: .  its - iteration number
302: .  fgnorm - 2-norm of residual
303: -  vf - viewer and format structure

305:    Notes:
306:    This routine prints the largest value in each row of the Jacobian

308:    Level: intermediate

310: .seealso: SNESMonitorSet(), SNESMonitorSolution()
311: @*/
312: PetscErrorCode  SNESMonitorScaling(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
313: {
315:   PetscViewer    viewer = vf->viewer;
316:   KSP            ksp;
317:   Mat            J;
318:   Vec            v;

322:   SNESGetKSP(snes,&ksp);
323:   KSPGetOperators(ksp,&J,NULL);
324:   MatCreateVecs(J,&v,NULL);
325:   MatGetRowMaxAbs(J,v,NULL);
326:   PetscViewerPushFormat(viewer,vf->format);
327:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
328:   PetscViewerASCIIPrintf(viewer,"%3D SNES Jacobian maximum row entries \n");
329:   VecView(v,viewer);
330:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
331:   PetscViewerPopFormat(viewer);
332:   VecDestroy(&v);
333:   return(0);
334: }

336: PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,PetscViewerAndFormat *vf)
337: {
338: #if defined(PETSC_MISSING_LAPACK_GEEV)
339:   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
340: #elif defined(PETSC_HAVE_ESSL)
341:   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - No support for ESSL Lapack Routines");
342: #else
343:   Vec            X;
344:   Mat            J,dJ,dJdense;
346:   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
347:   PetscInt       n,i;
348:   PetscBLASInt   nb,lwork;
349:   PetscReal      *eigr,*eigi;
350:   PetscScalar    *work;
351:   PetscScalar    *a;

354:   if (it == 0) return(0);
355:   /* create the difference between the current update and the current jacobian */
356:   SNESGetSolution(snes,&X);
357:   SNESGetJacobian(snes,NULL,&J,&func,NULL);
358:   MatDuplicate(J,MAT_COPY_VALUES,&dJ);
359:   SNESComputeJacobian(snes,X,dJ,dJ);
360:   MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);

362:   /* compute the spectrum directly */
363:   MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);
364:   MatGetSize(dJ,&n,NULL);
365:   PetscBLASIntCast(n,&nb);
366:   lwork = 3*nb;
367:   PetscMalloc1(n,&eigr);
368:   PetscMalloc1(n,&eigi);
369:   PetscMalloc1(lwork,&work);
370:   MatDenseGetArray(dJdense,&a);
371: #if !defined(PETSC_USE_COMPLEX)
372:   {
373:     PetscBLASInt lierr;
374:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
375:     PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,NULL,&nb,NULL,&nb,work,&lwork,&lierr));
376:     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr);
377:     PetscFPTrapPop();
378:   }
379: #else
380:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
381: #endif
382:   PetscPrintf(PetscObjectComm((PetscObject)snes),"Eigenvalues of J_%d - J_%d:\n",it,it-1);
383:   for (i=0;i<n;i++) {
384:     PetscPrintf(PetscObjectComm((PetscObject)snes),"%5d: %20.5g + %20.5gi\n",i,(double)eigr[i],(double)eigi[i]);
385:   }
386:   MatDenseRestoreArray(dJdense,&a);
387:   MatDestroy(&dJ);
388:   MatDestroy(&dJdense);
389:   PetscFree(eigr);
390:   PetscFree(eigi);
391:   PetscFree(work);
392:   return(0);
393: #endif
394: }

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

398: PetscErrorCode  SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
399: {
401:   Vec            resid;
402:   PetscReal      rmax,pwork;
403:   PetscInt       i,n,N;
404:   PetscScalar    *r;

407:   SNESGetFunction(snes,&resid,0,0);
408:   VecNorm(resid,NORM_INFINITY,&rmax);
409:   VecGetLocalSize(resid,&n);
410:   VecGetSize(resid,&N);
411:   VecGetArray(resid,&r);
412:   pwork = 0.0;
413:   for (i=0; i<n; i++) {
414:     pwork += (PetscAbsScalar(r[i]) > .20*rmax);
415:   }
416:   MPIU_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));
417:   VecRestoreArray(resid,&r);
418:   *per = *per/N;
419:   return(0);
420: }

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

425:    Collective on SNES

427:    Input Parameters:
428: +  snes   - iterative context
429: .  it    - iteration number
430: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
431: -  dummy - unused monitor context

433:    Options Database Key:
434: .  -snes_monitor_range - Activates SNESMonitorRange()

436:    Level: intermediate

438: .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
439: @*/
440: PetscErrorCode  SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,PetscViewerAndFormat *vf)
441: {
443:   PetscReal      perc,rel;
444:   PetscViewer    viewer = vf->viewer;
445:   /* should be in a MonitorRangeContext */
446:   static PetscReal prev;

450:   if (!it) prev = rnorm;
451:   SNESMonitorRange_Private(snes,it,&perc);

453:   rel  = (prev - rnorm)/prev;
454:   prev = rnorm;
455:   PetscViewerPushFormat(viewer,vf->format);
456:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
457:   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));
458:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
459:   PetscViewerPopFormat(viewer);
460:   return(0);
461: }

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

467:    Collective on SNES

469:    Input Parameters:
470: +  snes - the SNES context
471: .  its - iteration number
472: .  fgnorm - 2-norm of residual (or gradient)
473: -  dummy -  context of monitor

475:    Level: intermediate

477:    Notes:
478:     Insure that SNESMonitorRatio() is called when you set this monitor
479: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorRatio()
480: @*/
481: PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
482: {
483:   PetscErrorCode          ierr;
484:   PetscInt                len;
485:   PetscReal               *history;
486:   PetscViewer             viewer = vf->viewer;

489:   SNESGetConvergenceHistory(snes,&history,NULL,&len);
490:   PetscViewerPushFormat(viewer,vf->format);
491:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
492:   if (!its || !history || its > len) {
493:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
494:   } else {
495:     PetscReal ratio = fgnorm/history[its-1];
496:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);
497:   }
498:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
499:   PetscViewerPopFormat(viewer);
500:   return(0);
501: }

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

506:    Collective on SNES

508:    Input Parameters:
509: +   snes - the SNES context
510: -   viewer - the PetscViewer object (ignored)

512:    Level: intermediate

514: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorRatio()
515: @*/
516: PetscErrorCode  SNESMonitorRatioSetUp(SNES snes,PetscViewerAndFormat *vf)
517: {
518:   PetscErrorCode          ierr;
519:   PetscReal               *history;

522:   SNESGetConvergenceHistory(snes,&history,NULL,NULL);
523:   if (!history) {
524:     SNESSetConvergenceHistory(snes,NULL,NULL,100,PETSC_TRUE);
525:   }
526:   return(0);
527: }

529: /* ---------------------------------------------------------------- */
530: /*
531:      Default (short) SNES Monitor, same as SNESMonitorDefault() except
532:   it prints fewer digits of the residual as the residual gets smaller.
533:   This is because the later digits are meaningless and are often
534:   different on different machines; by using this routine different
535:   machines will usually generate the same output.

537:   Deprecated: Intentionally has no manual page
538: */
539: PetscErrorCode  SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
540: {
542:   PetscViewer    viewer = vf->viewer;

546:   PetscViewerPushFormat(viewer,vf->format);
547:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
548:   if (fgnorm > 1.e-9) {
549:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %g \n",its,(double)fgnorm);
550:   } else if (fgnorm > 1.e-11) {
551:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,(double)fgnorm);
552:   } else {
553:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);
554:   }
555:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
556:   PetscViewerPopFormat(viewer);
557:   return(0);
558: }

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

563:   Collective on SNES

565:   Input Parameters:
566: + snes   - the SNES context
567: . its    - iteration number
568: . fgnorm - 2-norm of residual
569: - ctx    - the PetscViewer

571:   Notes:
572:   This routine uses the DM attached to the residual vector

574:   Level: intermediate

576: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault()
577: @*/
578: PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
579: {
580:   PetscViewer    viewer = vf->viewer;
581:   Vec            r;
582:   DM             dm;
583:   PetscReal      res[256];
584:   PetscInt       tablevel;

589:   SNESGetFunction(snes, &r, NULL, NULL);
590:   VecGetDM(r, &dm);
591:   if (!dm) {SNESMonitorDefault(snes, its, fgnorm, vf);}
592:   else {
593:     PetscSection s, gs;
594:     PetscInt     Nf, f;

596:     DMGetLocalSection(dm, &s);
597:     DMGetGlobalSection(dm, &gs);
598:     if (!s || !gs) {SNESMonitorDefault(snes, its, fgnorm, vf);}
599:     PetscSectionGetNumFields(s, &Nf);
600:     if (Nf > 256) SETERRQ1(PetscObjectComm((PetscObject) snes), PETSC_ERR_SUP, "Do not support %d fields > 256", Nf);
601:     PetscSectionVecNorm(s, gs, r, NORM_2, res);
602:     PetscObjectGetTabLevel((PetscObject) snes, &tablevel);
603:     PetscViewerPushFormat(viewer,vf->format);
604:     PetscViewerASCIIAddTab(viewer, tablevel);
605:     PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);
606:     for (f = 0; f < Nf; ++f) {
607:       if (f) {PetscViewerASCIIPrintf(viewer, ", ");}
608:       PetscViewerASCIIPrintf(viewer, "%14.12e", res[f]);
609:     }
610:     PetscViewerASCIIPrintf(viewer, "] \n");
611:     PetscViewerASCIISubtractTab(viewer, tablevel);
612:     PetscViewerPopFormat(viewer);
613:   }
614:   return(0);
615: }
616: /* ---------------------------------------------------------------- */
617: /*@C
618:    SNESConvergedDefault - Convergence test of the solvers for
619:    systems of nonlinear equations (default).

621:    Collective on SNES

623:    Input Parameters:
624: +  snes - the SNES context
625: .  it - the iteration (0 indicates before any Newton steps)
626: .  xnorm - 2-norm of current iterate
627: .  snorm - 2-norm of current step
628: .  fnorm - 2-norm of function at current iterate
629: -  dummy - unused context

631:    Output Parameter:
632: .   reason  - one of
633: $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
634: $  SNES_CONVERGED_SNORM_RELATIVE  - (snorm < stol*xnorm),
635: $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
636: $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
637: $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
638: $  SNES_CONVERGED_ITERATING       - (otherwise),

640:    where
641: +    maxf - maximum number of function evaluations,
642:             set with SNESSetTolerances()
643: .    nfct - number of function evaluations,
644: .    abstol - absolute function norm tolerance,
645:             set with SNESSetTolerances()
646: -    rtol - relative function norm tolerance, set with SNESSetTolerances()

648:    Level: intermediate

650: .seealso: SNESSetConvergenceTest()
651: @*/
652: PetscErrorCode  SNESConvergedDefault(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
653: {


660:   *reason = SNES_CONVERGED_ITERATING;

662:   if (!it) {
663:     /* set parameter for default relative tolerance convergence test */
664:     snes->ttol   = fnorm*snes->rtol;
665:     snes->rnorm0 = fnorm;
666:   }
667:   if (PetscIsInfOrNanReal(fnorm)) {
668:     PetscInfo(snes,"Failed to converged, function norm is NaN\n");
669:     *reason = SNES_DIVERGED_FNORM_NAN;
670:   } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
671:     PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);
672:     *reason = SNES_CONVERGED_FNORM_ABS;
673:   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
674:     PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);
675:     *reason = SNES_DIVERGED_FUNCTION_COUNT;
676:   }

678:   if (it && !*reason) {
679:     if (fnorm <= snes->ttol) {
680:       PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);
681:       *reason = SNES_CONVERGED_FNORM_RELATIVE;
682:     } else if (snorm < snes->stol*xnorm) {
683:       PetscInfo3(snes,"Converged due to small update length: %14.12e < %14.12e * %14.12e\n",(double)snorm,(double)snes->stol,(double)xnorm);
684:       *reason = SNES_CONVERGED_SNORM_RELATIVE;
685:     } else if (snes->divtol > 0 && (fnorm > snes->divtol*snes->rnorm0)) {
686:       PetscInfo3(snes,"Diverged due to increase in function norm: %14.12e > %14.12e * %14.12e\n",(double)fnorm,(double)snes->divtol,(double)snes->rnorm0);
687:       *reason = SNES_DIVERGED_DTOL;
688:     }

690:   }
691:   return(0);
692: }

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

698:    Logically Collective on SNES

700:    Input Parameters:
701: +  snes - the SNES context
702: .  it - the iteration (0 indicates before any Newton steps)
703: .  xnorm - 2-norm of current iterate
704: .  snorm - 2-norm of current step
705: .  fnorm - 2-norm of function at current iterate
706: -  dummy - unused context

708:    Output Parameter:
709: .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN

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

714:    Level: advanced

716: .seealso: SNESSetConvergenceTest()
717: @*/
718: PetscErrorCode  SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
719: {


726:   *reason = SNES_CONVERGED_ITERATING;

728:   if (fnorm != fnorm) {
729:     PetscInfo(snes,"Failed to converged, function norm is NaN\n");
730:     *reason = SNES_DIVERGED_FNORM_NAN;
731:   } else if (it == snes->max_its) {
732:     *reason = SNES_CONVERGED_ITS;
733:   }
734:   return(0);
735: }

737: /*@C
738:   SNESSetWorkVecs - Gets a number of work vectors.

740:   Input Parameters:
741: + snes  - the SNES context
742: - nw - number of work vectors to allocate

744:   Level: developer

746: @*/
747: PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw)
748: {
749:   DM             dm;
750:   Vec            v;

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

757:   SNESGetDM(snes, &dm);
758:   DMGetGlobalVector(dm, &v);
759:   VecDuplicateVecs(v,snes->nwork,&snes->work);
760:   DMRestoreGlobalVector(dm, &v);
761:   PetscLogObjectParents(snes,nw,snes->work);
762:   return(0);
763: }