xref: /petsc/src/snes/interface/snes.c (revision 81540f2fa6ed413bc384182845a7ebe823db191f)
1 #define PETSCSNES_DLL
2 
3 #include "include/private/snesimpl.h"      /*I "petscsnes.h"  I*/
4 
5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
6 PetscFList SNESList              = PETSC_NULL;
7 
8 /* Logging support */
9 PetscCookie PETSCSNES_DLLEXPORT SNES_COOKIE;
10 PetscLogEvent  SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval;
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "SNESSetFunctionDomainError"
14 /*@
15    SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not
16      in the functions domain. For example, negative pressure.
17 
18    Collective on SNES
19 
20    Input Parameters:
21 .  SNES - the SNES context
22 
23    Level: advanced
24 
25 .keywords: SNES, view
26 
27 .seealso: SNESCreate(), SNESSetFunction()
28 @*/
29 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunctionDomainError(SNES snes)
30 {
31   PetscFunctionBegin;
32   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
33   snes->domainerror = PETSC_TRUE;
34   PetscFunctionReturn(0);
35 }
36 
37 #undef __FUNCT__
38 #define __FUNCT__ "SNESView"
39 /*@C
40    SNESView - Prints the SNES data structure.
41 
42    Collective on SNES
43 
44    Input Parameters:
45 +  SNES - the SNES context
46 -  viewer - visualization context
47 
48    Options Database Key:
49 .  -snes_view - Calls SNESView() at end of SNESSolve()
50 
51    Notes:
52    The available visualization contexts include
53 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
54 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
55          output where only the first processor opens
56          the file.  All other processors send their
57          data to the first processor to print.
58 
59    The user can open an alternative visualization context with
60    PetscViewerASCIIOpen() - output to a specified file.
61 
62    Level: beginner
63 
64 .keywords: SNES, view
65 
66 .seealso: PetscViewerASCIIOpen()
67 @*/
68 PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer)
69 {
70   SNESKSPEW           *kctx;
71   PetscErrorCode      ierr;
72   KSP                 ksp;
73   const SNESType      type;
74   PetscTruth          iascii,isstring;
75 
76   PetscFunctionBegin;
77   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
78   if (!viewer) {
79     ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr);
80   }
81   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
82   PetscCheckSameComm(snes,1,viewer,2);
83 
84   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
85   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
86   if (iascii) {
87     if (((PetscObject)snes)->prefix) {
88       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",((PetscObject)snes)->prefix);CHKERRQ(ierr);
89     } else {
90       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
91     }
92     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
93     if (type) {
94       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
95     } else {
96       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
97     }
98     if (snes->ops->view) {
99       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
100       ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr);
101       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
102     }
103     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
104     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, solution=%G\n",
105                  snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr);
106     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr);
107     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr);
108     if (snes->ksp_ewconv) {
109       kctx = (SNESKSPEW *)snes->kspconvctx;
110       if (kctx) {
111         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr);
112         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
113         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
114       }
115     }
116   } else if (isstring) {
117     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
118     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
119   }
120   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
121   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
122   ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
123   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
124   PetscFunctionReturn(0);
125 }
126 
127 /*
128   We retain a list of functions that also take SNES command
129   line options. These are called at the end SNESSetFromOptions()
130 */
131 #define MAXSETFROMOPTIONS 5
132 static PetscInt numberofsetfromoptions;
133 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
134 
135 #undef __FUNCT__
136 #define __FUNCT__ "SNESAddOptionsChecker"
137 /*@C
138   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
139 
140   Not Collective
141 
142   Input Parameter:
143 . snescheck - function that checks for options
144 
145   Level: developer
146 
147 .seealso: SNESSetFromOptions()
148 @*/
149 PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
150 {
151   PetscFunctionBegin;
152   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
153     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
154   }
155   othersetfromoptions[numberofsetfromoptions++] = snescheck;
156   PetscFunctionReturn(0);
157 }
158 
159 #undef __FUNCT__
160 #define __FUNCT__ "SNESSetFromOptions"
161 /*@
162    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
163 
164    Collective on SNES
165 
166    Input Parameter:
167 .  snes - the SNES context
168 
169    Options Database Keys:
170 +  -snes_type <type> - ls, tr, umls, umtr, test
171 .  -snes_stol - convergence tolerance in terms of the norm
172                 of the change in the solution between steps
173 .  -snes_atol <abstol> - absolute tolerance of residual norm
174 .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
175 .  -snes_max_it <max_it> - maximum number of iterations
176 .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
177 .  -snes_max_fail <max_fail> - maximum number of failures
178 .  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
179 .  -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
180 .  -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
181 .  -snes_trtol <trtol> - trust region tolerance
182 .  -snes_no_convergence_test - skip convergence test in nonlinear
183                                solver; hence iterations will continue until max_it
184                                or some other criterion is reached. Saves expense
185                                of convergence test
186 .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
187                                        filename given prints to stdout
188 .  -snes_monitor_solution - plots solution at each iteration
189 .  -snes_monitor_residual - plots residual (not its norm) at each iteration
190 .  -snes_monitor_solution_update - plots update to solution at each iteration
191 .  -snes_monitor_draw - plots residual norm at each iteration
192 .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
193 .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
194 -  -snes_converged_reason - print the reason for convergence/divergence after each solve
195 
196     Options Database for Eisenstat-Walker method:
197 +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
198 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
199 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
200 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
201 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
202 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
203 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
204 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
205 
206    Notes:
207    To see all options, run your program with the -help option or consult
208    the users manual.
209 
210    Level: beginner
211 
212 .keywords: SNES, nonlinear, set, options, database
213 
214 .seealso: SNESSetOptionsPrefix()
215 @*/
216 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes)
217 {
218   PetscTruth              flg;
219   PetscInt                i,indx,lag;
220   const char              *deft = SNESLS;
221   const char              *convtests[] = {"default","skip"};
222   SNESKSPEW               *kctx = NULL;
223   char                    type[256], monfilename[PETSC_MAX_PATH_LEN];
224   PetscViewerASCIIMonitor monviewer;
225   PetscErrorCode          ierr;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
229 
230   ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
231     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
232     if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; }
233     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
234     if (flg) {
235       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
236     } else if (!((PetscObject)snes)->type_name) {
237       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
238     }
239     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
240 
241     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
242     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr);
243 
244     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
245     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
246     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
247     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
248     ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr);
249 
250     ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr);
251     if (flg) {
252       ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr);
253     }
254     ierr = PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);CHKERRQ(ierr);
255     if (flg) {
256       ierr = SNESSetLagJacobian(snes,lag);CHKERRQ(ierr);
257     }
258 
259     ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr);
260     if (flg) {
261       switch (indx) {
262       case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break;
263       case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);    break;
264       }
265     }
266 
267     ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr);
268     if (flg) {
269       snes->printreason = PETSC_TRUE;
270     }
271 
272     kctx = (SNESKSPEW *)snes->kspconvctx;
273 
274     ierr = PetscOptionsTruth("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr);
275 
276     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
277     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
278     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
279     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
280     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
281     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
282     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
283 
284     ierr = PetscOptionsName("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",&flg);CHKERRQ(ierr);
285     if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);}
286 
287     ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
288     if (flg) {
289       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr);
290       ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
291     }
292 
293     ierr = PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
294     if (flg) {
295       ierr = SNESMonitorSet(snes,SNESMonitorRange,0,0);CHKERRQ(ierr);
296     }
297 
298     ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
299     if (flg) {
300       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr);
301       ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr);
302     }
303 
304     ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
305     if (flg) {
306       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr);
307       ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
308     }
309 
310     ierr = PetscOptionsName("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",&flg);CHKERRQ(ierr);
311     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);}
312     ierr = PetscOptionsName("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",&flg);CHKERRQ(ierr);
313     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);}
314     ierr = PetscOptionsName("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",&flg);CHKERRQ(ierr);
315     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);}
316     ierr = PetscOptionsName("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",&flg);CHKERRQ(ierr);
317     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
318     ierr = PetscOptionsName("-snes_monitor_range_draw","Plot function range at each iteration","SNESMonitorLG",&flg);CHKERRQ(ierr);
319     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLGRange,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
320 
321     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
322     if (flg) {
323       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
324       ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr);
325     }
326 
327     for(i = 0; i < numberofsetfromoptions; i++) {
328       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
329     }
330 
331     if (snes->ops->setfromoptions) {
332       ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr);
333     }
334   ierr = PetscOptionsEnd();CHKERRQ(ierr);
335 
336 
337   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
338   ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr);
339 
340   PetscFunctionReturn(0);
341 }
342 
343 
344 #undef __FUNCT__
345 #define __FUNCT__ "SNESSetApplicationContext"
346 /*@
347    SNESSetApplicationContext - Sets the optional user-defined context for
348    the nonlinear solvers.
349 
350    Collective on SNES
351 
352    Input Parameters:
353 +  snes - the SNES context
354 -  usrP - optional user context
355 
356    Level: intermediate
357 
358 .keywords: SNES, nonlinear, set, application, context
359 
360 .seealso: SNESGetApplicationContext()
361 @*/
362 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP)
363 {
364   PetscFunctionBegin;
365   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
366   snes->user		= usrP;
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "SNESGetApplicationContext"
372 /*@C
373    SNESGetApplicationContext - Gets the user-defined context for the
374    nonlinear solvers.
375 
376    Not Collective
377 
378    Input Parameter:
379 .  snes - SNES context
380 
381    Output Parameter:
382 .  usrP - user context
383 
384    Level: intermediate
385 
386 .keywords: SNES, nonlinear, get, application, context
387 
388 .seealso: SNESSetApplicationContext()
389 @*/
390 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP)
391 {
392   PetscFunctionBegin;
393   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
394   *usrP = snes->user;
395   PetscFunctionReturn(0);
396 }
397 
398 #undef __FUNCT__
399 #define __FUNCT__ "SNESGetIterationNumber"
400 /*@
401    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
402    at this time.
403 
404    Not Collective
405 
406    Input Parameter:
407 .  snes - SNES context
408 
409    Output Parameter:
410 .  iter - iteration number
411 
412    Notes:
413    For example, during the computation of iteration 2 this would return 1.
414 
415    This is useful for using lagged Jacobians (where one does not recompute the
416    Jacobian at each SNES iteration). For example, the code
417 .vb
418       ierr = SNESGetIterationNumber(snes,&it);
419       if (!(it % 2)) {
420         [compute Jacobian here]
421       }
422 .ve
423    can be used in your ComputeJacobian() function to cause the Jacobian to be
424    recomputed every second SNES iteration.
425 
426    Level: intermediate
427 
428 .keywords: SNES, nonlinear, get, iteration, number,
429 
430 .seealso:   SNESGetFunctionNorm(), SNESGetLinearSolveIterations()
431 @*/
432 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter)
433 {
434   PetscFunctionBegin;
435   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
436   PetscValidIntPointer(iter,2);
437   *iter = snes->iter;
438   PetscFunctionReturn(0);
439 }
440 
441 #undef __FUNCT__
442 #define __FUNCT__ "SNESGetFunctionNorm"
443 /*@
444    SNESGetFunctionNorm - Gets the norm of the current function that was set
445    with SNESSSetFunction().
446 
447    Collective on SNES
448 
449    Input Parameter:
450 .  snes - SNES context
451 
452    Output Parameter:
453 .  fnorm - 2-norm of function
454 
455    Level: intermediate
456 
457 .keywords: SNES, nonlinear, get, function, norm
458 
459 .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations()
460 @*/
461 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscReal *fnorm)
462 {
463   PetscFunctionBegin;
464   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
465   PetscValidScalarPointer(fnorm,2);
466   *fnorm = snes->norm;
467   PetscFunctionReturn(0);
468 }
469 
470 #undef __FUNCT__
471 #define __FUNCT__ "SNESGetNonlinearStepFailures"
472 /*@
473    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
474    attempted by the nonlinear solver.
475 
476    Not Collective
477 
478    Input Parameter:
479 .  snes - SNES context
480 
481    Output Parameter:
482 .  nfails - number of unsuccessful steps attempted
483 
484    Notes:
485    This counter is reset to zero for each successive call to SNESSolve().
486 
487    Level: intermediate
488 
489 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
490 
491 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
492           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
493 @*/
494 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails)
495 {
496   PetscFunctionBegin;
497   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
498   PetscValidIntPointer(nfails,2);
499   *nfails = snes->numFailures;
500   PetscFunctionReturn(0);
501 }
502 
503 #undef __FUNCT__
504 #define __FUNCT__ "SNESSetMaxNonlinearStepFailures"
505 /*@
506    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
507    attempted by the nonlinear solver before it gives up.
508 
509    Not Collective
510 
511    Input Parameters:
512 +  snes     - SNES context
513 -  maxFails - maximum of unsuccessful steps
514 
515    Level: intermediate
516 
517 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
518 
519 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
520           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
521 @*/
522 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
523 {
524   PetscFunctionBegin;
525   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
526   snes->maxFailures = maxFails;
527   PetscFunctionReturn(0);
528 }
529 
530 #undef __FUNCT__
531 #define __FUNCT__ "SNESGetMaxNonlinearStepFailures"
532 /*@
533    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
534    attempted by the nonlinear solver before it gives up.
535 
536    Not Collective
537 
538    Input Parameter:
539 .  snes     - SNES context
540 
541    Output Parameter:
542 .  maxFails - maximum of unsuccessful steps
543 
544    Level: intermediate
545 
546 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
547 
548 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
549           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
550 
551 @*/
552 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
553 {
554   PetscFunctionBegin;
555   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
556   PetscValidIntPointer(maxFails,2);
557   *maxFails = snes->maxFailures;
558   PetscFunctionReturn(0);
559 }
560 
561 #undef __FUNCT__
562 #define __FUNCT__ "SNESGetNumberFunctionEvals"
563 /*@
564    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
565      done by SNES.
566 
567    Not Collective
568 
569    Input Parameter:
570 .  snes     - SNES context
571 
572    Output Parameter:
573 .  nfuncs - number of evaluations
574 
575    Level: intermediate
576 
577 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
578 
579 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures()
580 @*/
581 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
582 {
583   PetscFunctionBegin;
584   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
585   PetscValidIntPointer(nfuncs,2);
586   *nfuncs = snes->nfuncs;
587   PetscFunctionReturn(0);
588 }
589 
590 #undef __FUNCT__
591 #define __FUNCT__ "SNESGetLinearSolveFailures"
592 /*@
593    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
594    linear solvers.
595 
596    Not Collective
597 
598    Input Parameter:
599 .  snes - SNES context
600 
601    Output Parameter:
602 .  nfails - number of failed solves
603 
604    Notes:
605    This counter is reset to zero for each successive call to SNESSolve().
606 
607    Level: intermediate
608 
609 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
610 
611 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
612 @*/
613 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails)
614 {
615   PetscFunctionBegin;
616   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
617   PetscValidIntPointer(nfails,2);
618   *nfails = snes->numLinearSolveFailures;
619   PetscFunctionReturn(0);
620 }
621 
622 #undef __FUNCT__
623 #define __FUNCT__ "SNESSetMaxLinearSolveFailures"
624 /*@
625    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
626    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
627 
628    Collective on SNES
629 
630    Input Parameters:
631 +  snes     - SNES context
632 -  maxFails - maximum allowed linear solve failures
633 
634    Level: intermediate
635 
636    Notes: By default this is 0; that is SNES returns on the first failed linear solve
637 
638 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
639 
640 .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
641 @*/
642 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
643 {
644   PetscFunctionBegin;
645   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
646   snes->maxLinearSolveFailures = maxFails;
647   PetscFunctionReturn(0);
648 }
649 
650 #undef __FUNCT__
651 #define __FUNCT__ "SNESGetMaxLinearSolveFailures"
652 /*@
653    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
654      are allowed before SNES terminates
655 
656    Not Collective
657 
658    Input Parameter:
659 .  snes     - SNES context
660 
661    Output Parameter:
662 .  maxFails - maximum of unsuccessful solves allowed
663 
664    Level: intermediate
665 
666    Notes: By default this is 1; that is SNES returns on the first failed linear solve
667 
668 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
669 
670 .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
671 @*/
672 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
673 {
674   PetscFunctionBegin;
675   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
676   PetscValidIntPointer(maxFails,2);
677   *maxFails = snes->maxLinearSolveFailures;
678   PetscFunctionReturn(0);
679 }
680 
681 #undef __FUNCT__
682 #define __FUNCT__ "SNESGetLinearSolveIterations"
683 /*@
684    SNESGetLinearSolveIterations - Gets the total number of linear iterations
685    used by the nonlinear solver.
686 
687    Not Collective
688 
689    Input Parameter:
690 .  snes - SNES context
691 
692    Output Parameter:
693 .  lits - number of linear iterations
694 
695    Notes:
696    This counter is reset to zero for each successive call to SNESSolve().
697 
698    Level: intermediate
699 
700 .keywords: SNES, nonlinear, get, number, linear, iterations
701 
702 .seealso:  SNESGetIterationNumber(), SNESGetFunctionNorm()S, NESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures()
703 @*/
704 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveIterations(SNES snes,PetscInt* lits)
705 {
706   PetscFunctionBegin;
707   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
708   PetscValidIntPointer(lits,2);
709   *lits = snes->linear_its;
710   PetscFunctionReturn(0);
711 }
712 
713 #undef __FUNCT__
714 #define __FUNCT__ "SNESGetKSP"
715 /*@
716    SNESGetKSP - Returns the KSP context for a SNES solver.
717 
718    Not Collective, but if SNES object is parallel, then KSP object is parallel
719 
720    Input Parameter:
721 .  snes - the SNES context
722 
723    Output Parameter:
724 .  ksp - the KSP context
725 
726    Notes:
727    The user can then directly manipulate the KSP context to set various
728    options, etc.  Likewise, the user can then extract and manipulate the
729    PC contexts as well.
730 
731    Level: beginner
732 
733 .keywords: SNES, nonlinear, get, KSP, context
734 
735 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
736 @*/
737 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp)
738 {
739   PetscErrorCode ierr;
740 
741   PetscFunctionBegin;
742   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
743   PetscValidPointer(ksp,2);
744 
745   if (!snes->ksp) {
746     ierr = KSPCreate(((PetscObject)snes)->comm,&snes->ksp);CHKERRQ(ierr);
747     ierr = PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);CHKERRQ(ierr);
748     ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr);
749   }
750   *ksp = snes->ksp;
751   PetscFunctionReturn(0);
752 }
753 
754 #undef __FUNCT__
755 #define __FUNCT__ "SNESSetKSP"
756 /*@
757    SNESSetKSP - Sets a KSP context for the SNES object to use
758 
759    Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
760 
761    Input Parameters:
762 +  snes - the SNES context
763 -  ksp - the KSP context
764 
765    Notes:
766    The SNES object already has its KSP object, you can obtain with SNESGetKSP()
767    so this routine is rarely needed.
768 
769    The KSP object that is already in the SNES object has its reference count
770    decreased by one.
771 
772    Level: developer
773 
774 .keywords: SNES, nonlinear, get, KSP, context
775 
776 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
777 @*/
778 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetKSP(SNES snes,KSP ksp)
779 {
780   PetscErrorCode ierr;
781 
782   PetscFunctionBegin;
783   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
784   PetscValidHeaderSpecific(ksp,KSP_COOKIE,2);
785   PetscCheckSameComm(snes,1,ksp,2);
786   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
787   if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);}
788   snes->ksp = ksp;
789   PetscFunctionReturn(0);
790 }
791 
792 #if 0
793 #undef __FUNCT__
794 #define __FUNCT__ "SNESPublish_Petsc"
795 static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
796 {
797   PetscFunctionBegin;
798   PetscFunctionReturn(0);
799 }
800 #endif
801 
802 /* -----------------------------------------------------------*/
803 #undef __FUNCT__
804 #define __FUNCT__ "SNESCreate"
805 /*@
806    SNESCreate - Creates a nonlinear solver context.
807 
808    Collective on MPI_Comm
809 
810    Input Parameters:
811 .  comm - MPI communicator
812 
813    Output Parameter:
814 .  outsnes - the new SNES context
815 
816    Options Database Keys:
817 +   -snes_mf - Activates default matrix-free Jacobian-vector products,
818                and no preconditioning matrix
819 .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
820                products, and a user-provided preconditioning matrix
821                as set by SNESSetJacobian()
822 -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
823 
824    Level: beginner
825 
826 .keywords: SNES, nonlinear, create, context
827 
828 .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()
829 
830 @*/
831 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes)
832 {
833   PetscErrorCode      ierr;
834   SNES                snes;
835   SNESKSPEW           *kctx;
836 
837   PetscFunctionBegin;
838   PetscValidPointer(outsnes,2);
839   *outsnes = PETSC_NULL;
840 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
841   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
842 #endif
843 
844   ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr);
845 
846   snes->ops->converged    = SNESDefaultConverged;
847   snes->max_its           = 50;
848   snes->max_funcs	  = 10000;
849   snes->norm		  = 0.0;
850   snes->rtol		  = 1.e-8;
851   snes->ttol              = 0.0;
852   snes->abstol		  = 1.e-50;
853   snes->xtol		  = 1.e-8;
854   snes->deltatol	  = 1.e-12;
855   snes->nfuncs            = 0;
856   snes->numFailures       = 0;
857   snes->maxFailures       = 1;
858   snes->linear_its        = 0;
859   snes->lagjacobian       = 1;
860   snes->lagpreconditioner = 1;
861   snes->numbermonitors    = 0;
862   snes->data              = 0;
863   snes->setupcalled       = PETSC_FALSE;
864   snes->ksp_ewconv        = PETSC_FALSE;
865   snes->vwork             = 0;
866   snes->nwork             = 0;
867   snes->conv_hist_len     = 0;
868   snes->conv_hist_max     = 0;
869   snes->conv_hist         = PETSC_NULL;
870   snes->conv_hist_its     = PETSC_NULL;
871   snes->conv_hist_reset   = PETSC_TRUE;
872   snes->reason            = SNES_CONVERGED_ITERATING;
873 
874   snes->numLinearSolveFailures = 0;
875   snes->maxLinearSolveFailures = 1;
876 
877   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
878   ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr);
879   snes->kspconvctx  = (void*)kctx;
880   kctx->version     = 2;
881   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
882                              this was too large for some test cases */
883   kctx->rtol_last   = 0;
884   kctx->rtol_max    = .9;
885   kctx->gamma       = 1.0;
886   kctx->alpha       = .5*(1.0 + sqrt(5.0));
887   kctx->alpha2      = kctx->alpha;
888   kctx->threshold   = .1;
889   kctx->lresid_last = 0;
890   kctx->norm_last   = 0;
891 
892   *outsnes = snes;
893   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
894   PetscFunctionReturn(0);
895 }
896 
897 #undef __FUNCT__
898 #define __FUNCT__ "SNESSetFunction"
899 /*@C
900    SNESSetFunction - Sets the function evaluation routine and function
901    vector for use by the SNES routines in solving systems of nonlinear
902    equations.
903 
904    Collective on SNES
905 
906    Input Parameters:
907 +  snes - the SNES context
908 .  r - vector to store function value
909 .  func - function evaluation routine
910 -  ctx - [optional] user-defined context for private data for the
911          function evaluation routine (may be PETSC_NULL)
912 
913    Calling sequence of func:
914 $    func (SNES snes,Vec x,Vec f,void *ctx);
915 
916 .  f - function vector
917 -  ctx - optional user-defined function context
918 
919    Notes:
920    The Newton-like methods typically solve linear systems of the form
921 $      f'(x) x = -f(x),
922    where f'(x) denotes the Jacobian matrix and f(x) is the function.
923 
924    Level: beginner
925 
926 .keywords: SNES, nonlinear, set, function
927 
928 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
929 @*/
930 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
931 {
932   PetscErrorCode ierr;
933   PetscFunctionBegin;
934   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
935   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
936   PetscCheckSameComm(snes,1,r,2);
937   ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr);
938   if (snes->vec_func) { ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr); }
939   snes->ops->computefunction = func;
940   snes->vec_func             = r;
941   snes->funP                 = ctx;
942   PetscFunctionReturn(0);
943 }
944 
945 /* --------------------------------------------------------------- */
946 #undef __FUNCT__
947 #define __FUNCT__ "SNESGetRhs"
948 /*@C
949    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
950    it assumes a zero right hand side.
951 
952    Collective on SNES
953 
954    Input Parameter:
955 .  snes - the SNES context
956 
957    Output Parameter:
958 .  rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null
959 
960    Level: intermediate
961 
962 .keywords: SNES, nonlinear, get, function, right hand side
963 
964 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
965 @*/
966 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs)
967 {
968   PetscFunctionBegin;
969   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
970   PetscValidPointer(rhs,2);
971   *rhs = snes->vec_rhs;
972   PetscFunctionReturn(0);
973 }
974 
975 #undef __FUNCT__
976 #define __FUNCT__ "SNESComputeFunction"
977 /*@
978    SNESComputeFunction - Calls the function that has been set with
979                          SNESSetFunction().
980 
981    Collective on SNES
982 
983    Input Parameters:
984 +  snes - the SNES context
985 -  x - input vector
986 
987    Output Parameter:
988 .  y - function vector, as set by SNESSetFunction()
989 
990    Notes:
991    SNESComputeFunction() is typically used within nonlinear solvers
992    implementations, so most users would not generally call this routine
993    themselves.
994 
995    Level: developer
996 
997 .keywords: SNES, nonlinear, compute, function
998 
999 .seealso: SNESSetFunction(), SNESGetFunction()
1000 @*/
1001 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y)
1002 {
1003   PetscErrorCode ierr;
1004 
1005   PetscFunctionBegin;
1006   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1007   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1008   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
1009   PetscCheckSameComm(snes,1,x,2);
1010   PetscCheckSameComm(snes,1,y,3);
1011 
1012   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1013   if (snes->ops->computefunction) {
1014     PetscStackPush("SNES user function");
1015     CHKMEMQ;
1016     ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP);
1017     CHKMEMQ;
1018     PetscStackPop;
1019     if (PetscExceptionValue(ierr)) {
1020       PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr);
1021     }
1022     CHKERRQ(ierr);
1023   } else if (snes->vec_rhs) {
1024     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
1025   } else {
1026     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
1027   }
1028   if (snes->vec_rhs) {
1029     ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr);
1030   }
1031   snes->nfuncs++;
1032   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1033   PetscFunctionReturn(0);
1034 }
1035 
1036 #undef __FUNCT__
1037 #define __FUNCT__ "SNESComputeJacobian"
1038 /*@
1039    SNESComputeJacobian - Computes the Jacobian matrix that has been
1040    set with SNESSetJacobian().
1041 
1042    Collective on SNES and Mat
1043 
1044    Input Parameters:
1045 +  snes - the SNES context
1046 -  x - input vector
1047 
1048    Output Parameters:
1049 +  A - Jacobian matrix
1050 .  B - optional preconditioning matrix
1051 -  flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER)
1052 
1053   Options Database Keys:
1054 +    -snes_lag_preconditioner <lag>
1055 -    -snes_lag_jacobian <lag>
1056 
1057    Notes:
1058    Most users should not need to explicitly call this routine, as it
1059    is used internally within the nonlinear solvers.
1060 
1061    See KSPSetOperators() for important information about setting the
1062    flag parameter.
1063 
1064    Level: developer
1065 
1066 .keywords: SNES, compute, Jacobian, matrix
1067 
1068 .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
1069 @*/
1070 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
1071 {
1072   PetscErrorCode ierr;
1073   PetscTruth     flag;
1074 
1075   PetscFunctionBegin;
1076   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1077   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
1078   PetscValidPointer(flg,5);
1079   PetscCheckSameComm(snes,1,X,2);
1080   if (!snes->ops->computejacobian) PetscFunctionReturn(0);
1081 
1082   /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
1083 
1084   if (snes->lagjacobian == -2) {
1085     snes->lagjacobian = -1;
1086     ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr);
1087   } else if (snes->lagjacobian == -1) {
1088     *flg = SAME_PRECONDITIONER;
1089     ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr);
1090     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1091     if (flag) {
1092       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1093       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1094     }
1095     PetscFunctionReturn(0);
1096   } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) {
1097     *flg = SAME_PRECONDITIONER;
1098     ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr);
1099     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1100     if (flag) {
1101       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1102       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1103     }
1104     PetscFunctionReturn(0);
1105   }
1106 
1107   *flg = DIFFERENT_NONZERO_PATTERN;
1108   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1109   PetscStackPush("SNES user Jacobian function");
1110   CHKMEMQ;
1111   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1112   CHKMEMQ;
1113   PetscStackPop;
1114   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1115 
1116   if (snes->lagpreconditioner == -2) {
1117     ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr);
1118     snes->lagpreconditioner = -1;
1119   } else if (snes->lagpreconditioner == -1) {
1120     *flg = SAME_PRECONDITIONER;
1121     ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr);
1122   } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) {
1123     *flg = SAME_PRECONDITIONER;
1124     ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr);
1125   }
1126 
1127   /* make sure user returned a correct Jacobian and preconditioner */
1128   /* PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
1129     PetscValidHeaderSpecific(*B,MAT_COOKIE,4);   */
1130   PetscFunctionReturn(0);
1131 }
1132 
1133 #undef __FUNCT__
1134 #define __FUNCT__ "SNESSetJacobian"
1135 /*@C
1136    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1137    location to store the matrix.
1138 
1139    Collective on SNES and Mat
1140 
1141    Input Parameters:
1142 +  snes - the SNES context
1143 .  A - Jacobian matrix
1144 .  B - preconditioner matrix (usually same as the Jacobian)
1145 .  func - Jacobian evaluation routine
1146 -  ctx - [optional] user-defined context for private data for the
1147          Jacobian evaluation routine (may be PETSC_NULL)
1148 
1149    Calling sequence of func:
1150 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1151 
1152 +  x - input vector
1153 .  A - Jacobian matrix
1154 .  B - preconditioner matrix, usually the same as A
1155 .  flag - flag indicating information about the preconditioner matrix
1156    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1157 -  ctx - [optional] user-defined Jacobian context
1158 
1159    Notes:
1160    See KSPSetOperators() for important information about setting the flag
1161    output parameter in the routine func().  Be sure to read this information!
1162 
1163    The routine func() takes Mat * as the matrix arguments rather than Mat.
1164    This allows the Jacobian evaluation routine to replace A and/or B with a
1165    completely new new matrix structure (not just different matrix elements)
1166    when appropriate, for instance, if the nonzero structure is changing
1167    throughout the global iterations.
1168 
1169    If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on
1170    each matrix.
1171 
1172    Level: beginner
1173 
1174 .keywords: SNES, nonlinear, set, Jacobian, matrix
1175 
1176 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
1177 @*/
1178 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1179 {
1180   PetscErrorCode ierr;
1181 
1182   PetscFunctionBegin;
1183   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1184   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
1185   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
1186   if (A) PetscCheckSameComm(snes,1,A,2);
1187   if (B) PetscCheckSameComm(snes,1,B,2);
1188   if (func) snes->ops->computejacobian = func;
1189   if (ctx)  snes->jacP                 = ctx;
1190   if (A) {
1191     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1192     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1193     snes->jacobian = A;
1194   }
1195   if (B) {
1196     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
1197     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1198     snes->jacobian_pre = B;
1199   }
1200   PetscFunctionReturn(0);
1201 }
1202 
1203 #undef __FUNCT__
1204 #define __FUNCT__ "SNESGetJacobian"
1205 /*@C
1206    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1207    provided context for evaluating the Jacobian.
1208 
1209    Not Collective, but Mat object will be parallel if SNES object is
1210 
1211    Input Parameter:
1212 .  snes - the nonlinear solver context
1213 
1214    Output Parameters:
1215 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1216 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1217 .  func - location to put Jacobian function (or PETSC_NULL)
1218 -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1219 
1220    Level: advanced
1221 
1222 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1223 @*/
1224 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1225 {
1226   PetscFunctionBegin;
1227   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1228   if (A)    *A    = snes->jacobian;
1229   if (B)    *B    = snes->jacobian_pre;
1230   if (func) *func = snes->ops->computejacobian;
1231   if (ctx)  *ctx  = snes->jacP;
1232   PetscFunctionReturn(0);
1233 }
1234 
1235 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1236 EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
1237 
1238 #undef __FUNCT__
1239 #define __FUNCT__ "SNESSetUp"
1240 /*@
1241    SNESSetUp - Sets up the internal data structures for the later use
1242    of a nonlinear solver.
1243 
1244    Collective on SNES
1245 
1246    Input Parameters:
1247 .  snes - the SNES context
1248 
1249    Notes:
1250    For basic use of the SNES solvers the user need not explicitly call
1251    SNESSetUp(), since these actions will automatically occur during
1252    the call to SNESSolve().  However, if one wishes to control this
1253    phase separately, SNESSetUp() should be called after SNESCreate()
1254    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1255 
1256    Level: advanced
1257 
1258 .keywords: SNES, nonlinear, setup
1259 
1260 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1261 @*/
1262 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
1263 {
1264   PetscErrorCode ierr;
1265   PetscTruth     flg;
1266 
1267   PetscFunctionBegin;
1268   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1269   if (snes->setupcalled) PetscFunctionReturn(0);
1270 
1271   if (!((PetscObject)snes)->type_name) {
1272     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
1273   }
1274 
1275   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1276   /*
1277       This version replaces the user provided Jacobian matrix with a
1278       matrix-free version but still employs the user-provided preconditioner matrix
1279   */
1280   if (flg) {
1281     Mat J;
1282     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
1283     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1284     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
1285     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
1286     ierr = MatDestroy(J);CHKERRQ(ierr);
1287   }
1288 
1289 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT)
1290   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
1291   if (flg) {
1292     Mat J;
1293     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1294     ierr = PetscInfo(snes,"Setting default matrix-free operator routines (version 2)\n");CHKERRQ(ierr);
1295     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
1296     ierr = MatDestroy(J);CHKERRQ(ierr);
1297   }
1298 #endif
1299 
1300   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1301   /*
1302       This version replaces both the user-provided Jacobian and the user-
1303       provided preconditioner matrix with the default matrix free version.
1304    */
1305   if (flg) {
1306     Mat  J;
1307     KSP ksp;
1308     PC   pc;
1309     /* create and set matrix-free operator */
1310     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
1311     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1312     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
1313     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
1314     ierr = MatDestroy(J);CHKERRQ(ierr);
1315     /* force no preconditioner */
1316     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1317     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1318     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1319     if (!flg) {
1320       ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines;\nThat is no preconditioner is being used\n");CHKERRQ(ierr);
1321       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1322     }
1323   }
1324 
1325   if (!snes->vec_func && !snes->vec_rhs) {
1326     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1327   }
1328   if (!snes->ops->computefunction && !snes->vec_rhs) {
1329     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1330   }
1331   if (!snes->jacobian) {
1332     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1333   }
1334   if (snes->vec_func == snes->vec_sol) {
1335     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1336   }
1337 
1338   if (snes->ops->setup) {
1339     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
1340   }
1341   snes->setupcalled = PETSC_TRUE;
1342   PetscFunctionReturn(0);
1343 }
1344 
1345 #undef __FUNCT__
1346 #define __FUNCT__ "SNESDestroy"
1347 /*@
1348    SNESDestroy - Destroys the nonlinear solver context that was created
1349    with SNESCreate().
1350 
1351    Collective on SNES
1352 
1353    Input Parameter:
1354 .  snes - the SNES context
1355 
1356    Level: beginner
1357 
1358 .keywords: SNES, nonlinear, destroy
1359 
1360 .seealso: SNESCreate(), SNESSolve()
1361 @*/
1362 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
1363 {
1364   PetscErrorCode ierr;
1365 
1366   PetscFunctionBegin;
1367   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1368   if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0);
1369 
1370   /* if memory was published with AMS then destroy it */
1371   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1372   if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);}
1373 
1374   if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);}
1375   if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);}
1376   if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);}
1377   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1378   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1379   if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);}
1380   ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);
1381   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1382   ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
1383   if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);}
1384   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
1385   PetscFunctionReturn(0);
1386 }
1387 
1388 /* ----------- Routines to set solver parameters ---------- */
1389 
1390 #undef __FUNCT__
1391 #define __FUNCT__ "SNESSetLagPreconditioner"
1392 /*@
1393    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
1394 
1395    Collective on SNES
1396 
1397    Input Parameters:
1398 +  snes - the SNES context
1399 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1400          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1401 
1402    Options Database Keys:
1403 .    -snes_lag_preconditioner <lag>
1404 
1405    Notes:
1406    The default is 1
1407    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1408    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
1409 
1410    Level: intermediate
1411 
1412 .keywords: SNES, nonlinear, set, convergence, tolerances
1413 
1414 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
1415 
1416 @*/
1417 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag)
1418 {
1419   PetscFunctionBegin;
1420   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1421   if (lag < -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1422   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1423   snes->lagpreconditioner = lag;
1424   PetscFunctionReturn(0);
1425 }
1426 
1427 #undef __FUNCT__
1428 #define __FUNCT__ "SNESGetLagPreconditioner"
1429 /*@
1430    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
1431 
1432    Collective on SNES
1433 
1434    Input Parameter:
1435 .  snes - the SNES context
1436 
1437    Output Parameter:
1438 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1439          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1440 
1441    Options Database Keys:
1442 .    -snes_lag_preconditioner <lag>
1443 
1444    Notes:
1445    The default is 1
1446    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1447 
1448    Level: intermediate
1449 
1450 .keywords: SNES, nonlinear, set, convergence, tolerances
1451 
1452 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
1453 
1454 @*/
1455 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
1456 {
1457   PetscFunctionBegin;
1458   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1459   *lag = snes->lagpreconditioner;
1460   PetscFunctionReturn(0);
1461 }
1462 
1463 #undef __FUNCT__
1464 #define __FUNCT__ "SNESSetLagJacobian"
1465 /*@
1466    SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
1467      often the preconditioner is rebuilt.
1468 
1469    Collective on SNES
1470 
1471    Input Parameters:
1472 +  snes - the SNES context
1473 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1474          the Jacobian is built etc. -2 means rebuild at next chance but then never again
1475 
1476    Options Database Keys:
1477 .    -snes_lag_jacobian <lag>
1478 
1479    Notes:
1480    The default is 1
1481    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1482    If  -1 is used before the very first nonlinear solve the CODE WILL FAIL! because no Jacobian is used, use -2 to indicate you want it recomputed
1483    at the next Newton step but never again (unless it is reset to another value)
1484 
1485    Level: intermediate
1486 
1487 .keywords: SNES, nonlinear, set, convergence, tolerances
1488 
1489 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
1490 
1491 @*/
1492 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagJacobian(SNES snes,PetscInt lag)
1493 {
1494   PetscFunctionBegin;
1495   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1496   if (lag < -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1497   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1498   snes->lagjacobian = lag;
1499   PetscFunctionReturn(0);
1500 }
1501 
1502 #undef __FUNCT__
1503 #define __FUNCT__ "SNESGetLagJacobian"
1504 /*@
1505    SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
1506 
1507    Collective on SNES
1508 
1509    Input Parameter:
1510 .  snes - the SNES context
1511 
1512    Output Parameter:
1513 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1514          the Jacobian is built etc.
1515 
1516    Options Database Keys:
1517 .    -snes_lag_jacobian <lag>
1518 
1519    Notes:
1520    The default is 1
1521    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1522 
1523    Level: intermediate
1524 
1525 .keywords: SNES, nonlinear, set, convergence, tolerances
1526 
1527 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
1528 
1529 @*/
1530 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagJacobian(SNES snes,PetscInt *lag)
1531 {
1532   PetscFunctionBegin;
1533   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1534   *lag = snes->lagjacobian;
1535   PetscFunctionReturn(0);
1536 }
1537 
1538 #undef __FUNCT__
1539 #define __FUNCT__ "SNESSetTolerances"
1540 /*@
1541    SNESSetTolerances - Sets various parameters used in convergence tests.
1542 
1543    Collective on SNES
1544 
1545    Input Parameters:
1546 +  snes - the SNES context
1547 .  abstol - absolute convergence tolerance
1548 .  rtol - relative convergence tolerance
1549 .  stol -  convergence tolerance in terms of the norm
1550            of the change in the solution between steps
1551 .  maxit - maximum number of iterations
1552 -  maxf - maximum number of function evaluations
1553 
1554    Options Database Keys:
1555 +    -snes_atol <abstol> - Sets abstol
1556 .    -snes_rtol <rtol> - Sets rtol
1557 .    -snes_stol <stol> - Sets stol
1558 .    -snes_max_it <maxit> - Sets maxit
1559 -    -snes_max_funcs <maxf> - Sets maxf
1560 
1561    Notes:
1562    The default maximum number of iterations is 50.
1563    The default maximum number of function evaluations is 1000.
1564 
1565    Level: intermediate
1566 
1567 .keywords: SNES, nonlinear, set, convergence, tolerances
1568 
1569 .seealso: SNESSetTrustRegionTolerance()
1570 @*/
1571 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
1572 {
1573   PetscFunctionBegin;
1574   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1575   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1576   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1577   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1578   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1579   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1580   PetscFunctionReturn(0);
1581 }
1582 
1583 #undef __FUNCT__
1584 #define __FUNCT__ "SNESGetTolerances"
1585 /*@
1586    SNESGetTolerances - Gets various parameters used in convergence tests.
1587 
1588    Not Collective
1589 
1590    Input Parameters:
1591 +  snes - the SNES context
1592 .  atol - absolute convergence tolerance
1593 .  rtol - relative convergence tolerance
1594 .  stol -  convergence tolerance in terms of the norm
1595            of the change in the solution between steps
1596 .  maxit - maximum number of iterations
1597 -  maxf - maximum number of function evaluations
1598 
1599    Notes:
1600    The user can specify PETSC_NULL for any parameter that is not needed.
1601 
1602    Level: intermediate
1603 
1604 .keywords: SNES, nonlinear, get, convergence, tolerances
1605 
1606 .seealso: SNESSetTolerances()
1607 @*/
1608 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
1609 {
1610   PetscFunctionBegin;
1611   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1612   if (atol)  *atol  = snes->abstol;
1613   if (rtol)  *rtol  = snes->rtol;
1614   if (stol)  *stol  = snes->xtol;
1615   if (maxit) *maxit = snes->max_its;
1616   if (maxf)  *maxf  = snes->max_funcs;
1617   PetscFunctionReturn(0);
1618 }
1619 
1620 #undef __FUNCT__
1621 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1622 /*@
1623    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1624 
1625    Collective on SNES
1626 
1627    Input Parameters:
1628 +  snes - the SNES context
1629 -  tol - tolerance
1630 
1631    Options Database Key:
1632 .  -snes_trtol <tol> - Sets tol
1633 
1634    Level: intermediate
1635 
1636 .keywords: SNES, nonlinear, set, trust region, tolerance
1637 
1638 .seealso: SNESSetTolerances()
1639 @*/
1640 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1641 {
1642   PetscFunctionBegin;
1643   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1644   snes->deltatol = tol;
1645   PetscFunctionReturn(0);
1646 }
1647 
1648 /*
1649    Duplicate the lg monitors for SNES from KSP; for some reason with
1650    dynamic libraries things don't work under Sun4 if we just use
1651    macros instead of functions
1652 */
1653 #undef __FUNCT__
1654 #define __FUNCT__ "SNESMonitorLG"
1655 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1656 {
1657   PetscErrorCode ierr;
1658 
1659   PetscFunctionBegin;
1660   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1661   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1662   PetscFunctionReturn(0);
1663 }
1664 
1665 #undef __FUNCT__
1666 #define __FUNCT__ "SNESMonitorLGCreate"
1667 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1668 {
1669   PetscErrorCode ierr;
1670 
1671   PetscFunctionBegin;
1672   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1673   PetscFunctionReturn(0);
1674 }
1675 
1676 #undef __FUNCT__
1677 #define __FUNCT__ "SNESMonitorLGDestroy"
1678 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw)
1679 {
1680   PetscErrorCode ierr;
1681 
1682   PetscFunctionBegin;
1683   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1684   PetscFunctionReturn(0);
1685 }
1686 
1687 extern PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
1688 #undef __FUNCT__
1689 #define __FUNCT__ "SNESMonitorLGRange"
1690 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
1691 {
1692   PetscDrawLG      lg;
1693   PetscErrorCode   ierr;
1694   PetscReal        x,y,per;
1695   PetscViewer      v = (PetscViewer)monctx;
1696   static PetscReal prev; /* should be in the context */
1697   PetscDraw        draw;
1698   PetscFunctionBegin;
1699   if (!monctx) {
1700     MPI_Comm    comm;
1701 
1702     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1703     v      = PETSC_VIEWER_DRAW_(comm);
1704   }
1705   ierr   = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr);
1706   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1707   ierr   = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1708   ierr   = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr);
1709   x = (PetscReal) n;
1710   if (rnorm > 0.0) y = log10(rnorm); else y = -15.0;
1711   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1712   if (n < 20 || !(n % 5)) {
1713     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1714   }
1715 
1716   ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr);
1717   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1718   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1719   ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr);
1720   ierr =  SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr);
1721   x = (PetscReal) n;
1722   y = 100.0*per;
1723   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1724   if (n < 20 || !(n % 5)) {
1725     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1726   }
1727 
1728   ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr);
1729   if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1730   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1731   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr);
1732   x = (PetscReal) n;
1733   y = (prev - rnorm)/prev;
1734   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1735   if (n < 20 || !(n % 5)) {
1736     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1737   }
1738 
1739   ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr);
1740   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1741   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1742   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
1743   x = (PetscReal) n;
1744   y = (prev - rnorm)/(prev*per);
1745   if (n > 2) { /*skip initial crazy value */
1746     ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1747   }
1748   if (n < 20 || !(n % 5)) {
1749     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1750   }
1751   prev = rnorm;
1752   PetscFunctionReturn(0);
1753 }
1754 
1755 #undef __FUNCT__
1756 #define __FUNCT__ "SNESMonitorLGRangeCreate"
1757 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1758 {
1759   PetscErrorCode ierr;
1760 
1761   PetscFunctionBegin;
1762   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1763   PetscFunctionReturn(0);
1764 }
1765 
1766 #undef __FUNCT__
1767 #define __FUNCT__ "SNESMonitorLGRangeDestroy"
1768 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeDestroy(PetscDrawLG draw)
1769 {
1770   PetscErrorCode ierr;
1771 
1772   PetscFunctionBegin;
1773   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1774   PetscFunctionReturn(0);
1775 }
1776 
1777 /* ------------ Routines to set performance monitoring options ----------- */
1778 
1779 #undef __FUNCT__
1780 #define __FUNCT__ "SNESMonitorSet"
1781 /*@C
1782    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
1783    iteration of the nonlinear solver to display the iteration's
1784    progress.
1785 
1786    Collective on SNES
1787 
1788    Input Parameters:
1789 +  snes - the SNES context
1790 .  func - monitoring routine
1791 .  mctx - [optional] user-defined context for private data for the
1792           monitor routine (use PETSC_NULL if no context is desired)
1793 -  monitordestroy - [optional] routine that frees monitor context
1794           (may be PETSC_NULL)
1795 
1796    Calling sequence of func:
1797 $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1798 
1799 +    snes - the SNES context
1800 .    its - iteration number
1801 .    norm - 2-norm function value (may be estimated)
1802 -    mctx - [optional] monitoring context
1803 
1804    Options Database Keys:
1805 +    -snes_monitor        - sets SNESMonitorDefault()
1806 .    -snes_monitor_draw    - sets line graph monitor,
1807                             uses SNESMonitorLGCreate()
1808 _    -snes_monitor_cancel - cancels all monitors that have
1809                             been hardwired into a code by
1810                             calls to SNESMonitorSet(), but
1811                             does not cancel those set via
1812                             the options database.
1813 
1814    Notes:
1815    Several different monitoring routines may be set by calling
1816    SNESMonitorSet() multiple times; all will be called in the
1817    order in which they were set.
1818 
1819    Fortran notes: Only a single monitor function can be set for each SNES object
1820 
1821    Level: intermediate
1822 
1823 .keywords: SNES, nonlinear, set, monitor
1824 
1825 .seealso: SNESMonitorDefault(), SNESMonitorCancel()
1826 @*/
1827 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
1828 {
1829   PetscInt i;
1830 
1831   PetscFunctionBegin;
1832   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1833   if (snes->numbermonitors >= MAXSNESMONITORS) {
1834     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1835   }
1836   for (i=0; i<snes->numbermonitors;i++) {
1837     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0);
1838 
1839     /* check if both default monitors that share common ASCII viewer */
1840     if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) {
1841       if (mctx && snes->monitorcontext[i]) {
1842         PetscErrorCode          ierr;
1843         PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx;
1844         PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i];
1845         if (viewer1->viewer == viewer2->viewer) {
1846           ierr = (*monitordestroy)(mctx);CHKERRQ(ierr);
1847           PetscFunctionReturn(0);
1848         }
1849       }
1850     }
1851   }
1852   snes->monitor[snes->numbermonitors]           = monitor;
1853   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1854   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1855   PetscFunctionReturn(0);
1856 }
1857 
1858 #undef __FUNCT__
1859 #define __FUNCT__ "SNESMonitorCancel"
1860 /*@C
1861    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
1862 
1863    Collective on SNES
1864 
1865    Input Parameters:
1866 .  snes - the SNES context
1867 
1868    Options Database Key:
1869 .  -snes_monitor_cancel - cancels all monitors that have been hardwired
1870     into a code by calls to SNESMonitorSet(), but does not cancel those
1871     set via the options database
1872 
1873    Notes:
1874    There is no way to clear one specific monitor from a SNES object.
1875 
1876    Level: intermediate
1877 
1878 .keywords: SNES, nonlinear, set, monitor
1879 
1880 .seealso: SNESMonitorDefault(), SNESMonitorSet()
1881 @*/
1882 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes)
1883 {
1884   PetscErrorCode ierr;
1885   PetscInt       i;
1886 
1887   PetscFunctionBegin;
1888   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1889   for (i=0; i<snes->numbermonitors; i++) {
1890     if (snes->monitordestroy[i]) {
1891       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1892     }
1893   }
1894   snes->numbermonitors = 0;
1895   PetscFunctionReturn(0);
1896 }
1897 
1898 #undef __FUNCT__
1899 #define __FUNCT__ "SNESSetConvergenceTest"
1900 /*@C
1901    SNESSetConvergenceTest - Sets the function that is to be used
1902    to test for convergence of the nonlinear iterative solution.
1903 
1904    Collective on SNES
1905 
1906    Input Parameters:
1907 +  snes - the SNES context
1908 .  func - routine to test for convergence
1909 .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
1910 -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
1911 
1912    Calling sequence of func:
1913 $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1914 
1915 +    snes - the SNES context
1916 .    it - current iteration (0 is the first and is before any Newton step)
1917 .    cctx - [optional] convergence context
1918 .    reason - reason for convergence/divergence
1919 .    xnorm - 2-norm of current iterate
1920 .    gnorm - 2-norm of current step
1921 -    f - 2-norm of function
1922 
1923    Level: advanced
1924 
1925 .keywords: SNES, nonlinear, set, convergence, test
1926 
1927 .seealso: SNESDefaultConverged(), SNESSkipConverged()
1928 @*/
1929 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1930 {
1931   PetscErrorCode ierr;
1932 
1933   PetscFunctionBegin;
1934   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1935   if (!func) func = SNESSkipConverged;
1936   if (snes->ops->convergeddestroy) {
1937     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
1938   }
1939   snes->ops->converged        = func;
1940   snes->ops->convergeddestroy = destroy;
1941   snes->cnvP                  = cctx;
1942   PetscFunctionReturn(0);
1943 }
1944 
1945 #undef __FUNCT__
1946 #define __FUNCT__ "SNESGetConvergedReason"
1947 /*@
1948    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1949 
1950    Not Collective
1951 
1952    Input Parameter:
1953 .  snes - the SNES context
1954 
1955    Output Parameter:
1956 .  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
1957             manual pages for the individual convergence tests for complete lists
1958 
1959    Level: intermediate
1960 
1961    Notes: Can only be called after the call the SNESSolve() is complete.
1962 
1963 .keywords: SNES, nonlinear, set, convergence, test
1964 
1965 .seealso: SNESSetConvergenceTest(), SNESConvergedReason
1966 @*/
1967 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1968 {
1969   PetscFunctionBegin;
1970   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1971   PetscValidPointer(reason,2);
1972   *reason = snes->reason;
1973   PetscFunctionReturn(0);
1974 }
1975 
1976 #undef __FUNCT__
1977 #define __FUNCT__ "SNESSetConvergenceHistory"
1978 /*@
1979    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1980 
1981    Collective on SNES
1982 
1983    Input Parameters:
1984 +  snes - iterative context obtained from SNESCreate()
1985 .  a   - array to hold history
1986 .  its - integer array holds the number of linear iterations for each solve.
1987 .  na  - size of a and its
1988 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1989            else it continues storing new values for new nonlinear solves after the old ones
1990 
1991    Notes:
1992    If set, this array will contain the function norms computed
1993    at each step.
1994 
1995    This routine is useful, e.g., when running a code for purposes
1996    of accurate performance monitoring, when no I/O should be done
1997    during the section of code that is being timed.
1998 
1999    Level: intermediate
2000 
2001 .keywords: SNES, set, convergence, history
2002 
2003 .seealso: SNESGetConvergenceHistory()
2004 
2005 @*/
2006 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset)
2007 {
2008   PetscFunctionBegin;
2009   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2010   if (na)  PetscValidScalarPointer(a,2);
2011   if (its) PetscValidIntPointer(its,3);
2012   snes->conv_hist       = a;
2013   snes->conv_hist_its   = its;
2014   snes->conv_hist_max   = na;
2015   snes->conv_hist_len   = 0;
2016   snes->conv_hist_reset = reset;
2017   PetscFunctionReturn(0);
2018 }
2019 
2020 #undef __FUNCT__
2021 #define __FUNCT__ "SNESGetConvergenceHistory"
2022 /*@C
2023    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
2024 
2025    Collective on SNES
2026 
2027    Input Parameter:
2028 .  snes - iterative context obtained from SNESCreate()
2029 
2030    Output Parameters:
2031 .  a   - array to hold history
2032 .  its - integer array holds the number of linear iterations (or
2033          negative if not converged) for each solve.
2034 -  na  - size of a and its
2035 
2036    Notes:
2037     The calling sequence for this routine in Fortran is
2038 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
2039 
2040    This routine is useful, e.g., when running a code for purposes
2041    of accurate performance monitoring, when no I/O should be done
2042    during the section of code that is being timed.
2043 
2044    Level: intermediate
2045 
2046 .keywords: SNES, get, convergence, history
2047 
2048 .seealso: SNESSetConvergencHistory()
2049 
2050 @*/
2051 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
2052 {
2053   PetscFunctionBegin;
2054   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2055   if (a)   *a   = snes->conv_hist;
2056   if (its) *its = snes->conv_hist_its;
2057   if (na)  *na  = snes->conv_hist_len;
2058   PetscFunctionReturn(0);
2059 }
2060 
2061 #undef __FUNCT__
2062 #define __FUNCT__ "SNESSetUpdate"
2063 /*@C
2064   SNESSetUpdate - Sets the general-purpose update function called
2065   at the beginning o every iteration of the nonlinear solve. Specifically
2066   it is called just before the Jacobian is "evaluated".
2067 
2068   Collective on SNES
2069 
2070   Input Parameters:
2071 . snes - The nonlinear solver context
2072 . func - The function
2073 
2074   Calling sequence of func:
2075 . func (SNES snes, PetscInt step);
2076 
2077 . step - The current step of the iteration
2078 
2079   Level: intermediate
2080 
2081 .keywords: SNES, update
2082 
2083 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
2084 @*/
2085 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
2086 {
2087   PetscFunctionBegin;
2088   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
2089   snes->ops->update = func;
2090   PetscFunctionReturn(0);
2091 }
2092 
2093 #undef __FUNCT__
2094 #define __FUNCT__ "SNESDefaultUpdate"
2095 /*@
2096   SNESDefaultUpdate - The default update function which does nothing.
2097 
2098   Not collective
2099 
2100   Input Parameters:
2101 . snes - The nonlinear solver context
2102 . step - The current step of the iteration
2103 
2104   Level: intermediate
2105 
2106 .keywords: SNES, update
2107 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
2108 @*/
2109 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
2110 {
2111   PetscFunctionBegin;
2112   PetscFunctionReturn(0);
2113 }
2114 
2115 #undef __FUNCT__
2116 #define __FUNCT__ "SNESScaleStep_Private"
2117 /*
2118    SNESScaleStep_Private - Scales a step so that its length is less than the
2119    positive parameter delta.
2120 
2121     Input Parameters:
2122 +   snes - the SNES context
2123 .   y - approximate solution of linear system
2124 .   fnorm - 2-norm of current function
2125 -   delta - trust region size
2126 
2127     Output Parameters:
2128 +   gpnorm - predicted function norm at the new point, assuming local
2129     linearization.  The value is zero if the step lies within the trust
2130     region, and exceeds zero otherwise.
2131 -   ynorm - 2-norm of the step
2132 
2133     Note:
2134     For non-trust region methods such as SNESLS, the parameter delta
2135     is set to be the maximum allowable step size.
2136 
2137 .keywords: SNES, nonlinear, scale, step
2138 */
2139 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
2140 {
2141   PetscReal      nrm;
2142   PetscScalar    cnorm;
2143   PetscErrorCode ierr;
2144 
2145   PetscFunctionBegin;
2146   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2147   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
2148   PetscCheckSameComm(snes,1,y,2);
2149 
2150   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
2151   if (nrm > *delta) {
2152      nrm = *delta/nrm;
2153      *gpnorm = (1.0 - nrm)*(*fnorm);
2154      cnorm = nrm;
2155      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
2156      *ynorm = *delta;
2157   } else {
2158      *gpnorm = 0.0;
2159      *ynorm = nrm;
2160   }
2161   PetscFunctionReturn(0);
2162 }
2163 
2164 #undef __FUNCT__
2165 #define __FUNCT__ "SNESSolve"
2166 /*@C
2167    SNESSolve - Solves a nonlinear system F(x) = b.
2168    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
2169 
2170    Collective on SNES
2171 
2172    Input Parameters:
2173 +  snes - the SNES context
2174 .  b - the constant part of the equation, or PETSC_NULL to use zero.
2175 -  x - the solution vector.
2176 
2177    Notes:
2178    The user should initialize the vector,x, with the initial guess
2179    for the nonlinear solve prior to calling SNESSolve.  In particular,
2180    to employ an initial guess of zero, the user should explicitly set
2181    this vector to zero by calling VecSet().
2182 
2183    Level: beginner
2184 
2185 .keywords: SNES, nonlinear, solve
2186 
2187 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian()
2188 @*/
2189 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
2190 {
2191   PetscErrorCode ierr;
2192   PetscTruth     flg;
2193   char           filename[PETSC_MAX_PATH_LEN];
2194   PetscViewer    viewer;
2195 
2196   PetscFunctionBegin;
2197   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2198   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
2199   PetscCheckSameComm(snes,1,x,3);
2200   if (b) PetscValidHeaderSpecific(b,VEC_COOKIE,2);
2201   if (b) PetscCheckSameComm(snes,1,b,2);
2202 
2203   /* set solution vector */
2204   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
2205   if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); }
2206   snes->vec_sol = x;
2207   /* set afine vector if provided */
2208   if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
2209   if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); }
2210   snes->vec_rhs = b;
2211 
2212   if (!snes->vec_func && snes->vec_rhs) {
2213     ierr = VecDuplicate(b, &snes->vec_func); CHKERRQ(ierr);
2214   }
2215 
2216   ierr = SNESSetUp(snes);CHKERRQ(ierr);
2217 
2218   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
2219   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
2220 
2221   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2222   ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
2223   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2224   if (snes->domainerror){
2225     snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
2226     snes->domainerror = PETSC_FALSE;
2227   }
2228 
2229   if (!snes->reason) {
2230     SETERRQ(PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
2231   }
2232 
2233   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
2234   if (flg && !PetscPreLoadingOn) {
2235     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
2236     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
2237     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
2238   }
2239 
2240   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
2241   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
2242   if (snes->printreason) {
2243     if (snes->reason > 0) {
2244       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2245     } else {
2246       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2247     }
2248   }
2249 
2250   PetscFunctionReturn(0);
2251 }
2252 
2253 /* --------- Internal routines for SNES Package --------- */
2254 
2255 #undef __FUNCT__
2256 #define __FUNCT__ "SNESSetType"
2257 /*@C
2258    SNESSetType - Sets the method for the nonlinear solver.
2259 
2260    Collective on SNES
2261 
2262    Input Parameters:
2263 +  snes - the SNES context
2264 -  type - a known method
2265 
2266    Options Database Key:
2267 .  -snes_type <type> - Sets the method; use -help for a list
2268    of available methods (for instance, ls or tr)
2269 
2270    Notes:
2271    See "petsc/include/petscsnes.h" for available methods (for instance)
2272 +    SNESLS - Newton's method with line search
2273      (systems of nonlinear equations)
2274 .    SNESTR - Newton's method with trust region
2275      (systems of nonlinear equations)
2276 
2277   Normally, it is best to use the SNESSetFromOptions() command and then
2278   set the SNES solver type from the options database rather than by using
2279   this routine.  Using the options database provides the user with
2280   maximum flexibility in evaluating the many nonlinear solvers.
2281   The SNESSetType() routine is provided for those situations where it
2282   is necessary to set the nonlinear solver independently of the command
2283   line or options database.  This might be the case, for example, when
2284   the choice of solver changes during the execution of the program,
2285   and the user's application is taking responsibility for choosing the
2286   appropriate method.
2287 
2288   Level: intermediate
2289 
2290 .keywords: SNES, set, type
2291 
2292 .seealso: SNESType, SNESCreate()
2293 
2294 @*/
2295 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type)
2296 {
2297   PetscErrorCode ierr,(*r)(SNES);
2298   PetscTruth     match;
2299 
2300   PetscFunctionBegin;
2301   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2302   PetscValidCharPointer(type,2);
2303 
2304   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
2305   if (match) PetscFunctionReturn(0);
2306 
2307   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
2308   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
2309   /* Destroy the previous private SNES context */
2310   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
2311   /* Reinitialize function pointers in SNESOps structure */
2312   snes->ops->setup          = 0;
2313   snes->ops->solve          = 0;
2314   snes->ops->view           = 0;
2315   snes->ops->setfromoptions = 0;
2316   snes->ops->destroy        = 0;
2317   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
2318   snes->setupcalled = PETSC_FALSE;
2319   ierr = (*r)(snes);CHKERRQ(ierr);
2320   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
2321   PetscFunctionReturn(0);
2322 }
2323 
2324 
2325 /* --------------------------------------------------------------------- */
2326 #undef __FUNCT__
2327 #define __FUNCT__ "SNESRegisterDestroy"
2328 /*@
2329    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2330    registered by SNESRegisterDynamic().
2331 
2332    Not Collective
2333 
2334    Level: advanced
2335 
2336 .keywords: SNES, nonlinear, register, destroy
2337 
2338 .seealso: SNESRegisterAll(), SNESRegisterAll()
2339 @*/
2340 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
2341 {
2342   PetscErrorCode ierr;
2343 
2344   PetscFunctionBegin;
2345   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
2346   SNESRegisterAllCalled = PETSC_FALSE;
2347   PetscFunctionReturn(0);
2348 }
2349 
2350 #undef __FUNCT__
2351 #define __FUNCT__ "SNESGetType"
2352 /*@C
2353    SNESGetType - Gets the SNES method type and name (as a string).
2354 
2355    Not Collective
2356 
2357    Input Parameter:
2358 .  snes - nonlinear solver context
2359 
2360    Output Parameter:
2361 .  type - SNES method (a character string)
2362 
2363    Level: intermediate
2364 
2365 .keywords: SNES, nonlinear, get, type, name
2366 @*/
2367 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type)
2368 {
2369   PetscFunctionBegin;
2370   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2371   PetscValidPointer(type,2);
2372   *type = ((PetscObject)snes)->type_name;
2373   PetscFunctionReturn(0);
2374 }
2375 
2376 #undef __FUNCT__
2377 #define __FUNCT__ "SNESGetSolution"
2378 /*@
2379    SNESGetSolution - Returns the vector where the approximate solution is
2380    stored.
2381 
2382    Not Collective, but Vec is parallel if SNES is parallel
2383 
2384    Input Parameter:
2385 .  snes - the SNES context
2386 
2387    Output Parameter:
2388 .  x - the solution
2389 
2390    Level: intermediate
2391 
2392 .keywords: SNES, nonlinear, get, solution
2393 
2394 .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
2395 @*/
2396 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
2397 {
2398   PetscFunctionBegin;
2399   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2400   PetscValidPointer(x,2);
2401   *x = snes->vec_sol;
2402   PetscFunctionReturn(0);
2403 }
2404 
2405 #undef __FUNCT__
2406 #define __FUNCT__ "SNESGetSolutionUpdate"
2407 /*@
2408    SNESGetSolutionUpdate - Returns the vector where the solution update is
2409    stored.
2410 
2411    Not Collective, but Vec is parallel if SNES is parallel
2412 
2413    Input Parameter:
2414 .  snes - the SNES context
2415 
2416    Output Parameter:
2417 .  x - the solution update
2418 
2419    Level: advanced
2420 
2421 .keywords: SNES, nonlinear, get, solution, update
2422 
2423 .seealso: SNESGetSolution(), SNESGetFunction()
2424 @*/
2425 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
2426 {
2427   PetscFunctionBegin;
2428   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2429   PetscValidPointer(x,2);
2430   *x = snes->vec_sol_update;
2431   PetscFunctionReturn(0);
2432 }
2433 
2434 #undef __FUNCT__
2435 #define __FUNCT__ "SNESGetFunction"
2436 /*@C
2437    SNESGetFunction - Returns the vector where the function is stored.
2438 
2439    Not Collective, but Vec is parallel if SNES is parallel
2440 
2441    Input Parameter:
2442 .  snes - the SNES context
2443 
2444    Output Parameter:
2445 +  r - the function (or PETSC_NULL)
2446 .  func - the function (or PETSC_NULL)
2447 -  ctx - the function context (or PETSC_NULL)
2448 
2449    Level: advanced
2450 
2451 .keywords: SNES, nonlinear, get, function
2452 
2453 .seealso: SNESSetFunction(), SNESGetSolution()
2454 @*/
2455 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
2456 {
2457   PetscFunctionBegin;
2458   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2459   if (r)    *r    = snes->vec_func;
2460   if (func) *func = snes->ops->computefunction;
2461   if (ctx)  *ctx  = snes->funP;
2462   PetscFunctionReturn(0);
2463 }
2464 
2465 #undef __FUNCT__
2466 #define __FUNCT__ "SNESSetOptionsPrefix"
2467 /*@C
2468    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2469    SNES options in the database.
2470 
2471    Collective on SNES
2472 
2473    Input Parameter:
2474 +  snes - the SNES context
2475 -  prefix - the prefix to prepend to all option names
2476 
2477    Notes:
2478    A hyphen (-) must NOT be given at the beginning of the prefix name.
2479    The first character of all runtime options is AUTOMATICALLY the hyphen.
2480 
2481    Level: advanced
2482 
2483 .keywords: SNES, set, options, prefix, database
2484 
2485 .seealso: SNESSetFromOptions()
2486 @*/
2487 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
2488 {
2489   PetscErrorCode ierr;
2490 
2491   PetscFunctionBegin;
2492   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2493   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2494   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2495   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2496   PetscFunctionReturn(0);
2497 }
2498 
2499 #undef __FUNCT__
2500 #define __FUNCT__ "SNESAppendOptionsPrefix"
2501 /*@C
2502    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2503    SNES options in the database.
2504 
2505    Collective on SNES
2506 
2507    Input Parameters:
2508 +  snes - the SNES context
2509 -  prefix - the prefix to prepend to all option names
2510 
2511    Notes:
2512    A hyphen (-) must NOT be given at the beginning of the prefix name.
2513    The first character of all runtime options is AUTOMATICALLY the hyphen.
2514 
2515    Level: advanced
2516 
2517 .keywords: SNES, append, options, prefix, database
2518 
2519 .seealso: SNESGetOptionsPrefix()
2520 @*/
2521 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
2522 {
2523   PetscErrorCode ierr;
2524 
2525   PetscFunctionBegin;
2526   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2527   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2528   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2529   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2530   PetscFunctionReturn(0);
2531 }
2532 
2533 #undef __FUNCT__
2534 #define __FUNCT__ "SNESGetOptionsPrefix"
2535 /*@C
2536    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2537    SNES options in the database.
2538 
2539    Not Collective
2540 
2541    Input Parameter:
2542 .  snes - the SNES context
2543 
2544    Output Parameter:
2545 .  prefix - pointer to the prefix string used
2546 
2547    Notes: On the fortran side, the user should pass in a string 'prefix' of
2548    sufficient length to hold the prefix.
2549 
2550    Level: advanced
2551 
2552 .keywords: SNES, get, options, prefix, database
2553 
2554 .seealso: SNESAppendOptionsPrefix()
2555 @*/
2556 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
2557 {
2558   PetscErrorCode ierr;
2559 
2560   PetscFunctionBegin;
2561   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2562   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2563   PetscFunctionReturn(0);
2564 }
2565 
2566 
2567 #undef __FUNCT__
2568 #define __FUNCT__ "SNESRegister"
2569 /*@C
2570   SNESRegister - See SNESRegisterDynamic()
2571 
2572   Level: advanced
2573 @*/
2574 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2575 {
2576   char           fullname[PETSC_MAX_PATH_LEN];
2577   PetscErrorCode ierr;
2578 
2579   PetscFunctionBegin;
2580   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2581   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2582   PetscFunctionReturn(0);
2583 }
2584 
2585 #undef __FUNCT__
2586 #define __FUNCT__ "SNESTestLocalMin"
2587 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2588 {
2589   PetscErrorCode ierr;
2590   PetscInt       N,i,j;
2591   Vec            u,uh,fh;
2592   PetscScalar    value;
2593   PetscReal      norm;
2594 
2595   PetscFunctionBegin;
2596   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2597   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2598   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2599 
2600   /* currently only works for sequential */
2601   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2602   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2603   for (i=0; i<N; i++) {
2604     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2605     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2606     for (j=-10; j<11; j++) {
2607       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2608       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2609       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2610       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2611       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2612       value = -value;
2613       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2614     }
2615   }
2616   ierr = VecDestroy(uh);CHKERRQ(ierr);
2617   ierr = VecDestroy(fh);CHKERRQ(ierr);
2618   PetscFunctionReturn(0);
2619 }
2620 
2621 #undef __FUNCT__
2622 #define __FUNCT__ "SNESKSPSetUseEW"
2623 /*@
2624    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
2625    computing relative tolerance for linear solvers within an inexact
2626    Newton method.
2627 
2628    Collective on SNES
2629 
2630    Input Parameters:
2631 +  snes - SNES context
2632 -  flag - PETSC_TRUE or PETSC_FALSE
2633 
2634    Notes:
2635    Currently, the default is to use a constant relative tolerance for
2636    the inner linear solvers.  Alternatively, one can use the
2637    Eisenstat-Walker method, where the relative convergence tolerance
2638    is reset at each Newton iteration according progress of the nonlinear
2639    solver.
2640 
2641    Level: advanced
2642 
2643    Reference:
2644    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2645    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2646 
2647 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2648 
2649 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2650 @*/
2651 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag)
2652 {
2653   PetscFunctionBegin;
2654   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2655   snes->ksp_ewconv = flag;
2656   PetscFunctionReturn(0);
2657 }
2658 
2659 #undef __FUNCT__
2660 #define __FUNCT__ "SNESKSPGetUseEW"
2661 /*@
2662    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
2663    for computing relative tolerance for linear solvers within an
2664    inexact Newton method.
2665 
2666    Not Collective
2667 
2668    Input Parameter:
2669 .  snes - SNES context
2670 
2671    Output Parameter:
2672 .  flag - PETSC_TRUE or PETSC_FALSE
2673 
2674    Notes:
2675    Currently, the default is to use a constant relative tolerance for
2676    the inner linear solvers.  Alternatively, one can use the
2677    Eisenstat-Walker method, where the relative convergence tolerance
2678    is reset at each Newton iteration according progress of the nonlinear
2679    solver.
2680 
2681    Level: advanced
2682 
2683    Reference:
2684    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2685    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2686 
2687 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2688 
2689 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2690 @*/
2691 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag)
2692 {
2693   PetscFunctionBegin;
2694   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2695   PetscValidPointer(flag,2);
2696   *flag = snes->ksp_ewconv;
2697   PetscFunctionReturn(0);
2698 }
2699 
2700 #undef __FUNCT__
2701 #define __FUNCT__ "SNESKSPSetParametersEW"
2702 /*@
2703    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
2704    convergence criteria for the linear solvers within an inexact
2705    Newton method.
2706 
2707    Collective on SNES
2708 
2709    Input Parameters:
2710 +    snes - SNES context
2711 .    version - version 1, 2 (default is 2) or 3
2712 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2713 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2714 .    gamma - multiplicative factor for version 2 rtol computation
2715              (0 <= gamma2 <= 1)
2716 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2717 .    alpha2 - power for safeguard
2718 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2719 
2720    Note:
2721    Version 3 was contributed by Luis Chacon, June 2006.
2722 
2723    Use PETSC_DEFAULT to retain the default for any of the parameters.
2724 
2725    Level: advanced
2726 
2727    Reference:
2728    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2729    inexact Newton method", Utah State University Math. Stat. Dept. Res.
2730    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
2731 
2732 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
2733 
2734 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
2735 @*/
2736 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
2737 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
2738 {
2739   SNESKSPEW *kctx;
2740   PetscFunctionBegin;
2741   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2742   kctx = (SNESKSPEW*)snes->kspconvctx;
2743   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2744 
2745   if (version != PETSC_DEFAULT)   kctx->version   = version;
2746   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
2747   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
2748   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
2749   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
2750   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
2751   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
2752 
2753   if (kctx->version < 1 || kctx->version > 3) {
2754     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
2755   }
2756   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
2757     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
2758   }
2759   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
2760     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
2761   }
2762   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
2763     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
2764   }
2765   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
2766     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
2767   }
2768   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
2769     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
2770   }
2771   PetscFunctionReturn(0);
2772 }
2773 
2774 #undef __FUNCT__
2775 #define __FUNCT__ "SNESKSPGetParametersEW"
2776 /*@
2777    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
2778    convergence criteria for the linear solvers within an inexact
2779    Newton method.
2780 
2781    Not Collective
2782 
2783    Input Parameters:
2784      snes - SNES context
2785 
2786    Output Parameters:
2787 +    version - version 1, 2 (default is 2) or 3
2788 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2789 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2790 .    gamma - multiplicative factor for version 2 rtol computation
2791              (0 <= gamma2 <= 1)
2792 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2793 .    alpha2 - power for safeguard
2794 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2795 
2796    Level: advanced
2797 
2798 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
2799 
2800 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
2801 @*/
2802 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
2803 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
2804 {
2805   SNESKSPEW *kctx;
2806   PetscFunctionBegin;
2807   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2808   kctx = (SNESKSPEW*)snes->kspconvctx;
2809   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2810   if(version)   *version   = kctx->version;
2811   if(rtol_0)    *rtol_0    = kctx->rtol_0;
2812   if(rtol_max)  *rtol_max  = kctx->rtol_max;
2813   if(gamma)     *gamma     = kctx->gamma;
2814   if(alpha)     *alpha     = kctx->alpha;
2815   if(alpha2)    *alpha2    = kctx->alpha2;
2816   if(threshold) *threshold = kctx->threshold;
2817   PetscFunctionReturn(0);
2818 }
2819 
2820 #undef __FUNCT__
2821 #define __FUNCT__ "SNESKSPEW_PreSolve"
2822 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
2823 {
2824   PetscErrorCode ierr;
2825   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2826   PetscReal      rtol=PETSC_DEFAULT,stol;
2827 
2828   PetscFunctionBegin;
2829   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2830   if (!snes->iter) { /* first time in, so use the original user rtol */
2831     rtol = kctx->rtol_0;
2832   } else {
2833     if (kctx->version == 1) {
2834       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
2835       if (rtol < 0.0) rtol = -rtol;
2836       stol = pow(kctx->rtol_last,kctx->alpha2);
2837       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2838     } else if (kctx->version == 2) {
2839       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2840       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
2841       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2842     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
2843       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2844       /* safeguard: avoid sharp decrease of rtol */
2845       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
2846       stol = PetscMax(rtol,stol);
2847       rtol = PetscMin(kctx->rtol_0,stol);
2848       /* safeguard: avoid oversolving */
2849       stol = kctx->gamma*(snes->ttol)/snes->norm;
2850       stol = PetscMax(rtol,stol);
2851       rtol = PetscMin(kctx->rtol_0,stol);
2852     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
2853   }
2854   /* safeguard: avoid rtol greater than one */
2855   rtol = PetscMin(rtol,kctx->rtol_max);
2856   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2857   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
2858   PetscFunctionReturn(0);
2859 }
2860 
2861 #undef __FUNCT__
2862 #define __FUNCT__ "SNESKSPEW_PostSolve"
2863 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
2864 {
2865   PetscErrorCode ierr;
2866   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2867   PCSide         pcside;
2868   Vec            lres;
2869 
2870   PetscFunctionBegin;
2871   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2872   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
2873   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
2874   if (kctx->version == 1) {
2875     ierr = KSPGetPreconditionerSide(ksp,&pcside);CHKERRQ(ierr);
2876     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
2877       /* KSP residual is true linear residual */
2878       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
2879     } else {
2880       /* KSP residual is preconditioned residual */
2881       /* compute true linear residual norm */
2882       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
2883       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
2884       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
2885       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
2886       ierr = VecDestroy(lres);CHKERRQ(ierr);
2887     }
2888   }
2889   PetscFunctionReturn(0);
2890 }
2891 
2892 #undef __FUNCT__
2893 #define __FUNCT__ "SNES_KSPSolve"
2894 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
2895 {
2896   PetscErrorCode ierr;
2897 
2898   PetscFunctionBegin;
2899   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
2900   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
2901   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
2902   PetscFunctionReturn(0);
2903 }
2904