Actual source code: snesut.c

petsc-3.7.5 2017-01-01
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 -  a viewer

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

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

 47:    Collective on SNES

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

 55:    Level: intermediate

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

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

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

 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:    Level: intermediate

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

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

104:   SNESGetSolutionUpdate(snes,&x);
105:   PetscViewerPushFormat(viewer,vf->format);
106:   VecView(x,viewer);
107:   PetscViewerPopFormat(viewer);
108:   return(0);
109: }

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

116:    Collective on KSP

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

124:    Level: intermediate

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

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

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

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

159: #include <petscdraw.h>

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

167:    Collective on KSP

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

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

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

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

186:    Level: intermediate

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

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

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

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

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

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

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

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

251: /*@
252:    KSPMonitorSNESLGResidualNormDestroy - Destroys a line graph context that was created
253:    with KSPMonitorSNESLGResidualNormCreate().

255:    Collective on KSP

257:    Input Parameter:
258: .  draw - the drawing context

260:    Level: intermediate

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

264: .seealso: KSPMonitorSNESLGResidualNormCreate(), KSPMonitorSNESLGTrueResidualDestroy(), KSPMonitorSet()
265: @*/
266: PetscErrorCode  KSPMonitorSNESLGResidualNormDestroy(PetscObject **objs)
267: {
269:   PetscDrawLG    lg = (PetscDrawLG) (*objs)[1];

272:   PetscDrawLGDestroy(&lg);
273:   PetscFree(*objs);
274:   return(0);
275: }

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

282:    Collective on SNES

284:    Input Parameters:
285: +  snes - the SNES context
286: .  its - iteration number
287: .  fgnorm - 2-norm of residual
288: -  vf - viewer and format structure

290:    Notes:
291:    This routine prints the residual norm at each iteration.

293:    Level: intermediate

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

297: .seealso: SNESMonitorSet(), SNESMonitorSolution()
298: @*/
299: PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
300: {
302:   PetscViewer    viewer = vf->viewer;

306:   PetscViewerPushFormat(viewer,vf->format);
307:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
308:   PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
309:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
310:   PetscViewerPopFormat(viewer);CHKERRQ(ierr) ;
311:   return(0);
312: }

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

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

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

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

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

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

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

409:    Collective on SNES

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

417:    Options Database Key:
418: .  -snes_monitor_range - Activates SNESMonitorRange()

420:    Level: intermediate

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

424: .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
425: @*/
426: PetscErrorCode  SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,PetscViewerAndFormat *vf)
427: {
429:   PetscReal      perc,rel;
430:   PetscViewer    viewer = vf->viewer;
431:   /* should be in a MonitorRangeContext */
432:   static PetscReal prev;

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

439:   rel  = (prev - rnorm)/prev;
440:   prev = rnorm;
441:   PetscViewerPushFormat(viewer,vf->format);
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:   PetscViewerPopFormat(viewer);
446:   return(0);
447: }

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

455:    Collective on SNES

457:    Input Parameters:
458: +  snes - the SNES context
459: .  its - iteration number
460: .  fgnorm - 2-norm of residual (or gradient)
461: -  dummy -  context of monitor

463:    Level: intermediate

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

468: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorRatio()
469: @*/
470: PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
471: {
472:   PetscErrorCode          ierr;
473:   PetscInt                len;
474:   PetscReal               *history;
475:   PetscViewer             viewer = vf->viewer;
476: 
478:   SNESGetConvergenceHistory(snes,&history,NULL,&len);
479:   PetscViewerPushFormat(viewer,vf->format);
480:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
481:   if (!its || !history || its > len) {
482:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
483:   } else {
484:     PetscReal ratio = fgnorm/history[its-1];
485:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);
486:   }
487:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
488:   PetscViewerPopFormat(viewer);
489:   return(0);
490: }

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

497:    Collective on SNES

499:    Input Parameters:
500: +   snes - the SNES context
501: -   viewer - the PetscViewer object (ignored)

503:    Level: intermediate

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

507: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorRatio()
508: @*/
509: PetscErrorCode  SNESMonitorRatioSetUp(SNES snes,PetscViewerAndFormat *vf)
510: {
511:   PetscErrorCode          ierr;
512:   PetscReal               *history;

515:   SNESGetConvergenceHistory(snes,&history,NULL,NULL);
516:   if (!history) {
517:     SNESSetConvergenceHistory(snes,NULL,NULL,100,PETSC_TRUE);
518:   }
519:   return(0);
520: }

522: /* ---------------------------------------------------------------- */
525: /*
526:      Default (short) SNES Monitor, same as SNESMonitorDefault() except
527:   it prints fewer digits of the residual as the residual gets smaller.
528:   This is because the later digits are meaningless and are often
529:   different on different machines; by using this routine different
530:   machines will usually generate the same output.
531: */
532: PetscErrorCode  SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *vf)
533: {
535:   PetscViewer    viewer = vf->viewer;

539:   PetscViewerPushFormat(viewer,vf->format);
540:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
541:   if (fgnorm > 1.e-9) {
542:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %g \n",its,(double)fgnorm);
543:   } else if (fgnorm > 1.e-11) {
544:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,(double)fgnorm);
545:   } else {
546:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);
547:   }
548:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
549:   PetscViewerPopFormat(viewer);
550:   return(0);
551: }

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

558:   Collective on SNES

560:   Input Parameters:
561: + snes   - the SNES context
562: . its    - iteration number
563: . fgnorm - 2-norm of residual
564: - ctx    - the PetscViewer

566:   Notes:
567:   This routine uses the DM attached to the residual vector

569:   Level: intermediate

571: .keywords: SNES, nonlinear, field, monitor, norm
572: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorDefaultShort()
573: @*/
574: PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
575: {
576:   PetscViewer    viewer = vf->viewer;
577:   Vec            r;
578:   DM             dm;
579:   PetscReal      res[256];
580:   PetscInt       tablevel;

585:   SNESGetFunction(snes, &r, NULL, NULL);
586:   VecGetDM(r, &dm);
587:   if (!dm) {SNESMonitorDefault(snes, its, fgnorm, vf);}
588:   else {
589:     PetscSection s, gs;
590:     PetscInt     Nf, f;

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

619:    Collective on SNES

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

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

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

646:    Level: intermediate

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

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:   }
666:   if (PetscIsInfOrNanReal(fnorm)) {
667:     PetscInfo(snes,"Failed to converged, function norm is NaN\n");
668:     *reason = SNES_DIVERGED_FNORM_NAN;
669:   } else if (fnorm < snes->abstol) {
670:     PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);
671:     *reason = SNES_CONVERGED_FNORM_ABS;
672:   } else if (snes->nfuncs >= snes->max_funcs) {
673:     PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);
674:     *reason = SNES_DIVERGED_FUNCTION_COUNT;
675:   }

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

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

695:    Logically Collective on SNES

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

705:    Output Parameter:
706: .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN

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

711:    Level: advanced

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

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


725:   *reason = SNES_CONVERGED_ITERATING;

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

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

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

745:    Level: developer

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

749: @*/
750: PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw)
751: {
752:   DM             dm;
753:   Vec            v;

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

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