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