xref: /petsc/src/snes/interface/snesut.c (revision 5ed2d596fd9914de76abbaed56c65f4792ae57ae)
1 /*$Id: snesut.c,v 1.66 2001/08/06 21:17:07 bsmith Exp $*/
2 
3 #include "src/snes/snesimpl.h"       /*I   "petscsnes.h"   I*/
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "SNESVecViewMonitor"
7 /*@C
8    SNESVecViewMonitor - Monitors progress of the SNES solvers by calling
9    VecView() for the approximate solution at each iteration.
10 
11    Collective on SNES
12 
13    Input Parameters:
14 +  snes - the SNES context
15 .  its - iteration number
16 .  fgnorm - 2-norm of residual (or gradient)
17 -  dummy - either a viewer or PETSC_NULL
18 
19    Level: intermediate
20 
21 .keywords: SNES, nonlinear, vector, monitor, view
22 
23 .seealso: SNESSetMonitor(), SNESDefaultMonitor(), VecView()
24 @*/
25 int SNESVecViewMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
26 {
27   int         ierr;
28   Vec         x;
29   PetscViewer viewer = (PetscViewer) dummy;
30 
31   PetscFunctionBegin;
32   ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr);
33   if (!viewer) {
34     MPI_Comm comm;
35     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
36     viewer = PETSC_VIEWER_DRAW_(comm);
37   }
38   ierr = VecView(x,viewer);CHKERRQ(ierr);
39 
40   PetscFunctionReturn(0);
41 }
42 
43 #undef __FUNCT__
44 #define __FUNCT__ "SNESVecViewResidualMonitor"
45 /*@C
46    SNESVecViewResidualMonitor - Monitors progress of the SNES solvers by calling
47    VecView() for the residual at each iteration.
48 
49    Collective on SNES
50 
51    Input Parameters:
52 +  snes - the SNES context
53 .  its - iteration number
54 .  fgnorm - 2-norm of residual (or gradient)
55 -  dummy - either a viewer or PETSC_NULL
56 
57    Level: intermediate
58 
59 .keywords: SNES, nonlinear, vector, monitor, view
60 
61 .seealso: SNESSetMonitor(), SNESDefaultMonitor(), VecView()
62 @*/
63 int SNESVecViewResidualMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
64 {
65   int         ierr;
66   Vec         x;
67   PetscViewer viewer = (PetscViewer) dummy;
68 
69   PetscFunctionBegin;
70   ierr = SNESGetFunction(snes,&x,0,0);CHKERRQ(ierr);
71   if (!viewer) {
72     MPI_Comm comm;
73     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
74     viewer = PETSC_VIEWER_DRAW_(comm);
75   }
76   ierr = VecView(x,viewer);CHKERRQ(ierr);
77 
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "SNESVecViewUpdateMonitor"
83 /*@C
84    SNESVecViewUpdateMonitor - Monitors progress of the SNES solvers by calling
85    VecView() for the UPDATE to the solution at each iteration.
86 
87    Collective on SNES
88 
89    Input Parameters:
90 +  snes - the SNES context
91 .  its - iteration number
92 .  fgnorm - 2-norm of residual (or gradient)
93 -  dummy - either a viewer or PETSC_NULL
94 
95    Level: intermediate
96 
97 .keywords: SNES, nonlinear, vector, monitor, view
98 
99 .seealso: SNESSetMonitor(), SNESDefaultMonitor(), VecView()
100 @*/
101 int SNESVecViewUpdateMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
102 {
103   int         ierr;
104   Vec         x;
105   PetscViewer viewer = (PetscViewer) dummy;
106 
107   PetscFunctionBegin;
108   ierr = SNESGetSolutionUpdate(snes,&x);CHKERRQ(ierr);
109   if (!viewer) {
110     MPI_Comm comm;
111     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
112     viewer = PETSC_VIEWER_DRAW_(comm);
113   }
114   ierr = VecView(x,viewer);CHKERRQ(ierr);
115 
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "SNESDefaultMonitor"
121 /*@C
122    SNESDefaultMonitor - Monitoring progress of the SNES solvers (default).
123 
124    Collective on SNES
125 
126    Input Parameters:
127 +  snes - the SNES context
128 .  its - iteration number
129 .  fgnorm - 2-norm of residual (or gradient)
130 -  dummy - unused context
131 
132    Notes:
133    For SNES_NONLINEAR_EQUATIONS methods the routine prints the
134    residual norm at each iteration.
135 
136    For SNES_UNCONSTRAINED_MINIMIZATION methods the routine prints the
137    function value and gradient norm at each iteration.
138 
139    Level: intermediate
140 
141 .keywords: SNES, nonlinear, default, monitor, norm
142 
143 .seealso: SNESSetMonitor(), SNESVecViewMonitor()
144 @*/
145 int SNESDefaultMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
146 {
147   int         ierr;
148   PetscViewer viewer = (PetscViewer) dummy;
149 
150   PetscFunctionBegin;
151   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm);
152 
153   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
154     ierr = PetscViewerASCIIPrintf(viewer,"%3d SNES Function norm %14.12e \n",its,fgnorm);CHKERRQ(ierr);
155   } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
156     ierr = PetscViewerASCIIPrintf(viewer,"%3d SNES Function value %14.12e, Gradient norm %14.12e \n",its,snes->fc,fgnorm);CHKERRQ(ierr);
157   } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class");
158   PetscFunctionReturn(0);
159 }
160 
161 #undef __FUNCT__
162 #define __FUNCT__ "SNESRatioMonitor"
163 /*@C
164    SNESRatioMonitor - Monitoring progress of the SNES solvers, prints ratio
165       of residual norm at each iteration to previous
166 
167    Collective on SNES
168 
169    Input Parameters:
170 +  snes - the SNES context
171 .  its - iteration number
172 .  fgnorm - 2-norm of residual (or gradient)
173 -  dummy - unused context
174 
175    Level: intermediate
176 
177 .keywords: SNES, nonlinear, monitor, norm
178 
179 .seealso: SNESSetMonitor(), SNESVecViewMonitor()
180 @*/
181 int SNESRatioMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
182 {
183   int         ierr,len;
184   PetscReal   *history;
185   PetscViewer viewer;
186 
187   PetscFunctionBegin;
188   viewer = PETSC_VIEWER_STDOUT_(snes->comm);
189 
190   ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,&len);CHKERRQ(ierr);
191   if (its == 0 || !history || its > len) {
192     ierr = PetscViewerASCIIPrintf(viewer,"%3d SNES Function norm %14.12e \n",its,fgnorm);CHKERRQ(ierr);
193   } else {
194     PetscReal ratio = fgnorm/history[its-1];
195     ierr = PetscViewerASCIIPrintf(viewer,"%3d SNES Function norm %14.12e %g \n",its,fgnorm,ratio);CHKERRQ(ierr);
196   }
197   PetscFunctionReturn(0);
198 }
199 
200 /*
201    If the we set the history monitor space then we must destroy it
202 */
203 #undef __FUNCT__
204 #define __FUNCT__ "SNESRatioMonitorDestroy"
205 int SNESRatioMonitorDestroy(void *history)
206 {
207   int         ierr;
208 
209   PetscFunctionBegin;
210   ierr = PetscFree(history);CHKERRQ(ierr);
211   PetscFunctionReturn(0);
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "SNESSetRatioMonitor"
216 /*@C
217    SNESSetRatioMonitor - Sets SNES to use a monitor that prints the
218      ratio of the function norm at each iteration
219 
220    Collective on SNES
221 
222    Input Parameters:
223 .   snes - the SNES context
224 
225    Level: intermediate
226 
227 .keywords: SNES, nonlinear, monitor, norm
228 
229 .seealso: SNESSetMonitor(), SNESVecViewMonitor(), SNESDefaultMonitor()
230 @*/
231 int SNESSetRatioMonitor(SNES snes)
232 {
233   int       ierr;
234   PetscReal *history;
235 
236   PetscFunctionBegin;
237 
238   ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
239   if (!history) {
240     ierr = PetscMalloc(100*sizeof(double),&history);CHKERRQ(ierr);
241     ierr = SNESSetConvergenceHistory(snes,history,0,100,PETSC_TRUE);CHKERRQ(ierr);
242     ierr = SNESSetMonitor(snes,SNESRatioMonitor,history,SNESRatioMonitorDestroy);CHKERRQ(ierr);
243   } else {
244     ierr = SNESSetMonitor(snes,SNESRatioMonitor,0,0);CHKERRQ(ierr);
245   }
246   PetscFunctionReturn(0);
247 }
248 
249 /* ---------------------------------------------------------------- */
250 #undef __FUNCT__
251 #define __FUNCT__ "SNESDefaultSMonitor"
252 /*
253      Default (short) SNES Monitor, same as SNESDefaultMonitor() except
254   it prints fewer digits of the residual as the residual gets smaller.
255   This is because the later digits are meaningless and are often
256   different on different machines; by using this routine different
257   machines will usually generate the same output.
258 */
259 int SNESDefaultSMonitor(SNES snes,int its,PetscReal fgnorm,void *dummy)
260 {
261   int ierr;
262 
263   PetscFunctionBegin;
264   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
265     if (fgnorm > 1.e-9) {
266       ierr = PetscPrintf(snes->comm,"%3d SNES Function norm %g \n",its,fgnorm);CHKERRQ(ierr);
267     } else if (fgnorm > 1.e-11){
268       ierr = PetscPrintf(snes->comm,"%3d SNES Function norm %5.3e \n",its,fgnorm);CHKERRQ(ierr);
269     } else {
270       ierr = PetscPrintf(snes->comm,"%3d SNES Function norm < 1.e-11\n",its);CHKERRQ(ierr);
271     }
272   } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
273     if (fgnorm > 1.e-9) {
274       ierr = PetscPrintf(snes->comm,"%3d SNES Function value %g, Gradient norm %g \n",its,snes->fc,fgnorm);CHKERRQ(ierr);
275     } else if (fgnorm > 1.e-11) {
276       ierr = PetscPrintf(snes->comm,"%3d SNES Function value %g, Gradient norm %5.3e \n",its,snes->fc,fgnorm);CHKERRQ(ierr);
277     } else {
278       ierr = PetscPrintf(snes->comm,"%3d SNES Function value %g, Gradient norm < 1.e-11\n",its,snes->fc);CHKERRQ(ierr);
279     }
280   } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class");
281   PetscFunctionReturn(0);
282 }
283 /* ---------------------------------------------------------------- */
284 #undef __FUNCT__
285 #define __FUNCT__ "SNESConverged_EQ_LS"
286 /*@C
287    SNESConverged_EQ_LS - Monitors the convergence of the solvers for
288    systems of nonlinear equations (default).
289 
290    Collective on SNES
291 
292    Input Parameters:
293 +  snes - the SNES context
294 .  xnorm - 2-norm of current iterate
295 .  pnorm - 2-norm of current step
296 .  fnorm - 2-norm of function
297 -  dummy - unused context
298 
299    Output Parameter:
300 .   reason  - one of
301 $  SNES_CONVERGED_FNORM_ABS       - (fnorm < atol),
302 $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
303 $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
304 $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
305 $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
306 $  SNES_CONVERGED_ITERATING       - (otherwise),
307 
308    where
309 +    maxf - maximum number of function evaluations,
310             set with SNESSetTolerances()
311 .    nfct - number of function evaluations,
312 .    atol - absolute function norm tolerance,
313             set with SNESSetTolerances()
314 -    rtol - relative function norm tolerance, set with SNESSetTolerances()
315 
316    Level: intermediate
317 
318 .keywords: SNES, nonlinear, default, converged, convergence
319 
320 .seealso: SNESSetConvergenceTest(), SNESEisenstatWalkerConverged()
321 @*/
322 int SNESConverged_EQ_LS(SNES snes,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
323 {
324   PetscFunctionBegin;
325   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
326      SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
327   }
328 
329   if (fnorm != fnorm) {
330     PetscLogInfo(snes,"SNESConverged_EQ_LS:Failed to converged, function norm is NaN\n");
331     *reason = SNES_DIVERGED_FNORM_NAN;
332   } else if (fnorm <= snes->ttol) {
333     PetscLogInfo(snes,"SNESConverged_EQ_LS:Converged due to function norm %g < %g (relative tolerance)\n",fnorm,snes->ttol);
334     *reason = SNES_CONVERGED_FNORM_RELATIVE;
335   } else if (fnorm < snes->atol) {
336     PetscLogInfo(snes,"SNESConverged_EQ_LS:Converged due to function norm %g < %g\n",fnorm,snes->atol);
337     *reason = SNES_CONVERGED_FNORM_ABS;
338   } else if (pnorm < snes->xtol*xnorm) {
339     PetscLogInfo(snes,"SNESConverged_EQ_LS:Converged due to small update length: %g < %g * %g\n",pnorm,snes->xtol,xnorm);
340     *reason = SNES_CONVERGED_PNORM_RELATIVE;
341   } else if (snes->nfuncs > snes->max_funcs) {
342     PetscLogInfo(snes,"SNESConverged_EQ_LS:Exceeded maximum number of function evaluations: %d > %d\n",snes->nfuncs,snes->max_funcs);
343     *reason = SNES_DIVERGED_FUNCTION_COUNT ;
344   } else {
345     *reason = SNES_CONVERGED_ITERATING;
346   }
347   PetscFunctionReturn(0);
348 }
349 /* ------------------------------------------------------------ */
350 #undef __FUNCT__
351 #define __FUNCT__ "SNES_KSP_SetConvergenceTestEW"
352 /*@
353    SNES_KSP_SetConvergenceTestEW - Sets alternative convergence test
354    for the linear solvers within an inexact Newton method.
355 
356    Collective on SNES
357 
358    Input Parameter:
359 .  snes - SNES context
360 
361    Notes:
362    Currently, the default is to use a constant relative tolerance for
363    the inner linear solvers.  Alternatively, one can use the
364    Eisenstat-Walker method, where the relative convergence tolerance
365    is reset at each Newton iteration according progress of the nonlinear
366    solver.
367 
368    Level: advanced
369 
370    Reference:
371    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
372    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
373 
374 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
375 @*/
376 int SNES_KSP_SetConvergenceTestEW(SNES snes)
377 {
378   PetscFunctionBegin;
379   snes->ksp_ewconv = PETSC_TRUE;
380   PetscFunctionReturn(0);
381 }
382 
383 #undef __FUNCT__
384 #define __FUNCT__ "SNES_KSP_SetParametersEW"
385 /*@
386    SNES_KSP_SetParametersEW - Sets parameters for Eisenstat-Walker
387    convergence criteria for the linear solvers within an inexact
388    Newton method.
389 
390    Collective on SNES
391 
392    Input Parameters:
393 +    snes - SNES context
394 .    version - version 1 or 2 (default is 2)
395 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
396 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
397 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
398 .    alpha2 - power for safeguard
399 .    gamma2 - multiplicative factor for version 2 rtol computation
400               (0 <= gamma2 <= 1)
401 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
402 
403    Note:
404    Use PETSC_DEFAULT to retain the default for any of the parameters.
405 
406    Level: advanced
407 
408    Reference:
409    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
410    inexact Newton method", Utah State University Math. Stat. Dept. Res.
411    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
412 
413 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
414 
415 .seealso: SNES_KSP_SetConvergenceTestEW()
416 @*/
417 int SNES_KSP_SetParametersEW(SNES snes,int version,PetscReal rtol_0,
418                              PetscReal rtol_max,PetscReal gamma2,PetscReal alpha,
419                              PetscReal alpha2,PetscReal threshold)
420 {
421   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
422 
423   PetscFunctionBegin;
424   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
425   if (version != PETSC_DEFAULT)   kctx->version   = version;
426   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
427   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
428   if (gamma2 != PETSC_DEFAULT)    kctx->gamma     = gamma2;
429   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
430   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
431   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
432   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
433     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %g",kctx->rtol_0);
434   }
435   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
436     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%g) < 1.0\n",kctx->rtol_max);
437   }
438   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
439     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%g) < 1.0\n",kctx->threshold);
440   }
441   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
442     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%g) <= 1.0\n",kctx->gamma);
443   }
444   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
445     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%g) <= 2.0\n",kctx->alpha);
446   }
447   if (kctx->version != 1 && kctx->version !=2) {
448     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1 and 2 are supported: %d",kctx->version);
449   }
450   PetscFunctionReturn(0);
451 }
452 
453 #undef __FUNCT__
454 #define __FUNCT__ "SNES_KSP_EW_ComputeRelativeTolerance_Private"
455 int SNES_KSP_EW_ComputeRelativeTolerance_Private(SNES snes,KSP ksp)
456 {
457   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
458   PetscReal           rtol = 0.0,stol;
459   int                 ierr;
460 
461   PetscFunctionBegin;
462   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
463   if (!snes->iter) { /* first time in, so use the original user rtol */
464     rtol = kctx->rtol_0;
465   } else {
466     if (kctx->version == 1) {
467       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
468       if (rtol < 0.0) rtol = -rtol;
469       stol = pow(kctx->rtol_last,kctx->alpha2);
470       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
471     } else if (kctx->version == 2) {
472       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
473       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
474       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
475     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1 or 2 are supported: %d",kctx->version);
476   }
477   rtol = PetscMin(rtol,kctx->rtol_max);
478   kctx->rtol_last = rtol;
479   PetscLogInfo(snes,"SNES_KSP_EW_ComputeRelativeTolerance_Private: iter %d, Eisenstat-Walker (version %d) KSP rtol = %g\n",snes->iter,kctx->version,rtol);
480   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
481   kctx->norm_last = snes->norm;
482   PetscFunctionReturn(0);
483 }
484 
485 #undef __FUNCT__
486 #define __FUNCT__ "SNES_KSP_EW_Converged_Private"
487 int SNES_KSP_EW_Converged_Private(KSP ksp,int n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
488 {
489   SNES                snes = (SNES)ctx;
490   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx*)snes->kspconvctx;
491   int                 ierr;
492 
493   PetscFunctionBegin;
494   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context set");
495   if (n == 0) {ierr = SNES_KSP_EW_ComputeRelativeTolerance_Private(snes,ksp);CHKERRQ(ierr);}
496   ierr = KSPDefaultConverged(ksp,n,rnorm,reason,ctx);CHKERRQ(ierr);
497   kctx->lresid_last = rnorm;
498   if (*reason) {
499     PetscLogInfo(snes,"SNES_KSP_EW_Converged_Private: KSP iterations=%d, rnorm=%g\n",n,rnorm);
500   }
501   PetscFunctionReturn(0);
502 }
503 
504 
505