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