Actual source code: snesut.c

petsc-3.3-p7 2013-05-11
  2: #include <petsc-private/snesimpl.h>       /*I   "petscsnes.h"   I*/

  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 - either a viewer or PETSC_NULL

 18:    Level: intermediate

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

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

 31:   SNESGetSolution(snes,&x);
 32:   if (!viewer) {
 33:     MPI_Comm comm;
 34:     PetscObjectGetComm((PetscObject)snes,&comm);
 35:     viewer = PETSC_VIEWER_DRAW_(comm);
 36:   }
 37:   VecView(x,viewer);

 39:   return(0);
 40: }

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

 48:    Collective on SNES

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

 56:    Level: intermediate

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

 60: .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
 61: @*/
 62: PetscErrorCode  SNESMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
 63: {
 65:   Vec            x;
 66:   PetscViewer    viewer = (PetscViewer) dummy;

 69:   SNESGetFunction(snes,&x,0,0);
 70:   if (!viewer) {
 71:     MPI_Comm comm;
 72:     PetscObjectGetComm((PetscObject)snes,&comm);
 73:     viewer = PETSC_VIEWER_DRAW_(comm);
 74:   }
 75:   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 PETSC_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);

115:   return(0);
116: }

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

123:    Collective on SNES

125:    Input Parameters:
126: +  snes - the SNES context
127: .  its - iteration number
128: .  fgnorm - 2-norm of residual
129: -  dummy - unused context

131:    Notes:
132:    This routine prints the residual norm at each iteration.

134:    Level: intermediate

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

138: .seealso: SNESMonitorSet(), SNESMonitorSolution()
139: @*/
140: PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
141: {
143:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);

146:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
147:   PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
148:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
149:   return(0);
150: }

154: PetscErrorCode  SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
155: {
156:   PetscErrorCode          ierr;
157:   Vec                     resid;
158:   PetscReal               rmax,pwork;
159:   PetscInt                i,n,N;
160:   PetscScalar             *r;

163:   SNESGetFunction(snes,&resid,0,0);
164:   VecNorm(resid,NORM_INFINITY,&rmax);
165:   VecGetLocalSize(resid,&n);
166:   VecGetSize(resid,&N);
167:   VecGetArray(resid,&r);
168:   pwork = 0.0;
169:   for (i=0; i<n; i++) {
170:     pwork += (PetscAbsScalar(r[i]) > .20*rmax);
171:   }
172:   MPI_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);
173:   VecRestoreArray(resid,&r);
174:   *per  = *per/N;
175:   return(0);
176: }

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

183:    Collective on SNES

185:    Input Parameters:
186: +  snes   - iterative context
187: .  it    - iteration number
188: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).  
189: -  dummy - unused monitor context 

191:    Options Database Key:
192: .  -snes_monitor_range - Activates SNESMonitorRange()

194:    Level: intermediate

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

198: .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
199: @*/
200: PetscErrorCode  SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,void *dummy)
201: {
202:   PetscErrorCode   ierr;
203:   PetscReal        perc,rel;
204:   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
205:   /* should be in a MonitorRangeContext */
206:   static PetscReal prev;

209:   if (!it) prev = rnorm;
210:   SNESMonitorRange_Private(snes,it,&perc);

212:   rel  = (prev - rnorm)/prev;
213:   prev = rnorm;
214:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
215:   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);
216:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
217:   return(0);
218: }

220: typedef struct {
221:   PetscViewer viewer;
222:   PetscReal   *history;
223: } SNESMonitorRatioContext;

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

231:    Collective on SNES

233:    Input Parameters:
234: +  snes - the SNES context
235: .  its - iteration number
236: .  fgnorm - 2-norm of residual (or gradient)
237: -  dummy -  context of monitor

239:    Level: intermediate

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

243: .seealso: SNESMonitorSet(), SNESMonitorSolution()
244: @*/
245: PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
246: {
247:   PetscErrorCode          ierr;
248:   PetscInt                len;
249:   PetscReal               *history;
250:   SNESMonitorRatioContext *ctx = (SNESMonitorRatioContext*)dummy;

253:   SNESGetConvergenceHistory(snes,&history,PETSC_NULL,&len);
254:   PetscViewerASCIIAddTab(ctx->viewer,((PetscObject)snes)->tablevel);
255:   if (!its || !history || its > len) {
256:     PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);
257:   } else {
258:     PetscReal ratio = fgnorm/history[its-1];
259:     PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);
260:   }
261:   PetscViewerASCIISubtractTab(ctx->viewer,((PetscObject)snes)->tablevel);
262:   return(0);
263: }

265: /*
266:    If the we set the history monitor space then we must destroy it
267: */
270: PetscErrorCode SNESMonitorRatioDestroy(void **ct)
271: {
272:   PetscErrorCode          ierr;
273:   SNESMonitorRatioContext *ctx = *(SNESMonitorRatioContext**)ct;

276:   PetscFree(ctx->history);
277:   PetscViewerDestroy(&ctx->viewer);
278:   PetscFree(ctx);
279:   return(0);
280: }

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

288:    Collective on SNES

290:    Input Parameters:
291: +   snes - the SNES context
292: -   viewer - ASCII viewer to print output

294:    Level: intermediate

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

298: .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault()
299: @*/
300: PetscErrorCode  SNESMonitorSetRatio(SNES snes,PetscViewer viewer)
301: {
302:   PetscErrorCode          ierr;
303:   SNESMonitorRatioContext *ctx;
304:   PetscReal               *history;

307:   if (!viewer) {
308:     PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&viewer);
309:     PetscObjectReference((PetscObject)viewer);
310:   }
311:   PetscNewLog(snes,SNESMonitorRatioContext,&ctx);
312:   SNESGetConvergenceHistory(snes,&history,PETSC_NULL,PETSC_NULL);
313:   if (!history) {
314:     PetscMalloc(100*sizeof(PetscReal),&ctx->history);
315:     SNESSetConvergenceHistory(snes,ctx->history,0,100,PETSC_TRUE);
316:   }
317:   ctx->viewer = viewer;
318:   SNESMonitorSet(snes,SNESMonitorRatio,ctx,SNESMonitorRatioDestroy);
319:   return(0);
320: }

322: /* ---------------------------------------------------------------- */
325: /*
326:      Default (short) SNES Monitor, same as SNESMonitorDefault() except
327:   it prints fewer digits of the residual as the residual gets smaller.
328:   This is because the later digits are meaningless and are often 
329:   different on different machines; by using this routine different 
330:   machines will usually generate the same output.
331: */
332: PetscErrorCode  SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
333: {
335:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);

338:   PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);
339:   if (fgnorm > 1.e-9) {
340:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %G \n",its,fgnorm);
341:   } else if (fgnorm > 1.e-11){
342:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,fgnorm);
343:   } else {
344:     PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);
345:   }
346:   PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);
347:   return(0);
348: }
349: /* ---------------------------------------------------------------- */
352: /*@C 
353:    SNESDefaultConverged - Convergence test of the solvers for
354:    systems of nonlinear equations (default).

356:    Collective on SNES

358:    Input Parameters:
359: +  snes - the SNES context
360: .  it - the iteration (0 indicates before any Newton steps)
361: .  xnorm - 2-norm of current iterate
362: .  snorm - 2-norm of current step 
363: .  fnorm - 2-norm of function at current iterate
364: -  dummy - unused context

366:    Output Parameter:
367: .   reason  - one of
368: $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
369: $  SNES_CONVERGED_SNORM_RELATIVE  - (snorm < stol*xnorm),
370: $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
371: $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
372: $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
373: $  SNES_CONVERGED_ITERATING       - (otherwise),

375:    where
376: +    maxf - maximum number of function evaluations,
377:             set with SNESSetTolerances()
378: .    nfct - number of function evaluations,
379: .    abstol - absolute function norm tolerance,
380:             set with SNESSetTolerances()
381: -    rtol - relative function norm tolerance, set with SNESSetTolerances()

383:    Level: intermediate

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

387: .seealso: SNESSetConvergenceTest()
388: @*/
389: PetscErrorCode  SNESDefaultConverged(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
390: {

396: 
397:   *reason = SNES_CONVERGED_ITERATING;

399:   if (!it) {
400:     /* set parameter for default relative tolerance convergence test */
401:     snes->ttol = fnorm*snes->rtol;
402:   }
403:   if (fnorm != fnorm) {
404:     PetscInfo(snes,"Failed to converged, function norm is NaN\n");
405:     *reason = SNES_DIVERGED_FNORM_NAN;
406:   } else if (fnorm < snes->abstol) {
407:     PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);
408:     *reason = SNES_CONVERGED_FNORM_ABS;
409:   } else if (snes->nfuncs >= snes->max_funcs) {
410:     PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);
411:     *reason = SNES_DIVERGED_FUNCTION_COUNT;
412:   }

414:   if (it && !*reason) {
415:     if (fnorm <= snes->ttol) {
416:       PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);
417:       *reason = SNES_CONVERGED_FNORM_RELATIVE;
418:     } else if (snorm < snes->stol*xnorm) {
419:       PetscInfo3(snes,"Converged due to small update length: %14.12e < %14.12e * %14.12e\n",(double)snorm,(double)snes->stol,(double)xnorm);
420:       *reason = SNES_CONVERGED_SNORM_RELATIVE;
421:     }
422:   }
423:   return(0);
424: }

428: /*@C 
429:    SNESSkipConverged - Convergence test for SNES that NEVER returns as
430:    converged, UNLESS the maximum number of iteration have been reached.

432:    Logically Collective on SNES

434:    Input Parameters:
435: +  snes - the SNES context
436: .  it - the iteration (0 indicates before any Newton steps)
437: .  xnorm - 2-norm of current iterate
438: .  snorm - 2-norm of current step 
439: .  fnorm - 2-norm of function at current iterate
440: -  dummy - unused context

442:    Output Parameter:
443: .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN

445:    Notes:
446:    Convergence is then declared after a fixed number of iterations have been used.
447:                     
448:    Level: advanced

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

452: .seealso: SNESSetConvergenceTest()
453: @*/
454: PetscErrorCode  SNESSkipConverged(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
455: {


462:   *reason = SNES_CONVERGED_ITERATING;

464:   if (fnorm != fnorm) {
465:     PetscInfo(snes,"Failed to converged, function norm is NaN\n");
466:     *reason = SNES_DIVERGED_FNORM_NAN;
467:   } else if(it == snes->max_its) {
468:     *reason = SNES_CONVERGED_ITS;
469:   }
470:   return(0);
471: }

475: /*@
476:   SNESDefaultGetWork - Gets a number of work vectors.

478:   Input Parameters:
479: . snes  - the SNES context
480: . nw - number of work vectors to allocate

482:    Level: developer

484:   Notes:
485:   Call this only if no work vectors have been allocated
486: @*/
487: PetscErrorCode SNESDefaultGetWork(SNES snes,PetscInt nw)
488: {

492:   if (snes->work) {VecDestroyVecs(snes->nwork,&snes->work);}
493:   snes->nwork = nw;
494:   VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);
495:   PetscLogObjectParents(snes,nw,snes->work);
496:   return(0);
497: }