Actual source code: snesut.c

  1: /*$Id: snesut.c,v 1.64 2001/04/10 19:36:48 bsmith Exp $*/

  3: #include "src/snes/snesimpl.h"       /*I   "petscsnes.h"   I*/

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

  9:    Collective on SNES

 11:    Input Parameters:
 12: +  snes - the SNES context
 13: .  its - iteration number
 14: .  fgnorm - 2-norm of residual (or gradient)
 15: -  dummy - either a viewer or PETSC_NULL

 17:    Level: intermediate

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

 21: .seealso: SNESSetMonitor(), SNESDefaultMonitor(), VecView()
 22: @*/
 23: int SNESVecViewMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
 24: {
 25:   int         ierr;
 26:   Vec         x;
 27:   PetscViewer viewer = (PetscViewer) dummy;

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

 38:   return(0);
 39: }

 41: /*@C
 42:    SNESVecViewUpdateMonitor - Monitors progress of the SNES solvers by calling 
 43:    VecView() for the UPDATE to the solution 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 (or gradient)
 51: -  dummy - either a viewer or PETSC_NULL

 53:    Level: intermediate

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

 57: .seealso: SNESSetMonitor(), SNESDefaultMonitor(), VecView()
 58: @*/
 59: int SNESVecViewUpdateMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
 60: {
 61:   int         ierr;
 62:   Vec         x;
 63:   PetscViewer viewer = (PetscViewer) dummy;

 66:   SNESGetSolutionUpdate(snes,&x);
 67:   if (!viewer) {
 68:     MPI_Comm comm;
 69:     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);
 70:     viewer = PETSC_VIEWER_DRAW_(comm);
 71:   }
 72:   VecView(x,viewer);

 74:   return(0);
 75: }

 77: /*@C
 78:    SNESDefaultMonitor - Monitoring progress of the SNES solvers (default).

 80:    Collective on SNES

 82:    Input Parameters:
 83: +  snes - the SNES context
 84: .  its - iteration number
 85: .  fgnorm - 2-norm of residual (or gradient)
 86: -  dummy - unused context

 88:    Notes:
 89:    For SNES_NONLINEAR_EQUATIONS methods the routine prints the 
 90:    residual norm at each iteration.

 92:    For SNES_UNCONSTRAINED_MINIMIZATION methods the routine prints the
 93:    function value and gradient norm at each iteration.

 95:    Level: intermediate

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

 99: .seealso: SNESSetMonitor(), SNESVecViewMonitor()
100: @*/
101: int SNESDefaultMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
102: {
103:   int         ierr;
104:   PetscViewer viewer = (PetscViewer) dummy;

107:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm);

109:   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
110:     PetscViewerASCIIPrintf(viewer,"%3d SNES Function norm %14.12e n",its,fgnorm);
111:   } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
112:     PetscViewerASCIIPrintf(viewer,"%3d SNES Function value %14.12e, Gradient norm %14.12e n",its,snes->fc,fgnorm);
113:   } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class");
114:   return(0);
115: }

117: /* ---------------------------------------------------------------- */
118: /*
119:      Default (short) SNES Monitor, same as SNESDefaultMonitor() except
120:   it prints fewer digits of the residual as the residual gets smaller.
121:   This is because the later digits are meaningless and are often 
122:   different on different machines; by using this routine different 
123:   machines will usually generate the same output.
124: */
125: int SNESDefaultSMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
126: {

130:   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
131:     if (fgnorm > 1.e-9) {
132:       PetscPrintf(snes->comm,"%3d SNES Function norm %g n",its,fgnorm);
133:     } else if (fgnorm > 1.e-11){
134:       PetscPrintf(snes->comm,"%3d SNES Function norm %5.3e n",its,fgnorm);
135:     } else {
136:       PetscPrintf(snes->comm,"%3d SNES Function norm < 1.e-11n",its);
137:     }
138:   } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
139:     if (fgnorm > 1.e-9) {
140:       PetscPrintf(snes->comm,"%3d SNES Function value %g, Gradient norm %g n",its,snes->fc,fgnorm);
141:     } else if (fgnorm > 1.e-11) {
142:       PetscPrintf(snes->comm,"%3d SNES Function value %g, Gradient norm %5.3e n",its,snes->fc,fgnorm);
143:     } else {
144:       PetscPrintf(snes->comm,"%3d SNES Function value %g, Gradient norm < 1.e-11n",its,snes->fc);
145:     }
146:   } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class");
147:   return(0);
148: }
149: /* ---------------------------------------------------------------- */
150: /*@C 
151:    SNESConverged_EQ_LS - Monitors the convergence of the solvers for
152:    systems of nonlinear equations (default).

154:    Collective on SNES

156:    Input Parameters:
157: +  snes - the SNES context
158: .  xnorm - 2-norm of current iterate
159: .  pnorm - 2-norm of current step 
160: .  fnorm - 2-norm of function
161: -  dummy - unused context

163:    Output Parameter:
164: .   reason  - one of
165: $  SNES_CONVERGED_FNORM_ABS       - (fnorm < atol),
166: $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
167: $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
168: $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
169: $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
170: $  SNES_CONVERGED_ITERATING       - (otherwise),

172:    where
173: +    maxf - maximum number of function evaluations,
174:             set with SNESSetTolerances()
175: .    nfct - number of function evaluations,
176: .    atol - absolute function norm tolerance,
177:             set with SNESSetTolerances()
178: -    rtol - relative function norm tolerance, set with SNESSetTolerances()

180:    Level: intermediate

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

184: .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
185: @*/
186: int SNESConverged_EQ_LS(SNES snes,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
187: {
189:   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
190:      SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
191:   }

193:   if (fnorm != fnorm) {
194:     PetscLogInfo(snes,"SNESConverged_EQ_LS:Failed to converged, function norm is NaNn");
195:     *reason = SNES_DIVERGED_FNORM_NAN;
196:   } else if (fnorm <= snes->ttol) {
197:     PetscLogInfo(snes,"SNESConverged_EQ_LS:Converged due to function norm %g < %g (relative tolerance)n",fnorm,snes->ttol);
198:     *reason = SNES_CONVERGED_FNORM_RELATIVE;
199:   } else if (fnorm < snes->atol) {
200:     PetscLogInfo(snes,"SNESConverged_EQ_LS:Converged due to function norm %g < %gn",fnorm,snes->atol);
201:     *reason = SNES_CONVERGED_FNORM_ABS;
202:   } else if (pnorm < snes->xtol*(xnorm)) {
203:     PetscLogInfo(snes,"SNESConverged_EQ_LS:Converged due to small update length: %g < %g * %gn",pnorm,snes->xtol,xnorm);
204:     *reason = SNES_CONVERGED_PNORM_RELATIVE;
205:   } else if (snes->nfuncs > snes->max_funcs) {
206:     PetscLogInfo(snes,"SNESConverged_EQ_LS:Exceeded maximum number of function evaluations: %d > %dn",snes->nfuncs,snes->max_funcs);
207:     *reason = SNES_DIVERGED_FUNCTION_COUNT ;
208:   } else {
209:     *reason = SNES_CONVERGED_ITERATING;
210:   }
211:   return(0);
212: }
213: /* ------------------------------------------------------------ */
214: /*@
215:    SNES_KSP_SetConvergenceTestEW - Sets alternative convergence test
216:    for the linear solvers within an inexact Newton method.  

218:    Collective on SNES

220:    Input Parameter:
221: .  snes - SNES context

223:    Notes:
224:    Currently, the default is to use a constant relative tolerance for 
225:    the inner linear solvers.  Alternatively, one can use the 
226:    Eisenstat-Walker method, where the relative convergence tolerance 
227:    is reset at each Newton iteration according progress of the nonlinear 
228:    solver. 

230:    Level: advanced

232:    Reference:
233:    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 
234:    inexact Newton method", SISC 17 (1), pp.16-32, 1996.

236: .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
237: @*/
238: int SNES_KSP_SetConvergenceTestEW(SNES snes)
239: {
241:   snes->ksp_ewconv = PETSC_TRUE;
242:   return(0);
243: }

245: /*@
246:    SNES_KSP_SetParametersEW - Sets parameters for Eisenstat-Walker
247:    convergence criteria for the linear solvers within an inexact
248:    Newton method.

250:    Collective on SNES
251:  
252:    Input Parameters:
253: +    snes - SNES context
254: .    version - version 1 or 2 (default is 2)
255: .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
256: .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
257: .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
258: .    alpha2 - power for safeguard
259: .    gamma2 - multiplicative factor for version 2 rtol computation
260:               (0 <= gamma2 <= 1)
261: -    threshold - threshold for imposing safeguard (0 < threshold < 1)

263:    Note:
264:    Use PETSC_DEFAULT to retain the default for any of the parameters.

266:    Level: advanced

268:    Reference:
269:    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 
270:    inexact Newton method", Utah State University Math. Stat. Dept. Res. 
271:    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 

273: .keywords: SNES, KSP, Eisenstat, Walker, set, parameters

275: .seealso: SNES_KSP_SetConvergenceTestEW()
276: @*/
277: int SNES_KSP_SetParametersEW(SNES snes,int version,PetscReal rtol_0,
278:                              PetscReal rtol_max,PetscReal gamma2,PetscReal alpha,
279:                              PetscReal alpha2,PetscReal threshold)
280: {
281:   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;

284:   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
285:   if (version != PETSC_DEFAULT)   kctx->version   = version;
286:   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
287:   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
288:   if (gamma2 != PETSC_DEFAULT)    kctx->gamma     = gamma2;
289:   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
290:   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
291:   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
292:   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
293:     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %g",kctx->rtol_0);
294:   }
295:   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
296:     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max < 1.0n",kctx->rtol_max);
297:   }
298:   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
299:     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold < 1.0n",kctx->threshold);
300:   }
301:   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
302:     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= alpha <= 1.0n",kctx->gamma);
303:   }
304:   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
305:     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha <= 2.0n",kctx->alpha);
306:   }
307:   if (kctx->version != 1 && kctx->version !=2) {
308:     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1 and 2 are supported: %d",kctx->version);
309:   }
310:   return(0);
311: }

313: int SNES_KSP_EW_ComputeRelativeTolerance_Private(SNES snes,KSP ksp)
314: {
315:   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
316:   PetscReal           rtol = 0.0,stol;
317:   int                 ierr;

320:   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
321:   if (!snes->iter) { /* first time in, so use the original user rtol */
322:     rtol = kctx->rtol_0;
323:   } else {
324:     if (kctx->version == 1) {
325:       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
326:       if (rtol < 0.0) rtol = -rtol;
327:       stol = pow(kctx->rtol_last,kctx->alpha2);
328:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
329:     } else if (kctx->version == 2) {
330:       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
331:       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
332:       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
333:     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1 or 2 are supported: %d",kctx->version);
334:   }
335:   rtol = PetscMin(rtol,kctx->rtol_max);
336:   kctx->rtol_last = rtol;
337:   PetscLogInfo(snes,"SNES_KSP_EW_ComputeRelativeTolerance_Private: iter %d, Eisenstat-Walker (version %d) KSP rtol = %gn",snes->iter,kctx->version,rtol);
338:   KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
339:   kctx->norm_last = snes->norm;
340:   return(0);
341: }

343: int SNES_KSP_EW_Converged_Private(KSP ksp,int n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
344: {
345:   SNES                snes = (SNES)ctx;
346:   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
347:   int                 ierr;

350:   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context set");
351:   if (n == 0) {SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);}
352:   KSPDefaultConverged(ksp,n,rnorm,reason,ctx);
353:   kctx->lresid_last = rnorm;
354:   if (*reason) {
355:     PetscLogInfo(snes,"SNES_KSP_EW_Converged_Private: KSP iterations=%d, rnorm=%gn",n,rnorm);
356:   }
357:   return(0);
358: }