xref: /petsc/src/snes/interface/snesut.c (revision b59f628eff0cd0541d0457de20841eec6c6bc428)
1 
2 #include <petsc-private/snesimpl.h>       /*I   "petsc-private/snesimpl.h"   I*/
3 #include <petscdm.h>
4 #include <petscblaslapack.h>
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "SNESMonitorSolution"
8 /*@C
9    SNESMonitorSolution - Monitors progress of the SNES solvers by calling
10    VecView() for the approximate solution at each iteration.
11 
12    Collective on SNES
13 
14    Input Parameters:
15 +  snes - the SNES context
16 .  its - iteration number
17 .  fgnorm - 2-norm of residual
18 -  dummy - either a viewer or NULL
19 
20    Level: intermediate
21 
22 .keywords: SNES, nonlinear, vector, monitor, view
23 
24 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
25 @*/
26 PetscErrorCode  SNESMonitorSolution(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
27 {
28   PetscErrorCode ierr;
29   Vec            x;
30   PetscViewer    viewer = (PetscViewer) dummy;
31 
32   PetscFunctionBegin;
33   ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr);
34   if (!viewer) {
35     MPI_Comm comm;
36     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
37     viewer = PETSC_VIEWER_DRAW_(comm);
38   }
39   ierr = VecView(x,viewer);CHKERRQ(ierr);
40   PetscFunctionReturn(0);
41 }
42 
43 #undef __FUNCT__
44 #define __FUNCT__ "SNESMonitorResidual"
45 /*@C
46    SNESMonitorResidual - 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
55 -  dummy - either a viewer or NULL
56 
57    Level: intermediate
58 
59 .keywords: SNES, nonlinear, vector, monitor, view
60 
61 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
62 @*/
63 PetscErrorCode  SNESMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
64 {
65   PetscErrorCode 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   PetscFunctionReturn(0);
78 }
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "SNESMonitorSolutionUpdate"
82 /*@C
83    SNESMonitorSolutionUpdate - Monitors progress of the SNES solvers by calling
84    VecView() for the UPDATE to the solution at each iteration.
85 
86    Collective on SNES
87 
88    Input Parameters:
89 +  snes - the SNES context
90 .  its - iteration number
91 .  fgnorm - 2-norm of residual
92 -  dummy - either a viewer or NULL
93 
94    Level: intermediate
95 
96 .keywords: SNES, nonlinear, vector, monitor, view
97 
98 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
99 @*/
100 PetscErrorCode  SNESMonitorSolutionUpdate(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
101 {
102   PetscErrorCode ierr;
103   Vec            x;
104   PetscViewer    viewer = (PetscViewer) dummy;
105 
106   PetscFunctionBegin;
107   ierr = SNESGetSolutionUpdate(snes,&x);CHKERRQ(ierr);
108   if (!viewer) {
109     MPI_Comm comm;
110     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
111     viewer = PETSC_VIEWER_DRAW_(comm);
112   }
113   ierr = VecView(x,viewer);CHKERRQ(ierr);
114   PetscFunctionReturn(0);
115 }
116 
117 #undef __FUNCT__
118 #define __FUNCT__ "KSPMonitorSNES"
119 /*@C
120    KSPMonitorSNES - Print the residual norm of the nonlinear function at each iteration of the linear iterative solver.
121 
122    Collective on KSP
123 
124    Input Parameters:
125 +  ksp   - iterative context
126 .  n     - iteration number
127 .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
128 -  dummy - unused monitor context
129 
130    Level: intermediate
131 
132 .keywords: KSP, default, monitor, residual
133 
134 .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGResidualNormCreate()
135 @*/
136 PetscErrorCode  KSPMonitorSNES(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
137 {
138   PetscErrorCode ierr;
139   PetscViewer    viewer;
140   SNES           snes = (SNES) dummy;
141   Vec            snes_solution,work1,work2;
142   PetscReal      snorm;
143 
144   PetscFunctionBegin;
145   ierr = SNESGetSolution(snes,&snes_solution);CHKERRQ(ierr);
146   ierr = VecDuplicate(snes_solution,&work1);CHKERRQ(ierr);
147   ierr = VecDuplicate(snes_solution,&work2);CHKERRQ(ierr);
148   ierr = KSPBuildSolution(ksp,work1,NULL);CHKERRQ(ierr);
149   ierr = VecAYPX(work1,-1.0,snes_solution);CHKERRQ(ierr);
150   ierr = SNESComputeFunction(snes,work1,work2);CHKERRQ(ierr);
151   ierr = VecNorm(work2,NORM_2,&snorm);CHKERRQ(ierr);
152   ierr = VecDestroy(&work1);CHKERRQ(ierr);
153   ierr = VecDestroy(&work2);CHKERRQ(ierr);
154 
155   ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp),&viewer);CHKERRQ(ierr);
156   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);CHKERRQ(ierr);
157   if (n == 0 && ((PetscObject)ksp)->prefix) {
158     ierr = PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);CHKERRQ(ierr);
159   }
160   ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Residual norm %5.3e KSP Residual norm %5.3e \n",n,(double)snorm,(double)rnorm);CHKERRQ(ierr);
161   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);CHKERRQ(ierr);
162   PetscFunctionReturn(0);
163 }
164 
165 #include <petscdraw.h>
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "KSPMonitorSNESLGResidualNormCreate"
169 /*@C
170    KSPMonitorSNESLGResidualNormCreate - Creates a line graph context for use with
171    KSP to monitor convergence of preconditioned residual norms.
172 
173    Collective on KSP
174 
175    Input Parameters:
176 +  host - the X display to open, or null for the local machine
177 .  label - the title to put in the title bar
178 .  x, y - the screen coordinates of the upper left coordinate of
179           the window
180 -  m, n - the screen width and height in pixels
181 
182    Output Parameter:
183 .  draw - the drawing context
184 
185    Options Database Key:
186 .  -ksp_monitor_lg_residualnorm - Sets line graph monitor
187 
188    Notes:
189    Use KSPMonitorSNESLGResidualNormDestroy() to destroy this line graph; do not use PetscDrawLGDestroy().
190 
191    Level: intermediate
192 
193 .keywords: KSP, monitor, line graph, residual, create
194 
195 .seealso: KSPMonitorSNESLGResidualNormDestroy(), KSPMonitorSet(), KSPMonitorSNESLGTrueResidualCreate()
196 @*/
197 PetscErrorCode  KSPMonitorSNESLGResidualNormCreate(const char host[],const char label[],int x,int y,int m,int n,PetscObject **objs)
198 {
199   PetscDraw      draw;
200   PetscErrorCode ierr;
201   PetscDrawAxis  axis;
202   PetscDrawLG    drawlg;
203   const char     *names[] = {"Linear residual","Nonlinear residual"};
204 
205   PetscFunctionBegin;
206   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&draw);CHKERRQ(ierr);
207   ierr = PetscDrawSetFromOptions(draw);CHKERRQ(ierr);
208   ierr = PetscDrawLGCreate(draw,2,&drawlg);CHKERRQ(ierr);
209   ierr = PetscDrawLGGetAxis(drawlg,&axis);CHKERRQ(ierr);
210   ierr = PetscDrawAxisSetLabels(axis,"Convergence of Residual Norm","Iteration","Residual Norm");CHKERRQ(ierr);
211   ierr = PetscDrawLGSetLegend(drawlg,names);CHKERRQ(ierr);
212   ierr = PetscLogObjectParent((PetscObject)drawlg,(PetscObject)draw);CHKERRQ(ierr);
213 
214   ierr = PetscMalloc(2,objs);CHKERRQ(ierr);
215   (*objs)[1] = (PetscObject)drawlg;
216   PetscFunctionReturn(0);
217 }
218 
219 #undef __FUNCT__
220 #define __FUNCT__ "KSPMonitorSNESLGResidualNorm"
221 PetscErrorCode  KSPMonitorSNESLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscObject *objs)
222 {
223   PetscDrawLG    lg = (PetscDrawLG) objs[1];
224   PetscErrorCode ierr;
225   PetscReal      y[2];
226   SNES           snes = (SNES) objs[0];
227   Vec            snes_solution,work1,work2;
228 
229   PetscFunctionBegin;
230   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
231   else y[0] = -15.0;
232 
233   ierr = SNESGetSolution(snes,&snes_solution);CHKERRQ(ierr);
234   ierr = VecDuplicate(snes_solution,&work1);CHKERRQ(ierr);
235   ierr = VecDuplicate(snes_solution,&work2);CHKERRQ(ierr);
236   ierr = KSPBuildSolution(ksp,work1,NULL);CHKERRQ(ierr);
237   ierr = VecAYPX(work1,-1.0,snes_solution);CHKERRQ(ierr);
238   ierr = SNESComputeFunction(snes,work1,work2);CHKERRQ(ierr);
239   ierr = VecNorm(work2,NORM_2,y+1);CHKERRQ(ierr);
240   if (y[1] > 0.0) y[1] = PetscLog10Real(y[1]);
241   else y[1] = -15.0;
242   ierr = VecDestroy(&work1);CHKERRQ(ierr);
243   ierr = VecDestroy(&work2);CHKERRQ(ierr);
244 
245   ierr = PetscDrawLGAddPoint(lg,NULL,y);CHKERRQ(ierr);
246   if (n < 20 || !(n % 5)) {
247     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
248   }
249   PetscFunctionReturn(0);
250 }
251 
252 #undef __FUNCT__
253 #define __FUNCT__ "KSPMonitorSNESLGResidualNormDestroy"
254 /*@
255    KSPMonitorSNESLGResidualNormDestroy - Destroys a line graph context that was created
256    with KSPMonitorSNESLGResidualNormCreate().
257 
258    Collective on KSP
259 
260    Input Parameter:
261 .  draw - the drawing context
262 
263    Level: intermediate
264 
265 .keywords: KSP, monitor, line graph, destroy
266 
267 .seealso: KSPMonitorSNESLGResidualNormCreate(), KSPMonitorSNESLGTrueResidualDestroy(), KSPMonitorSet()
268 @*/
269 PetscErrorCode  KSPMonitorSNESLGResidualNormDestroy(PetscObject **objs)
270 {
271   PetscDraw      draw;
272   PetscErrorCode ierr;
273   PetscDrawLG    drawlg = (PetscDrawLG) (*objs)[1];
274 
275   PetscFunctionBegin;
276   ierr = PetscDrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr);
277   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
278   ierr = PetscDrawLGDestroy(&drawlg);CHKERRQ(ierr);
279   ierr = PetscFree(*objs);CHKERRQ(ierr);
280   PetscFunctionReturn(0);
281 }
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "SNESMonitorDefault"
285 /*@C
286    SNESMonitorDefault - Monitors progress of the SNES solvers (default).
287 
288    Collective on SNES
289 
290    Input Parameters:
291 +  snes - the SNES context
292 .  its - iteration number
293 .  fgnorm - 2-norm of residual
294 -  dummy - unused context
295 
296    Notes:
297    This routine prints the residual norm at each iteration.
298 
299    Level: intermediate
300 
301 .keywords: SNES, nonlinear, default, monitor, norm
302 
303 .seealso: SNESMonitorSet(), SNESMonitorSolution()
304 @*/
305 PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
306 {
307   PetscErrorCode ierr;
308   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
309 
310   PetscFunctionBegin;
311   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
312   ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr);
313   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
314   PetscFunctionReturn(0);
315 }
316 
317 #undef __FUNCT__
318 #define __FUNCT__ "SNESMonitorJacUpdateSpectrum"
319 PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,void *ctx)
320 {
321 #if defined(PETSC_MISSING_LAPACK_GEEV)
322   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
323 #elif defined(PETSC_HAVE_ESSL)
324   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - No support for ESSL Lapack Routines");
325 #else
326   Vec            X;
327   Mat            J,dJ,dJdense;
328   PetscErrorCode ierr;
329   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
330   PetscInt       n,i;
331   PetscBLASInt   nb,lwork;
332   PetscReal      *eigr,*eigi;
333   PetscScalar    *work;
334   PetscScalar    *a;
335 
336   PetscFunctionBegin;
337   if (it == 0) PetscFunctionReturn(0);
338   /* create the difference between the current update and the current jacobian */
339   ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr);
340   ierr = SNESGetJacobian(snes,NULL,&J,&func,NULL);CHKERRQ(ierr);
341   ierr = MatDuplicate(J,MAT_COPY_VALUES,&dJ);CHKERRQ(ierr);
342   ierr = SNESComputeJacobian(snes,X,dJ,dJ);CHKERRQ(ierr);
343   ierr = MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
344 
345   /* compute the spectrum directly */
346   ierr  = MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);CHKERRQ(ierr);
347   ierr  = MatGetSize(dJ,&n,NULL);CHKERRQ(ierr);
348   ierr  = PetscBLASIntCast(n,&nb);CHKERRQ(ierr);
349   lwork = 3*nb;
350   ierr  = PetscMalloc1(n,&eigr);CHKERRQ(ierr);
351   ierr  = PetscMalloc1(n,&eigi);CHKERRQ(ierr);
352   ierr  = PetscMalloc1(lwork,&work);CHKERRQ(ierr);
353   ierr  = MatDenseGetArray(dJdense,&a);CHKERRQ(ierr);
354 #if !defined(PETSC_USE_COMPLEX)
355   {
356     PetscBLASInt lierr;
357     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
358     PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,NULL,&nb,NULL,&nb,work,&lwork,&lierr));
359     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr);
360     ierr = PetscFPTrapPop();CHKERRQ(ierr);
361   }
362 #else
363   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
364 #endif
365   PetscPrintf(PetscObjectComm((PetscObject)snes),"Eigenvalues of J_%d - J_%d:\n",it,it-1);CHKERRQ(ierr);
366   for (i=0;i<n;i++) {
367     PetscPrintf(PetscObjectComm((PetscObject)snes),"%5d: %20.5g + %20.5gi\n",i,(double)eigr[i],(double)eigi[i]);CHKERRQ(ierr);
368   }
369   ierr = MatDenseRestoreArray(dJdense,&a);CHKERRQ(ierr);
370   ierr = MatDestroy(&dJ);CHKERRQ(ierr);
371   ierr = MatDestroy(&dJdense);CHKERRQ(ierr);
372   ierr = PetscFree(eigr);CHKERRQ(ierr);
373   ierr = PetscFree(eigi);CHKERRQ(ierr);
374   ierr = PetscFree(work);CHKERRQ(ierr);
375   PetscFunctionReturn(0);
376 #endif
377 }
378 
379 #undef __FUNCT__
380 #define __FUNCT__ "SNESMonitorRange_Private"
381 PetscErrorCode  SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
382 {
383   PetscErrorCode ierr;
384   Vec            resid;
385   PetscReal      rmax,pwork;
386   PetscInt       i,n,N;
387   PetscScalar    *r;
388 
389   PetscFunctionBegin;
390   ierr  = SNESGetFunction(snes,&resid,0,0);CHKERRQ(ierr);
391   ierr  = VecNorm(resid,NORM_INFINITY,&rmax);CHKERRQ(ierr);
392   ierr  = VecGetLocalSize(resid,&n);CHKERRQ(ierr);
393   ierr  = VecGetSize(resid,&N);CHKERRQ(ierr);
394   ierr  = VecGetArray(resid,&r);CHKERRQ(ierr);
395   pwork = 0.0;
396   for (i=0; i<n; i++) {
397     pwork += (PetscAbsScalar(r[i]) > .20*rmax);
398   }
399   ierr = MPI_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
400   ierr = VecRestoreArray(resid,&r);CHKERRQ(ierr);
401   *per = *per/N;
402   PetscFunctionReturn(0);
403 }
404 
405 #undef __FUNCT__
406 #define __FUNCT__ "SNESMonitorRange"
407 /*@C
408    SNESMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value.
409 
410    Collective on SNES
411 
412    Input Parameters:
413 +  snes   - iterative context
414 .  it    - iteration number
415 .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
416 -  dummy - unused monitor context
417 
418    Options Database Key:
419 .  -snes_monitor_range - Activates SNESMonitorRange()
420 
421    Level: intermediate
422 
423 .keywords: SNES, default, monitor, residual
424 
425 .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
426 @*/
427 PetscErrorCode  SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,void *dummy)
428 {
429   PetscErrorCode ierr;
430   PetscReal      perc,rel;
431   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
432   /* should be in a MonitorRangeContext */
433   static PetscReal prev;
434 
435   PetscFunctionBegin;
436   if (!it) prev = rnorm;
437   ierr = SNESMonitorRange_Private(snes,it,&perc);CHKERRQ(ierr);
438 
439   rel  = (prev - rnorm)/prev;
440   prev = rnorm;
441   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
442   ierr = 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));CHKERRQ(ierr);
443   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
444   PetscFunctionReturn(0);
445 }
446 
447 typedef struct {
448   PetscViewer viewer;
449   PetscReal   *history;
450 } SNESMonitorRatioContext;
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "SNESMonitorRatio"
454 /*@C
455    SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio
456    of residual norm at each iteration to the previous.
457 
458    Collective on SNES
459 
460    Input Parameters:
461 +  snes - the SNES context
462 .  its - iteration number
463 .  fgnorm - 2-norm of residual (or gradient)
464 -  dummy -  context of monitor
465 
466    Level: intermediate
467 
468 .keywords: SNES, nonlinear, monitor, norm
469 
470 .seealso: SNESMonitorSet(), SNESMonitorSolution()
471 @*/
472 PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
473 {
474   PetscErrorCode          ierr;
475   PetscInt                len;
476   PetscReal               *history;
477   SNESMonitorRatioContext *ctx = (SNESMonitorRatioContext*)dummy;
478 
479   PetscFunctionBegin;
480   ierr = SNESGetConvergenceHistory(snes,&history,NULL,&len);CHKERRQ(ierr);
481   ierr = PetscViewerASCIIAddTab(ctx->viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
482   if (!its || !history || its > len) {
483     ierr = PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr);
484   } else {
485     PetscReal ratio = fgnorm/history[its-1];
486     ierr = PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);CHKERRQ(ierr);
487   }
488   ierr = PetscViewerASCIISubtractTab(ctx->viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
489   PetscFunctionReturn(0);
490 }
491 
492 /*
493    If the we set the history monitor space then we must destroy it
494 */
495 #undef __FUNCT__
496 #define __FUNCT__ "SNESMonitorRatioDestroy"
497 PetscErrorCode SNESMonitorRatioDestroy(void **ct)
498 {
499   PetscErrorCode          ierr;
500   SNESMonitorRatioContext *ctx = *(SNESMonitorRatioContext**)ct;
501 
502   PetscFunctionBegin;
503   ierr = PetscFree(ctx->history);CHKERRQ(ierr);
504   ierr = PetscViewerDestroy(&ctx->viewer);CHKERRQ(ierr);
505   ierr = PetscFree(ctx);CHKERRQ(ierr);
506   PetscFunctionReturn(0);
507 }
508 
509 #undef __FUNCT__
510 #define __FUNCT__ "SNESMonitorSetRatio"
511 /*@C
512    SNESMonitorSetRatio - Sets SNES to use a monitor that prints the
513    ratio of the function norm at each iteration.
514 
515    Collective on SNES
516 
517    Input Parameters:
518 +   snes - the SNES context
519 -   viewer - ASCII viewer to print output
520 
521    Level: intermediate
522 
523 .keywords: SNES, nonlinear, monitor, norm
524 
525 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault()
526 @*/
527 PetscErrorCode  SNESMonitorSetRatio(SNES snes,PetscViewer viewer)
528 {
529   PetscErrorCode          ierr;
530   SNESMonitorRatioContext *ctx;
531   PetscReal               *history;
532 
533   PetscFunctionBegin;
534   if (!viewer) {
535     ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)snes),"stdout",&viewer);CHKERRQ(ierr);
536     ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);
537   }
538   ierr = PetscNewLog(snes,&ctx);CHKERRQ(ierr);
539   ierr = SNESGetConvergenceHistory(snes,&history,NULL,NULL);CHKERRQ(ierr);
540   if (!history) {
541     ierr = PetscMalloc1(100,&ctx->history);CHKERRQ(ierr);
542     ierr = SNESSetConvergenceHistory(snes,ctx->history,0,100,PETSC_TRUE);CHKERRQ(ierr);
543   }
544   ctx->viewer = viewer;
545 
546   ierr = SNESMonitorSet(snes,SNESMonitorRatio,ctx,SNESMonitorRatioDestroy);CHKERRQ(ierr);
547   PetscFunctionReturn(0);
548 }
549 
550 /* ---------------------------------------------------------------- */
551 #undef __FUNCT__
552 #define __FUNCT__ "SNESMonitorDefaultShort"
553 /*
554      Default (short) SNES Monitor, same as SNESMonitorDefault() except
555   it prints fewer digits of the residual as the residual gets smaller.
556   This is because the later digits are meaningless and are often
557   different on different machines; by using this routine different
558   machines will usually generate the same output.
559 */
560 PetscErrorCode  SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
561 {
562   PetscErrorCode ierr;
563   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
564 
565   PetscFunctionBegin;
566   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
567   if (fgnorm > 1.e-9) {
568     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %g \n",its,(double)fgnorm);CHKERRQ(ierr);
569   } else if (fgnorm > 1.e-11) {
570     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,(double)fgnorm);CHKERRQ(ierr);
571   } else {
572     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);CHKERRQ(ierr);
573   }
574   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
575   PetscFunctionReturn(0);
576 }
577 /* ---------------------------------------------------------------- */
578 #undef __FUNCT__
579 #define __FUNCT__ "SNESConvergedDefault"
580 /*@C
581    SNESConvergedDefault - Convergence test of the solvers for
582    systems of nonlinear equations (default).
583 
584    Collective on SNES
585 
586    Input Parameters:
587 +  snes - the SNES context
588 .  it - the iteration (0 indicates before any Newton steps)
589 .  xnorm - 2-norm of current iterate
590 .  snorm - 2-norm of current step
591 .  fnorm - 2-norm of function at current iterate
592 -  dummy - unused context
593 
594    Output Parameter:
595 .   reason  - one of
596 $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
597 $  SNES_CONVERGED_SNORM_RELATIVE  - (snorm < stol*xnorm),
598 $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
599 $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
600 $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
601 $  SNES_CONVERGED_ITERATING       - (otherwise),
602 
603    where
604 +    maxf - maximum number of function evaluations,
605             set with SNESSetTolerances()
606 .    nfct - number of function evaluations,
607 .    abstol - absolute function norm tolerance,
608             set with SNESSetTolerances()
609 -    rtol - relative function norm tolerance, set with SNESSetTolerances()
610 
611    Level: intermediate
612 
613 .keywords: SNES, nonlinear, default, converged, convergence
614 
615 .seealso: SNESSetConvergenceTest()
616 @*/
617 PetscErrorCode  SNESConvergedDefault(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
618 {
619   PetscErrorCode ierr;
620 
621   PetscFunctionBegin;
622   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
623   PetscValidPointer(reason,6);
624 
625   *reason = SNES_CONVERGED_ITERATING;
626 
627   if (!it) {
628     /* set parameter for default relative tolerance convergence test */
629     snes->ttol = fnorm*snes->rtol;
630   }
631   if (PetscIsInfOrNanReal(fnorm)) {
632     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
633     *reason = SNES_DIVERGED_FNORM_NAN;
634   } else if (fnorm < snes->abstol) {
635     ierr    = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
636     *reason = SNES_CONVERGED_FNORM_ABS;
637   } else if (snes->nfuncs >= snes->max_funcs) {
638     ierr    = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
639     *reason = SNES_DIVERGED_FUNCTION_COUNT;
640   }
641 
642   if (it && !*reason) {
643     if (fnorm <= snes->ttol) {
644       ierr    = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
645       *reason = SNES_CONVERGED_FNORM_RELATIVE;
646     } else if (snorm < snes->stol*xnorm) {
647       ierr    = PetscInfo3(snes,"Converged due to small update length: %14.12e < %14.12e * %14.12e\n",(double)snorm,(double)snes->stol,(double)xnorm);CHKERRQ(ierr);
648       *reason = SNES_CONVERGED_SNORM_RELATIVE;
649     }
650   }
651   PetscFunctionReturn(0);
652 }
653 
654 #undef __FUNCT__
655 #define __FUNCT__ "SNESConvergedSkip"
656 /*@C
657    SNESConvergedSkip - Convergence test for SNES that NEVER returns as
658    converged, UNLESS the maximum number of iteration have been reached.
659 
660    Logically Collective on SNES
661 
662    Input Parameters:
663 +  snes - the SNES context
664 .  it - the iteration (0 indicates before any Newton steps)
665 .  xnorm - 2-norm of current iterate
666 .  snorm - 2-norm of current step
667 .  fnorm - 2-norm of function at current iterate
668 -  dummy - unused context
669 
670    Output Parameter:
671 .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
672 
673    Notes:
674    Convergence is then declared after a fixed number of iterations have been used.
675 
676    Level: advanced
677 
678 .keywords: SNES, nonlinear, skip, converged, convergence
679 
680 .seealso: SNESSetConvergenceTest()
681 @*/
682 PetscErrorCode  SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
683 {
684   PetscErrorCode ierr;
685 
686   PetscFunctionBegin;
687   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
688   PetscValidPointer(reason,6);
689 
690   *reason = SNES_CONVERGED_ITERATING;
691 
692   if (fnorm != fnorm) {
693     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
694     *reason = SNES_DIVERGED_FNORM_NAN;
695   } else if (it == snes->max_its) {
696     *reason = SNES_CONVERGED_ITS;
697   }
698   PetscFunctionReturn(0);
699 }
700 
701 #undef __FUNCT__
702 #define __FUNCT__ "SNESSetWorkVecs"
703 /*@C
704   SNESSetWorkVecs - Gets a number of work vectors.
705 
706   Input Parameters:
707 . snes  - the SNES context
708 . nw - number of work vectors to allocate
709 
710    Level: developer
711 
712    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations
713 
714 @*/
715 PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw)
716 {
717   DM             dm;
718   Vec            v;
719   PetscErrorCode ierr;
720 
721   PetscFunctionBegin;
722   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
723   snes->nwork = nw;
724 
725   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
726   ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr);
727   ierr = VecDuplicateVecs(v,snes->nwork,&snes->work);CHKERRQ(ierr);
728   ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr);
729   ierr = PetscLogObjectParents(snes,nw,snes->work);CHKERRQ(ierr);
730   PetscFunctionReturn(0);
731 }
732