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