xref: /petsc/src/snes/interface/snesut.c (revision 36851e7ff8b7b138dacde63871c4a6f4b13dee77)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*36851e7fSLois Curfman McInnes static char vcid[] = "$Id: snesut.c,v 1.43 1999/01/13 21:46:21 bsmith Exp curfman $";
3e7e93795SLois Curfman McInnes #endif
4e7e93795SLois Curfman McInnes 
570f55243SBarry Smith #include "src/snes/snesimpl.h"       /*I   "snes.h"   I*/
6e7e93795SLois Curfman McInnes 
75615d1e5SSatish Balay #undef __FUNC__
83f1db9ecSBarry Smith #define __FUNC__ "SNESVecViewMonitor"
93f1db9ecSBarry Smith /*@C
10*36851e7fSLois Curfman McInnes    SNESVecViewMonitor - Monitors progress of the SNES solvers by calling
11*36851e7fSLois Curfman McInnes    VecView() for the approximate solution at each iteration.
123f1db9ecSBarry Smith 
133f1db9ecSBarry Smith    Collective on SNES
143f1db9ecSBarry Smith 
153f1db9ecSBarry Smith    Input Parameters:
163f1db9ecSBarry Smith +  snes - the SNES context
173f1db9ecSBarry Smith .  its - iteration number
183f1db9ecSBarry Smith .  fgnorm - 2-norm of residual (or gradient)
193f1db9ecSBarry Smith -  dummy - either a viewer or PETSC_NULL
203f1db9ecSBarry Smith 
21*36851e7fSLois Curfman McInnes    Level: intermediate
223f1db9ecSBarry Smith 
23*36851e7fSLois Curfman McInnes .keywords: SNES, nonlinear, vector, monitor, view
243f1db9ecSBarry Smith 
25*36851e7fSLois Curfman McInnes .seealso: SNESSetMonitor(), SNESDefaultMonitor(), VecView()
263f1db9ecSBarry Smith @*/
273f1db9ecSBarry Smith int SNESVecViewMonitor(SNES snes,int its,double fgnorm,void *dummy)
283f1db9ecSBarry Smith {
293f1db9ecSBarry Smith   int    ierr;
303f1db9ecSBarry Smith   Vec    x;
313f1db9ecSBarry Smith   Viewer viewer = (Viewer) dummy;
323f1db9ecSBarry Smith 
333f1db9ecSBarry Smith   PetscFunctionBegin;
343f1db9ecSBarry Smith   ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr);
353f1db9ecSBarry Smith   if (!viewer) {
363f1db9ecSBarry Smith     MPI_Comm comm;
373f1db9ecSBarry Smith     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
38c655490fSBarry Smith     viewer = VIEWER_DRAW_(comm);
393f1db9ecSBarry Smith   }
403f1db9ecSBarry Smith   ierr = VecView(x,viewer);CHKERRQ(ierr);
413f1db9ecSBarry Smith 
423f1db9ecSBarry Smith   PetscFunctionReturn(0);
433f1db9ecSBarry Smith }
443f1db9ecSBarry Smith 
453f1db9ecSBarry Smith #undef __FUNC__
46d4bb536fSBarry Smith #define __FUNC__ "SNESDefaultMonitor"
474b828684SBarry Smith /*@C
48f525115eSLois Curfman McInnes    SNESDefaultMonitor - Monitoring progress of the SNES solvers (default).
49e7e93795SLois Curfman McInnes 
50c7afd0dbSLois Curfman McInnes    Collective on SNES
51c7afd0dbSLois Curfman McInnes 
52e7e93795SLois Curfman McInnes    Input Parameters:
53c7afd0dbSLois Curfman McInnes +  snes - the SNES context
54e7e93795SLois Curfman McInnes .  its - iteration number
55e7e93795SLois Curfman McInnes .  fgnorm - 2-norm of residual (or gradient)
56c7afd0dbSLois Curfman McInnes -  dummy - unused context
57fee21e36SBarry Smith 
58e7e93795SLois Curfman McInnes    Notes:
59e7e93795SLois Curfman McInnes    For SNES_NONLINEAR_EQUATIONS methods the routine prints the
60e7e93795SLois Curfman McInnes    residual norm at each iteration.
61e7e93795SLois Curfman McInnes 
62e7e93795SLois Curfman McInnes    For SNES_UNCONSTRAINED_MINIMIZATION methods the routine prints the
63e7e93795SLois Curfman McInnes    function value and gradient norm at each iteration.
64e7e93795SLois Curfman McInnes 
65*36851e7fSLois Curfman McInnes    Level: intermediate
66*36851e7fSLois Curfman McInnes 
67e7e93795SLois Curfman McInnes .keywords: SNES, nonlinear, default, monitor, norm
68e7e93795SLois Curfman McInnes 
69*36851e7fSLois Curfman McInnes .seealso: SNESSetMonitor(), SNESVecViewMonitor()
70e7e93795SLois Curfman McInnes @*/
71e7e93795SLois Curfman McInnes int SNESDefaultMonitor(SNES snes,int its,double fgnorm,void *dummy)
72e7e93795SLois Curfman McInnes {
733a40ed3dSBarry Smith   PetscFunctionBegin;
7476be9ce4SBarry Smith   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
7577c4ece6SBarry Smith     PetscPrintf(snes->comm, "iter = %d, SNES Function norm %g \n",its,fgnorm);
7676be9ce4SBarry Smith   } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
773a40ed3dSBarry Smith     PetscPrintf(snes->comm,"iter = %d, SNES Function value %g, Gradient norm %g \n",its,snes->fc,fgnorm);
7876be9ce4SBarry Smith   } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
793a40ed3dSBarry Smith   PetscFunctionReturn(0);
80e7e93795SLois Curfman McInnes }
813f1db9ecSBarry Smith 
82e7e93795SLois Curfman McInnes /* ---------------------------------------------------------------- */
835615d1e5SSatish Balay #undef __FUNC__
84d4bb536fSBarry Smith #define __FUNC__ "SNESDefaultSMonitor"
85be1f7002SBarry Smith /*
86be1f7002SBarry Smith      Default (short) SNES Monitor, same as SNESDefaultMonitor() except
87be1f7002SBarry Smith   it prints fewer digits of the residual as the residual gets smaller.
88be1f7002SBarry Smith   This is because the later digits are meaningless and are often
89be1f7002SBarry Smith   different on different machines; by using this routine different
90be1f7002SBarry Smith   machines will usually generate the same output.
91be1f7002SBarry Smith */
92e7e93795SLois Curfman McInnes int SNESDefaultSMonitor(SNES snes,int its, double fgnorm,void *dummy)
93e7e93795SLois Curfman McInnes {
943a40ed3dSBarry Smith   PetscFunctionBegin;
95e7e93795SLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
968f240d10SBarry Smith     if (fgnorm > 1.e-9) {
97c7ab52efSLois Curfman McInnes       PetscPrintf(snes->comm, "iter = %d, SNES Function norm %g \n",its,fgnorm);
983a40ed3dSBarry Smith     } else if (fgnorm > 1.e-11){
99c7ab52efSLois Curfman McInnes       PetscPrintf(snes->comm, "iter = %d, SNES Function norm %5.3e \n",its,fgnorm);
1003a40ed3dSBarry Smith     } else {
101c7ab52efSLois Curfman McInnes       PetscPrintf(snes->comm, "iter = %d, SNES Function norm < 1.e-11\n",its);
102e7e93795SLois Curfman McInnes     }
103e7e93795SLois Curfman McInnes   } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1048f240d10SBarry Smith     if (fgnorm > 1.e-9) {
10577c4ece6SBarry Smith       PetscPrintf(snes->comm,
1063a40ed3dSBarry Smith        "iter = %d, SNES Function value %g, Gradient norm %g \n",its,snes->fc,fgnorm);
1073a40ed3dSBarry Smith     } else if (fgnorm > 1.e-11) {
10877c4ece6SBarry Smith       PetscPrintf(snes->comm,
1093a40ed3dSBarry Smith         "iter = %d, SNES Function value %g, Gradient norm %5.3e \n",its,snes->fc,fgnorm);
1103a40ed3dSBarry Smith     } else {
11177c4ece6SBarry Smith       PetscPrintf(snes->comm,
1123a40ed3dSBarry Smith         "iter = %d, SNES Function value %g, Gradient norm < 1.e-11\n",its,snes->fc);
113e7e93795SLois Curfman McInnes     }
114a8c6a408SBarry Smith   } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
1153a40ed3dSBarry Smith   PetscFunctionReturn(0);
116e7e93795SLois Curfman McInnes }
117e7e93795SLois Curfman McInnes /* ---------------------------------------------------------------- */
1185615d1e5SSatish Balay #undef __FUNC__
1195615d1e5SSatish Balay #define __FUNC__ "SNESConverged_EQ_LS"
1204b828684SBarry Smith /*@C
121f525115eSLois Curfman McInnes    SNESConverged_EQ_LS - Monitors the convergence of the solvers for
122f525115eSLois Curfman McInnes    systems of nonlinear equations (default).
123e7e93795SLois Curfman McInnes 
124c7afd0dbSLois Curfman McInnes    Collective on SNES
125c7afd0dbSLois Curfman McInnes 
126e7e93795SLois Curfman McInnes    Input Parameters:
127c7afd0dbSLois Curfman McInnes +  snes - the SNES context
128e7e93795SLois Curfman McInnes .  xnorm - 2-norm of current iterate
129e7e93795SLois Curfman McInnes .  pnorm - 2-norm of current step
130e7e93795SLois Curfman McInnes .  fnorm - 2-norm of function
131c7afd0dbSLois Curfman McInnes -  dummy - unused context
132e7e93795SLois Curfman McInnes 
133e7e93795SLois Curfman McInnes    Returns:
134c7afd0dbSLois Curfman McInnes +  2  - if  ( fnorm < atol ),
135c7afd0dbSLois Curfman McInnes .  3  - if  ( pnorm < xtol*xnorm ),
136c7afd0dbSLois Curfman McInnes .  4  - if  ( fnorm < rtol*fnorm0 ),
137c7afd0dbSLois Curfman McInnes . -2  - if  ( nfct > maxf ),
138c7afd0dbSLois Curfman McInnes -  0  - otherwise,
139e7e93795SLois Curfman McInnes 
140e7e93795SLois Curfman McInnes    where
141c7afd0dbSLois Curfman McInnes +    maxf - maximum number of function evaluations,
142c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
143c7afd0dbSLois Curfman McInnes .    nfct - number of function evaluations,
144c7afd0dbSLois Curfman McInnes .    atol - absolute function norm tolerance,
145c7afd0dbSLois Curfman McInnes             set with SNESSetTolerances()
146c7afd0dbSLois Curfman McInnes -    rtol - relative function norm tolerance, set with SNESSetTolerances()
147fee21e36SBarry Smith 
148*36851e7fSLois Curfman McInnes    Level: intermediate
149*36851e7fSLois Curfman McInnes 
150e7e93795SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence
151e7e93795SLois Curfman McInnes 
152e7e93795SLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
153e7e93795SLois Curfman McInnes @*/
15440191667SLois Curfman McInnes int SNESConverged_EQ_LS(SNES snes,double xnorm,double pnorm,double fnorm,void *dummy)
155e7e93795SLois Curfman McInnes {
1563a40ed3dSBarry Smith   PetscFunctionBegin;
157d252947aSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
158a8c6a408SBarry Smith      SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
159d252947aSBarry Smith   }
160082acdaeSLois Curfman McInnes   /* Note:  Reserve return code 1, -1 for compatibility with SNESConverged_EQ_TR */
161d252947aSBarry Smith   if (fnorm != fnorm) {
162981c4779SBarry Smith     PLogInfo(snes,"SNESConverged_EQ_LS:Failed to converged, function norm is NaN\n");
1633a40ed3dSBarry Smith     PetscFunctionReturn(-3);
164d252947aSBarry Smith   }
1655d2e0e51SBarry Smith   if (fnorm <= snes->ttol) {
16694a424c1SBarry Smith     PLogInfo(snes,
167981c4779SBarry Smith     "SNESConverged_EQ_LS:Converged due to function norm %g < %g (relative tolerance)\n",fnorm,snes->ttol);
1683a40ed3dSBarry Smith     PetscFunctionReturn(4);
1695d2e0e51SBarry Smith   }
1705d2e0e51SBarry Smith 
171e7e93795SLois Curfman McInnes   if (fnorm < snes->atol) {
17294a424c1SBarry Smith     PLogInfo(snes,
173981c4779SBarry Smith       "SNESConverged_EQ_LS: Converged due to function norm %g < %g\n",fnorm,snes->atol);
1743a40ed3dSBarry Smith     PetscFunctionReturn(2);
175e7e93795SLois Curfman McInnes   }
176e7e93795SLois Curfman McInnes   if (pnorm < snes->xtol*(xnorm)) {
17794a424c1SBarry Smith     PLogInfo(snes,
178981c4779SBarry Smith       "SNESConverged_EQ_LS: Converged due to small update length: %g < %g * %g\n",
179e7e93795SLois Curfman McInnes        pnorm,snes->xtol,xnorm);
1803a40ed3dSBarry Smith     PetscFunctionReturn(3);
181e7e93795SLois Curfman McInnes   }
182e7e93795SLois Curfman McInnes   if (snes->nfuncs > snes->max_funcs) {
183981c4779SBarry Smith     PLogInfo(snes,"SNESConverged_EQ_LS: Exceeded maximum number of function evaluations: %d > %d\n",
184e7e93795SLois Curfman McInnes       snes->nfuncs, snes->max_funcs );
1853a40ed3dSBarry Smith     PetscFunctionReturn(-2);
186e7e93795SLois Curfman McInnes   }
1873a40ed3dSBarry Smith   PetscFunctionReturn(0);
188e7e93795SLois Curfman McInnes }
189e7e93795SLois Curfman McInnes /* ------------------------------------------------------------ */
1905615d1e5SSatish Balay #undef __FUNC__
1915615d1e5SSatish Balay #define __FUNC__ "SNES_KSP_SetConvergenceTestEW"
192e7e93795SLois Curfman McInnes /*@
193f525115eSLois Curfman McInnes    SNES_KSP_SetConvergenceTestEW - Sets alternative convergence test
194e7e93795SLois Curfman McInnes    for the linear solvers within an inexact Newton method.
195e7e93795SLois Curfman McInnes 
196c7afd0dbSLois Curfman McInnes    Collective on SNES
197c7afd0dbSLois Curfman McInnes 
198e7e93795SLois Curfman McInnes    Input Parameter:
199e7e93795SLois Curfman McInnes .  snes - SNES context
200e7e93795SLois Curfman McInnes 
201e7e93795SLois Curfman McInnes    Notes:
202e7e93795SLois Curfman McInnes    Currently, the default is to use a constant relative tolerance for
203e7e93795SLois Curfman McInnes    the inner linear solvers.  Alternatively, one can use the
204e7e93795SLois Curfman McInnes    Eisenstat-Walker method, where the relative convergence tolerance
205e7e93795SLois Curfman McInnes    is reset at each Newton iteration according progress of the nonlinear
206e7e93795SLois Curfman McInnes    solver.
207e7e93795SLois Curfman McInnes 
208*36851e7fSLois Curfman McInnes    Level: advanced
209*36851e7fSLois Curfman McInnes 
210e7e93795SLois Curfman McInnes    Reference:
211e7e93795SLois Curfman McInnes    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
212e30ad881SLois Curfman McInnes    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
213e7e93795SLois Curfman McInnes 
214e7e93795SLois Curfman McInnes .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
215e7e93795SLois Curfman McInnes @*/
216e7e93795SLois Curfman McInnes int SNES_KSP_SetConvergenceTestEW(SNES snes)
217e7e93795SLois Curfman McInnes {
2183a40ed3dSBarry Smith   PetscFunctionBegin;
219e7e93795SLois Curfman McInnes   snes->ksp_ewconv = 1;
2203a40ed3dSBarry Smith   PetscFunctionReturn(0);
221e7e93795SLois Curfman McInnes }
222e7e93795SLois Curfman McInnes 
2235615d1e5SSatish Balay #undef __FUNC__
2245615d1e5SSatish Balay #define __FUNC__ "SNES_KSP_SetParametersEW"
225e7e93795SLois Curfman McInnes /*@
226e7e93795SLois Curfman McInnes    SNES_KSP_SetParametersEW - Sets parameters for Eisenstat-Walker
227e7e93795SLois Curfman McInnes    convergence criteria for the linear solvers within an inexact
228e7e93795SLois Curfman McInnes    Newton method.
229e7e93795SLois Curfman McInnes 
230c7afd0dbSLois Curfman McInnes    Collective on SNES
231c7afd0dbSLois Curfman McInnes 
232e7e93795SLois Curfman McInnes    Input Parameters:
233c7afd0dbSLois Curfman McInnes +    snes - SNES context
234e7e93795SLois Curfman McInnes .    version - version 1 or 2 (default is 2)
235c7afd0dbSLois Curfman McInnes .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
236c7afd0dbSLois Curfman McInnes .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
237c7afd0dbSLois Curfman McInnes .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
238e7e93795SLois Curfman McInnes .    alpha2 - power for safeguard
239e7e93795SLois Curfman McInnes .    gamma2 - multiplicative factor for version 2 rtol computation
240c7afd0dbSLois Curfman McInnes               (0 <= gamma2 <= 1)
241c7afd0dbSLois Curfman McInnes -    threshold - threshold for imposing safeguard (0 < threshold < 1)
242fee21e36SBarry Smith 
243e7e93795SLois Curfman McInnes    Note:
244e7e93795SLois Curfman McInnes    Use PETSC_DEFAULT to retain the default for any of the parameters.
245e7e93795SLois Curfman McInnes 
246*36851e7fSLois Curfman McInnes    Level: advanced
247*36851e7fSLois Curfman McInnes 
248e7e93795SLois Curfman McInnes    Reference:
249e7e93795SLois Curfman McInnes    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
250e7e93795SLois Curfman McInnes    inexact Newton method", Utah State University Math. Stat. Dept. Res.
251e7e93795SLois Curfman McInnes    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
252e7e93795SLois Curfman McInnes 
253e7e93795SLois Curfman McInnes .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
254e7e93795SLois Curfman McInnes 
255e7e93795SLois Curfman McInnes .seealso: SNES_KSP_SetConvergenceTestEW()
256e7e93795SLois Curfman McInnes @*/
257e7e93795SLois Curfman McInnes int SNES_KSP_SetParametersEW(SNES snes,int version,double rtol_0,
258e7e93795SLois Curfman McInnes                              double rtol_max,double gamma2,double alpha,
259e7e93795SLois Curfman McInnes                              double alpha2,double threshold)
260e7e93795SLois Curfman McInnes {
261e7e93795SLois Curfman McInnes   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
2623a40ed3dSBarry Smith 
2633a40ed3dSBarry Smith   PetscFunctionBegin;
264a8c6a408SBarry Smith   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"No Eisenstat-Walker context existing");
265e7e93795SLois Curfman McInnes   if (version != PETSC_DEFAULT)   kctx->version   = version;
266e7e93795SLois Curfman McInnes   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
267e7e93795SLois Curfman McInnes   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
268e7e93795SLois Curfman McInnes   if (gamma2 != PETSC_DEFAULT)    kctx->gamma     = gamma2;
269e7e93795SLois Curfman McInnes   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
270e7e93795SLois Curfman McInnes   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
271e7e93795SLois Curfman McInnes   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
272a8c6a408SBarry Smith   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
273596552b5SBarry Smith     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"0.0 <= rtol_0 < 1.0: %g",kctx->rtol_0);
274a8c6a408SBarry Smith   }
275a8c6a408SBarry Smith   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
276596552b5SBarry Smith     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"0.0 <= rtol_max < 1.0\n",kctx->rtol_max);
277a8c6a408SBarry Smith   }
278a8c6a408SBarry Smith   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
279596552b5SBarry Smith     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"0.0 < threshold < 1.0\n",kctx->threshold);
280a8c6a408SBarry Smith   }
281a8c6a408SBarry Smith   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
282596552b5SBarry Smith     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"0.0 <= alpha <= 1.0\n",kctx->gamma);
283a8c6a408SBarry Smith   }
284a8c6a408SBarry Smith   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
285596552b5SBarry Smith     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"1.0 < alpha <= 2.0\n",kctx->alpha);
286a8c6a408SBarry Smith   }
287a8c6a408SBarry Smith   if (kctx->version != 1 && kctx->version !=2) {
288596552b5SBarry Smith     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"Only versions 1 and 2 are supported: %d",kctx->version);
289a8c6a408SBarry Smith   }
2903a40ed3dSBarry Smith   PetscFunctionReturn(0);
291e7e93795SLois Curfman McInnes }
292e7e93795SLois Curfman McInnes 
2935615d1e5SSatish Balay #undef __FUNC__
2945615d1e5SSatish Balay #define __FUNC__ "SNES_KSP_EW_ComputeRelativeTolerance_Private"
295e7e93795SLois Curfman McInnes int SNES_KSP_EW_ComputeRelativeTolerance_Private(SNES snes,KSP ksp)
296e7e93795SLois Curfman McInnes {
297e7e93795SLois Curfman McInnes   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
2983f1db9ecSBarry Smith   double              rtol = 0.0, stol;
299e7e93795SLois Curfman McInnes   int                 ierr;
3003a40ed3dSBarry Smith 
3013a40ed3dSBarry Smith   PetscFunctionBegin;
302a8c6a408SBarry Smith   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"No Eisenstat-Walker context exists");
303e7e93795SLois Curfman McInnes   if (snes->iter == 1) {
304e7e93795SLois Curfman McInnes     rtol = kctx->rtol_0;
305e7e93795SLois Curfman McInnes   } else {
306e7e93795SLois Curfman McInnes     if (kctx->version == 1) {
307e7e93795SLois Curfman McInnes       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
308e7e93795SLois Curfman McInnes       if (rtol < 0.0) rtol = -rtol;
309e7e93795SLois Curfman McInnes       stol = pow(kctx->rtol_last,kctx->alpha2);
3100452661fSBarry Smith       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
311e7e93795SLois Curfman McInnes     } else if (kctx->version == 2) {
312e7e93795SLois Curfman McInnes       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
313e7e93795SLois Curfman McInnes       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
3140452661fSBarry Smith       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
315596552b5SBarry Smith     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"Only versions 1 or 2 are supported: %d",kctx->version);
316e7e93795SLois Curfman McInnes   }
3170452661fSBarry Smith   rtol = PetscMin(rtol,kctx->rtol_max);
318e7e93795SLois Curfman McInnes   kctx->rtol_last = rtol;
319596552b5SBarry Smith   PLogInfo(snes,"SNESConverged_EQ_LS: iter %d, Eisenstat-Walker (version %d) KSP rtol = %g\n",snes->iter,kctx->version,rtol);
3203131a8b6SLois Curfman McInnes   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); CHKERRQ(ierr);
321e7e93795SLois Curfman McInnes   kctx->norm_last = snes->norm;
3223a40ed3dSBarry Smith   PetscFunctionReturn(0);
323e7e93795SLois Curfman McInnes }
324e7e93795SLois Curfman McInnes 
3255615d1e5SSatish Balay #undef __FUNC__
3265615d1e5SSatish Balay #define __FUNC__ "SNES_KSP_EW_Converged_Private"
327e7e93795SLois Curfman McInnes int SNES_KSP_EW_Converged_Private(KSP ksp,int n,double rnorm,void *ctx)
328e7e93795SLois Curfman McInnes {
329e7e93795SLois Curfman McInnes   SNES                snes = (SNES)ctx;
330e7e93795SLois Curfman McInnes   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
331e7e93795SLois Curfman McInnes   int                 convinfo;
332e7e93795SLois Curfman McInnes 
3333a40ed3dSBarry Smith   PetscFunctionBegin;
334a8c6a408SBarry Smith   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"No Eisenstat-Walker context set");
335e7e93795SLois Curfman McInnes   if (n == 0) SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);
336e7e93795SLois Curfman McInnes   convinfo = KSPDefaultConverged(ksp,n,rnorm,ctx);
337e7e93795SLois Curfman McInnes   kctx->lresid_last = rnorm;
3383a40ed3dSBarry Smith   if (convinfo) {
339981c4779SBarry Smith     PLogInfo(snes,"SNES_KSP_EW_Converged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm);
3403a40ed3dSBarry Smith   }
3413a40ed3dSBarry Smith   PetscFunctionReturn(convinfo);
342e7e93795SLois Curfman McInnes }
343e7e93795SLois Curfman McInnes 
344e7e93795SLois Curfman McInnes 
345