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