xref: /petsc/src/snes/interface/snesut.c (revision 00d931fe9835bef04c3bcd2a9a1bf118d64cc4c2)
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    lg;
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,&lg);CHKERRQ(ierr);
198   ierr = PetscDrawLGSetLegend(lg,names);CHKERRQ(ierr);
199   ierr = PetscDrawLGSetFromOptions(lg);CHKERRQ(ierr);
200   ierr = PetscDrawLGGetAxis(lg,&axis);CHKERRQ(ierr);
201   ierr = PetscDrawAxisSetLabels(axis,"Convergence of Residual Norm","Iteration","Residual Norm");CHKERRQ(ierr);
202   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
203 
204   ierr = PetscMalloc1(2,objs);CHKERRQ(ierr);
205   (*objs)[1] = (PetscObject)lg;
206   PetscFunctionReturn(0);
207 }
208 
209 #undef __FUNCT__
210 #define __FUNCT__ "KSPMonitorSNESLGResidualNorm"
211 PetscErrorCode  KSPMonitorSNESLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscObject *objs)
212 {
213   SNES           snes = (SNES) objs[0];
214   PetscDrawLG    lg   = (PetscDrawLG) objs[1];
215   PetscErrorCode ierr;
216   PetscReal      y[2];
217   Vec            snes_solution,work1,work2;
218 
219   PetscFunctionBegin;
220   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
221   else y[0] = -15.0;
222 
223   ierr = SNESGetSolution(snes,&snes_solution);CHKERRQ(ierr);
224   ierr = VecDuplicate(snes_solution,&work1);CHKERRQ(ierr);
225   ierr = VecDuplicate(snes_solution,&work2);CHKERRQ(ierr);
226   ierr = KSPBuildSolution(ksp,work1,NULL);CHKERRQ(ierr);
227   ierr = VecAYPX(work1,-1.0,snes_solution);CHKERRQ(ierr);
228   ierr = SNESComputeFunction(snes,work1,work2);CHKERRQ(ierr);
229   ierr = VecNorm(work2,NORM_2,y+1);CHKERRQ(ierr);
230   if (y[1] > 0.0) y[1] = PetscLog10Real(y[1]);
231   else y[1] = -15.0;
232   ierr = VecDestroy(&work1);CHKERRQ(ierr);
233   ierr = VecDestroy(&work2);CHKERRQ(ierr);
234 
235   ierr = PetscDrawLGAddPoint(lg,NULL,y);CHKERRQ(ierr);
236   if (n < 20 || !(n % 5)) {
237     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
238   }
239   PetscFunctionReturn(0);
240 }
241 
242 #undef __FUNCT__
243 #define __FUNCT__ "KSPMonitorSNESLGResidualNormDestroy"
244 /*@
245    KSPMonitorSNESLGResidualNormDestroy - Destroys a line graph context that was created
246    with KSPMonitorSNESLGResidualNormCreate().
247 
248    Collective on KSP
249 
250    Input Parameter:
251 .  draw - the drawing context
252 
253    Level: intermediate
254 
255 .keywords: KSP, monitor, line graph, destroy
256 
257 .seealso: KSPMonitorSNESLGResidualNormCreate(), KSPMonitorSNESLGTrueResidualDestroy(), KSPMonitorSet()
258 @*/
259 PetscErrorCode  KSPMonitorSNESLGResidualNormDestroy(PetscObject **objs)
260 {
261   PetscErrorCode ierr;
262   PetscDrawLG    lg = (PetscDrawLG) (*objs)[1];
263 
264   PetscFunctionBegin;
265   ierr = PetscDrawLGDestroy(&lg);CHKERRQ(ierr);
266   ierr = PetscFree(*objs);CHKERRQ(ierr);
267   PetscFunctionReturn(0);
268 }
269 
270 #undef __FUNCT__
271 #define __FUNCT__ "SNESMonitorDefault"
272 /*@C
273    SNESMonitorDefault - Monitors progress of the SNES solvers (default).
274 
275    Collective on SNES
276 
277    Input Parameters:
278 +  snes - the SNES context
279 .  its - iteration number
280 .  fgnorm - 2-norm of residual
281 -  dummy - unused context
282 
283    Notes:
284    This routine prints the residual norm at each iteration.
285 
286    Level: intermediate
287 
288 .keywords: SNES, nonlinear, default, monitor, norm
289 
290 .seealso: SNESMonitorSet(), SNESMonitorSolution()
291 @*/
292 PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
293 {
294   PetscErrorCode ierr;
295   PetscViewer    viewer = (PetscViewer) dummy;
296 
297   PetscFunctionBegin;
298   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
299   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
300   ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr);
301   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
302   PetscFunctionReturn(0);
303 }
304 
305 #undef __FUNCT__
306 #define __FUNCT__ "SNESMonitorJacUpdateSpectrum"
307 PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,void *ctx)
308 {
309 #if defined(PETSC_MISSING_LAPACK_GEEV)
310   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
311 #elif defined(PETSC_HAVE_ESSL)
312   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - No support for ESSL Lapack Routines");
313 #else
314   Vec            X;
315   Mat            J,dJ,dJdense;
316   PetscErrorCode ierr;
317   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
318   PetscInt       n,i;
319   PetscBLASInt   nb,lwork;
320   PetscReal      *eigr,*eigi;
321   PetscScalar    *work;
322   PetscScalar    *a;
323 
324   PetscFunctionBegin;
325   if (it == 0) PetscFunctionReturn(0);
326   /* create the difference between the current update and the current jacobian */
327   ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr);
328   ierr = SNESGetJacobian(snes,NULL,&J,&func,NULL);CHKERRQ(ierr);
329   ierr = MatDuplicate(J,MAT_COPY_VALUES,&dJ);CHKERRQ(ierr);
330   ierr = SNESComputeJacobian(snes,X,dJ,dJ);CHKERRQ(ierr);
331   ierr = MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
332 
333   /* compute the spectrum directly */
334   ierr  = MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);CHKERRQ(ierr);
335   ierr  = MatGetSize(dJ,&n,NULL);CHKERRQ(ierr);
336   ierr  = PetscBLASIntCast(n,&nb);CHKERRQ(ierr);
337   lwork = 3*nb;
338   ierr  = PetscMalloc1(n,&eigr);CHKERRQ(ierr);
339   ierr  = PetscMalloc1(n,&eigi);CHKERRQ(ierr);
340   ierr  = PetscMalloc1(lwork,&work);CHKERRQ(ierr);
341   ierr  = MatDenseGetArray(dJdense,&a);CHKERRQ(ierr);
342 #if !defined(PETSC_USE_COMPLEX)
343   {
344     PetscBLASInt lierr;
345     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
346     PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,NULL,&nb,NULL,&nb,work,&lwork,&lierr));
347     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr);
348     ierr = PetscFPTrapPop();CHKERRQ(ierr);
349   }
350 #else
351   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
352 #endif
353   ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Eigenvalues of J_%d - J_%d:\n",it,it-1);CHKERRQ(ierr);
354   for (i=0;i<n;i++) {
355     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"%5d: %20.5g + %20.5gi\n",i,(double)eigr[i],(double)eigi[i]);CHKERRQ(ierr);
356   }
357   ierr = MatDenseRestoreArray(dJdense,&a);CHKERRQ(ierr);
358   ierr = MatDestroy(&dJ);CHKERRQ(ierr);
359   ierr = MatDestroy(&dJdense);CHKERRQ(ierr);
360   ierr = PetscFree(eigr);CHKERRQ(ierr);
361   ierr = PetscFree(eigi);CHKERRQ(ierr);
362   ierr = PetscFree(work);CHKERRQ(ierr);
363   PetscFunctionReturn(0);
364 #endif
365 }
366 
367 #undef __FUNCT__
368 #define __FUNCT__ "SNESMonitorRange_Private"
369 PetscErrorCode  SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
370 {
371   PetscErrorCode ierr;
372   Vec            resid;
373   PetscReal      rmax,pwork;
374   PetscInt       i,n,N;
375   PetscScalar    *r;
376 
377   PetscFunctionBegin;
378   ierr  = SNESGetFunction(snes,&resid,0,0);CHKERRQ(ierr);
379   ierr  = VecNorm(resid,NORM_INFINITY,&rmax);CHKERRQ(ierr);
380   ierr  = VecGetLocalSize(resid,&n);CHKERRQ(ierr);
381   ierr  = VecGetSize(resid,&N);CHKERRQ(ierr);
382   ierr  = VecGetArray(resid,&r);CHKERRQ(ierr);
383   pwork = 0.0;
384   for (i=0; i<n; i++) {
385     pwork += (PetscAbsScalar(r[i]) > .20*rmax);
386   }
387   ierr = MPIU_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
388   ierr = VecRestoreArray(resid,&r);CHKERRQ(ierr);
389   *per = *per/N;
390   PetscFunctionReturn(0);
391 }
392 
393 #undef __FUNCT__
394 #define __FUNCT__ "SNESMonitorRange"
395 /*@C
396    SNESMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value.
397 
398    Collective on SNES
399 
400    Input Parameters:
401 +  snes   - iterative context
402 .  it    - iteration number
403 .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
404 -  dummy - unused monitor context
405 
406    Options Database Key:
407 .  -snes_monitor_range - Activates SNESMonitorRange()
408 
409    Level: intermediate
410 
411 .keywords: SNES, default, monitor, residual
412 
413 .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
414 @*/
415 PetscErrorCode  SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,void *dummy)
416 {
417   PetscErrorCode ierr;
418   PetscReal      perc,rel;
419   PetscViewer    viewer = (PetscViewer) dummy;
420   /* should be in a MonitorRangeContext */
421   static PetscReal prev;
422 
423   PetscFunctionBegin;
424   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
425   if (!it) prev = rnorm;
426   ierr = SNESMonitorRange_Private(snes,it,&perc);CHKERRQ(ierr);
427 
428   rel  = (prev - rnorm)/prev;
429   prev = rnorm;
430   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
431   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);
432   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
433   PetscFunctionReturn(0);
434 }
435 
436 #undef __FUNCT__
437 #define __FUNCT__ "SNESMonitorRatio"
438 /*@C
439    SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio
440    of residual norm at each iteration to the previous.
441 
442    Collective on SNES
443 
444    Input Parameters:
445 +  snes - the SNES context
446 .  its - iteration number
447 .  fgnorm - 2-norm of residual (or gradient)
448 -  dummy -  context of monitor
449 
450    Level: intermediate
451 
452    Notes: Insure that SNESMonitorRatio() is called when you set this monitor
453 .keywords: SNES, nonlinear, monitor, norm
454 
455 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorRatio()
456 @*/
457 PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
458 {
459   PetscErrorCode          ierr;
460   PetscInt                len;
461   PetscReal               *history;
462   PetscViewer             viewer = (PetscViewer) dummy;
463 
464   PetscFunctionBegin;
465   ierr = SNESGetConvergenceHistory(snes,&history,NULL,&len);CHKERRQ(ierr);
466   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
467   if (!its || !history || its > len) {
468     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr);
469   } else {
470     PetscReal ratio = fgnorm/history[its-1];
471     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);CHKERRQ(ierr);
472   }
473   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
474   PetscFunctionReturn(0);
475 }
476 
477 #undef __FUNCT__
478 #define __FUNCT__ "SNESMonitorRatioSetUp"
479 /*@C
480    SNESMonitorRatioSetUp - Insures the SNES object is saving its history since this monitor needs access to it
481 
482    Collective on SNES
483 
484    Input Parameters:
485 +   snes - the SNES context
486 -   viewer - the PetscViewer object (ignored)
487 
488    Level: intermediate
489 
490 .keywords: SNES, nonlinear, monitor, norm
491 
492 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorRatio()
493 @*/
494 PetscErrorCode  SNESMonitorRatioSetUp(SNES snes,PetscViewer viewer)
495 {
496   PetscErrorCode          ierr;
497   PetscReal               *history;
498 
499   PetscFunctionBegin;
500   ierr = SNESGetConvergenceHistory(snes,&history,NULL,NULL);CHKERRQ(ierr);
501   if (!history) {
502     ierr = SNESSetConvergenceHistory(snes,NULL,NULL,100,PETSC_TRUE);CHKERRQ(ierr);
503   }
504   PetscFunctionReturn(0);
505 }
506 
507 /* ---------------------------------------------------------------- */
508 #undef __FUNCT__
509 #define __FUNCT__ "SNESMonitorDefaultShort"
510 /*
511      Default (short) SNES Monitor, same as SNESMonitorDefault() except
512   it prints fewer digits of the residual as the residual gets smaller.
513   This is because the later digits are meaningless and are often
514   different on different machines; by using this routine different
515   machines will usually generate the same output.
516 */
517 PetscErrorCode  SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
518 {
519   PetscErrorCode ierr;
520   PetscViewer    viewer = (PetscViewer) dummy;
521 
522   PetscFunctionBegin;
523   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
524   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
525   if (fgnorm > 1.e-9) {
526     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %g \n",its,(double)fgnorm);CHKERRQ(ierr);
527   } else if (fgnorm > 1.e-11) {
528     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,(double)fgnorm);CHKERRQ(ierr);
529   } else {
530     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);CHKERRQ(ierr);
531   }
532   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
533   PetscFunctionReturn(0);
534 }
535 
536 #undef __FUNCT__
537 #define __FUNCT__ "SNESMonitorDefaultField"
538 /*@C
539   SNESMonitorDefaultField - Monitors progress of the SNES solvers, separated into fields.
540 
541   Collective on SNES
542 
543   Input Parameters:
544 + snes   - the SNES context
545 . its    - iteration number
546 . fgnorm - 2-norm of residual
547 - ctx    - the PetscViewer
548 
549   Notes:
550   This routine uses the DM attached to the residual vector
551 
552   Level: intermediate
553 
554 .keywords: SNES, nonlinear, field, monitor, norm
555 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorDefaultShort()
556 @*/
557 PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, void *ctx)
558 {
559   PetscViewer    viewer = (PetscViewer) ctx;
560   Vec            r;
561   DM             dm;
562   PetscReal      res[256];
563   PetscInt       tablevel;
564   PetscErrorCode ierr;
565 
566   PetscFunctionBegin;
567   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
568   ierr = SNESGetFunction(snes, &r, NULL, NULL);CHKERRQ(ierr);
569   ierr = VecGetDM(r, &dm);CHKERRQ(ierr);
570   if (!dm) {ierr = SNESMonitorDefault(snes, its, fgnorm, ctx);CHKERRQ(ierr);}
571   else {
572     PetscSection s, gs;
573     PetscInt     Nf, f;
574 
575     ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
576     ierr = DMGetDefaultGlobalSection(dm, &gs);CHKERRQ(ierr);
577     if (!s || !gs) {ierr = SNESMonitorDefault(snes, its, fgnorm, ctx);CHKERRQ(ierr);}
578     ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr);
579     if (Nf > 256) SETERRQ1(PetscObjectComm((PetscObject) snes), PETSC_ERR_SUP, "Do not support %d fields > 256", Nf);
580     ierr = PetscSectionVecNorm(s, gs, r, NORM_2, res);CHKERRQ(ierr);
581     ierr = PetscObjectGetTabLevel((PetscObject) snes, &tablevel);CHKERRQ(ierr);
582     ierr = PetscViewerASCIIAddTab(viewer, tablevel);CHKERRQ(ierr);
583     ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
584     for (f = 0; f < Nf; ++f) {
585       if (f) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
586       ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", res[f]);CHKERRQ(ierr);
587     }
588     ierr = PetscViewerASCIIPrintf(viewer, "] \n");CHKERRQ(ierr);
589     ierr = PetscViewerASCIISubtractTab(viewer, tablevel);CHKERRQ(ierr);
590   }
591   PetscFunctionReturn(0);
592 }
593 /* ---------------------------------------------------------------- */
594 #undef __FUNCT__
595 #define __FUNCT__ "SNESConvergedDefault"
596 /*@C
597    SNESConvergedDefault - Convergence test of the solvers for
598    systems of nonlinear equations (default).
599 
600    Collective on SNES
601 
602    Input Parameters:
603 +  snes - the SNES context
604 .  it - the iteration (0 indicates before any Newton steps)
605 .  xnorm - 2-norm of current iterate
606 .  snorm - 2-norm of current step
607 .  fnorm - 2-norm of function at current iterate
608 -  dummy - unused context
609 
610    Output Parameter:
611 .   reason  - one of
612 $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
613 $  SNES_CONVERGED_SNORM_RELATIVE  - (snorm < stol*xnorm),
614 $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
615 $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
616 $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
617 $  SNES_CONVERGED_ITERATING       - (otherwise),
618 
619    where
620 +    maxf - maximum number of function evaluations,
621             set with SNESSetTolerances()
622 .    nfct - number of function evaluations,
623 .    abstol - absolute function norm tolerance,
624             set with SNESSetTolerances()
625 -    rtol - relative function norm tolerance, set with SNESSetTolerances()
626 
627    Level: intermediate
628 
629 .keywords: SNES, nonlinear, default, converged, convergence
630 
631 .seealso: SNESSetConvergenceTest()
632 @*/
633 PetscErrorCode  SNESConvergedDefault(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
634 {
635   PetscErrorCode ierr;
636 
637   PetscFunctionBegin;
638   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
639   PetscValidPointer(reason,6);
640 
641   *reason = SNES_CONVERGED_ITERATING;
642 
643   if (!it) {
644     /* set parameter for default relative tolerance convergence test */
645     snes->ttol = fnorm*snes->rtol;
646   }
647   if (PetscIsInfOrNanReal(fnorm)) {
648     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
649     *reason = SNES_DIVERGED_FNORM_NAN;
650   } else if (fnorm < snes->abstol) {
651     ierr    = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
652     *reason = SNES_CONVERGED_FNORM_ABS;
653   } else if (snes->nfuncs >= snes->max_funcs) {
654     ierr    = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
655     *reason = SNES_DIVERGED_FUNCTION_COUNT;
656   }
657 
658   if (it && !*reason) {
659     if (fnorm <= snes->ttol) {
660       ierr    = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
661       *reason = SNES_CONVERGED_FNORM_RELATIVE;
662     } else if (snorm < snes->stol*xnorm) {
663       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);
664       *reason = SNES_CONVERGED_SNORM_RELATIVE;
665     }
666   }
667   PetscFunctionReturn(0);
668 }
669 
670 #undef __FUNCT__
671 #define __FUNCT__ "SNESConvergedSkip"
672 /*@C
673    SNESConvergedSkip - Convergence test for SNES that NEVER returns as
674    converged, UNLESS the maximum number of iteration have been reached.
675 
676    Logically Collective on SNES
677 
678    Input Parameters:
679 +  snes - the SNES context
680 .  it - the iteration (0 indicates before any Newton steps)
681 .  xnorm - 2-norm of current iterate
682 .  snorm - 2-norm of current step
683 .  fnorm - 2-norm of function at current iterate
684 -  dummy - unused context
685 
686    Output Parameter:
687 .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
688 
689    Notes:
690    Convergence is then declared after a fixed number of iterations have been used.
691 
692    Level: advanced
693 
694 .keywords: SNES, nonlinear, skip, converged, convergence
695 
696 .seealso: SNESSetConvergenceTest()
697 @*/
698 PetscErrorCode  SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
699 {
700   PetscErrorCode ierr;
701 
702   PetscFunctionBegin;
703   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
704   PetscValidPointer(reason,6);
705 
706   *reason = SNES_CONVERGED_ITERATING;
707 
708   if (fnorm != fnorm) {
709     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
710     *reason = SNES_DIVERGED_FNORM_NAN;
711   } else if (it == snes->max_its) {
712     *reason = SNES_CONVERGED_ITS;
713   }
714   PetscFunctionReturn(0);
715 }
716 
717 #undef __FUNCT__
718 #define __FUNCT__ "SNESSetWorkVecs"
719 /*@C
720   SNESSetWorkVecs - Gets a number of work vectors.
721 
722   Input Parameters:
723 . snes  - the SNES context
724 . nw - number of work vectors to allocate
725 
726    Level: developer
727 
728    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations
729 
730 @*/
731 PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw)
732 {
733   DM             dm;
734   Vec            v;
735   PetscErrorCode ierr;
736 
737   PetscFunctionBegin;
738   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
739   snes->nwork = nw;
740 
741   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
742   ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr);
743   ierr = VecDuplicateVecs(v,snes->nwork,&snes->work);CHKERRQ(ierr);
744   ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr);
745   ierr = PetscLogObjectParents(snes,nw,snes->work);CHKERRQ(ierr);
746   PetscFunctionReturn(0);
747 }
748