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