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