Actual source code: snesut.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: #include <petsc/private/snesimpl.h>       /*I   "petsc/private/snesimpl.h"   I*/
  3: #include <petscdm.h>
  4: #include <petscblaslapack.h>

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

 12:    Collective on SNES

 14:    Input Parameters:
 15: +  snes - the SNES context
 16: .  its - iteration number
 17: .  fgnorm - 2-norm of residual
 18: -  dummy - either a viewer or NULL

 20:    Level: intermediate

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

 24: .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
 25: @*/
 26: PetscErrorCode  SNESMonitorSolution(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
 27: {
 29:   Vec            x;
 30:   PetscViewer    viewer = (PetscViewer) dummy;

 33:   SNESGetSolution(snes,&x);
 34:   if (!viewer) {
 35:     MPI_Comm comm;
 36:     PetscObjectGetComm((PetscObject)snes,&comm);
 37:     viewer = PETSC_VIEWER_DRAW_(comm);
 38:   }
 39:   VecView(x,viewer);
 40:   return(0);
 41: }

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

 49:    Collective on SNES

 51:    Input Parameters:
 52: +  snes - the SNES context
 53: .  its - iteration number
 54: .  fgnorm - 2-norm of residual
 55: -  dummy - either a viewer or NULL

 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,void *dummy)
 64: {
 66:   Vec            x;
 67:   PetscViewer    viewer = (PetscViewer) dummy;

 70:   SNESGetFunction(snes,&x,0,0);
 71:   if (!viewer) {
 72:     MPI_Comm comm;
 73:     PetscObjectGetComm((PetscObject)snes,&comm);
 74:     viewer = PETSC_VIEWER_DRAW_(comm);
 75:   }
 76:   VecView(x,viewer);
 77:   return(0);
 78: }

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

 86:    Collective on SNES

 88:    Input Parameters:
 89: +  snes - the SNES context
 90: .  its - iteration number
 91: .  fgnorm - 2-norm of residual
 92: -  dummy - either a viewer or NULL

 94:    Level: intermediate

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

 98: .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
 99: @*/
100: PetscErrorCode  SNESMonitorSolutionUpdate(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
101: {
103:   Vec            x;
104:   PetscViewer    viewer = (PetscViewer) dummy;

107:   SNESGetSolutionUpdate(snes,&x);
108:   if (!viewer) {
109:     MPI_Comm comm;
110:     PetscObjectGetComm((PetscObject)snes,&comm);
111:     viewer = PETSC_VIEWER_DRAW_(comm);
112:   }
113:   VecView(x,viewer);
114:   return(0);
115: }

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

122:    Collective on KSP

124:    Input Parameters:
125: +  ksp   - iterative context
126: .  n     - iteration number
127: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
128: -  dummy - unused monitor context

130:    Level: intermediate

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

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

145:   SNESGetSolution(snes,&snes_solution);
146:   VecDuplicate(snes_solution,&work1);
147:   VecDuplicate(snes_solution,&work2);
148:   KSPBuildSolution(ksp,work1,NULL);
149:   VecAYPX(work1,-1.0,snes_solution);
150:   SNESComputeFunction(snes,work1,work2);
151:   VecNorm(work2,NORM_2,&snorm);
152:   VecDestroy(&work1);
153:   VecDestroy(&work2);

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

165: #include <petscdraw.h>

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

173:    Collective on KSP

175:    Input Parameters:
176: +  host - the X display to open, or null for the local machine
177: .  label - the title to put in the title bar
178: .  x, y - the screen coordinates of the upper left coordinate of
179:           the window
180: -  m, n - the screen width and height in pixels

182:    Output Parameter:
183: .  draw - the drawing context

185:    Options Database Key:
186: .  -ksp_monitor_lg_residualnorm - Sets line graph monitor

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

191:    Level: intermediate

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

195: .seealso: KSPMonitorSNESLGResidualNormDestroy(), KSPMonitorSet(), KSPMonitorSNESLGTrueResidualCreate()
196: @*/
197: PetscErrorCode  KSPMonitorSNESLGResidualNormCreate(const char host[],const char label[],int x,int y,int m,int n,PetscObject **objs)
198: {
199:   PetscDraw      draw;
201:   PetscDrawAxis  axis;
202:   PetscDrawLG    drawlg;
203:   const char     *names[] = {"Linear residual","Nonlinear residual"};

206:   PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&draw);
207:   PetscDrawSetFromOptions(draw);
208:   PetscDrawLGCreate(draw,2,&drawlg);
209:   PetscDrawLGSetFromOptions(drawlg);
210:   PetscDrawLGGetAxis(drawlg,&axis);
211:   PetscDrawAxisSetLabels(axis,"Convergence of Residual Norm","Iteration","Residual Norm");
212:   PetscDrawLGSetLegend(drawlg,names);
213:   PetscLogObjectParent((PetscObject)drawlg,(PetscObject)draw);

215:   PetscMalloc1(3,objs);
216:   (*objs)[1] = (PetscObject)drawlg;
217:   (*objs)[2] = (PetscObject)draw;
218:   return(0);
219: }

223: PetscErrorCode  KSPMonitorSNESLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscObject *objs)
224: {
225:   PetscDrawLG    lg = (PetscDrawLG) objs[1];
227:   PetscReal      y[2];
228:   SNES           snes = (SNES) objs[0];
229:   Vec            snes_solution,work1,work2;

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

235:   SNESGetSolution(snes,&snes_solution);
236:   VecDuplicate(snes_solution,&work1);
237:   VecDuplicate(snes_solution,&work2);
238:   KSPBuildSolution(ksp,work1,NULL);
239:   VecAYPX(work1,-1.0,snes_solution);
240:   SNESComputeFunction(snes,work1,work2);
241:   VecNorm(work2,NORM_2,y+1);
242:   if (y[1] > 0.0) y[1] = PetscLog10Real(y[1]);
243:   else y[1] = -15.0;
244:   VecDestroy(&work1);
245:   VecDestroy(&work2);

247:   PetscDrawLGAddPoint(lg,NULL,y);
248:   if (n < 20 || !(n % 5)) {
249:     PetscDrawLGDraw(lg);
250:   }
251:   return(0);
252: }

256: /*@
257:    KSPMonitorSNESLGResidualNormDestroy - Destroys a line graph context that was created
258:    with KSPMonitorSNESLGResidualNormCreate().

260:    Collective on KSP

262:    Input Parameter:
263: .  draw - the drawing context

265:    Level: intermediate

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

269: .seealso: KSPMonitorSNESLGResidualNormCreate(), KSPMonitorSNESLGTrueResidualDestroy(), KSPMonitorSet()
270: @*/
271: PetscErrorCode  KSPMonitorSNESLGResidualNormDestroy(PetscObject **objs)
272: {
274:   PetscDrawLG    drawlg = (PetscDrawLG) (*objs)[1];
275:   PetscDraw      draw = (PetscDraw) (*objs)[2];

278:   PetscDrawDestroy(&draw);
279:   PetscDrawLGDestroy(&drawlg);
280:   PetscFree(*objs);
281:   return(0);
282: }

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

289:    Collective on SNES

291:    Input Parameters:
292: +  snes - the SNES context
293: .  its - iteration number
294: .  fgnorm - 2-norm of residual
295: -  dummy - unused context

297:    Notes:
298:    This routine prints the residual norm at each iteration.

300:    Level: intermediate

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

304: .seealso: SNESMonitorSet(), SNESMonitorSolution()
305: @*/
306: PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
307: {
309:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));

312:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
313:   PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
314:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
315:   return(0);
316: }

320: PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,void *ctx)
321: {
322: #if defined(PETSC_MISSING_LAPACK_GEEV)
323:   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
324: #elif defined(PETSC_HAVE_ESSL)
325:   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - No support for ESSL Lapack Routines");
326: #else
327:   Vec            X;
328:   Mat            J,dJ,dJdense;
330:   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
331:   PetscInt       n,i;
332:   PetscBLASInt   nb,lwork;
333:   PetscReal      *eigr,*eigi;
334:   PetscScalar    *work;
335:   PetscScalar    *a;

338:   if (it == 0) return(0);
339:   /* create the difference between the current update and the current jacobian */
340:   SNESGetSolution(snes,&X);
341:   SNESGetJacobian(snes,NULL,&J,&func,NULL);
342:   MatDuplicate(J,MAT_COPY_VALUES,&dJ);
343:   SNESComputeJacobian(snes,X,dJ,dJ);
344:   MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);

346:   /* compute the spectrum directly */
347:   MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);
348:   MatGetSize(dJ,&n,NULL);
349:   PetscBLASIntCast(n,&nb);
350:   lwork = 3*nb;
351:   PetscMalloc1(n,&eigr);
352:   PetscMalloc1(n,&eigi);
353:   PetscMalloc1(lwork,&work);
354:   MatDenseGetArray(dJdense,&a);
355: #if !defined(PETSC_USE_COMPLEX)
356:   {
357:     PetscBLASInt lierr;
358:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
359:     PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,NULL,&nb,NULL,&nb,work,&lwork,&lierr));
360:     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr);
361:     PetscFPTrapPop();
362:   }
363: #else
364:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
365: #endif
366:   PetscPrintf(PetscObjectComm((PetscObject)snes),"Eigenvalues of J_%d - J_%d:\n",it,it-1);
367:   for (i=0;i<n;i++) {
368:     PetscPrintf(PetscObjectComm((PetscObject)snes),"%5d: %20.5g + %20.5gi\n",i,(double)eigr[i],(double)eigi[i]);
369:   }
370:   MatDenseRestoreArray(dJdense,&a);
371:   MatDestroy(&dJ);
372:   MatDestroy(&dJdense);
373:   PetscFree(eigr);
374:   PetscFree(eigi);
375:   PetscFree(work);
376:   return(0);
377: #endif
378: }

382: PetscErrorCode  SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
383: {
385:   Vec            resid;
386:   PetscReal      rmax,pwork;
387:   PetscInt       i,n,N;
388:   PetscScalar    *r;

391:   SNESGetFunction(snes,&resid,0,0);
392:   VecNorm(resid,NORM_INFINITY,&rmax);
393:   VecGetLocalSize(resid,&n);
394:   VecGetSize(resid,&N);
395:   VecGetArray(resid,&r);
396:   pwork = 0.0;
397:   for (i=0; i<n; i++) {
398:     pwork += (PetscAbsScalar(r[i]) > .20*rmax);
399:   }
400:   MPI_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));
401:   VecRestoreArray(resid,&r);
402:   *per = *per/N;
403:   return(0);
404: }

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

411:    Collective on SNES

413:    Input Parameters:
414: +  snes   - iterative context
415: .  it    - iteration number
416: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
417: -  dummy - unused monitor context

419:    Options Database Key:
420: .  -snes_monitor_range - Activates SNESMonitorRange()

422:    Level: intermediate

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

426: .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
427: @*/
428: PetscErrorCode  SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,void *dummy)
429: {
431:   PetscReal      perc,rel;
432:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
433:   /* should be in a MonitorRangeContext */
434:   static PetscReal prev;

437:   if (!it) prev = rnorm;
438:   SNESMonitorRange_Private(snes,it,&perc);

440:   rel  = (prev - rnorm)/prev;
441:   prev = rnorm;
442:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
443:   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));
444:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
445:   return(0);
446: }

448: typedef struct {
449:   PetscViewer viewer;
450:   PetscReal   *history;
451: } SNESMonitorRatioContext;

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

459:    Collective on SNES

461:    Input Parameters:
462: +  snes - the SNES context
463: .  its - iteration number
464: .  fgnorm - 2-norm of residual (or gradient)
465: -  dummy -  context of monitor

467:    Level: intermediate

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

471: .seealso: SNESMonitorSet(), SNESMonitorSolution()
472: @*/
473: PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
474: {
475:   PetscErrorCode          ierr;
476:   PetscInt                len;
477:   PetscReal               *history;
478:   SNESMonitorRatioContext *ctx = (SNESMonitorRatioContext*)dummy;

481:   SNESGetConvergenceHistory(snes,&history,NULL,&len);
482:   PetscViewerASCIIAddTab(ctx->viewer,((PetscObject)snes)->tablevel);
483:   if (!its || !history || its > len) {
484:     PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
485:   } else {
486:     PetscReal ratio = fgnorm/history[its-1];
487:     PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);
488:   }
489:   PetscViewerASCIISubtractTab(ctx->viewer,((PetscObject)snes)->tablevel);
490:   return(0);
491: }

493: /*
494:    If the we set the history monitor space then we must destroy it
495: */
498: PetscErrorCode SNESMonitorRatioDestroy(void **ct)
499: {
500:   PetscErrorCode          ierr;
501:   SNESMonitorRatioContext *ctx = *(SNESMonitorRatioContext**)ct;

504:   PetscFree(ctx->history);
505:   PetscViewerDestroy(&ctx->viewer);
506:   PetscFree(ctx);
507:   return(0);
508: }

512: /*@C
513:    SNESMonitorSetRatio - Sets SNES to use a monitor that prints the
514:    ratio of the function norm at each iteration.

516:    Collective on SNES

518:    Input Parameters:
519: +   snes - the SNES context
520: -   viewer - ASCII viewer to print output

522:    Level: intermediate

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

526: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault()
527: @*/
528: PetscErrorCode  SNESMonitorSetRatio(SNES snes,PetscViewer viewer)
529: {
530:   PetscErrorCode          ierr;
531:   SNESMonitorRatioContext *ctx;
532:   PetscReal               *history;

535:   if (!viewer) {
536:     PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),"stdout",&viewer);
537:     PetscObjectReference((PetscObject)viewer);
538:   }
539:   PetscNewLog(snes,&ctx);
540:   SNESGetConvergenceHistory(snes,&history,NULL,NULL);
541:   if (!history) {
542:     PetscMalloc1(100,&ctx->history);
543:     SNESSetConvergenceHistory(snes,ctx->history,0,100,PETSC_TRUE);
544:   }
545:   ctx->viewer = viewer;

547:   SNESMonitorSet(snes,SNESMonitorRatio,ctx,SNESMonitorRatioDestroy);
548:   return(0);
549: }

551: /* ---------------------------------------------------------------- */
554: /*
555:      Default (short) SNES Monitor, same as SNESMonitorDefault() except
556:   it prints fewer digits of the residual as the residual gets smaller.
557:   This is because the later digits are meaningless and are often
558:   different on different machines; by using this routine different
559:   machines will usually generate the same output.
560: */
561: PetscErrorCode  SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
562: {
564:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));

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:   return(0);
577: }

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(), SNESMonitorDefaultShort()
599: @*/
600: PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, void *ctx)
601: {
602:   PetscViewer    viewer = ctx ? (PetscViewer) ctx : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject) snes));
603:   Vec            r;
604:   DM             dm;
605:   PetscReal      res[256];
606:   PetscInt       tablevel;

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

617:     DMGetDefaultSection(dm, &s);
618:     DMGetDefaultGlobalSection(dm, &gs);
619:     if (!s || !gs) {SNESMonitorDefault(snes, its, fgnorm, ctx);}
620:     PetscSectionGetNumFields(s, &Nf);
621:     if (Nf > 256) SETERRQ1(PetscObjectComm((PetscObject) snes), PETSC_ERR_SUP, "Do not support %d fields > 256", Nf);
622:     PetscSectionVecNorm(s, gs, r, NORM_2, res);
623:     PetscObjectGetTabLevel((PetscObject) snes, &tablevel);
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:   }
633:   return(0);
634: }
635: /* ---------------------------------------------------------------- */
638: /*@C
639:    SNESConvergedDefault - Convergence test of the solvers for
640:    systems of nonlinear equations (default).

642:    Collective on SNES

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

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

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

669:    Level: intermediate

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

673: .seealso: SNESSetConvergenceTest()
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:   }
689:   if (PetscIsInfOrNanReal(fnorm)) {
690:     PetscInfo(snes,"Failed to converged, function norm is NaN\n");
691:     *reason = SNES_DIVERGED_FNORM_NAN;
692:   } else if (fnorm < snes->abstol) {
693:     PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);
694:     *reason = SNES_CONVERGED_FNORM_ABS;
695:   } else if (snes->nfuncs >= snes->max_funcs) {
696:     PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);
697:     *reason = SNES_DIVERGED_FUNCTION_COUNT;
698:   }

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

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

718:    Logically Collective on SNES

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

728:    Output Parameter:
729: .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN

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

734:    Level: advanced

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

738: .seealso: SNESSetConvergenceTest()
739: @*/
740: PetscErrorCode  SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
741: {


748:   *reason = SNES_CONVERGED_ITERATING;

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

761: /*@C
762:   SNESSetWorkVecs - Gets a number of work vectors.

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

768:    Level: developer

770:    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations

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: }