xref: /petsc/src/snes/interface/snes.c (revision e35cf81ddd3d80ce338419283d20072b93e92362)
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 == -1) {
1074     *flg = SAME_PRECONDITIONER;
1075     ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr);
1076     PetscFunctionReturn(0);
1077   } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) {
1078     *flg = SAME_PRECONDITIONER;
1079     ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr);
1080     PetscFunctionReturn(0);
1081   }
1082 
1083   *flg = DIFFERENT_NONZERO_PATTERN;
1084   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1085   PetscStackPush("SNES user Jacobian function");
1086   CHKMEMQ;
1087   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1088   CHKMEMQ;
1089   PetscStackPop;
1090   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1091 
1092   if (snes->lagpreconditioner == -1) {
1093     *flg = SAME_PRECONDITIONER;
1094     ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr);
1095   } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) {
1096     *flg = SAME_PRECONDITIONER;
1097     ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr);
1098   }
1099 
1100   /* make sure user returned a correct Jacobian and preconditioner */
1101   /* PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
1102     PetscValidHeaderSpecific(*B,MAT_COOKIE,4);   */
1103   PetscFunctionReturn(0);
1104 }
1105 
1106 #undef __FUNCT__
1107 #define __FUNCT__ "SNESSetJacobian"
1108 /*@C
1109    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1110    location to store the matrix.
1111 
1112    Collective on SNES and Mat
1113 
1114    Input Parameters:
1115 +  snes - the SNES context
1116 .  A - Jacobian matrix
1117 .  B - preconditioner matrix (usually same as the Jacobian)
1118 .  func - Jacobian evaluation routine
1119 -  ctx - [optional] user-defined context for private data for the
1120          Jacobian evaluation routine (may be PETSC_NULL)
1121 
1122    Calling sequence of func:
1123 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1124 
1125 +  x - input vector
1126 .  A - Jacobian matrix
1127 .  B - preconditioner matrix, usually the same as A
1128 .  flag - flag indicating information about the preconditioner matrix
1129    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1130 -  ctx - [optional] user-defined Jacobian context
1131 
1132    Notes:
1133    See KSPSetOperators() for important information about setting the flag
1134    output parameter in the routine func().  Be sure to read this information!
1135 
1136    The routine func() takes Mat * as the matrix arguments rather than Mat.
1137    This allows the Jacobian evaluation routine to replace A and/or B with a
1138    completely new new matrix structure (not just different matrix elements)
1139    when appropriate, for instance, if the nonzero structure is changing
1140    throughout the global iterations.
1141 
1142    If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on
1143    each matrix.
1144 
1145    Level: beginner
1146 
1147 .keywords: SNES, nonlinear, set, Jacobian, matrix
1148 
1149 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
1150 @*/
1151 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1152 {
1153   PetscErrorCode ierr;
1154 
1155   PetscFunctionBegin;
1156   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1157   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
1158   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
1159   if (A) PetscCheckSameComm(snes,1,A,2);
1160   if (B) PetscCheckSameComm(snes,1,B,2);
1161   if (func) snes->ops->computejacobian = func;
1162   if (ctx)  snes->jacP                 = ctx;
1163   if (A) {
1164     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1165     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1166     snes->jacobian = A;
1167   }
1168   if (B) {
1169     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
1170     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1171     snes->jacobian_pre = B;
1172   }
1173   PetscFunctionReturn(0);
1174 }
1175 
1176 #undef __FUNCT__
1177 #define __FUNCT__ "SNESGetJacobian"
1178 /*@C
1179    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1180    provided context for evaluating the Jacobian.
1181 
1182    Not Collective, but Mat object will be parallel if SNES object is
1183 
1184    Input Parameter:
1185 .  snes - the nonlinear solver context
1186 
1187    Output Parameters:
1188 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1189 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1190 .  func - location to put Jacobian function (or PETSC_NULL)
1191 -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1192 
1193    Level: advanced
1194 
1195 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1196 @*/
1197 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1198 {
1199   PetscFunctionBegin;
1200   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1201   if (A)    *A    = snes->jacobian;
1202   if (B)    *B    = snes->jacobian_pre;
1203   if (func) *func = snes->ops->computejacobian;
1204   if (ctx)  *ctx  = snes->jacP;
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1209 EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
1210 
1211 #undef __FUNCT__
1212 #define __FUNCT__ "SNESSetUp"
1213 /*@
1214    SNESSetUp - Sets up the internal data structures for the later use
1215    of a nonlinear solver.
1216 
1217    Collective on SNES
1218 
1219    Input Parameters:
1220 .  snes - the SNES context
1221 
1222    Notes:
1223    For basic use of the SNES solvers the user need not explicitly call
1224    SNESSetUp(), since these actions will automatically occur during
1225    the call to SNESSolve().  However, if one wishes to control this
1226    phase separately, SNESSetUp() should be called after SNESCreate()
1227    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1228 
1229    Level: advanced
1230 
1231 .keywords: SNES, nonlinear, setup
1232 
1233 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1234 @*/
1235 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
1236 {
1237   PetscErrorCode ierr;
1238   PetscTruth     flg;
1239 
1240   PetscFunctionBegin;
1241   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1242   if (snes->setupcalled) PetscFunctionReturn(0);
1243 
1244   if (!((PetscObject)snes)->type_name) {
1245     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
1246   }
1247 
1248   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1249   /*
1250       This version replaces the user provided Jacobian matrix with a
1251       matrix-free version but still employs the user-provided preconditioner matrix
1252   */
1253   if (flg) {
1254     Mat J;
1255     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
1256     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1257     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
1258     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
1259     ierr = MatDestroy(J);CHKERRQ(ierr);
1260   }
1261 
1262 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT)
1263   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
1264   if (flg) {
1265     Mat J;
1266     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1267     ierr = PetscInfo(snes,"Setting default matrix-free operator routines (version 2)\n");CHKERRQ(ierr);
1268     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
1269     ierr = MatDestroy(J);CHKERRQ(ierr);
1270   }
1271 #endif
1272 
1273   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1274   /*
1275       This version replaces both the user-provided Jacobian and the user-
1276       provided preconditioner matrix with the default matrix free version.
1277    */
1278   if (flg) {
1279     Mat  J;
1280     KSP ksp;
1281     PC   pc;
1282     /* create and set matrix-free operator */
1283     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
1284     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1285     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
1286     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
1287     ierr = MatDestroy(J);CHKERRQ(ierr);
1288     /* force no preconditioner */
1289     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1290     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1291     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1292     if (!flg) {
1293       ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines;\nThat is no preconditioner is being used\n");CHKERRQ(ierr);
1294       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1295     }
1296   }
1297 
1298   if (!snes->vec_func && !snes->vec_rhs) {
1299     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1300   }
1301   if (!snes->ops->computefunction && !snes->vec_rhs) {
1302     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1303   }
1304   if (!snes->jacobian) {
1305     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1306   }
1307   if (snes->vec_func == snes->vec_sol) {
1308     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1309   }
1310 
1311   if (snes->ops->setup) {
1312     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
1313   }
1314   snes->setupcalled = PETSC_TRUE;
1315   PetscFunctionReturn(0);
1316 }
1317 
1318 #undef __FUNCT__
1319 #define __FUNCT__ "SNESDestroy"
1320 /*@
1321    SNESDestroy - Destroys the nonlinear solver context that was created
1322    with SNESCreate().
1323 
1324    Collective on SNES
1325 
1326    Input Parameter:
1327 .  snes - the SNES context
1328 
1329    Level: beginner
1330 
1331 .keywords: SNES, nonlinear, destroy
1332 
1333 .seealso: SNESCreate(), SNESSolve()
1334 @*/
1335 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
1336 {
1337   PetscErrorCode ierr;
1338 
1339   PetscFunctionBegin;
1340   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1341   if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0);
1342 
1343   /* if memory was published with AMS then destroy it */
1344   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1345   if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);}
1346 
1347   if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);}
1348   if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);}
1349   if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);}
1350   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1351   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1352   if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);}
1353   ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);
1354   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1355   ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
1356   if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);}
1357   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
1358   PetscFunctionReturn(0);
1359 }
1360 
1361 /* ----------- Routines to set solver parameters ---------- */
1362 
1363 #undef __FUNCT__
1364 #define __FUNCT__ "SNESSetLagPreconditioner"
1365 /*@
1366    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
1367 
1368    Collective on SNES
1369 
1370    Input Parameters:
1371 +  snes - the SNES context
1372 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1373          the Jacobian is built etc.
1374 
1375    Options Database Keys:
1376 .    -snes_lag_preconditioner <lag>
1377 
1378    Notes:
1379    The default is 1
1380    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1381    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
1382 
1383    Level: intermediate
1384 
1385 .keywords: SNES, nonlinear, set, convergence, tolerances
1386 
1387 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
1388 
1389 @*/
1390 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag)
1391 {
1392   PetscFunctionBegin;
1393   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1394   if (lag < -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -1, 1 or greater");
1395   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1396   snes->lagpreconditioner = lag;
1397   PetscFunctionReturn(0);
1398 }
1399 
1400 #undef __FUNCT__
1401 #define __FUNCT__ "SNESGetLagPreconditioner"
1402 /*@
1403    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
1404 
1405    Collective on SNES
1406 
1407    Input Parameter:
1408 .  snes - the SNES context
1409 
1410    Output Parameter:
1411 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1412          the Jacobian is built etc.
1413 
1414    Options Database Keys:
1415 .    -snes_lag_preconditioner <lag>
1416 
1417    Notes:
1418    The default is 1
1419    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1420 
1421    Level: intermediate
1422 
1423 .keywords: SNES, nonlinear, set, convergence, tolerances
1424 
1425 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
1426 
1427 @*/
1428 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
1429 {
1430   PetscFunctionBegin;
1431   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1432   *lag = snes->lagpreconditioner;
1433   PetscFunctionReturn(0);
1434 }
1435 
1436 #undef __FUNCT__
1437 #define __FUNCT__ "SNESSetLagJacobian"
1438 /*@
1439    SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
1440      often the preconditioner is rebuilt.
1441 
1442    Collective on SNES
1443 
1444    Input Parameters:
1445 +  snes - the SNES context
1446 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1447          the Jacobian is built etc.
1448 
1449    Options Database Keys:
1450 .    -snes_lag_jacobian <lag>
1451 
1452    Notes:
1453    The default is 1
1454    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1455    If  -1 is used before the very first nonlinear solve the CODE WILL FAIL! because no Jacobian is used
1456 
1457    Level: intermediate
1458 
1459 .keywords: SNES, nonlinear, set, convergence, tolerances
1460 
1461 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
1462 
1463 @*/
1464 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagJacobian(SNES snes,PetscInt lag)
1465 {
1466   PetscFunctionBegin;
1467   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1468   if (lag < -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -1, 1 or greater");
1469   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1470   snes->lagjacobian = lag;
1471   PetscFunctionReturn(0);
1472 }
1473 
1474 #undef __FUNCT__
1475 #define __FUNCT__ "SNESGetLagJacobian"
1476 /*@
1477    SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
1478 
1479    Collective on SNES
1480 
1481    Input Parameter:
1482 .  snes - the SNES context
1483 
1484    Output Parameter:
1485 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1486          the Jacobian is built etc.
1487 
1488    Options Database Keys:
1489 .    -snes_lag_jacobian <lag>
1490 
1491    Notes:
1492    The default is 1
1493    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1494 
1495    Level: intermediate
1496 
1497 .keywords: SNES, nonlinear, set, convergence, tolerances
1498 
1499 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
1500 
1501 @*/
1502 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagJacobian(SNES snes,PetscInt *lag)
1503 {
1504   PetscFunctionBegin;
1505   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1506   *lag = snes->lagjacobian;
1507   PetscFunctionReturn(0);
1508 }
1509 
1510 #undef __FUNCT__
1511 #define __FUNCT__ "SNESSetTolerances"
1512 /*@
1513    SNESSetTolerances - Sets various parameters used in convergence tests.
1514 
1515    Collective on SNES
1516 
1517    Input Parameters:
1518 +  snes - the SNES context
1519 .  abstol - absolute convergence tolerance
1520 .  rtol - relative convergence tolerance
1521 .  stol -  convergence tolerance in terms of the norm
1522            of the change in the solution between steps
1523 .  maxit - maximum number of iterations
1524 -  maxf - maximum number of function evaluations
1525 
1526    Options Database Keys:
1527 +    -snes_atol <abstol> - Sets abstol
1528 .    -snes_rtol <rtol> - Sets rtol
1529 .    -snes_stol <stol> - Sets stol
1530 .    -snes_max_it <maxit> - Sets maxit
1531 -    -snes_max_funcs <maxf> - Sets maxf
1532 
1533    Notes:
1534    The default maximum number of iterations is 50.
1535    The default maximum number of function evaluations is 1000.
1536 
1537    Level: intermediate
1538 
1539 .keywords: SNES, nonlinear, set, convergence, tolerances
1540 
1541 .seealso: SNESSetTrustRegionTolerance()
1542 @*/
1543 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
1544 {
1545   PetscFunctionBegin;
1546   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1547   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1548   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1549   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1550   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1551   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1552   PetscFunctionReturn(0);
1553 }
1554 
1555 #undef __FUNCT__
1556 #define __FUNCT__ "SNESGetTolerances"
1557 /*@
1558    SNESGetTolerances - Gets various parameters used in convergence tests.
1559 
1560    Not Collective
1561 
1562    Input Parameters:
1563 +  snes - the SNES context
1564 .  atol - absolute convergence tolerance
1565 .  rtol - relative convergence tolerance
1566 .  stol -  convergence tolerance in terms of the norm
1567            of the change in the solution between steps
1568 .  maxit - maximum number of iterations
1569 -  maxf - maximum number of function evaluations
1570 
1571    Notes:
1572    The user can specify PETSC_NULL for any parameter that is not needed.
1573 
1574    Level: intermediate
1575 
1576 .keywords: SNES, nonlinear, get, convergence, tolerances
1577 
1578 .seealso: SNESSetTolerances()
1579 @*/
1580 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
1581 {
1582   PetscFunctionBegin;
1583   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1584   if (atol)  *atol  = snes->abstol;
1585   if (rtol)  *rtol  = snes->rtol;
1586   if (stol)  *stol  = snes->xtol;
1587   if (maxit) *maxit = snes->max_its;
1588   if (maxf)  *maxf  = snes->max_funcs;
1589   PetscFunctionReturn(0);
1590 }
1591 
1592 #undef __FUNCT__
1593 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1594 /*@
1595    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1596 
1597    Collective on SNES
1598 
1599    Input Parameters:
1600 +  snes - the SNES context
1601 -  tol - tolerance
1602 
1603    Options Database Key:
1604 .  -snes_trtol <tol> - Sets tol
1605 
1606    Level: intermediate
1607 
1608 .keywords: SNES, nonlinear, set, trust region, tolerance
1609 
1610 .seealso: SNESSetTolerances()
1611 @*/
1612 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1613 {
1614   PetscFunctionBegin;
1615   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1616   snes->deltatol = tol;
1617   PetscFunctionReturn(0);
1618 }
1619 
1620 /*
1621    Duplicate the lg monitors for SNES from KSP; for some reason with
1622    dynamic libraries things don't work under Sun4 if we just use
1623    macros instead of functions
1624 */
1625 #undef __FUNCT__
1626 #define __FUNCT__ "SNESMonitorLG"
1627 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1628 {
1629   PetscErrorCode ierr;
1630 
1631   PetscFunctionBegin;
1632   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1633   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1634   PetscFunctionReturn(0);
1635 }
1636 
1637 #undef __FUNCT__
1638 #define __FUNCT__ "SNESMonitorLGCreate"
1639 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1640 {
1641   PetscErrorCode ierr;
1642 
1643   PetscFunctionBegin;
1644   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1645   PetscFunctionReturn(0);
1646 }
1647 
1648 #undef __FUNCT__
1649 #define __FUNCT__ "SNESMonitorLGDestroy"
1650 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw)
1651 {
1652   PetscErrorCode ierr;
1653 
1654   PetscFunctionBegin;
1655   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1656   PetscFunctionReturn(0);
1657 }
1658 
1659 /* ------------ Routines to set performance monitoring options ----------- */
1660 
1661 #undef __FUNCT__
1662 #define __FUNCT__ "SNESMonitorSet"
1663 /*@C
1664    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
1665    iteration of the nonlinear solver to display the iteration's
1666    progress.
1667 
1668    Collective on SNES
1669 
1670    Input Parameters:
1671 +  snes - the SNES context
1672 .  func - monitoring routine
1673 .  mctx - [optional] user-defined context for private data for the
1674           monitor routine (use PETSC_NULL if no context is desired)
1675 -  monitordestroy - [optional] routine that frees monitor context
1676           (may be PETSC_NULL)
1677 
1678    Calling sequence of func:
1679 $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1680 
1681 +    snes - the SNES context
1682 .    its - iteration number
1683 .    norm - 2-norm function value (may be estimated)
1684 -    mctx - [optional] monitoring context
1685 
1686    Options Database Keys:
1687 +    -snes_monitor        - sets SNESMonitorDefault()
1688 .    -snes_monitor_draw    - sets line graph monitor,
1689                             uses SNESMonitorLGCreate()
1690 _    -snes_monitor_cancel - cancels all monitors that have
1691                             been hardwired into a code by
1692                             calls to SNESMonitorSet(), but
1693                             does not cancel those set via
1694                             the options database.
1695 
1696    Notes:
1697    Several different monitoring routines may be set by calling
1698    SNESMonitorSet() multiple times; all will be called in the
1699    order in which they were set.
1700 
1701    Fortran notes: Only a single monitor function can be set for each SNES object
1702 
1703    Level: intermediate
1704 
1705 .keywords: SNES, nonlinear, set, monitor
1706 
1707 .seealso: SNESMonitorDefault(), SNESMonitorCancel()
1708 @*/
1709 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
1710 {
1711   PetscInt i;
1712 
1713   PetscFunctionBegin;
1714   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1715   if (snes->numbermonitors >= MAXSNESMONITORS) {
1716     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1717   }
1718   for (i=0; i<snes->numbermonitors;i++) {
1719     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0);
1720 
1721     /* check if both default monitors that share common ASCII viewer */
1722     if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) {
1723       if (mctx && snes->monitorcontext[i]) {
1724         PetscErrorCode          ierr;
1725         PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx;
1726         PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i];
1727         if (viewer1->viewer == viewer2->viewer) {
1728           ierr = (*monitordestroy)(mctx);CHKERRQ(ierr);
1729           PetscFunctionReturn(0);
1730         }
1731       }
1732     }
1733   }
1734   snes->monitor[snes->numbermonitors]           = monitor;
1735   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1736   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1737   PetscFunctionReturn(0);
1738 }
1739 
1740 #undef __FUNCT__
1741 #define __FUNCT__ "SNESMonitorCancel"
1742 /*@C
1743    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
1744 
1745    Collective on SNES
1746 
1747    Input Parameters:
1748 .  snes - the SNES context
1749 
1750    Options Database Key:
1751 .  -snes_monitor_cancel - cancels all monitors that have been hardwired
1752     into a code by calls to SNESMonitorSet(), but does not cancel those
1753     set via the options database
1754 
1755    Notes:
1756    There is no way to clear one specific monitor from a SNES object.
1757 
1758    Level: intermediate
1759 
1760 .keywords: SNES, nonlinear, set, monitor
1761 
1762 .seealso: SNESMonitorDefault(), SNESMonitorSet()
1763 @*/
1764 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes)
1765 {
1766   PetscErrorCode ierr;
1767   PetscInt       i;
1768 
1769   PetscFunctionBegin;
1770   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1771   for (i=0; i<snes->numbermonitors; i++) {
1772     if (snes->monitordestroy[i]) {
1773       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1774     }
1775   }
1776   snes->numbermonitors = 0;
1777   PetscFunctionReturn(0);
1778 }
1779 
1780 #undef __FUNCT__
1781 #define __FUNCT__ "SNESSetConvergenceTest"
1782 /*@C
1783    SNESSetConvergenceTest - Sets the function that is to be used
1784    to test for convergence of the nonlinear iterative solution.
1785 
1786    Collective on SNES
1787 
1788    Input Parameters:
1789 +  snes - the SNES context
1790 .  func - routine to test for convergence
1791 .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
1792 -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
1793 
1794    Calling sequence of func:
1795 $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1796 
1797 +    snes - the SNES context
1798 .    it - current iteration (0 is the first and is before any Newton step)
1799 .    cctx - [optional] convergence context
1800 .    reason - reason for convergence/divergence
1801 .    xnorm - 2-norm of current iterate
1802 .    gnorm - 2-norm of current step
1803 -    f - 2-norm of function
1804 
1805    Level: advanced
1806 
1807 .keywords: SNES, nonlinear, set, convergence, test
1808 
1809 .seealso: SNESDefaultConverged(), SNESSkipConverged()
1810 @*/
1811 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1812 {
1813   PetscErrorCode ierr;
1814 
1815   PetscFunctionBegin;
1816   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1817   if (!func) func = SNESSkipConverged;
1818   if (snes->ops->convergeddestroy) {
1819     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
1820   }
1821   snes->ops->converged        = func;
1822   snes->ops->convergeddestroy = destroy;
1823   snes->cnvP                  = cctx;
1824   PetscFunctionReturn(0);
1825 }
1826 
1827 #undef __FUNCT__
1828 #define __FUNCT__ "SNESGetConvergedReason"
1829 /*@
1830    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1831 
1832    Not Collective
1833 
1834    Input Parameter:
1835 .  snes - the SNES context
1836 
1837    Output Parameter:
1838 .  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
1839             manual pages for the individual convergence tests for complete lists
1840 
1841    Level: intermediate
1842 
1843    Notes: Can only be called after the call the SNESSolve() is complete.
1844 
1845 .keywords: SNES, nonlinear, set, convergence, test
1846 
1847 .seealso: SNESSetConvergenceTest(), SNESConvergedReason
1848 @*/
1849 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1850 {
1851   PetscFunctionBegin;
1852   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1853   PetscValidPointer(reason,2);
1854   *reason = snes->reason;
1855   PetscFunctionReturn(0);
1856 }
1857 
1858 #undef __FUNCT__
1859 #define __FUNCT__ "SNESSetConvergenceHistory"
1860 /*@
1861    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1862 
1863    Collective on SNES
1864 
1865    Input Parameters:
1866 +  snes - iterative context obtained from SNESCreate()
1867 .  a   - array to hold history
1868 .  its - integer array holds the number of linear iterations for each solve.
1869 .  na  - size of a and its
1870 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1871            else it continues storing new values for new nonlinear solves after the old ones
1872 
1873    Notes:
1874    If set, this array will contain the function norms computed
1875    at each step.
1876 
1877    This routine is useful, e.g., when running a code for purposes
1878    of accurate performance monitoring, when no I/O should be done
1879    during the section of code that is being timed.
1880 
1881    Level: intermediate
1882 
1883 .keywords: SNES, set, convergence, history
1884 
1885 .seealso: SNESGetConvergenceHistory()
1886 
1887 @*/
1888 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset)
1889 {
1890   PetscFunctionBegin;
1891   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1892   if (na)  PetscValidScalarPointer(a,2);
1893   if (its) PetscValidIntPointer(its,3);
1894   snes->conv_hist       = a;
1895   snes->conv_hist_its   = its;
1896   snes->conv_hist_max   = na;
1897   snes->conv_hist_len   = 0;
1898   snes->conv_hist_reset = reset;
1899   PetscFunctionReturn(0);
1900 }
1901 
1902 #undef __FUNCT__
1903 #define __FUNCT__ "SNESGetConvergenceHistory"
1904 /*@C
1905    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1906 
1907    Collective on SNES
1908 
1909    Input Parameter:
1910 .  snes - iterative context obtained from SNESCreate()
1911 
1912    Output Parameters:
1913 .  a   - array to hold history
1914 .  its - integer array holds the number of linear iterations (or
1915          negative if not converged) for each solve.
1916 -  na  - size of a and its
1917 
1918    Notes:
1919     The calling sequence for this routine in Fortran is
1920 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1921 
1922    This routine is useful, e.g., when running a code for purposes
1923    of accurate performance monitoring, when no I/O should be done
1924    during the section of code that is being timed.
1925 
1926    Level: intermediate
1927 
1928 .keywords: SNES, get, convergence, history
1929 
1930 .seealso: SNESSetConvergencHistory()
1931 
1932 @*/
1933 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
1934 {
1935   PetscFunctionBegin;
1936   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1937   if (a)   *a   = snes->conv_hist;
1938   if (its) *its = snes->conv_hist_its;
1939   if (na)  *na  = snes->conv_hist_len;
1940   PetscFunctionReturn(0);
1941 }
1942 
1943 #undef __FUNCT__
1944 #define __FUNCT__ "SNESSetUpdate"
1945 /*@C
1946   SNESSetUpdate - Sets the general-purpose update function called
1947   at the beginning o every iteration of the nonlinear solve. Specifically
1948   it is called just before the Jacobian is "evaluated".
1949 
1950   Collective on SNES
1951 
1952   Input Parameters:
1953 . snes - The nonlinear solver context
1954 . func - The function
1955 
1956   Calling sequence of func:
1957 . func (SNES snes, PetscInt step);
1958 
1959 . step - The current step of the iteration
1960 
1961   Level: intermediate
1962 
1963 .keywords: SNES, update
1964 
1965 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
1966 @*/
1967 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
1968 {
1969   PetscFunctionBegin;
1970   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
1971   snes->ops->update = func;
1972   PetscFunctionReturn(0);
1973 }
1974 
1975 #undef __FUNCT__
1976 #define __FUNCT__ "SNESDefaultUpdate"
1977 /*@
1978   SNESDefaultUpdate - The default update function which does nothing.
1979 
1980   Not collective
1981 
1982   Input Parameters:
1983 . snes - The nonlinear solver context
1984 . step - The current step of the iteration
1985 
1986   Level: intermediate
1987 
1988 .keywords: SNES, update
1989 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
1990 @*/
1991 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
1992 {
1993   PetscFunctionBegin;
1994   PetscFunctionReturn(0);
1995 }
1996 
1997 #undef __FUNCT__
1998 #define __FUNCT__ "SNESScaleStep_Private"
1999 /*
2000    SNESScaleStep_Private - Scales a step so that its length is less than the
2001    positive parameter delta.
2002 
2003     Input Parameters:
2004 +   snes - the SNES context
2005 .   y - approximate solution of linear system
2006 .   fnorm - 2-norm of current function
2007 -   delta - trust region size
2008 
2009     Output Parameters:
2010 +   gpnorm - predicted function norm at the new point, assuming local
2011     linearization.  The value is zero if the step lies within the trust
2012     region, and exceeds zero otherwise.
2013 -   ynorm - 2-norm of the step
2014 
2015     Note:
2016     For non-trust region methods such as SNESLS, the parameter delta
2017     is set to be the maximum allowable step size.
2018 
2019 .keywords: SNES, nonlinear, scale, step
2020 */
2021 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
2022 {
2023   PetscReal      nrm;
2024   PetscScalar    cnorm;
2025   PetscErrorCode ierr;
2026 
2027   PetscFunctionBegin;
2028   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2029   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
2030   PetscCheckSameComm(snes,1,y,2);
2031 
2032   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
2033   if (nrm > *delta) {
2034      nrm = *delta/nrm;
2035      *gpnorm = (1.0 - nrm)*(*fnorm);
2036      cnorm = nrm;
2037      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
2038      *ynorm = *delta;
2039   } else {
2040      *gpnorm = 0.0;
2041      *ynorm = nrm;
2042   }
2043   PetscFunctionReturn(0);
2044 }
2045 
2046 #undef __FUNCT__
2047 #define __FUNCT__ "SNESSolve"
2048 /*@C
2049    SNESSolve - Solves a nonlinear system F(x) = b.
2050    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
2051 
2052    Collective on SNES
2053 
2054    Input Parameters:
2055 +  snes - the SNES context
2056 .  b - the constant part of the equation, or PETSC_NULL to use zero.
2057 -  x - the solution vector.
2058 
2059    Notes:
2060    The user should initialize the vector,x, with the initial guess
2061    for the nonlinear solve prior to calling SNESSolve.  In particular,
2062    to employ an initial guess of zero, the user should explicitly set
2063    this vector to zero by calling VecSet().
2064 
2065    Level: beginner
2066 
2067 .keywords: SNES, nonlinear, solve
2068 
2069 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian()
2070 @*/
2071 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
2072 {
2073   PetscErrorCode ierr;
2074   PetscTruth     flg;
2075   char           filename[PETSC_MAX_PATH_LEN];
2076   PetscViewer    viewer;
2077 
2078   PetscFunctionBegin;
2079   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2080   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
2081   PetscCheckSameComm(snes,1,x,3);
2082   if (b) PetscValidHeaderSpecific(b,VEC_COOKIE,2);
2083   if (b) PetscCheckSameComm(snes,1,b,2);
2084 
2085   /* set solution vector */
2086   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
2087   if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); }
2088   snes->vec_sol = x;
2089   /* set afine vector if provided */
2090   if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
2091   if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); }
2092   snes->vec_rhs = b;
2093 
2094   if (!snes->vec_func && snes->vec_rhs) {
2095     ierr = VecDuplicate(b, &snes->vec_func); CHKERRQ(ierr);
2096   }
2097 
2098   ierr = SNESSetUp(snes);CHKERRQ(ierr);
2099 
2100   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
2101   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
2102 
2103   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2104   ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
2105   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2106   if (snes->domainerror){
2107     snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
2108     snes->domainerror = PETSC_FALSE;
2109   }
2110 
2111   if (!snes->reason) {
2112     SETERRQ(PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
2113   }
2114 
2115   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
2116   if (flg && !PetscPreLoadingOn) {
2117     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
2118     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
2119     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
2120   }
2121 
2122   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
2123   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
2124   if (snes->printreason) {
2125     if (snes->reason > 0) {
2126       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2127     } else {
2128       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2129     }
2130   }
2131 
2132   PetscFunctionReturn(0);
2133 }
2134 
2135 /* --------- Internal routines for SNES Package --------- */
2136 
2137 #undef __FUNCT__
2138 #define __FUNCT__ "SNESSetType"
2139 /*@C
2140    SNESSetType - Sets the method for the nonlinear solver.
2141 
2142    Collective on SNES
2143 
2144    Input Parameters:
2145 +  snes - the SNES context
2146 -  type - a known method
2147 
2148    Options Database Key:
2149 .  -snes_type <type> - Sets the method; use -help for a list
2150    of available methods (for instance, ls or tr)
2151 
2152    Notes:
2153    See "petsc/include/petscsnes.h" for available methods (for instance)
2154 +    SNESLS - Newton's method with line search
2155      (systems of nonlinear equations)
2156 .    SNESTR - Newton's method with trust region
2157      (systems of nonlinear equations)
2158 
2159   Normally, it is best to use the SNESSetFromOptions() command and then
2160   set the SNES solver type from the options database rather than by using
2161   this routine.  Using the options database provides the user with
2162   maximum flexibility in evaluating the many nonlinear solvers.
2163   The SNESSetType() routine is provided for those situations where it
2164   is necessary to set the nonlinear solver independently of the command
2165   line or options database.  This might be the case, for example, when
2166   the choice of solver changes during the execution of the program,
2167   and the user's application is taking responsibility for choosing the
2168   appropriate method.
2169 
2170   Level: intermediate
2171 
2172 .keywords: SNES, set, type
2173 
2174 .seealso: SNESType, SNESCreate()
2175 
2176 @*/
2177 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type)
2178 {
2179   PetscErrorCode ierr,(*r)(SNES);
2180   PetscTruth     match;
2181 
2182   PetscFunctionBegin;
2183   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2184   PetscValidCharPointer(type,2);
2185 
2186   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
2187   if (match) PetscFunctionReturn(0);
2188 
2189   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
2190   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
2191   /* Destroy the previous private SNES context */
2192   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
2193   /* Reinitialize function pointers in SNESOps structure */
2194   snes->ops->setup          = 0;
2195   snes->ops->solve          = 0;
2196   snes->ops->view           = 0;
2197   snes->ops->setfromoptions = 0;
2198   snes->ops->destroy        = 0;
2199   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
2200   snes->setupcalled = PETSC_FALSE;
2201   ierr = (*r)(snes);CHKERRQ(ierr);
2202   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
2203   PetscFunctionReturn(0);
2204 }
2205 
2206 
2207 /* --------------------------------------------------------------------- */
2208 #undef __FUNCT__
2209 #define __FUNCT__ "SNESRegisterDestroy"
2210 /*@
2211    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2212    registered by SNESRegisterDynamic().
2213 
2214    Not Collective
2215 
2216    Level: advanced
2217 
2218 .keywords: SNES, nonlinear, register, destroy
2219 
2220 .seealso: SNESRegisterAll(), SNESRegisterAll()
2221 @*/
2222 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
2223 {
2224   PetscErrorCode ierr;
2225 
2226   PetscFunctionBegin;
2227   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
2228   SNESRegisterAllCalled = PETSC_FALSE;
2229   PetscFunctionReturn(0);
2230 }
2231 
2232 #undef __FUNCT__
2233 #define __FUNCT__ "SNESGetType"
2234 /*@C
2235    SNESGetType - Gets the SNES method type and name (as a string).
2236 
2237    Not Collective
2238 
2239    Input Parameter:
2240 .  snes - nonlinear solver context
2241 
2242    Output Parameter:
2243 .  type - SNES method (a character string)
2244 
2245    Level: intermediate
2246 
2247 .keywords: SNES, nonlinear, get, type, name
2248 @*/
2249 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type)
2250 {
2251   PetscFunctionBegin;
2252   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2253   PetscValidPointer(type,2);
2254   *type = ((PetscObject)snes)->type_name;
2255   PetscFunctionReturn(0);
2256 }
2257 
2258 #undef __FUNCT__
2259 #define __FUNCT__ "SNESGetSolution"
2260 /*@
2261    SNESGetSolution - Returns the vector where the approximate solution is
2262    stored.
2263 
2264    Not Collective, but Vec is parallel if SNES is parallel
2265 
2266    Input Parameter:
2267 .  snes - the SNES context
2268 
2269    Output Parameter:
2270 .  x - the solution
2271 
2272    Level: intermediate
2273 
2274 .keywords: SNES, nonlinear, get, solution
2275 
2276 .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
2277 @*/
2278 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
2279 {
2280   PetscFunctionBegin;
2281   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2282   PetscValidPointer(x,2);
2283   *x = snes->vec_sol;
2284   PetscFunctionReturn(0);
2285 }
2286 
2287 #undef __FUNCT__
2288 #define __FUNCT__ "SNESGetSolutionUpdate"
2289 /*@
2290    SNESGetSolutionUpdate - Returns the vector where the solution update is
2291    stored.
2292 
2293    Not Collective, but Vec is parallel if SNES is parallel
2294 
2295    Input Parameter:
2296 .  snes - the SNES context
2297 
2298    Output Parameter:
2299 .  x - the solution update
2300 
2301    Level: advanced
2302 
2303 .keywords: SNES, nonlinear, get, solution, update
2304 
2305 .seealso: SNESGetSolution(), SNESGetFunction()
2306 @*/
2307 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
2308 {
2309   PetscFunctionBegin;
2310   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2311   PetscValidPointer(x,2);
2312   *x = snes->vec_sol_update;
2313   PetscFunctionReturn(0);
2314 }
2315 
2316 #undef __FUNCT__
2317 #define __FUNCT__ "SNESGetFunction"
2318 /*@C
2319    SNESGetFunction - Returns the vector where the function is stored.
2320 
2321    Not Collective, but Vec is parallel if SNES is parallel
2322 
2323    Input Parameter:
2324 .  snes - the SNES context
2325 
2326    Output Parameter:
2327 +  r - the function (or PETSC_NULL)
2328 .  func - the function (or PETSC_NULL)
2329 -  ctx - the function context (or PETSC_NULL)
2330 
2331    Level: advanced
2332 
2333 .keywords: SNES, nonlinear, get, function
2334 
2335 .seealso: SNESSetFunction(), SNESGetSolution()
2336 @*/
2337 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
2338 {
2339   PetscFunctionBegin;
2340   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2341   if (r)    *r    = snes->vec_func;
2342   if (func) *func = snes->ops->computefunction;
2343   if (ctx)  *ctx  = snes->funP;
2344   PetscFunctionReturn(0);
2345 }
2346 
2347 #undef __FUNCT__
2348 #define __FUNCT__ "SNESSetOptionsPrefix"
2349 /*@C
2350    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2351    SNES options in the database.
2352 
2353    Collective on SNES
2354 
2355    Input Parameter:
2356 +  snes - the SNES context
2357 -  prefix - the prefix to prepend to all option names
2358 
2359    Notes:
2360    A hyphen (-) must NOT be given at the beginning of the prefix name.
2361    The first character of all runtime options is AUTOMATICALLY the hyphen.
2362 
2363    Level: advanced
2364 
2365 .keywords: SNES, set, options, prefix, database
2366 
2367 .seealso: SNESSetFromOptions()
2368 @*/
2369 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
2370 {
2371   PetscErrorCode ierr;
2372 
2373   PetscFunctionBegin;
2374   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2375   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2376   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2377   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2378   PetscFunctionReturn(0);
2379 }
2380 
2381 #undef __FUNCT__
2382 #define __FUNCT__ "SNESAppendOptionsPrefix"
2383 /*@C
2384    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2385    SNES options in the database.
2386 
2387    Collective on SNES
2388 
2389    Input Parameters:
2390 +  snes - the SNES context
2391 -  prefix - the prefix to prepend to all option names
2392 
2393    Notes:
2394    A hyphen (-) must NOT be given at the beginning of the prefix name.
2395    The first character of all runtime options is AUTOMATICALLY the hyphen.
2396 
2397    Level: advanced
2398 
2399 .keywords: SNES, append, options, prefix, database
2400 
2401 .seealso: SNESGetOptionsPrefix()
2402 @*/
2403 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
2404 {
2405   PetscErrorCode ierr;
2406 
2407   PetscFunctionBegin;
2408   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2409   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2410   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2411   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2412   PetscFunctionReturn(0);
2413 }
2414 
2415 #undef __FUNCT__
2416 #define __FUNCT__ "SNESGetOptionsPrefix"
2417 /*@C
2418    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2419    SNES options in the database.
2420 
2421    Not Collective
2422 
2423    Input Parameter:
2424 .  snes - the SNES context
2425 
2426    Output Parameter:
2427 .  prefix - pointer to the prefix string used
2428 
2429    Notes: On the fortran side, the user should pass in a string 'prefix' of
2430    sufficient length to hold the prefix.
2431 
2432    Level: advanced
2433 
2434 .keywords: SNES, get, options, prefix, database
2435 
2436 .seealso: SNESAppendOptionsPrefix()
2437 @*/
2438 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
2439 {
2440   PetscErrorCode ierr;
2441 
2442   PetscFunctionBegin;
2443   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2444   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2445   PetscFunctionReturn(0);
2446 }
2447 
2448 
2449 #undef __FUNCT__
2450 #define __FUNCT__ "SNESRegister"
2451 /*@C
2452   SNESRegister - See SNESRegisterDynamic()
2453 
2454   Level: advanced
2455 @*/
2456 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2457 {
2458   char           fullname[PETSC_MAX_PATH_LEN];
2459   PetscErrorCode ierr;
2460 
2461   PetscFunctionBegin;
2462   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2463   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2464   PetscFunctionReturn(0);
2465 }
2466 
2467 #undef __FUNCT__
2468 #define __FUNCT__ "SNESTestLocalMin"
2469 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2470 {
2471   PetscErrorCode ierr;
2472   PetscInt       N,i,j;
2473   Vec            u,uh,fh;
2474   PetscScalar    value;
2475   PetscReal      norm;
2476 
2477   PetscFunctionBegin;
2478   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2479   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2480   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2481 
2482   /* currently only works for sequential */
2483   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2484   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2485   for (i=0; i<N; i++) {
2486     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2487     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2488     for (j=-10; j<11; j++) {
2489       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2490       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2491       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2492       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2493       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2494       value = -value;
2495       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2496     }
2497   }
2498   ierr = VecDestroy(uh);CHKERRQ(ierr);
2499   ierr = VecDestroy(fh);CHKERRQ(ierr);
2500   PetscFunctionReturn(0);
2501 }
2502 
2503 #undef __FUNCT__
2504 #define __FUNCT__ "SNESKSPSetUseEW"
2505 /*@
2506    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
2507    computing relative tolerance for linear solvers within an inexact
2508    Newton method.
2509 
2510    Collective on SNES
2511 
2512    Input Parameters:
2513 +  snes - SNES context
2514 -  flag - PETSC_TRUE or PETSC_FALSE
2515 
2516    Notes:
2517    Currently, the default is to use a constant relative tolerance for
2518    the inner linear solvers.  Alternatively, one can use the
2519    Eisenstat-Walker method, where the relative convergence tolerance
2520    is reset at each Newton iteration according progress of the nonlinear
2521    solver.
2522 
2523    Level: advanced
2524 
2525    Reference:
2526    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2527    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2528 
2529 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2530 
2531 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2532 @*/
2533 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag)
2534 {
2535   PetscFunctionBegin;
2536   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2537   snes->ksp_ewconv = flag;
2538   PetscFunctionReturn(0);
2539 }
2540 
2541 #undef __FUNCT__
2542 #define __FUNCT__ "SNESKSPGetUseEW"
2543 /*@
2544    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
2545    for computing relative tolerance for linear solvers within an
2546    inexact Newton method.
2547 
2548    Not Collective
2549 
2550    Input Parameter:
2551 .  snes - SNES context
2552 
2553    Output Parameter:
2554 .  flag - PETSC_TRUE or PETSC_FALSE
2555 
2556    Notes:
2557    Currently, the default is to use a constant relative tolerance for
2558    the inner linear solvers.  Alternatively, one can use the
2559    Eisenstat-Walker method, where the relative convergence tolerance
2560    is reset at each Newton iteration according progress of the nonlinear
2561    solver.
2562 
2563    Level: advanced
2564 
2565    Reference:
2566    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2567    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2568 
2569 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2570 
2571 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2572 @*/
2573 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag)
2574 {
2575   PetscFunctionBegin;
2576   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2577   PetscValidPointer(flag,2);
2578   *flag = snes->ksp_ewconv;
2579   PetscFunctionReturn(0);
2580 }
2581 
2582 #undef __FUNCT__
2583 #define __FUNCT__ "SNESKSPSetParametersEW"
2584 /*@
2585    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
2586    convergence criteria for the linear solvers within an inexact
2587    Newton method.
2588 
2589    Collective on SNES
2590 
2591    Input Parameters:
2592 +    snes - SNES context
2593 .    version - version 1, 2 (default is 2) or 3
2594 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2595 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2596 .    gamma - multiplicative factor for version 2 rtol computation
2597              (0 <= gamma2 <= 1)
2598 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2599 .    alpha2 - power for safeguard
2600 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2601 
2602    Note:
2603    Version 3 was contributed by Luis Chacon, June 2006.
2604 
2605    Use PETSC_DEFAULT to retain the default for any of the parameters.
2606 
2607    Level: advanced
2608 
2609    Reference:
2610    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2611    inexact Newton method", Utah State University Math. Stat. Dept. Res.
2612    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
2613 
2614 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
2615 
2616 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
2617 @*/
2618 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
2619 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
2620 {
2621   SNESKSPEW *kctx;
2622   PetscFunctionBegin;
2623   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2624   kctx = (SNESKSPEW*)snes->kspconvctx;
2625   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2626 
2627   if (version != PETSC_DEFAULT)   kctx->version   = version;
2628   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
2629   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
2630   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
2631   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
2632   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
2633   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
2634 
2635   if (kctx->version < 1 || kctx->version > 3) {
2636     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
2637   }
2638   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
2639     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
2640   }
2641   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
2642     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
2643   }
2644   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
2645     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
2646   }
2647   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
2648     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
2649   }
2650   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
2651     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
2652   }
2653   PetscFunctionReturn(0);
2654 }
2655 
2656 #undef __FUNCT__
2657 #define __FUNCT__ "SNESKSPGetParametersEW"
2658 /*@
2659    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
2660    convergence criteria for the linear solvers within an inexact
2661    Newton method.
2662 
2663    Not Collective
2664 
2665    Input Parameters:
2666      snes - SNES context
2667 
2668    Output Parameters:
2669 +    version - version 1, 2 (default is 2) or 3
2670 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2671 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2672 .    gamma - multiplicative factor for version 2 rtol computation
2673              (0 <= gamma2 <= 1)
2674 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2675 .    alpha2 - power for safeguard
2676 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2677 
2678    Level: advanced
2679 
2680 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
2681 
2682 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
2683 @*/
2684 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
2685 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
2686 {
2687   SNESKSPEW *kctx;
2688   PetscFunctionBegin;
2689   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2690   kctx = (SNESKSPEW*)snes->kspconvctx;
2691   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2692   if(version)   *version   = kctx->version;
2693   if(rtol_0)    *rtol_0    = kctx->rtol_0;
2694   if(rtol_max)  *rtol_max  = kctx->rtol_max;
2695   if(gamma)     *gamma     = kctx->gamma;
2696   if(alpha)     *alpha     = kctx->alpha;
2697   if(alpha2)    *alpha2    = kctx->alpha2;
2698   if(threshold) *threshold = kctx->threshold;
2699   PetscFunctionReturn(0);
2700 }
2701 
2702 #undef __FUNCT__
2703 #define __FUNCT__ "SNESKSPEW_PreSolve"
2704 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
2705 {
2706   PetscErrorCode ierr;
2707   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2708   PetscReal      rtol=PETSC_DEFAULT,stol;
2709 
2710   PetscFunctionBegin;
2711   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2712   if (!snes->iter) { /* first time in, so use the original user rtol */
2713     rtol = kctx->rtol_0;
2714   } else {
2715     if (kctx->version == 1) {
2716       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
2717       if (rtol < 0.0) rtol = -rtol;
2718       stol = pow(kctx->rtol_last,kctx->alpha2);
2719       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2720     } else if (kctx->version == 2) {
2721       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2722       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
2723       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2724     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
2725       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2726       /* safeguard: avoid sharp decrease of rtol */
2727       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
2728       stol = PetscMax(rtol,stol);
2729       rtol = PetscMin(kctx->rtol_0,stol);
2730       /* safeguard: avoid oversolving */
2731       stol = kctx->gamma*(snes->ttol)/snes->norm;
2732       stol = PetscMax(rtol,stol);
2733       rtol = PetscMin(kctx->rtol_0,stol);
2734     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
2735   }
2736   /* safeguard: avoid rtol greater than one */
2737   rtol = PetscMin(rtol,kctx->rtol_max);
2738   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2739   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
2740   PetscFunctionReturn(0);
2741 }
2742 
2743 #undef __FUNCT__
2744 #define __FUNCT__ "SNESKSPEW_PostSolve"
2745 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
2746 {
2747   PetscErrorCode ierr;
2748   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2749   PCSide         pcside;
2750   Vec            lres;
2751 
2752   PetscFunctionBegin;
2753   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2754   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
2755   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
2756   if (kctx->version == 1) {
2757     ierr = KSPGetPreconditionerSide(ksp,&pcside);CHKERRQ(ierr);
2758     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
2759       /* KSP residual is true linear residual */
2760       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
2761     } else {
2762       /* KSP residual is preconditioned residual */
2763       /* compute true linear residual norm */
2764       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
2765       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
2766       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
2767       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
2768       ierr = VecDestroy(lres);CHKERRQ(ierr);
2769     }
2770   }
2771   PetscFunctionReturn(0);
2772 }
2773 
2774 #undef __FUNCT__
2775 #define __FUNCT__ "SNES_KSPSolve"
2776 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
2777 {
2778   PetscErrorCode ierr;
2779 
2780   PetscFunctionBegin;
2781   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
2782   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
2783   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
2784   PetscFunctionReturn(0);
2785 }
2786