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