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