xref: /petsc/src/snes/interface/snesut.c (revision 009bbdc485cd9ad46be9940d3549e2dde9cdc322)
1 
2 #include <petsc-private/snesimpl.h>       /*I   "petscsnes.h"   I*/
3 #include <petscblaslapack.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "SNESMonitorSolution"
7 /*@C
8    SNESMonitorSolution - Monitors progress of the SNES solvers by calling
9    VecView() for the approximate solution at each iteration.
10 
11    Collective on SNES
12 
13    Input Parameters:
14 +  snes - the SNES context
15 .  its - iteration number
16 .  fgnorm - 2-norm of residual
17 -  dummy - either a viewer or PETSC_NULL
18 
19    Level: intermediate
20 
21 .keywords: SNES, nonlinear, vector, monitor, view
22 
23 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
24 @*/
25 PetscErrorCode  SNESMonitorSolution(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
26 {
27   PetscErrorCode ierr;
28   Vec            x;
29   PetscViewer    viewer = (PetscViewer) dummy;
30 
31   PetscFunctionBegin;
32   ierr = SNESGetSolution(snes,&x);CHKERRQ(ierr);
33   if (!viewer) {
34     MPI_Comm comm;
35     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
36     viewer = PETSC_VIEWER_DRAW_(comm);
37   }
38   ierr = VecView(x,viewer);CHKERRQ(ierr);
39   PetscFunctionReturn(0);
40 }
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "SNESMonitorResidual"
44 /*@C
45    SNESMonitorResidual - Monitors progress of the SNES solvers by calling
46    VecView() for the residual at each iteration.
47 
48    Collective on SNES
49 
50    Input Parameters:
51 +  snes - the SNES context
52 .  its - iteration number
53 .  fgnorm - 2-norm of residual
54 -  dummy - either a viewer or PETSC_NULL
55 
56    Level: intermediate
57 
58 .keywords: SNES, nonlinear, vector, monitor, view
59 
60 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
61 @*/
62 PetscErrorCode  SNESMonitorResidual(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
63 {
64   PetscErrorCode ierr;
65   Vec            x;
66   PetscViewer    viewer = (PetscViewer) dummy;
67 
68   PetscFunctionBegin;
69   ierr = SNESGetFunction(snes,&x,0,0);CHKERRQ(ierr);
70   if (!viewer) {
71     MPI_Comm comm;
72     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
73     viewer = PETSC_VIEWER_DRAW_(comm);
74   }
75   ierr = VecView(x,viewer);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "SNESMonitorSolutionUpdate"
81 /*@C
82    SNESMonitorSolutionUpdate - Monitors progress of the SNES solvers by calling
83    VecView() for the UPDATE to the solution at each iteration.
84 
85    Collective on SNES
86 
87    Input Parameters:
88 +  snes - the SNES context
89 .  its - iteration number
90 .  fgnorm - 2-norm of residual
91 -  dummy - either a viewer or PETSC_NULL
92 
93    Level: intermediate
94 
95 .keywords: SNES, nonlinear, vector, monitor, view
96 
97 .seealso: SNESMonitorSet(), SNESMonitorDefault(), VecView()
98 @*/
99 PetscErrorCode  SNESMonitorSolutionUpdate(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
100 {
101   PetscErrorCode ierr;
102   Vec            x;
103   PetscViewer    viewer = (PetscViewer) dummy;
104 
105   PetscFunctionBegin;
106   ierr = SNESGetSolutionUpdate(snes,&x);CHKERRQ(ierr);
107   if (!viewer) {
108     MPI_Comm comm;
109     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
110     viewer = PETSC_VIEWER_DRAW_(comm);
111   }
112   ierr = VecView(x,viewer);CHKERRQ(ierr);
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "SNESMonitorDefault"
118 /*@C
119    SNESMonitorDefault - Monitors progress of the SNES solvers (default).
120 
121    Collective on SNES
122 
123    Input Parameters:
124 +  snes - the SNES context
125 .  its - iteration number
126 .  fgnorm - 2-norm of residual
127 -  dummy - unused context
128 
129    Notes:
130    This routine prints the residual norm at each iteration.
131 
132    Level: intermediate
133 
134 .keywords: SNES, nonlinear, default, monitor, norm
135 
136 .seealso: SNESMonitorSet(), SNESMonitorSolution()
137 @*/
138 PetscErrorCode  SNESMonitorDefault(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
139 {
140   PetscErrorCode ierr;
141   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
142 
143   PetscFunctionBegin;
144   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
145   ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr);
146   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "SNESMonitorJacUpdateSpectrum"
152 PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES snes,PetscInt it,PetscReal fnorm,void *ctx)
153 {
154 #if defined(PETSC_MISSING_LAPACK_GEEV)
155   SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable\nNot able to provide eigen values.");
156 #elif defined(PETSC_HAVE_ESSL)
157   SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_SUP,"GEEV - No support for ESSL Lapack Routines");
158 #else
159   Vec            X;
160   Mat            J,dJ,dJdense;
161   PetscErrorCode ierr;
162   PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
163   PetscInt       n,i;
164   PetscBLASInt   nb,lwork;
165   PetscReal      *eigr,*eigi;
166   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
167   PetscScalar    *work;
168   PetscScalar    *a;
169 
170   PetscFunctionBegin;
171   if (it == 0) PetscFunctionReturn(0);
172   /* create the difference between the current update and the current jacobian */
173   ierr = SNESGetSolution(snes,&X);CHKERRQ(ierr);
174   ierr = SNESGetJacobian(snes,&J,PETSC_NULL,&func,PETSC_NULL);CHKERRQ(ierr);
175   ierr = MatDuplicate(J,MAT_COPY_VALUES,&dJ);CHKERRQ(ierr);
176   ierr = SNESComputeJacobian(snes,X,&dJ,&dJ,&flg);CHKERRQ(ierr);
177   ierr = MatAXPY(dJ,-1.0,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
178   /* compute the spectrum directly */
179   ierr = MatConvert(dJ,MATSEQDENSE,MAT_INITIAL_MATRIX,&dJdense);CHKERRQ(ierr);
180   ierr = MatGetSize(dJ,&n,PETSC_NULL);CHKERRQ(ierr);
181   ierr = PetscBLASIntCast(n,&nb);CHKERRQ(ierr);
182   lwork = 3*nb;
183   ierr = PetscMalloc(n*sizeof(PetscReal),&eigr);CHKERRQ(ierr);
184   ierr = PetscMalloc(n*sizeof(PetscReal),&eigi);CHKERRQ(ierr);
185   ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
186   ierr = MatDenseGetArray(dJdense,&a);CHKERRQ(ierr);
187 #if !defined(PETSC_USE_COMPLEX)
188   {
189     PetscBLASInt lierr;
190     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
191     LAPACKgeev_("N","N",&nb,a,&nb,eigr,eigi,PETSC_NULL,&nb,PETSC_NULL,&nb,work,&lwork,&lierr);
192     if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"geev() error %d",lierr);
193     ierr = PetscFPTrapPop();CHKERRQ(ierr);
194   }
195 #else
196   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not coded for complex");
197 #endif
198   PetscPrintf(((PetscObject)snes)->comm,"Eigenvalues of J_%d - J_%d:\n",it,it-1);CHKERRQ(ierr);
199   for (i=0;i<n;i++) {
200     PetscPrintf(((PetscObject)snes)->comm,"%5d: %20.5g + %20.5gi\n",i,eigr[i],eigi[i]);CHKERRQ(ierr);
201   }
202   ierr = MatDenseRestoreArray(dJdense,&a);CHKERRQ(ierr);
203   ierr = MatDestroy(&dJ);CHKERRQ(ierr);
204   ierr = MatDestroy(&dJdense);CHKERRQ(ierr);
205   ierr = PetscFree(eigr);CHKERRQ(ierr);
206   ierr = PetscFree(eigi);CHKERRQ(ierr);
207   ierr = PetscFree(work);CHKERRQ(ierr);
208   PetscFunctionReturn(0);
209 #endif
210 }
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "SNESMonitorRange_Private"
214 PetscErrorCode  SNESMonitorRange_Private(SNES snes,PetscInt it,PetscReal *per)
215 {
216   PetscErrorCode          ierr;
217   Vec                     resid;
218   PetscReal               rmax,pwork;
219   PetscInt                i,n,N;
220   PetscScalar             *r;
221 
222   PetscFunctionBegin;
223   ierr = SNESGetFunction(snes,&resid,0,0);CHKERRQ(ierr);
224   ierr = VecNorm(resid,NORM_INFINITY,&rmax);CHKERRQ(ierr);
225   ierr = VecGetLocalSize(resid,&n);CHKERRQ(ierr);
226   ierr = VecGetSize(resid,&N);CHKERRQ(ierr);
227   ierr = VecGetArray(resid,&r);CHKERRQ(ierr);
228   pwork = 0.0;
229   for (i=0; i<n; i++) {
230     pwork += (PetscAbsScalar(r[i]) > .20*rmax);
231   }
232   ierr = MPI_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,((PetscObject)snes)->comm);CHKERRQ(ierr);
233   ierr = VecRestoreArray(resid,&r);CHKERRQ(ierr);
234   *per  = *per/N;
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "SNESMonitorRange"
240 /*@C
241    SNESMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value.
242 
243    Collective on SNES
244 
245    Input Parameters:
246 +  snes   - iterative context
247 .  it    - iteration number
248 .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
249 -  dummy - unused monitor context
250 
251    Options Database Key:
252 .  -snes_monitor_range - Activates SNESMonitorRange()
253 
254    Level: intermediate
255 
256 .keywords: SNES, default, monitor, residual
257 
258 .seealso: SNESMonitorSet(), SNESMonitorDefault(), SNESMonitorLGCreate()
259 @*/
260 PetscErrorCode  SNESMonitorRange(SNES snes,PetscInt it,PetscReal rnorm,void *dummy)
261 {
262   PetscErrorCode   ierr;
263   PetscReal        perc,rel;
264   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
265   /* should be in a MonitorRangeContext */
266   static PetscReal prev;
267 
268   PetscFunctionBegin;
269   if (!it) prev = rnorm;
270   ierr = SNESMonitorRange_Private(snes,it,&perc);CHKERRQ(ierr);
271 
272   rel  = (prev - rnorm)/prev;
273   prev = rnorm;
274   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
275   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);
276   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 
280 typedef struct {
281   PetscViewer viewer;
282   PetscReal   *history;
283 } SNESMonitorRatioContext;
284 
285 #undef __FUNCT__
286 #define __FUNCT__ "SNESMonitorRatio"
287 /*@C
288    SNESMonitorRatio - Monitors progress of the SNES solvers by printing the ratio
289    of residual norm at each iteration to the previous.
290 
291    Collective on SNES
292 
293    Input Parameters:
294 +  snes - the SNES context
295 .  its - iteration number
296 .  fgnorm - 2-norm of residual (or gradient)
297 -  dummy -  context of monitor
298 
299    Level: intermediate
300 
301 .keywords: SNES, nonlinear, monitor, norm
302 
303 .seealso: SNESMonitorSet(), SNESMonitorSolution()
304 @*/
305 PetscErrorCode  SNESMonitorRatio(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
306 {
307   PetscErrorCode          ierr;
308   PetscInt                len;
309   PetscReal               *history;
310   SNESMonitorRatioContext *ctx = (SNESMonitorRatioContext*)dummy;
311 
312   PetscFunctionBegin;
313   ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,&len);CHKERRQ(ierr);
314   ierr = PetscViewerASCIIAddTab(ctx->viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
315   if (!its || !history || its > len) {
316     ierr = PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e \n",its,(double)fgnorm);CHKERRQ(ierr);
317   } else {
318     PetscReal ratio = fgnorm/history[its-1];
319     ierr = PetscViewerASCIIPrintf(ctx->viewer,"%3D SNES Function norm %14.12e %14.12e \n",its,(double)fgnorm,(double)ratio);CHKERRQ(ierr);
320   }
321   ierr = PetscViewerASCIISubtractTab(ctx->viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
322   PetscFunctionReturn(0);
323 }
324 
325 /*
326    If the we set the history monitor space then we must destroy it
327 */
328 #undef __FUNCT__
329 #define __FUNCT__ "SNESMonitorRatioDestroy"
330 PetscErrorCode SNESMonitorRatioDestroy(void **ct)
331 {
332   PetscErrorCode          ierr;
333   SNESMonitorRatioContext *ctx = *(SNESMonitorRatioContext**)ct;
334 
335   PetscFunctionBegin;
336   ierr = PetscFree(ctx->history);CHKERRQ(ierr);
337   ierr = PetscViewerDestroy(&ctx->viewer);CHKERRQ(ierr);
338   ierr = PetscFree(ctx);CHKERRQ(ierr);
339   PetscFunctionReturn(0);
340 }
341 
342 #undef __FUNCT__
343 #define __FUNCT__ "SNESMonitorSetRatio"
344 /*@C
345    SNESMonitorSetRatio - Sets SNES to use a monitor that prints the
346    ratio of the function norm at each iteration.
347 
348    Collective on SNES
349 
350    Input Parameters:
351 +   snes - the SNES context
352 -   viewer - ASCII viewer to print output
353 
354    Level: intermediate
355 
356 .keywords: SNES, nonlinear, monitor, norm
357 
358 .seealso: SNESMonitorSet(), SNESMonitorSolution(), SNESMonitorDefault()
359 @*/
360 PetscErrorCode  SNESMonitorSetRatio(SNES snes,PetscViewer viewer)
361 {
362   PetscErrorCode          ierr;
363   SNESMonitorRatioContext *ctx;
364   PetscReal               *history;
365 
366   PetscFunctionBegin;
367   if (!viewer) {
368     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&viewer);CHKERRQ(ierr);
369     ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);
370   }
371   ierr = PetscNewLog(snes,SNESMonitorRatioContext,&ctx);CHKERRQ(ierr);
372   ierr = SNESGetConvergenceHistory(snes,&history,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
373   if (!history) {
374     ierr = PetscMalloc(100*sizeof(PetscReal),&ctx->history);CHKERRQ(ierr);
375     ierr = SNESSetConvergenceHistory(snes,ctx->history,0,100,PETSC_TRUE);CHKERRQ(ierr);
376   }
377   ctx->viewer = viewer;
378   ierr = SNESMonitorSet(snes,SNESMonitorRatio,ctx,SNESMonitorRatioDestroy);CHKERRQ(ierr);
379   PetscFunctionReturn(0);
380 }
381 
382 /* ---------------------------------------------------------------- */
383 #undef __FUNCT__
384 #define __FUNCT__ "SNESMonitorDefaultShort"
385 /*
386      Default (short) SNES Monitor, same as SNESMonitorDefault() except
387   it prints fewer digits of the residual as the residual gets smaller.
388   This is because the later digits are meaningless and are often
389   different on different machines; by using this routine different
390   machines will usually generate the same output.
391 */
392 PetscErrorCode  SNESMonitorDefaultShort(SNES snes,PetscInt its,PetscReal fgnorm,void *dummy)
393 {
394   PetscErrorCode ierr;
395   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);
396 
397   PetscFunctionBegin;
398   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
399   if (fgnorm > 1.e-9) {
400     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %G \n",its,fgnorm);CHKERRQ(ierr);
401   } else if (fgnorm > 1.e-11) {
402     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm %5.3e \n",its,fgnorm);CHKERRQ(ierr);
403   } else {
404     ierr = PetscViewerASCIIPrintf(viewer,"%3D SNES Function norm < 1.e-11\n",its);CHKERRQ(ierr);
405   }
406   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
407   PetscFunctionReturn(0);
408 }
409 /* ---------------------------------------------------------------- */
410 #undef __FUNCT__
411 #define __FUNCT__ "SNESDefaultConverged"
412 /*@C
413    SNESDefaultConverged - Convergence test of the solvers for
414    systems of nonlinear equations (default).
415 
416    Collective on SNES
417 
418    Input Parameters:
419 +  snes - the SNES context
420 .  it - the iteration (0 indicates before any Newton steps)
421 .  xnorm - 2-norm of current iterate
422 .  snorm - 2-norm of current step
423 .  fnorm - 2-norm of function at current iterate
424 -  dummy - unused context
425 
426    Output Parameter:
427 .   reason  - one of
428 $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
429 $  SNES_CONVERGED_SNORM_RELATIVE  - (snorm < stol*xnorm),
430 $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
431 $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
432 $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
433 $  SNES_CONVERGED_ITERATING       - (otherwise),
434 
435    where
436 +    maxf - maximum number of function evaluations,
437             set with SNESSetTolerances()
438 .    nfct - number of function evaluations,
439 .    abstol - absolute function norm tolerance,
440             set with SNESSetTolerances()
441 -    rtol - relative function norm tolerance, set with SNESSetTolerances()
442 
443    Level: intermediate
444 
445 .keywords: SNES, nonlinear, default, converged, convergence
446 
447 .seealso: SNESSetConvergenceTest()
448 @*/
449 PetscErrorCode  SNESDefaultConverged(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
450 {
451   PetscErrorCode ierr;
452 
453   PetscFunctionBegin;
454   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
455   PetscValidPointer(reason,6);
456 
457   *reason = SNES_CONVERGED_ITERATING;
458 
459   if (!it) {
460     /* set parameter for default relative tolerance convergence test */
461     snes->ttol = fnorm*snes->rtol;
462   }
463   if (PetscIsInfOrNanReal(fnorm)) {
464     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
465     *reason = SNES_DIVERGED_FNORM_NAN;
466   } else if (fnorm < snes->abstol) {
467     ierr = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)snes->abstol);CHKERRQ(ierr);
468     *reason = SNES_CONVERGED_FNORM_ABS;
469   } else if (snes->nfuncs >= snes->max_funcs) {
470     ierr = PetscInfo2(snes,"Exceeded maximum number of function evaluations: %D > %D\n",snes->nfuncs,snes->max_funcs);CHKERRQ(ierr);
471     *reason = SNES_DIVERGED_FUNCTION_COUNT;
472   }
473 
474   if (it && !*reason) {
475     if (fnorm <= snes->ttol) {
476       ierr = PetscInfo2(snes,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)snes->ttol);CHKERRQ(ierr);
477       *reason = SNES_CONVERGED_FNORM_RELATIVE;
478     } else if (snorm < snes->stol*xnorm) {
479       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);
480       *reason = SNES_CONVERGED_SNORM_RELATIVE;
481     }
482   }
483   PetscFunctionReturn(0);
484 }
485 
486 #undef __FUNCT__
487 #define __FUNCT__ "SNESSkipConverged"
488 /*@C
489    SNESSkipConverged - Convergence test for SNES that NEVER returns as
490    converged, UNLESS the maximum number of iteration have been reached.
491 
492    Logically Collective on SNES
493 
494    Input Parameters:
495 +  snes - the SNES context
496 .  it - the iteration (0 indicates before any Newton steps)
497 .  xnorm - 2-norm of current iterate
498 .  snorm - 2-norm of current step
499 .  fnorm - 2-norm of function at current iterate
500 -  dummy - unused context
501 
502    Output Parameter:
503 .   reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
504 
505    Notes:
506    Convergence is then declared after a fixed number of iterations have been used.
507 
508    Level: advanced
509 
510 .keywords: SNES, nonlinear, skip, converged, convergence
511 
512 .seealso: SNESSetConvergenceTest()
513 @*/
514 PetscErrorCode  SNESSkipConverged(SNES snes,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
515 {
516   PetscErrorCode ierr;
517 
518   PetscFunctionBegin;
519   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
520   PetscValidPointer(reason,6);
521 
522   *reason = SNES_CONVERGED_ITERATING;
523 
524   if (fnorm != fnorm) {
525     ierr = PetscInfo(snes,"Failed to converged, function norm is NaN\n");CHKERRQ(ierr);
526     *reason = SNES_DIVERGED_FNORM_NAN;
527   } else if (it == snes->max_its) {
528     *reason = SNES_CONVERGED_ITS;
529   }
530   PetscFunctionReturn(0);
531 }
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "SNESDefaultGetWork"
535 /*@
536   SNESDefaultGetWork - Gets a number of work vectors.
537 
538   Input Parameters:
539 . snes  - the SNES context
540 . nw - number of work vectors to allocate
541 
542    Level: developer
543 
544   Notes:
545   Call this only if no work vectors have been allocated
546 @*/
547 PetscErrorCode SNESDefaultGetWork(SNES snes,PetscInt nw)
548 {
549   PetscErrorCode ierr;
550 
551   PetscFunctionBegin;
552   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
553   snes->nwork = nw;
554   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
555   ierr = PetscLogObjectParents(snes,nw,snes->work);CHKERRQ(ierr);
556   PetscFunctionReturn(0);
557 }
558