xref: /petsc/src/snes/interface/snes.c (revision 83b03829eb385a3cc4c5beebc2599dacdb04a507)
1 #define PETSCSNES_DLL
2 
3 #include "private/snesimpl.h"      /*I "petscsnes.h"  I*/
4 
5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
6 PetscFList SNESList              = PETSC_NULL;
7 
8 /* Logging support */
9 PetscCookie PETSCSNES_DLLEXPORT SNES_COOKIE;
10 PetscLogEvent  SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval;
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "SNESSetFunctionDomainError"
14 /*@
15    SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not
16      in the functions domain. For example, negative pressure.
17 
18    Collective on SNES
19 
20    Input Parameters:
21 .  SNES - the SNES context
22 
23    Level: advanced
24 
25 .keywords: SNES, view
26 
27 .seealso: SNESCreate(), SNESSetFunction()
28 @*/
29 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunctionDomainError(SNES snes)
30 {
31   PetscFunctionBegin;
32   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
33   snes->domainerror = PETSC_TRUE;
34   PetscFunctionReturn(0);
35 }
36 
37 #undef __FUNCT__
38 #define __FUNCT__ "SNESView"
39 /*@C
40    SNESView - Prints the SNES data structure.
41 
42    Collective on SNES
43 
44    Input Parameters:
45 +  SNES - the SNES context
46 -  viewer - visualization context
47 
48    Options Database Key:
49 .  -snes_view - Calls SNESView() at end of SNESSolve()
50 
51    Notes:
52    The available visualization contexts include
53 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
54 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
55          output where only the first processor opens
56          the file.  All other processors send their
57          data to the first processor to print.
58 
59    The user can open an alternative visualization context with
60    PetscViewerASCIIOpen() - output to a specified file.
61 
62    Level: beginner
63 
64 .keywords: SNES, view
65 
66 .seealso: PetscViewerASCIIOpen()
67 @*/
68 PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer)
69 {
70   SNESKSPEW           *kctx;
71   PetscErrorCode      ierr;
72   KSP                 ksp;
73   const SNESType      type;
74   PetscTruth          iascii,isstring;
75 
76   PetscFunctionBegin;
77   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
78   if (!viewer) {
79     ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr);
80   }
81   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
82   PetscCheckSameComm(snes,1,viewer,2);
83 
84   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
85   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
86   if (iascii) {
87     if (((PetscObject)snes)->prefix) {
88       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",((PetscObject)snes)->prefix);CHKERRQ(ierr);
89     } else {
90       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
91     }
92     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
93     if (type) {
94       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
95     } else {
96       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
97     }
98     if (snes->ops->view) {
99       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
100       ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr);
101       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
102     }
103     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
104     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, solution=%G\n",
105                  snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr);
106     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr);
107     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr);
108     if (snes->ksp_ewconv) {
109       kctx = (SNESKSPEW *)snes->kspconvctx;
110       if (kctx) {
111         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr);
112         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
113         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
114       }
115     }
116   } else if (isstring) {
117     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
118     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
119   }
120   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
121   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
122   ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
123   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
124   PetscFunctionReturn(0);
125 }
126 
127 /*
128   We retain a list of functions that also take SNES command
129   line options. These are called at the end SNESSetFromOptions()
130 */
131 #define MAXSETFROMOPTIONS 5
132 static PetscInt numberofsetfromoptions;
133 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
134 
135 #undef __FUNCT__
136 #define __FUNCT__ "SNESAddOptionsChecker"
137 /*@C
138   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
139 
140   Not Collective
141 
142   Input Parameter:
143 . snescheck - function that checks for options
144 
145   Level: developer
146 
147 .seealso: SNESSetFromOptions()
148 @*/
149 PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
150 {
151   PetscFunctionBegin;
152   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
153     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
154   }
155   othersetfromoptions[numberofsetfromoptions++] = snescheck;
156   PetscFunctionReturn(0);
157 }
158 
159 #undef __FUNCT__
160 #define __FUNCT__ "SNESSetFromOptions"
161 /*@
162    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
163 
164    Collective on SNES
165 
166    Input Parameter:
167 .  snes - the SNES context
168 
169    Options Database Keys:
170 +  -snes_type <type> - ls, tr, umls, umtr, test
171 .  -snes_stol - convergence tolerance in terms of the norm
172                 of the change in the solution between steps
173 .  -snes_atol <abstol> - absolute tolerance of residual norm
174 .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
175 .  -snes_max_it <max_it> - maximum number of iterations
176 .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
177 .  -snes_max_fail <max_fail> - maximum number of failures
178 .  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
179 .  -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
180 .  -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
181 .  -snes_trtol <trtol> - trust region tolerance
182 .  -snes_no_convergence_test - skip convergence test in nonlinear
183                                solver; hence iterations will continue until max_it
184                                or some other criterion is reached. Saves expense
185                                of convergence test
186 .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
187                                        filename given prints to stdout
188 .  -snes_monitor_solution - plots solution at each iteration
189 .  -snes_monitor_residual - plots residual (not its norm) at each iteration
190 .  -snes_monitor_solution_update - plots update to solution at each iteration
191 .  -snes_monitor_draw - plots residual norm at each iteration
192 .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
193 .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
194 -  -snes_converged_reason - print the reason for convergence/divergence after each solve
195 
196     Options Database for Eisenstat-Walker method:
197 +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
198 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
199 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
200 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
201 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
202 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
203 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
204 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
205 
206    Notes:
207    To see all options, run your program with the -help option or consult
208    the users manual.
209 
210    Level: beginner
211 
212 .keywords: SNES, nonlinear, set, options, database
213 
214 .seealso: SNESSetOptionsPrefix()
215 @*/
216 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes)
217 {
218   PetscTruth              flg;
219   PetscInt                i,indx,lag;
220   const char              *deft = SNESLS;
221   const char              *convtests[] = {"default","skip"};
222   SNESKSPEW               *kctx = NULL;
223   char                    type[256], monfilename[PETSC_MAX_PATH_LEN];
224   PetscViewerASCIIMonitor monviewer;
225   PetscErrorCode          ierr;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
229 
230   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
231   ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
232     if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; }
233     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
234     if (flg) {
235       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
236     } else if (!((PetscObject)snes)->type_name) {
237       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
238     }
239     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     mf, mfOp, mfOp2;
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 = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) auxiliary options","SNES_EX");CHKERRQ(ierr);
1276     ierr = PetscOptionsName("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","MatCreateSNESMF",&mf);CHKERRQ(ierr);
1277     ierr = PetscOptionsName("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",&mfOp);CHKERRQ(ierr);
1278     ierr = PetscOptionsName("-snes_mf_operator2","Use a Matrix-Free (version 2) Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",&mfOp2);CHKERRQ(ierr);
1279   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1280 
1281   /*
1282       This version replaces the user provided Jacobian matrix with a
1283       matrix-free version but still employs the user-provided preconditioner matrix
1284   */
1285   if (mfOp) {
1286     Mat J;
1287     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
1288     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1289     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
1290     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
1291     ierr = MatDestroy(J);CHKERRQ(ierr);
1292   }
1293 
1294 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT)
1295   if (mfOp2) {
1296     Mat J;
1297     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1298     ierr = PetscInfo(snes,"Setting default matrix-free operator routines (version 2)\n");CHKERRQ(ierr);
1299     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
1300     ierr = MatDestroy(J);CHKERRQ(ierr);
1301   }
1302 #endif
1303 
1304   /*
1305       This version replaces both the user-provided Jacobian and the user-
1306       provided preconditioner matrix with the default matrix free version.
1307    */
1308   if (mf) {
1309     Mat  J;
1310     KSP ksp;
1311     PC   pc;
1312     PetscTruth flg;
1313     /* create and set matrix-free operator */
1314     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
1315     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1316     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
1317     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
1318     ierr = MatDestroy(J);CHKERRQ(ierr);
1319     /* force no preconditioner */
1320     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1321     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1322     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1323     if (!flg) {
1324       ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines;\nThat is no preconditioner is being used\n");CHKERRQ(ierr);
1325       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1326     }
1327   }
1328 
1329   if (!snes->vec_func && !snes->vec_rhs) {
1330     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1331   }
1332   if (!snes->ops->computefunction && !snes->vec_rhs) {
1333     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1334   }
1335   if (snes->vec_func == snes->vec_sol) {
1336     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1337   }
1338 
1339   if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);}
1340 
1341   if (snes->ops->setup) {
1342     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
1343   }
1344   snes->setupcalled = PETSC_TRUE;
1345   PetscFunctionReturn(0);
1346 }
1347 
1348 #undef __FUNCT__
1349 #define __FUNCT__ "SNESDestroy"
1350 /*@
1351    SNESDestroy - Destroys the nonlinear solver context that was created
1352    with SNESCreate().
1353 
1354    Collective on SNES
1355 
1356    Input Parameter:
1357 .  snes - the SNES context
1358 
1359    Level: beginner
1360 
1361 .keywords: SNES, nonlinear, destroy
1362 
1363 .seealso: SNESCreate(), SNESSolve()
1364 @*/
1365 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
1366 {
1367   PetscErrorCode ierr;
1368 
1369   PetscFunctionBegin;
1370   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1371   if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0);
1372 
1373   /* if memory was published with AMS then destroy it */
1374   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1375   if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);}
1376 
1377   if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);}
1378   if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);}
1379   if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);}
1380   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1381   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1382   if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);}
1383   ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);
1384   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1385   ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
1386   if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);}
1387   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
1388   PetscFunctionReturn(0);
1389 }
1390 
1391 /* ----------- Routines to set solver parameters ---------- */
1392 
1393 #undef __FUNCT__
1394 #define __FUNCT__ "SNESSetLagPreconditioner"
1395 /*@
1396    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
1397 
1398    Collective on SNES
1399 
1400    Input Parameters:
1401 +  snes - the SNES context
1402 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1403          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1404 
1405    Options Database Keys:
1406 .    -snes_lag_preconditioner <lag>
1407 
1408    Notes:
1409    The default is 1
1410    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1411    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
1412 
1413    Level: intermediate
1414 
1415 .keywords: SNES, nonlinear, set, convergence, tolerances
1416 
1417 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
1418 
1419 @*/
1420 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag)
1421 {
1422   PetscFunctionBegin;
1423   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1424   if (lag < -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1425   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1426   snes->lagpreconditioner = lag;
1427   PetscFunctionReturn(0);
1428 }
1429 
1430 #undef __FUNCT__
1431 #define __FUNCT__ "SNESGetLagPreconditioner"
1432 /*@
1433    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
1434 
1435    Collective on SNES
1436 
1437    Input Parameter:
1438 .  snes - the SNES context
1439 
1440    Output Parameter:
1441 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1442          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1443 
1444    Options Database Keys:
1445 .    -snes_lag_preconditioner <lag>
1446 
1447    Notes:
1448    The default is 1
1449    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1450 
1451    Level: intermediate
1452 
1453 .keywords: SNES, nonlinear, set, convergence, tolerances
1454 
1455 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
1456 
1457 @*/
1458 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
1459 {
1460   PetscFunctionBegin;
1461   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1462   *lag = snes->lagpreconditioner;
1463   PetscFunctionReturn(0);
1464 }
1465 
1466 #undef __FUNCT__
1467 #define __FUNCT__ "SNESSetLagJacobian"
1468 /*@
1469    SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
1470      often the preconditioner is rebuilt.
1471 
1472    Collective on SNES
1473 
1474    Input Parameters:
1475 +  snes - the SNES context
1476 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1477          the Jacobian is built etc. -2 means rebuild at next chance but then never again
1478 
1479    Options Database Keys:
1480 .    -snes_lag_jacobian <lag>
1481 
1482    Notes:
1483    The default is 1
1484    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1485    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
1486    at the next Newton step but never again (unless it is reset to another value)
1487 
1488    Level: intermediate
1489 
1490 .keywords: SNES, nonlinear, set, convergence, tolerances
1491 
1492 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
1493 
1494 @*/
1495 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagJacobian(SNES snes,PetscInt lag)
1496 {
1497   PetscFunctionBegin;
1498   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1499   if (lag < -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1500   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1501   snes->lagjacobian = lag;
1502   PetscFunctionReturn(0);
1503 }
1504 
1505 #undef __FUNCT__
1506 #define __FUNCT__ "SNESGetLagJacobian"
1507 /*@
1508    SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
1509 
1510    Collective on SNES
1511 
1512    Input Parameter:
1513 .  snes - the SNES context
1514 
1515    Output Parameter:
1516 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1517          the Jacobian is built etc.
1518 
1519    Options Database Keys:
1520 .    -snes_lag_jacobian <lag>
1521 
1522    Notes:
1523    The default is 1
1524    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1525 
1526    Level: intermediate
1527 
1528 .keywords: SNES, nonlinear, set, convergence, tolerances
1529 
1530 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
1531 
1532 @*/
1533 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagJacobian(SNES snes,PetscInt *lag)
1534 {
1535   PetscFunctionBegin;
1536   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1537   *lag = snes->lagjacobian;
1538   PetscFunctionReturn(0);
1539 }
1540 
1541 #undef __FUNCT__
1542 #define __FUNCT__ "SNESSetTolerances"
1543 /*@
1544    SNESSetTolerances - Sets various parameters used in convergence tests.
1545 
1546    Collective on SNES
1547 
1548    Input Parameters:
1549 +  snes - the SNES context
1550 .  abstol - absolute convergence tolerance
1551 .  rtol - relative convergence tolerance
1552 .  stol -  convergence tolerance in terms of the norm
1553            of the change in the solution between steps
1554 .  maxit - maximum number of iterations
1555 -  maxf - maximum number of function evaluations
1556 
1557    Options Database Keys:
1558 +    -snes_atol <abstol> - Sets abstol
1559 .    -snes_rtol <rtol> - Sets rtol
1560 .    -snes_stol <stol> - Sets stol
1561 .    -snes_max_it <maxit> - Sets maxit
1562 -    -snes_max_funcs <maxf> - Sets maxf
1563 
1564    Notes:
1565    The default maximum number of iterations is 50.
1566    The default maximum number of function evaluations is 1000.
1567 
1568    Level: intermediate
1569 
1570 .keywords: SNES, nonlinear, set, convergence, tolerances
1571 
1572 .seealso: SNESSetTrustRegionTolerance()
1573 @*/
1574 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
1575 {
1576   PetscFunctionBegin;
1577   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1578   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1579   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1580   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1581   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1582   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1583   PetscFunctionReturn(0);
1584 }
1585 
1586 #undef __FUNCT__
1587 #define __FUNCT__ "SNESGetTolerances"
1588 /*@
1589    SNESGetTolerances - Gets various parameters used in convergence tests.
1590 
1591    Not Collective
1592 
1593    Input Parameters:
1594 +  snes - the SNES context
1595 .  atol - absolute convergence tolerance
1596 .  rtol - relative convergence tolerance
1597 .  stol -  convergence tolerance in terms of the norm
1598            of the change in the solution between steps
1599 .  maxit - maximum number of iterations
1600 -  maxf - maximum number of function evaluations
1601 
1602    Notes:
1603    The user can specify PETSC_NULL for any parameter that is not needed.
1604 
1605    Level: intermediate
1606 
1607 .keywords: SNES, nonlinear, get, convergence, tolerances
1608 
1609 .seealso: SNESSetTolerances()
1610 @*/
1611 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
1612 {
1613   PetscFunctionBegin;
1614   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1615   if (atol)  *atol  = snes->abstol;
1616   if (rtol)  *rtol  = snes->rtol;
1617   if (stol)  *stol  = snes->xtol;
1618   if (maxit) *maxit = snes->max_its;
1619   if (maxf)  *maxf  = snes->max_funcs;
1620   PetscFunctionReturn(0);
1621 }
1622 
1623 #undef __FUNCT__
1624 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1625 /*@
1626    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1627 
1628    Collective on SNES
1629 
1630    Input Parameters:
1631 +  snes - the SNES context
1632 -  tol - tolerance
1633 
1634    Options Database Key:
1635 .  -snes_trtol <tol> - Sets tol
1636 
1637    Level: intermediate
1638 
1639 .keywords: SNES, nonlinear, set, trust region, tolerance
1640 
1641 .seealso: SNESSetTolerances()
1642 @*/
1643 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1644 {
1645   PetscFunctionBegin;
1646   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1647   snes->deltatol = tol;
1648   PetscFunctionReturn(0);
1649 }
1650 
1651 /*
1652    Duplicate the lg monitors for SNES from KSP; for some reason with
1653    dynamic libraries things don't work under Sun4 if we just use
1654    macros instead of functions
1655 */
1656 #undef __FUNCT__
1657 #define __FUNCT__ "SNESMonitorLG"
1658 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1659 {
1660   PetscErrorCode ierr;
1661 
1662   PetscFunctionBegin;
1663   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1664   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1665   PetscFunctionReturn(0);
1666 }
1667 
1668 #undef __FUNCT__
1669 #define __FUNCT__ "SNESMonitorLGCreate"
1670 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1671 {
1672   PetscErrorCode ierr;
1673 
1674   PetscFunctionBegin;
1675   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1676   PetscFunctionReturn(0);
1677 }
1678 
1679 #undef __FUNCT__
1680 #define __FUNCT__ "SNESMonitorLGDestroy"
1681 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw)
1682 {
1683   PetscErrorCode ierr;
1684 
1685   PetscFunctionBegin;
1686   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1687   PetscFunctionReturn(0);
1688 }
1689 
1690 extern PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
1691 #undef __FUNCT__
1692 #define __FUNCT__ "SNESMonitorLGRange"
1693 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
1694 {
1695   PetscDrawLG      lg;
1696   PetscErrorCode   ierr;
1697   PetscReal        x,y,per;
1698   PetscViewer      v = (PetscViewer)monctx;
1699   static PetscReal prev; /* should be in the context */
1700   PetscDraw        draw;
1701   PetscFunctionBegin;
1702   if (!monctx) {
1703     MPI_Comm    comm;
1704 
1705     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1706     v      = PETSC_VIEWER_DRAW_(comm);
1707   }
1708   ierr   = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr);
1709   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1710   ierr   = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1711   ierr   = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr);
1712   x = (PetscReal) n;
1713   if (rnorm > 0.0) y = log10(rnorm); else y = -15.0;
1714   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1715   if (n < 20 || !(n % 5)) {
1716     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1717   }
1718 
1719   ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr);
1720   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1721   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1722   ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr);
1723   ierr =  SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr);
1724   x = (PetscReal) n;
1725   y = 100.0*per;
1726   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1727   if (n < 20 || !(n % 5)) {
1728     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1729   }
1730 
1731   ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr);
1732   if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1733   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1734   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr);
1735   x = (PetscReal) n;
1736   y = (prev - rnorm)/prev;
1737   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1738   if (n < 20 || !(n % 5)) {
1739     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1740   }
1741 
1742   ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr);
1743   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1744   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1745   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
1746   x = (PetscReal) n;
1747   y = (prev - rnorm)/(prev*per);
1748   if (n > 2) { /*skip initial crazy value */
1749     ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1750   }
1751   if (n < 20 || !(n % 5)) {
1752     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1753   }
1754   prev = rnorm;
1755   PetscFunctionReturn(0);
1756 }
1757 
1758 #undef __FUNCT__
1759 #define __FUNCT__ "SNESMonitorLGRangeCreate"
1760 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1761 {
1762   PetscErrorCode ierr;
1763 
1764   PetscFunctionBegin;
1765   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1766   PetscFunctionReturn(0);
1767 }
1768 
1769 #undef __FUNCT__
1770 #define __FUNCT__ "SNESMonitorLGRangeDestroy"
1771 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeDestroy(PetscDrawLG draw)
1772 {
1773   PetscErrorCode ierr;
1774 
1775   PetscFunctionBegin;
1776   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1777   PetscFunctionReturn(0);
1778 }
1779 
1780 /* ------------ Routines to set performance monitoring options ----------- */
1781 
1782 #undef __FUNCT__
1783 #define __FUNCT__ "SNESMonitorSet"
1784 /*@C
1785    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
1786    iteration of the nonlinear solver to display the iteration's
1787    progress.
1788 
1789    Collective on SNES
1790 
1791    Input Parameters:
1792 +  snes - the SNES context
1793 .  func - monitoring routine
1794 .  mctx - [optional] user-defined context for private data for the
1795           monitor routine (use PETSC_NULL if no context is desired)
1796 -  monitordestroy - [optional] routine that frees monitor context
1797           (may be PETSC_NULL)
1798 
1799    Calling sequence of func:
1800 $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1801 
1802 +    snes - the SNES context
1803 .    its - iteration number
1804 .    norm - 2-norm function value (may be estimated)
1805 -    mctx - [optional] monitoring context
1806 
1807    Options Database Keys:
1808 +    -snes_monitor        - sets SNESMonitorDefault()
1809 .    -snes_monitor_draw    - sets line graph monitor,
1810                             uses SNESMonitorLGCreate()
1811 _    -snes_monitor_cancel - cancels all monitors that have
1812                             been hardwired into a code by
1813                             calls to SNESMonitorSet(), but
1814                             does not cancel those set via
1815                             the options database.
1816 
1817    Notes:
1818    Several different monitoring routines may be set by calling
1819    SNESMonitorSet() multiple times; all will be called in the
1820    order in which they were set.
1821 
1822    Fortran notes: Only a single monitor function can be set for each SNES object
1823 
1824    Level: intermediate
1825 
1826 .keywords: SNES, nonlinear, set, monitor
1827 
1828 .seealso: SNESMonitorDefault(), SNESMonitorCancel()
1829 @*/
1830 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
1831 {
1832   PetscInt i;
1833 
1834   PetscFunctionBegin;
1835   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1836   if (snes->numbermonitors >= MAXSNESMONITORS) {
1837     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1838   }
1839   for (i=0; i<snes->numbermonitors;i++) {
1840     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0);
1841 
1842     /* check if both default monitors that share common ASCII viewer */
1843     if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) {
1844       if (mctx && snes->monitorcontext[i]) {
1845         PetscErrorCode          ierr;
1846         PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx;
1847         PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i];
1848         if (viewer1->viewer == viewer2->viewer) {
1849           ierr = (*monitordestroy)(mctx);CHKERRQ(ierr);
1850           PetscFunctionReturn(0);
1851         }
1852       }
1853     }
1854   }
1855   snes->monitor[snes->numbermonitors]           = monitor;
1856   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1857   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1858   PetscFunctionReturn(0);
1859 }
1860 
1861 #undef __FUNCT__
1862 #define __FUNCT__ "SNESMonitorCancel"
1863 /*@C
1864    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
1865 
1866    Collective on SNES
1867 
1868    Input Parameters:
1869 .  snes - the SNES context
1870 
1871    Options Database Key:
1872 .  -snes_monitor_cancel - cancels all monitors that have been hardwired
1873     into a code by calls to SNESMonitorSet(), but does not cancel those
1874     set via the options database
1875 
1876    Notes:
1877    There is no way to clear one specific monitor from a SNES object.
1878 
1879    Level: intermediate
1880 
1881 .keywords: SNES, nonlinear, set, monitor
1882 
1883 .seealso: SNESMonitorDefault(), SNESMonitorSet()
1884 @*/
1885 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes)
1886 {
1887   PetscErrorCode ierr;
1888   PetscInt       i;
1889 
1890   PetscFunctionBegin;
1891   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1892   for (i=0; i<snes->numbermonitors; i++) {
1893     if (snes->monitordestroy[i]) {
1894       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1895     }
1896   }
1897   snes->numbermonitors = 0;
1898   PetscFunctionReturn(0);
1899 }
1900 
1901 #undef __FUNCT__
1902 #define __FUNCT__ "SNESSetConvergenceTest"
1903 /*@C
1904    SNESSetConvergenceTest - Sets the function that is to be used
1905    to test for convergence of the nonlinear iterative solution.
1906 
1907    Collective on SNES
1908 
1909    Input Parameters:
1910 +  snes - the SNES context
1911 .  func - routine to test for convergence
1912 .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
1913 -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
1914 
1915    Calling sequence of func:
1916 $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1917 
1918 +    snes - the SNES context
1919 .    it - current iteration (0 is the first and is before any Newton step)
1920 .    cctx - [optional] convergence context
1921 .    reason - reason for convergence/divergence
1922 .    xnorm - 2-norm of current iterate
1923 .    gnorm - 2-norm of current step
1924 -    f - 2-norm of function
1925 
1926    Level: advanced
1927 
1928 .keywords: SNES, nonlinear, set, convergence, test
1929 
1930 .seealso: SNESDefaultConverged(), SNESSkipConverged()
1931 @*/
1932 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1933 {
1934   PetscErrorCode ierr;
1935 
1936   PetscFunctionBegin;
1937   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1938   if (!func) func = SNESSkipConverged;
1939   if (snes->ops->convergeddestroy) {
1940     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
1941   }
1942   snes->ops->converged        = func;
1943   snes->ops->convergeddestroy = destroy;
1944   snes->cnvP                  = cctx;
1945   PetscFunctionReturn(0);
1946 }
1947 
1948 #undef __FUNCT__
1949 #define __FUNCT__ "SNESGetConvergedReason"
1950 /*@
1951    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1952 
1953    Not Collective
1954 
1955    Input Parameter:
1956 .  snes - the SNES context
1957 
1958    Output Parameter:
1959 .  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
1960             manual pages for the individual convergence tests for complete lists
1961 
1962    Level: intermediate
1963 
1964    Notes: Can only be called after the call the SNESSolve() is complete.
1965 
1966 .keywords: SNES, nonlinear, set, convergence, test
1967 
1968 .seealso: SNESSetConvergenceTest(), SNESConvergedReason
1969 @*/
1970 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1971 {
1972   PetscFunctionBegin;
1973   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1974   PetscValidPointer(reason,2);
1975   *reason = snes->reason;
1976   PetscFunctionReturn(0);
1977 }
1978 
1979 #undef __FUNCT__
1980 #define __FUNCT__ "SNESSetConvergenceHistory"
1981 /*@
1982    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1983 
1984    Collective on SNES
1985 
1986    Input Parameters:
1987 +  snes - iterative context obtained from SNESCreate()
1988 .  a   - array to hold history
1989 .  its - integer array holds the number of linear iterations for each solve.
1990 .  na  - size of a and its
1991 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1992            else it continues storing new values for new nonlinear solves after the old ones
1993 
1994    Notes:
1995    If set, this array will contain the function norms computed
1996    at each step.
1997 
1998    This routine is useful, e.g., when running a code for purposes
1999    of accurate performance monitoring, when no I/O should be done
2000    during the section of code that is being timed.
2001 
2002    Level: intermediate
2003 
2004 .keywords: SNES, set, convergence, history
2005 
2006 .seealso: SNESGetConvergenceHistory()
2007 
2008 @*/
2009 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset)
2010 {
2011   PetscFunctionBegin;
2012   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2013   if (na)  PetscValidScalarPointer(a,2);
2014   if (its) PetscValidIntPointer(its,3);
2015   snes->conv_hist       = a;
2016   snes->conv_hist_its   = its;
2017   snes->conv_hist_max   = na;
2018   snes->conv_hist_len   = 0;
2019   snes->conv_hist_reset = reset;
2020   PetscFunctionReturn(0);
2021 }
2022 
2023 #undef __FUNCT__
2024 #define __FUNCT__ "SNESGetConvergenceHistory"
2025 /*@C
2026    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
2027 
2028    Collective on SNES
2029 
2030    Input Parameter:
2031 .  snes - iterative context obtained from SNESCreate()
2032 
2033    Output Parameters:
2034 .  a   - array to hold history
2035 .  its - integer array holds the number of linear iterations (or
2036          negative if not converged) for each solve.
2037 -  na  - size of a and its
2038 
2039    Notes:
2040     The calling sequence for this routine in Fortran is
2041 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
2042 
2043    This routine is useful, e.g., when running a code for purposes
2044    of accurate performance monitoring, when no I/O should be done
2045    during the section of code that is being timed.
2046 
2047    Level: intermediate
2048 
2049 .keywords: SNES, get, convergence, history
2050 
2051 .seealso: SNESSetConvergencHistory()
2052 
2053 @*/
2054 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
2055 {
2056   PetscFunctionBegin;
2057   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2058   if (a)   *a   = snes->conv_hist;
2059   if (its) *its = snes->conv_hist_its;
2060   if (na)  *na  = snes->conv_hist_len;
2061   PetscFunctionReturn(0);
2062 }
2063 
2064 #undef __FUNCT__
2065 #define __FUNCT__ "SNESSetUpdate"
2066 /*@C
2067   SNESSetUpdate - Sets the general-purpose update function called
2068   at the beginning o every iteration of the nonlinear solve. Specifically
2069   it is called just before the Jacobian is "evaluated".
2070 
2071   Collective on SNES
2072 
2073   Input Parameters:
2074 . snes - The nonlinear solver context
2075 . func - The function
2076 
2077   Calling sequence of func:
2078 . func (SNES snes, PetscInt step);
2079 
2080 . step - The current step of the iteration
2081 
2082   Level: intermediate
2083 
2084 .keywords: SNES, update
2085 
2086 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
2087 @*/
2088 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
2089 {
2090   PetscFunctionBegin;
2091   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
2092   snes->ops->update = func;
2093   PetscFunctionReturn(0);
2094 }
2095 
2096 #undef __FUNCT__
2097 #define __FUNCT__ "SNESDefaultUpdate"
2098 /*@
2099   SNESDefaultUpdate - The default update function which does nothing.
2100 
2101   Not collective
2102 
2103   Input Parameters:
2104 . snes - The nonlinear solver context
2105 . step - The current step of the iteration
2106 
2107   Level: intermediate
2108 
2109 .keywords: SNES, update
2110 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
2111 @*/
2112 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
2113 {
2114   PetscFunctionBegin;
2115   PetscFunctionReturn(0);
2116 }
2117 
2118 #undef __FUNCT__
2119 #define __FUNCT__ "SNESScaleStep_Private"
2120 /*
2121    SNESScaleStep_Private - Scales a step so that its length is less than the
2122    positive parameter delta.
2123 
2124     Input Parameters:
2125 +   snes - the SNES context
2126 .   y - approximate solution of linear system
2127 .   fnorm - 2-norm of current function
2128 -   delta - trust region size
2129 
2130     Output Parameters:
2131 +   gpnorm - predicted function norm at the new point, assuming local
2132     linearization.  The value is zero if the step lies within the trust
2133     region, and exceeds zero otherwise.
2134 -   ynorm - 2-norm of the step
2135 
2136     Note:
2137     For non-trust region methods such as SNESLS, the parameter delta
2138     is set to be the maximum allowable step size.
2139 
2140 .keywords: SNES, nonlinear, scale, step
2141 */
2142 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
2143 {
2144   PetscReal      nrm;
2145   PetscScalar    cnorm;
2146   PetscErrorCode ierr;
2147 
2148   PetscFunctionBegin;
2149   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2150   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
2151   PetscCheckSameComm(snes,1,y,2);
2152 
2153   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
2154   if (nrm > *delta) {
2155      nrm = *delta/nrm;
2156      *gpnorm = (1.0 - nrm)*(*fnorm);
2157      cnorm = nrm;
2158      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
2159      *ynorm = *delta;
2160   } else {
2161      *gpnorm = 0.0;
2162      *ynorm = nrm;
2163   }
2164   PetscFunctionReturn(0);
2165 }
2166 
2167 #undef __FUNCT__
2168 #define __FUNCT__ "SNESSolve"
2169 /*@C
2170    SNESSolve - Solves a nonlinear system F(x) = b.
2171    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
2172 
2173    Collective on SNES
2174 
2175    Input Parameters:
2176 +  snes - the SNES context
2177 .  b - the constant part of the equation, or PETSC_NULL to use zero.
2178 -  x - the solution vector.
2179 
2180    Notes:
2181    The user should initialize the vector,x, with the initial guess
2182    for the nonlinear solve prior to calling SNESSolve.  In particular,
2183    to employ an initial guess of zero, the user should explicitly set
2184    this vector to zero by calling VecSet().
2185 
2186    Level: beginner
2187 
2188 .keywords: SNES, nonlinear, solve
2189 
2190 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian()
2191 @*/
2192 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
2193 {
2194   PetscErrorCode ierr;
2195   PetscTruth     flg;
2196   char           filename[PETSC_MAX_PATH_LEN];
2197   PetscViewer    viewer;
2198 
2199   PetscFunctionBegin;
2200   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2201   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
2202   PetscCheckSameComm(snes,1,x,3);
2203   if (b) PetscValidHeaderSpecific(b,VEC_COOKIE,2);
2204   if (b) PetscCheckSameComm(snes,1,b,2);
2205 
2206   /* set solution vector */
2207   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
2208   if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); }
2209   snes->vec_sol = x;
2210   /* set afine vector if provided */
2211   if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
2212   if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); }
2213   snes->vec_rhs = b;
2214 
2215   if (!snes->vec_func && snes->vec_rhs) {
2216     ierr = VecDuplicate(b, &snes->vec_func);CHKERRQ(ierr);
2217   }
2218 
2219   ierr = SNESSetUp(snes);CHKERRQ(ierr);
2220 
2221   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
2222   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
2223 
2224   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2225   ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
2226   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2227   if (snes->domainerror){
2228     snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
2229     snes->domainerror = PETSC_FALSE;
2230   }
2231 
2232   if (!snes->reason) {
2233     SETERRQ(PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
2234   }
2235 
2236   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
2237   if (flg && !PetscPreLoadingOn) {
2238     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
2239     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
2240     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
2241   }
2242 
2243   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
2244   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
2245   if (snes->printreason) {
2246     if (snes->reason > 0) {
2247       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2248     } else {
2249       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2250     }
2251   }
2252 
2253   PetscFunctionReturn(0);
2254 }
2255 
2256 /* --------- Internal routines for SNES Package --------- */
2257 
2258 #undef __FUNCT__
2259 #define __FUNCT__ "SNESSetType"
2260 /*@C
2261    SNESSetType - Sets the method for the nonlinear solver.
2262 
2263    Collective on SNES
2264 
2265    Input Parameters:
2266 +  snes - the SNES context
2267 -  type - a known method
2268 
2269    Options Database Key:
2270 .  -snes_type <type> - Sets the method; use -help for a list
2271    of available methods (for instance, ls or tr)
2272 
2273    Notes:
2274    See "petsc/include/petscsnes.h" for available methods (for instance)
2275 +    SNESLS - Newton's method with line search
2276      (systems of nonlinear equations)
2277 .    SNESTR - Newton's method with trust region
2278      (systems of nonlinear equations)
2279 
2280   Normally, it is best to use the SNESSetFromOptions() command and then
2281   set the SNES solver type from the options database rather than by using
2282   this routine.  Using the options database provides the user with
2283   maximum flexibility in evaluating the many nonlinear solvers.
2284   The SNESSetType() routine is provided for those situations where it
2285   is necessary to set the nonlinear solver independently of the command
2286   line or options database.  This might be the case, for example, when
2287   the choice of solver changes during the execution of the program,
2288   and the user's application is taking responsibility for choosing the
2289   appropriate method.
2290 
2291   Level: intermediate
2292 
2293 .keywords: SNES, set, type
2294 
2295 .seealso: SNESType, SNESCreate()
2296 
2297 @*/
2298 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type)
2299 {
2300   PetscErrorCode ierr,(*r)(SNES);
2301   PetscTruth     match;
2302 
2303   PetscFunctionBegin;
2304   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2305   PetscValidCharPointer(type,2);
2306 
2307   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
2308   if (match) PetscFunctionReturn(0);
2309 
2310   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
2311   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
2312   /* Destroy the previous private SNES context */
2313   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
2314   /* Reinitialize function pointers in SNESOps structure */
2315   snes->ops->setup          = 0;
2316   snes->ops->solve          = 0;
2317   snes->ops->view           = 0;
2318   snes->ops->setfromoptions = 0;
2319   snes->ops->destroy        = 0;
2320   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
2321   snes->setupcalled = PETSC_FALSE;
2322   ierr = (*r)(snes);CHKERRQ(ierr);
2323   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
2324   PetscFunctionReturn(0);
2325 }
2326 
2327 
2328 /* --------------------------------------------------------------------- */
2329 #undef __FUNCT__
2330 #define __FUNCT__ "SNESRegisterDestroy"
2331 /*@
2332    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2333    registered by SNESRegisterDynamic().
2334 
2335    Not Collective
2336 
2337    Level: advanced
2338 
2339 .keywords: SNES, nonlinear, register, destroy
2340 
2341 .seealso: SNESRegisterAll(), SNESRegisterAll()
2342 @*/
2343 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
2344 {
2345   PetscErrorCode ierr;
2346 
2347   PetscFunctionBegin;
2348   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
2349   SNESRegisterAllCalled = PETSC_FALSE;
2350   PetscFunctionReturn(0);
2351 }
2352 
2353 #undef __FUNCT__
2354 #define __FUNCT__ "SNESGetType"
2355 /*@C
2356    SNESGetType - Gets the SNES method type and name (as a string).
2357 
2358    Not Collective
2359 
2360    Input Parameter:
2361 .  snes - nonlinear solver context
2362 
2363    Output Parameter:
2364 .  type - SNES method (a character string)
2365 
2366    Level: intermediate
2367 
2368 .keywords: SNES, nonlinear, get, type, name
2369 @*/
2370 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type)
2371 {
2372   PetscFunctionBegin;
2373   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2374   PetscValidPointer(type,2);
2375   *type = ((PetscObject)snes)->type_name;
2376   PetscFunctionReturn(0);
2377 }
2378 
2379 #undef __FUNCT__
2380 #define __FUNCT__ "SNESGetSolution"
2381 /*@
2382    SNESGetSolution - Returns the vector where the approximate solution is
2383    stored.
2384 
2385    Not Collective, but Vec is parallel if SNES is parallel
2386 
2387    Input Parameter:
2388 .  snes - the SNES context
2389 
2390    Output Parameter:
2391 .  x - the solution
2392 
2393    Level: intermediate
2394 
2395 .keywords: SNES, nonlinear, get, solution
2396 
2397 .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
2398 @*/
2399 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
2400 {
2401   PetscFunctionBegin;
2402   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2403   PetscValidPointer(x,2);
2404   *x = snes->vec_sol;
2405   PetscFunctionReturn(0);
2406 }
2407 
2408 #undef __FUNCT__
2409 #define __FUNCT__ "SNESGetSolutionUpdate"
2410 /*@
2411    SNESGetSolutionUpdate - Returns the vector where the solution update is
2412    stored.
2413 
2414    Not Collective, but Vec is parallel if SNES is parallel
2415 
2416    Input Parameter:
2417 .  snes - the SNES context
2418 
2419    Output Parameter:
2420 .  x - the solution update
2421 
2422    Level: advanced
2423 
2424 .keywords: SNES, nonlinear, get, solution, update
2425 
2426 .seealso: SNESGetSolution(), SNESGetFunction()
2427 @*/
2428 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
2429 {
2430   PetscFunctionBegin;
2431   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2432   PetscValidPointer(x,2);
2433   *x = snes->vec_sol_update;
2434   PetscFunctionReturn(0);
2435 }
2436 
2437 #undef __FUNCT__
2438 #define __FUNCT__ "SNESGetFunction"
2439 /*@C
2440    SNESGetFunction - Returns the vector where the function is stored.
2441 
2442    Not Collective, but Vec is parallel if SNES is parallel
2443 
2444    Input Parameter:
2445 .  snes - the SNES context
2446 
2447    Output Parameter:
2448 +  r - the function (or PETSC_NULL)
2449 .  func - the function (or PETSC_NULL)
2450 -  ctx - the function context (or PETSC_NULL)
2451 
2452    Level: advanced
2453 
2454 .keywords: SNES, nonlinear, get, function
2455 
2456 .seealso: SNESSetFunction(), SNESGetSolution()
2457 @*/
2458 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
2459 {
2460   PetscFunctionBegin;
2461   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2462   if (r)    *r    = snes->vec_func;
2463   if (func) *func = snes->ops->computefunction;
2464   if (ctx)  *ctx  = snes->funP;
2465   PetscFunctionReturn(0);
2466 }
2467 
2468 #undef __FUNCT__
2469 #define __FUNCT__ "SNESSetOptionsPrefix"
2470 /*@C
2471    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2472    SNES options in the database.
2473 
2474    Collective on SNES
2475 
2476    Input Parameter:
2477 +  snes - the SNES context
2478 -  prefix - the prefix to prepend to all option names
2479 
2480    Notes:
2481    A hyphen (-) must NOT be given at the beginning of the prefix name.
2482    The first character of all runtime options is AUTOMATICALLY the hyphen.
2483 
2484    Level: advanced
2485 
2486 .keywords: SNES, set, options, prefix, database
2487 
2488 .seealso: SNESSetFromOptions()
2489 @*/
2490 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
2491 {
2492   PetscErrorCode ierr;
2493 
2494   PetscFunctionBegin;
2495   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2496   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2497   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2498   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2499   PetscFunctionReturn(0);
2500 }
2501 
2502 #undef __FUNCT__
2503 #define __FUNCT__ "SNESAppendOptionsPrefix"
2504 /*@C
2505    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2506    SNES options in the database.
2507 
2508    Collective on SNES
2509 
2510    Input Parameters:
2511 +  snes - the SNES context
2512 -  prefix - the prefix to prepend to all option names
2513 
2514    Notes:
2515    A hyphen (-) must NOT be given at the beginning of the prefix name.
2516    The first character of all runtime options is AUTOMATICALLY the hyphen.
2517 
2518    Level: advanced
2519 
2520 .keywords: SNES, append, options, prefix, database
2521 
2522 .seealso: SNESGetOptionsPrefix()
2523 @*/
2524 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
2525 {
2526   PetscErrorCode ierr;
2527 
2528   PetscFunctionBegin;
2529   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2530   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2531   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2532   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2533   PetscFunctionReturn(0);
2534 }
2535 
2536 #undef __FUNCT__
2537 #define __FUNCT__ "SNESGetOptionsPrefix"
2538 /*@C
2539    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2540    SNES options in the database.
2541 
2542    Not Collective
2543 
2544    Input Parameter:
2545 .  snes - the SNES context
2546 
2547    Output Parameter:
2548 .  prefix - pointer to the prefix string used
2549 
2550    Notes: On the fortran side, the user should pass in a string 'prefix' of
2551    sufficient length to hold the prefix.
2552 
2553    Level: advanced
2554 
2555 .keywords: SNES, get, options, prefix, database
2556 
2557 .seealso: SNESAppendOptionsPrefix()
2558 @*/
2559 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
2560 {
2561   PetscErrorCode ierr;
2562 
2563   PetscFunctionBegin;
2564   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2565   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2566   PetscFunctionReturn(0);
2567 }
2568 
2569 
2570 #undef __FUNCT__
2571 #define __FUNCT__ "SNESRegister"
2572 /*@C
2573   SNESRegister - See SNESRegisterDynamic()
2574 
2575   Level: advanced
2576 @*/
2577 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2578 {
2579   char           fullname[PETSC_MAX_PATH_LEN];
2580   PetscErrorCode ierr;
2581 
2582   PetscFunctionBegin;
2583   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2584   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2585   PetscFunctionReturn(0);
2586 }
2587 
2588 #undef __FUNCT__
2589 #define __FUNCT__ "SNESTestLocalMin"
2590 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2591 {
2592   PetscErrorCode ierr;
2593   PetscInt       N,i,j;
2594   Vec            u,uh,fh;
2595   PetscScalar    value;
2596   PetscReal      norm;
2597 
2598   PetscFunctionBegin;
2599   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2600   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2601   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2602 
2603   /* currently only works for sequential */
2604   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2605   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2606   for (i=0; i<N; i++) {
2607     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2608     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2609     for (j=-10; j<11; j++) {
2610       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2611       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2612       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2613       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2614       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2615       value = -value;
2616       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2617     }
2618   }
2619   ierr = VecDestroy(uh);CHKERRQ(ierr);
2620   ierr = VecDestroy(fh);CHKERRQ(ierr);
2621   PetscFunctionReturn(0);
2622 }
2623 
2624 #undef __FUNCT__
2625 #define __FUNCT__ "SNESKSPSetUseEW"
2626 /*@
2627    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
2628    computing relative tolerance for linear solvers within an inexact
2629    Newton method.
2630 
2631    Collective on SNES
2632 
2633    Input Parameters:
2634 +  snes - SNES context
2635 -  flag - PETSC_TRUE or PETSC_FALSE
2636 
2637    Notes:
2638    Currently, the default is to use a constant relative tolerance for
2639    the inner linear solvers.  Alternatively, one can use the
2640    Eisenstat-Walker method, where the relative convergence tolerance
2641    is reset at each Newton iteration according progress of the nonlinear
2642    solver.
2643 
2644    Level: advanced
2645 
2646    Reference:
2647    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2648    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2649 
2650 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2651 
2652 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2653 @*/
2654 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag)
2655 {
2656   PetscFunctionBegin;
2657   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2658   snes->ksp_ewconv = flag;
2659   PetscFunctionReturn(0);
2660 }
2661 
2662 #undef __FUNCT__
2663 #define __FUNCT__ "SNESKSPGetUseEW"
2664 /*@
2665    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
2666    for computing relative tolerance for linear solvers within an
2667    inexact Newton method.
2668 
2669    Not Collective
2670 
2671    Input Parameter:
2672 .  snes - SNES context
2673 
2674    Output Parameter:
2675 .  flag - PETSC_TRUE or PETSC_FALSE
2676 
2677    Notes:
2678    Currently, the default is to use a constant relative tolerance for
2679    the inner linear solvers.  Alternatively, one can use the
2680    Eisenstat-Walker method, where the relative convergence tolerance
2681    is reset at each Newton iteration according progress of the nonlinear
2682    solver.
2683 
2684    Level: advanced
2685 
2686    Reference:
2687    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2688    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2689 
2690 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2691 
2692 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2693 @*/
2694 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag)
2695 {
2696   PetscFunctionBegin;
2697   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2698   PetscValidPointer(flag,2);
2699   *flag = snes->ksp_ewconv;
2700   PetscFunctionReturn(0);
2701 }
2702 
2703 #undef __FUNCT__
2704 #define __FUNCT__ "SNESKSPSetParametersEW"
2705 /*@
2706    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
2707    convergence criteria for the linear solvers within an inexact
2708    Newton method.
2709 
2710    Collective on SNES
2711 
2712    Input Parameters:
2713 +    snes - SNES context
2714 .    version - version 1, 2 (default is 2) or 3
2715 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2716 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2717 .    gamma - multiplicative factor for version 2 rtol computation
2718              (0 <= gamma2 <= 1)
2719 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2720 .    alpha2 - power for safeguard
2721 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2722 
2723    Note:
2724    Version 3 was contributed by Luis Chacon, June 2006.
2725 
2726    Use PETSC_DEFAULT to retain the default for any of the parameters.
2727 
2728    Level: advanced
2729 
2730    Reference:
2731    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2732    inexact Newton method", Utah State University Math. Stat. Dept. Res.
2733    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
2734 
2735 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
2736 
2737 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
2738 @*/
2739 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
2740 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
2741 {
2742   SNESKSPEW *kctx;
2743   PetscFunctionBegin;
2744   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2745   kctx = (SNESKSPEW*)snes->kspconvctx;
2746   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2747 
2748   if (version != PETSC_DEFAULT)   kctx->version   = version;
2749   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
2750   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
2751   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
2752   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
2753   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
2754   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
2755 
2756   if (kctx->version < 1 || kctx->version > 3) {
2757     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
2758   }
2759   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
2760     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
2761   }
2762   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
2763     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
2764   }
2765   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
2766     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
2767   }
2768   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
2769     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
2770   }
2771   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
2772     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
2773   }
2774   PetscFunctionReturn(0);
2775 }
2776 
2777 #undef __FUNCT__
2778 #define __FUNCT__ "SNESKSPGetParametersEW"
2779 /*@
2780    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
2781    convergence criteria for the linear solvers within an inexact
2782    Newton method.
2783 
2784    Not Collective
2785 
2786    Input Parameters:
2787      snes - SNES context
2788 
2789    Output Parameters:
2790 +    version - version 1, 2 (default is 2) or 3
2791 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2792 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2793 .    gamma - multiplicative factor for version 2 rtol computation
2794              (0 <= gamma2 <= 1)
2795 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2796 .    alpha2 - power for safeguard
2797 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2798 
2799    Level: advanced
2800 
2801 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
2802 
2803 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
2804 @*/
2805 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
2806 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
2807 {
2808   SNESKSPEW *kctx;
2809   PetscFunctionBegin;
2810   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2811   kctx = (SNESKSPEW*)snes->kspconvctx;
2812   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2813   if(version)   *version   = kctx->version;
2814   if(rtol_0)    *rtol_0    = kctx->rtol_0;
2815   if(rtol_max)  *rtol_max  = kctx->rtol_max;
2816   if(gamma)     *gamma     = kctx->gamma;
2817   if(alpha)     *alpha     = kctx->alpha;
2818   if(alpha2)    *alpha2    = kctx->alpha2;
2819   if(threshold) *threshold = kctx->threshold;
2820   PetscFunctionReturn(0);
2821 }
2822 
2823 #undef __FUNCT__
2824 #define __FUNCT__ "SNESKSPEW_PreSolve"
2825 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
2826 {
2827   PetscErrorCode ierr;
2828   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2829   PetscReal      rtol=PETSC_DEFAULT,stol;
2830 
2831   PetscFunctionBegin;
2832   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2833   if (!snes->iter) { /* first time in, so use the original user rtol */
2834     rtol = kctx->rtol_0;
2835   } else {
2836     if (kctx->version == 1) {
2837       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
2838       if (rtol < 0.0) rtol = -rtol;
2839       stol = pow(kctx->rtol_last,kctx->alpha2);
2840       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2841     } else if (kctx->version == 2) {
2842       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2843       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
2844       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2845     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
2846       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2847       /* safeguard: avoid sharp decrease of rtol */
2848       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
2849       stol = PetscMax(rtol,stol);
2850       rtol = PetscMin(kctx->rtol_0,stol);
2851       /* safeguard: avoid oversolving */
2852       stol = kctx->gamma*(snes->ttol)/snes->norm;
2853       stol = PetscMax(rtol,stol);
2854       rtol = PetscMin(kctx->rtol_0,stol);
2855     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
2856   }
2857   /* safeguard: avoid rtol greater than one */
2858   rtol = PetscMin(rtol,kctx->rtol_max);
2859   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2860   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
2861   PetscFunctionReturn(0);
2862 }
2863 
2864 #undef __FUNCT__
2865 #define __FUNCT__ "SNESKSPEW_PostSolve"
2866 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
2867 {
2868   PetscErrorCode ierr;
2869   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2870   PCSide         pcside;
2871   Vec            lres;
2872 
2873   PetscFunctionBegin;
2874   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2875   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
2876   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
2877   if (kctx->version == 1) {
2878     ierr = KSPGetPreconditionerSide(ksp,&pcside);CHKERRQ(ierr);
2879     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
2880       /* KSP residual is true linear residual */
2881       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
2882     } else {
2883       /* KSP residual is preconditioned residual */
2884       /* compute true linear residual norm */
2885       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
2886       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
2887       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
2888       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
2889       ierr = VecDestroy(lres);CHKERRQ(ierr);
2890     }
2891   }
2892   PetscFunctionReturn(0);
2893 }
2894 
2895 #undef __FUNCT__
2896 #define __FUNCT__ "SNES_KSPSolve"
2897 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
2898 {
2899   PetscErrorCode ierr;
2900 
2901   PetscFunctionBegin;
2902   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
2903   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
2904   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
2905   PetscFunctionReturn(0);
2906 }
2907