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