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