xref: /petsc/src/snes/interface/snesut.c (revision ca960df886e40c602abfa1f92aa06f03adc14bf5)
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) || snes->reason) {
237     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
238     ierr = PetscDrawLGSave(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    lg = (PetscDrawLG) (*objs)[1];
264 
265   PetscFunctionBegin;
266   ierr = PetscDrawLGDestroy(&lg);CHKERRQ(ierr);
267   ierr = PetscFree(*objs);CHKERRQ(ierr);
268   PetscFunctionReturn(0);
269 }
270 
271 #undef __FUNCT__
272 #define __FUNCT__ "SNESMonitorDefault"
273 /*@C
274    SNESMonitorDefault - Monitors progress of the SNES solvers (default).
275 
276    Collective on SNES
277 
278    Input Parameters:
279 +  snes - the SNES context
280 .  its - iteration number
281 .  fgnorm - 2-norm of residual
282 -  dummy - unused context
283 
284    Notes:
285    This routine prints the residual norm at each iteration.
286 
287    Level: intermediate
288 
289 .keywords: SNES, nonlinear, default, monitor, norm
290 
291 .seealso: SNESMonitorSet(), SNESMonitorSolution()
292 @*/
293 PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
294 {
295   PetscErrorCode ierr;
296   PetscViewer    viewer = (PetscViewer) dummy;
297 
298   PetscFunctionBegin;
299   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
300   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
301   ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr);
302   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
303   PetscFunctionReturn(0);
304 }
305 
306 #undef __FUNCT__
307 #define __FUNCT__ "SNESMonitorJacUpdateSpectrum"
308 PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,void *ctx)
309 {
310 #if defined(PETSC_MISSING_LAPACK_GEEV)
311   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
312 #elif defined(PETSC_HAVE_ESSL)
313   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"GEEV - No support for ESSL Lapack Routines");
314 #else
315   Vec            X;
316   Mat            J,dJ,dJdense;
317   PetscErrorCode ierr;
318   PetscErrorCode (*func)(SNES,Vec,Mat,Mat,void*);
319   PetscInt       n,i;
320   PetscBLASInt   nb,lwork;
321   PetscReal      *eigr,*eigi;
322   PetscScalar    *work;
323   PetscScalar    *a;
324 
325   PetscFunctionBegin;
326   if (it == 0) PetscFunctionReturn(0);
327   /* create the difference between the current update and the current jacobian */
328   ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr);
329   ierr = SNESGetJacobian(snes,NULL,&J,&func,NULL);CHKERRQ(ierr);
330   ierr = MatDuplicate(J,MAT_COPY_VALUES,&dJ);CHKERRQ(ierr);
331   ierr = SNESComputeJacobian(snes,X,dJ,dJ);CHKERRQ(ierr);
332   ierr = MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
333 
334   /* compute the spectrum directly */
335   ierr  = MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);CHKERRQ(ierr);
336   ierr  = MatGetSize(dJ,&n,NULL);CHKERRQ(ierr);
337   ierr  = PetscBLASIntCast(n,&nb);CHKERRQ(ierr);
338   lwork = 3*nb;
339   ierr  = PetscMalloc1(n,&eigr);CHKERRQ(ierr);
340   ierr  = PetscMalloc1(n,&eigi);CHKERRQ(ierr);
341   ierr  = PetscMalloc1(lwork,&work);CHKERRQ(ierr);
342   ierr  = MatDenseGetArray(dJdense,&a);CHKERRQ(ierr);
343 #if !defined(PETSC_USE_COMPLEX)
344   {
345     PetscBLASInt lierr;
346     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
347     PetscStackCallBLAS("LAPACKgeev",LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,NULL,&nb,NULL,&nb,work,&lwork,&lierr));
348     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr);
349     ierr = PetscFPTrapPop();CHKERRQ(ierr);
350   }
351 #else
352   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
353 #endif
354   ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Eigenvalues of J_%d - J_%d:\n",it,it-1);CHKERRQ(ierr);
355   for (i=0;i<n;i++) {
356     ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"%5d: %20.5g + %20.5gi\n",i,(double)eigr[i],(double)eigi[i]);CHKERRQ(ierr);
357   }
358   ierr = MatDenseRestoreArray(dJdense,&a);CHKERRQ(ierr);
359   ierr = MatDestroy(&dJ);CHKERRQ(ierr);
360   ierr = MatDestroy(&dJdense);CHKERRQ(ierr);
361   ierr = PetscFree(eigr);CHKERRQ(ierr);
362   ierr = PetscFree(eigi);CHKERRQ(ierr);
363   ierr = PetscFree(work);CHKERRQ(ierr);
364   PetscFunctionReturn(0);
365 #endif
366 }
367 
368 PETSC_INTERN PetscErrorCode  SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
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 = MPIU_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 #undef __FUNCT__
440 #define __FUNCT__ "SNESMonitorRatio"
441 /*@C
442    SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio
443    of residual norm at each iteration to the previous.
444 
445    Collective on SNES
446 
447    Input Parameters:
448 +  snes - the SNES context
449 .  its - iteration number
450 .  fgnorm - 2-norm of residual (or gradient)
451 -  dummy -  context of monitor
452 
453    Level: intermediate
454 
455    Notes: Insure that SNESMonitorRatio() is called when you set this monitor
456 .keywords: SNES, nonlinear, monitor, norm
457 
458 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorRatio()
459 @*/
460 PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
461 {
462   PetscErrorCode          ierr;
463   PetscInt                len;
464   PetscReal               *history;
465   PetscViewer             viewer = (PetscViewer) dummy;
466 
467   PetscFunctionBegin;
468   ierr = SNESGetConvergenceHistory(snes,&history,NULL,&len);CHKERRQ(ierr);
469   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
470   if (!its || !history || its > len) {
471     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr);
472   } else {
473     PetscReal ratio = fgnorm/history[its-1];
474     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);CHKERRQ(ierr);
475   }
476   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
477   PetscFunctionReturn(0);
478 }
479 
480 #undef __FUNCT__
481 #define __FUNCT__ "SNESMonitorRatioSetUp"
482 /*@C
483    SNESMonitorRatioSetUp - Insures the SNES object is saving its history since this monitor needs access to it
484 
485    Collective on SNES
486 
487    Input Parameters:
488 +   snes - the SNES context
489 -   viewer - the PetscViewer object (ignored)
490 
491    Level: intermediate
492 
493 .keywords: SNES, nonlinear, monitor, norm
494 
495 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorRatio()
496 @*/
497 PetscErrorCode  SNESMonitorRatioSetUp(SNES snes,PetscViewer viewer)
498 {
499   PetscErrorCode          ierr;
500   PetscReal               *history;
501 
502   PetscFunctionBegin;
503   ierr = SNESGetConvergenceHistory(snes,&history,NULL,NULL);CHKERRQ(ierr);
504   if (!history) {
505     ierr = SNESSetConvergenceHistory(snes,NULL,NULL,100,PETSC_TRUE);CHKERRQ(ierr);
506   }
507   PetscFunctionReturn(0);
508 }
509 
510 /* ---------------------------------------------------------------- */
511 #undef __FUNCT__
512 #define __FUNCT__ "SNESMonitorDefaultShort"
513 /*
514      Default (short) SNES Monitor, same as SNESMonitorDefault() except
515   it prints fewer digits of the residual as the residual gets smaller.
516   This is because the later digits are meaningless and are often
517   different on different machines; by using this routine different
518   machines will usually generate the same output.
519 */
520 PetscErrorCode  SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
521 {
522   PetscErrorCode ierr;
523   PetscViewer    viewer = (PetscViewer) dummy;
524 
525   PetscFunctionBegin;
526   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
527   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
528   if (fgnorm > 1.e-9) {
529     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %g \n",its,(double)fgnorm);CHKERRQ(ierr);
530   } else if (fgnorm > 1.e-11) {
531     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,(double)fgnorm);CHKERRQ(ierr);
532   } else {
533     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);CHKERRQ(ierr);
534   }
535   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
536   PetscFunctionReturn(0);
537 }
538 
539 #undef __FUNCT__
540 #define __FUNCT__ "SNESMonitorDefaultField"
541 /*@C
542   SNESMonitorDefaultField - Monitors progress of the SNES solvers, separated into fields.
543 
544   Collective on SNES
545 
546   Input Parameters:
547 + snes   - the SNES context
548 . its    - iteration number
549 . fgnorm - 2-norm of residual
550 - ctx    - the PetscViewer
551 
552   Notes:
553   This routine uses the DM attached to the residual vector
554 
555   Level: intermediate
556 
557 .keywords: SNES, nonlinear, field, monitor, norm
558 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault(), SNESMonitorDefaultShort()
559 @*/
560 PetscErrorCode SNESMonitorDefaultField(SNES snes, PetscInt its, PetscReal fgnorm, void *ctx)
561 {
562   PetscViewer    viewer = (PetscViewer) ctx;
563   Vec            r;
564   DM             dm;
565   PetscReal      res[256];
566   PetscInt       tablevel;
567   PetscErrorCode ierr;
568 
569   PetscFunctionBegin;
570   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
571   ierr = SNESGetFunction(snes, &r, NULL, NULL);CHKERRQ(ierr);
572   ierr = VecGetDM(r, &dm);CHKERRQ(ierr);
573   if (!dm) {ierr = SNESMonitorDefault(snes, its, fgnorm, ctx);CHKERRQ(ierr);}
574   else {
575     PetscSection s, gs;
576     PetscInt     Nf, f;
577 
578     ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
579     ierr = DMGetDefaultGlobalSection(dm, &gs);CHKERRQ(ierr);
580     if (!s || !gs) {ierr = SNESMonitorDefault(snes, its, fgnorm, ctx);CHKERRQ(ierr);}
581     ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr);
582     if (Nf > 256) SETERRQ1(PetscObjectComm((PetscObject) snes), PETSC_ERR_SUP, "Do not support %d fields > 256", Nf);
583     ierr = PetscSectionVecNorm(s, gs, r, NORM_2, res);CHKERRQ(ierr);
584     ierr = PetscObjectGetTabLevel((PetscObject) snes, &tablevel);CHKERRQ(ierr);
585     ierr = PetscViewerASCIIAddTab(viewer, tablevel);CHKERRQ(ierr);
586     ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
587     for (f = 0; f < Nf; ++f) {
588       if (f) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
589       ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", res[f]);CHKERRQ(ierr);
590     }
591     ierr = PetscViewerASCIIPrintf(viewer, "] \n");CHKERRQ(ierr);
592     ierr = PetscViewerASCIISubtractTab(viewer, tablevel);CHKERRQ(ierr);
593   }
594   PetscFunctionReturn(0);
595 }
596 /* ---------------------------------------------------------------- */
597 #undef __FUNCT__
598 #define __FUNCT__ "SNESConvergedDefault"
599 /*@C
600    SNESConvergedDefault - Convergence test of the solvers for
601    systems of nonlinear equations (default).
602 
603    Collective on SNES
604 
605    Input Parameters:
606 +  snes - the SNES context
607 .  it - the iteration (0 indicates before any Newton steps)
608 .  xnorm - 2-norm of current iterate
609 .  snorm - 2-norm of current step
610 .  fnorm - 2-norm of function at current iterate
611 -  dummy - unused context
612 
613    Output Parameter:
614 .   reason  - one of
615 $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
616 $  SNES_CONVERGED_SNORM_RELATIVE  - (snorm < stol*xnorm),
617 $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
618 $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
619 $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
620 $  SNES_CONVERGED_ITERATING       - (otherwise),
621 
622    where
623 +    maxf - maximum number of function evaluations,
624             set with SNESSetTolerances()
625 .    nfct - number of function evaluations,
626 .    abstol - absolute function norm tolerance,
627             set with SNESSetTolerances()
628 -    rtol - relative function norm tolerance, set with SNESSetTolerances()
629 
630    Level: intermediate
631 
632 .keywords: SNES, nonlinear, default, converged, convergence
633 
634 .seealso: SNESSetConvergenceTest()
635 @*/
636 PetscErrorCode  SNESConvergedDefault(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
637 {
638   PetscErrorCode ierr;
639 
640   PetscFunctionBegin;
641   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
642   PetscValidPointer(reason,6);
643 
644   *reason = SNES_CONVERGED_ITERATING;
645 
646   if (!it) {
647     /* set parameter for default relative tolerance convergence test */
648     snes->ttol = fnorm*snes->rtol;
649   }
650   if (PetscIsInfOrNanReal(fnorm)) {
651     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
652     *reason = SNES_DIVERGED_FNORM_NAN;
653   } else if (fnorm < snes->abstol) {
654     ierr    = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
655     *reason = SNES_CONVERGED_FNORM_ABS;
656   } else if (snes->nfuncs >= snes->max_funcs) {
657     ierr    = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
658     *reason = SNES_DIVERGED_FUNCTION_COUNT;
659   }
660 
661   if (it && !*reason) {
662     if (fnorm <= snes->ttol) {
663       ierr    = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
664       *reason = SNES_CONVERGED_FNORM_RELATIVE;
665     } else if (snorm < snes->stol*xnorm) {
666       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);
667       *reason = SNES_CONVERGED_SNORM_RELATIVE;
668     }
669   }
670   PetscFunctionReturn(0);
671 }
672 
673 #undef __FUNCT__
674 #define __FUNCT__ "SNESConvergedSkip"
675 /*@C
676    SNESConvergedSkip - Convergence test for SNES that NEVER returns as
677    converged, UNLESS the maximum number of iteration have been reached.
678 
679    Logically Collective on SNES
680 
681    Input Parameters:
682 +  snes - the SNES context
683 .  it - the iteration (0 indicates before any Newton steps)
684 .  xnorm - 2-norm of current iterate
685 .  snorm - 2-norm of current step
686 .  fnorm - 2-norm of function at current iterate
687 -  dummy - unused context
688 
689    Output Parameter:
690 .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
691 
692    Notes:
693    Convergence is then declared after a fixed number of iterations have been used.
694 
695    Level: advanced
696 
697 .keywords: SNES, nonlinear, skip, converged, convergence
698 
699 .seealso: SNESSetConvergenceTest()
700 @*/
701 PetscErrorCode  SNESConvergedSkip(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
702 {
703   PetscErrorCode ierr;
704 
705   PetscFunctionBegin;
706   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
707   PetscValidPointer(reason,6);
708 
709   *reason = SNES_CONVERGED_ITERATING;
710 
711   if (fnorm != fnorm) {
712     ierr    = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
713     *reason = SNES_DIVERGED_FNORM_NAN;
714   } else if (it == snes->max_its) {
715     *reason = SNES_CONVERGED_ITS;
716   }
717   PetscFunctionReturn(0);
718 }
719 
720 #undef __FUNCT__
721 #define __FUNCT__ "SNESSetWorkVecs"
722 /*@C
723   SNESSetWorkVecs - Gets a number of work vectors.
724 
725   Input Parameters:
726 . snes  - the SNES context
727 . nw - number of work vectors to allocate
728 
729    Level: developer
730 
731    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations
732 
733 @*/
734 PetscErrorCode SNESSetWorkVecs(SNES snes,PetscInt nw)
735 {
736   DM             dm;
737   Vec            v;
738   PetscErrorCode ierr;
739 
740   PetscFunctionBegin;
741   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
742   snes->nwork = nw;
743 
744   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
745   ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr);
746   ierr = VecDuplicateVecs(v,snes->nwork,&snes->work);CHKERRQ(ierr);
747   ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr);
748   ierr = PetscLogObjectParents(snes,nw,snes->work);CHKERRQ(ierr);
749   PetscFunctionReturn(0);
750 }
751