xref: /petsc/src/snes/interface/snes.c (revision f22f69f01ec1ca0bb29797ff3752b65f1e1c47db)
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_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_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_ERR_SUP, "matrix-free operator rutines (version 2)");
199 #endif
200   } else {
201     SETERRQ(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 {
1117     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
1118   }
1119   if (snes->vec_rhs) {
1120     ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr);
1121   }
1122   snes->nfuncs++;
1123   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1124   PetscFunctionReturn(0);
1125 }
1126 
1127 #undef __FUNCT__
1128 #define __FUNCT__ "SNESComputeJacobian"
1129 /*@
1130    SNESComputeJacobian - Computes the Jacobian matrix that has been
1131    set with SNESSetJacobian().
1132 
1133    Collective on SNES and Mat
1134 
1135    Input Parameters:
1136 +  snes - the SNES context
1137 -  x - input vector
1138 
1139    Output Parameters:
1140 +  A - Jacobian matrix
1141 .  B - optional preconditioning matrix
1142 -  flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER)
1143 
1144   Options Database Keys:
1145 +    -snes_lag_preconditioner <lag>
1146 -    -snes_lag_jacobian <lag>
1147 
1148    Notes:
1149    Most users should not need to explicitly call this routine, as it
1150    is used internally within the nonlinear solvers.
1151 
1152    See KSPSetOperators() for important information about setting the
1153    flag parameter.
1154 
1155    Level: developer
1156 
1157 .keywords: SNES, compute, Jacobian, matrix
1158 
1159 .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
1160 @*/
1161 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
1162 {
1163   PetscErrorCode ierr;
1164   PetscTruth     flag;
1165 
1166   PetscFunctionBegin;
1167   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1168   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1169   PetscValidPointer(flg,5);
1170   PetscCheckSameComm(snes,1,X,2);
1171   if (!snes->ops->computejacobian) PetscFunctionReturn(0);
1172 
1173   /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
1174 
1175   if (snes->lagjacobian == -2) {
1176     snes->lagjacobian = -1;
1177     ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr);
1178   } else if (snes->lagjacobian == -1) {
1179     *flg = SAME_PRECONDITIONER;
1180     ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr);
1181     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1182     if (flag) {
1183       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1184       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1185     }
1186     PetscFunctionReturn(0);
1187   } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) {
1188     *flg = SAME_PRECONDITIONER;
1189     ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr);
1190     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1191     if (flag) {
1192       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1193       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1194     }
1195     PetscFunctionReturn(0);
1196   }
1197 
1198   *flg = DIFFERENT_NONZERO_PATTERN;
1199   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1200   PetscStackPush("SNES user Jacobian function");
1201   CHKMEMQ;
1202   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1203   CHKMEMQ;
1204   PetscStackPop;
1205   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1206 
1207   if (snes->lagpreconditioner == -2) {
1208     ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr);
1209     snes->lagpreconditioner = -1;
1210   } else if (snes->lagpreconditioner == -1) {
1211     *flg = SAME_PRECONDITIONER;
1212     ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr);
1213   } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) {
1214     *flg = SAME_PRECONDITIONER;
1215     ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr);
1216   }
1217 
1218   /* make sure user returned a correct Jacobian and preconditioner */
1219   /* PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
1220     PetscValidHeaderSpecific(*B,MAT_CLASSID,4);   */
1221   PetscFunctionReturn(0);
1222 }
1223 
1224 #undef __FUNCT__
1225 #define __FUNCT__ "SNESSetJacobian"
1226 /*@C
1227    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1228    location to store the matrix.
1229 
1230    Collective on SNES and Mat
1231 
1232    Input Parameters:
1233 +  snes - the SNES context
1234 .  A - Jacobian matrix
1235 .  B - preconditioner matrix (usually same as the Jacobian)
1236 .  func - Jacobian evaluation routine
1237 -  ctx - [optional] user-defined context for private data for the
1238          Jacobian evaluation routine (may be PETSC_NULL)
1239 
1240    Calling sequence of func:
1241 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1242 
1243 +  x - input vector
1244 .  A - Jacobian matrix
1245 .  B - preconditioner matrix, usually the same as A
1246 .  flag - flag indicating information about the preconditioner matrix
1247    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1248 -  ctx - [optional] user-defined Jacobian context
1249 
1250    Notes:
1251    See KSPSetOperators() for important information about setting the flag
1252    output parameter in the routine func().  Be sure to read this information!
1253 
1254    The routine func() takes Mat * as the matrix arguments rather than Mat.
1255    This allows the Jacobian evaluation routine to replace A and/or B with a
1256    completely new new matrix structure (not just different matrix elements)
1257    when appropriate, for instance, if the nonzero structure is changing
1258    throughout the global iterations.
1259 
1260    If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on
1261    each matrix.
1262 
1263    Level: beginner
1264 
1265 .keywords: SNES, nonlinear, set, Jacobian, matrix
1266 
1267 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
1268 @*/
1269 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1270 {
1271   PetscErrorCode ierr;
1272 
1273   PetscFunctionBegin;
1274   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1275   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
1276   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
1277   if (A) PetscCheckSameComm(snes,1,A,2);
1278   if (B) PetscCheckSameComm(snes,1,B,3);
1279   if (func) snes->ops->computejacobian = func;
1280   if (ctx)  snes->jacP                 = ctx;
1281   if (A) {
1282     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1283     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1284     snes->jacobian = A;
1285   }
1286   if (B) {
1287     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
1288     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1289     snes->jacobian_pre = B;
1290   }
1291   PetscFunctionReturn(0);
1292 }
1293 
1294 #undef __FUNCT__
1295 #define __FUNCT__ "SNESGetJacobian"
1296 /*@C
1297    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1298    provided context for evaluating the Jacobian.
1299 
1300    Not Collective, but Mat object will be parallel if SNES object is
1301 
1302    Input Parameter:
1303 .  snes - the nonlinear solver context
1304 
1305    Output Parameters:
1306 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1307 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1308 .  func - location to put Jacobian function (or PETSC_NULL)
1309 -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1310 
1311    Level: advanced
1312 
1313 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1314 @*/
1315 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1316 {
1317   PetscFunctionBegin;
1318   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1319   if (A)    *A    = snes->jacobian;
1320   if (B)    *B    = snes->jacobian_pre;
1321   if (func) *func = snes->ops->computejacobian;
1322   if (ctx)  *ctx  = snes->jacP;
1323   PetscFunctionReturn(0);
1324 }
1325 
1326 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1327 
1328 #undef __FUNCT__
1329 #define __FUNCT__ "SNESSetUp"
1330 /*@
1331    SNESSetUp - Sets up the internal data structures for the later use
1332    of a nonlinear solver.
1333 
1334    Collective on SNES
1335 
1336    Input Parameters:
1337 .  snes - the SNES context
1338 
1339    Notes:
1340    For basic use of the SNES solvers the user need not explicitly call
1341    SNESSetUp(), since these actions will automatically occur during
1342    the call to SNESSolve().  However, if one wishes to control this
1343    phase separately, SNESSetUp() should be called after SNESCreate()
1344    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1345 
1346    Level: advanced
1347 
1348 .keywords: SNES, nonlinear, setup
1349 
1350 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1351 @*/
1352 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
1353 {
1354   PetscErrorCode ierr;
1355 
1356   PetscFunctionBegin;
1357   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1358   if (snes->setupcalled) PetscFunctionReturn(0);
1359 
1360   if (!((PetscObject)snes)->type_name) {
1361     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
1362   }
1363 
1364   if (!snes->vec_func && !snes->vec_rhs) {
1365     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1366   }
1367   if (!snes->ops->computefunction && !snes->vec_rhs) {
1368     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1369   }
1370   if (snes->vec_func == snes->vec_sol) {
1371     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1372   }
1373   if (!snes->ops->computejacobian && snes->dm) {
1374     Mat           J;
1375     ISColoring    coloring;
1376     MatFDColoring fd;
1377 
1378     ierr = DMGetMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr);
1379     ierr = DMGetColoring(snes->dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr);
1380     ierr = MatFDColoringCreate(J,coloring,&fd);CHKERRQ(ierr);
1381     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
1382     ierr = SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fd);CHKERRQ(ierr);
1383     ierr = ISColoringDestroy(coloring);CHKERRQ(ierr);
1384   }
1385   if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);}
1386 
1387   if (snes->ops->setup) {
1388     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
1389   }
1390   snes->setupcalled = PETSC_TRUE;
1391   PetscFunctionReturn(0);
1392 }
1393 
1394 #undef __FUNCT__
1395 #define __FUNCT__ "SNESDestroy"
1396 /*@
1397    SNESDestroy - Destroys the nonlinear solver context that was created
1398    with SNESCreate().
1399 
1400    Collective on SNES
1401 
1402    Input Parameter:
1403 .  snes - the SNES context
1404 
1405    Level: beginner
1406 
1407 .keywords: SNES, nonlinear, destroy
1408 
1409 .seealso: SNESCreate(), SNESSolve()
1410 @*/
1411 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
1412 {
1413   PetscErrorCode ierr;
1414 
1415   PetscFunctionBegin;
1416   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1417   if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0);
1418 
1419   /* if memory was published with AMS then destroy it */
1420   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1421 
1422   if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);}
1423   if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);}
1424 
1425   if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);}
1426   if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);}
1427   if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);}
1428   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1429   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1430   if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);}
1431   ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);
1432   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1433   ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
1434   if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);}
1435   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
1436   PetscFunctionReturn(0);
1437 }
1438 
1439 /* ----------- Routines to set solver parameters ---------- */
1440 
1441 #undef __FUNCT__
1442 #define __FUNCT__ "SNESSetLagPreconditioner"
1443 /*@
1444    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
1445 
1446    Collective on SNES
1447 
1448    Input Parameters:
1449 +  snes - the SNES context
1450 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1451          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1452 
1453    Options Database Keys:
1454 .    -snes_lag_preconditioner <lag>
1455 
1456    Notes:
1457    The default is 1
1458    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1459    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
1460 
1461    Level: intermediate
1462 
1463 .keywords: SNES, nonlinear, set, convergence, tolerances
1464 
1465 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
1466 
1467 @*/
1468 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag)
1469 {
1470   PetscFunctionBegin;
1471   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1472   if (lag < -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1473   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1474   snes->lagpreconditioner = lag;
1475   PetscFunctionReturn(0);
1476 }
1477 
1478 #undef __FUNCT__
1479 #define __FUNCT__ "SNESGetLagPreconditioner"
1480 /*@
1481    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
1482 
1483    Collective on SNES
1484 
1485    Input Parameter:
1486 .  snes - the SNES context
1487 
1488    Output Parameter:
1489 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1490          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1491 
1492    Options Database Keys:
1493 .    -snes_lag_preconditioner <lag>
1494 
1495    Notes:
1496    The default is 1
1497    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1498 
1499    Level: intermediate
1500 
1501 .keywords: SNES, nonlinear, set, convergence, tolerances
1502 
1503 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
1504 
1505 @*/
1506 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
1507 {
1508   PetscFunctionBegin;
1509   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1510   *lag = snes->lagpreconditioner;
1511   PetscFunctionReturn(0);
1512 }
1513 
1514 #undef __FUNCT__
1515 #define __FUNCT__ "SNESSetLagJacobian"
1516 /*@
1517    SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
1518      often the preconditioner is rebuilt.
1519 
1520    Collective on SNES
1521 
1522    Input Parameters:
1523 +  snes - the SNES context
1524 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1525          the Jacobian is built etc. -2 means rebuild at next chance but then never again
1526 
1527    Options Database Keys:
1528 .    -snes_lag_jacobian <lag>
1529 
1530    Notes:
1531    The default is 1
1532    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1533    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
1534    at the next Newton step but never again (unless it is reset to another value)
1535 
1536    Level: intermediate
1537 
1538 .keywords: SNES, nonlinear, set, convergence, tolerances
1539 
1540 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
1541 
1542 @*/
1543 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagJacobian(SNES snes,PetscInt lag)
1544 {
1545   PetscFunctionBegin;
1546   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1547   if (lag < -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1548   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1549   snes->lagjacobian = lag;
1550   PetscFunctionReturn(0);
1551 }
1552 
1553 #undef __FUNCT__
1554 #define __FUNCT__ "SNESGetLagJacobian"
1555 /*@
1556    SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
1557 
1558    Collective on SNES
1559 
1560    Input Parameter:
1561 .  snes - the SNES context
1562 
1563    Output Parameter:
1564 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1565          the Jacobian is built etc.
1566 
1567    Options Database Keys:
1568 .    -snes_lag_jacobian <lag>
1569 
1570    Notes:
1571    The default is 1
1572    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1573 
1574    Level: intermediate
1575 
1576 .keywords: SNES, nonlinear, set, convergence, tolerances
1577 
1578 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
1579 
1580 @*/
1581 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagJacobian(SNES snes,PetscInt *lag)
1582 {
1583   PetscFunctionBegin;
1584   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1585   *lag = snes->lagjacobian;
1586   PetscFunctionReturn(0);
1587 }
1588 
1589 #undef __FUNCT__
1590 #define __FUNCT__ "SNESSetTolerances"
1591 /*@
1592    SNESSetTolerances - Sets various parameters used in convergence tests.
1593 
1594    Collective on SNES
1595 
1596    Input Parameters:
1597 +  snes - the SNES context
1598 .  abstol - absolute convergence tolerance
1599 .  rtol - relative convergence tolerance
1600 .  stol -  convergence tolerance in terms of the norm
1601            of the change in the solution between steps
1602 .  maxit - maximum number of iterations
1603 -  maxf - maximum number of function evaluations
1604 
1605    Options Database Keys:
1606 +    -snes_atol <abstol> - Sets abstol
1607 .    -snes_rtol <rtol> - Sets rtol
1608 .    -snes_stol <stol> - Sets stol
1609 .    -snes_max_it <maxit> - Sets maxit
1610 -    -snes_max_funcs <maxf> - Sets maxf
1611 
1612    Notes:
1613    The default maximum number of iterations is 50.
1614    The default maximum number of function evaluations is 1000.
1615 
1616    Level: intermediate
1617 
1618 .keywords: SNES, nonlinear, set, convergence, tolerances
1619 
1620 .seealso: SNESSetTrustRegionTolerance()
1621 @*/
1622 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
1623 {
1624   PetscFunctionBegin;
1625   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1626   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1627   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1628   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1629   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1630   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1631   PetscFunctionReturn(0);
1632 }
1633 
1634 #undef __FUNCT__
1635 #define __FUNCT__ "SNESGetTolerances"
1636 /*@
1637    SNESGetTolerances - Gets various parameters used in convergence tests.
1638 
1639    Not Collective
1640 
1641    Input Parameters:
1642 +  snes - the SNES context
1643 .  atol - absolute convergence tolerance
1644 .  rtol - relative convergence tolerance
1645 .  stol -  convergence tolerance in terms of the norm
1646            of the change in the solution between steps
1647 .  maxit - maximum number of iterations
1648 -  maxf - maximum number of function evaluations
1649 
1650    Notes:
1651    The user can specify PETSC_NULL for any parameter that is not needed.
1652 
1653    Level: intermediate
1654 
1655 .keywords: SNES, nonlinear, get, convergence, tolerances
1656 
1657 .seealso: SNESSetTolerances()
1658 @*/
1659 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
1660 {
1661   PetscFunctionBegin;
1662   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1663   if (atol)  *atol  = snes->abstol;
1664   if (rtol)  *rtol  = snes->rtol;
1665   if (stol)  *stol  = snes->xtol;
1666   if (maxit) *maxit = snes->max_its;
1667   if (maxf)  *maxf  = snes->max_funcs;
1668   PetscFunctionReturn(0);
1669 }
1670 
1671 #undef __FUNCT__
1672 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1673 /*@
1674    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1675 
1676    Collective on SNES
1677 
1678    Input Parameters:
1679 +  snes - the SNES context
1680 -  tol - tolerance
1681 
1682    Options Database Key:
1683 .  -snes_trtol <tol> - Sets tol
1684 
1685    Level: intermediate
1686 
1687 .keywords: SNES, nonlinear, set, trust region, tolerance
1688 
1689 .seealso: SNESSetTolerances()
1690 @*/
1691 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1692 {
1693   PetscFunctionBegin;
1694   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1695   snes->deltatol = tol;
1696   PetscFunctionReturn(0);
1697 }
1698 
1699 /*
1700    Duplicate the lg monitors for SNES from KSP; for some reason with
1701    dynamic libraries things don't work under Sun4 if we just use
1702    macros instead of functions
1703 */
1704 #undef __FUNCT__
1705 #define __FUNCT__ "SNESMonitorLG"
1706 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1707 {
1708   PetscErrorCode ierr;
1709 
1710   PetscFunctionBegin;
1711   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1712   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1713   PetscFunctionReturn(0);
1714 }
1715 
1716 #undef __FUNCT__
1717 #define __FUNCT__ "SNESMonitorLGCreate"
1718 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1719 {
1720   PetscErrorCode ierr;
1721 
1722   PetscFunctionBegin;
1723   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1724   PetscFunctionReturn(0);
1725 }
1726 
1727 #undef __FUNCT__
1728 #define __FUNCT__ "SNESMonitorLGDestroy"
1729 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw)
1730 {
1731   PetscErrorCode ierr;
1732 
1733   PetscFunctionBegin;
1734   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1735   PetscFunctionReturn(0);
1736 }
1737 
1738 extern PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
1739 #undef __FUNCT__
1740 #define __FUNCT__ "SNESMonitorLGRange"
1741 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
1742 {
1743   PetscDrawLG      lg;
1744   PetscErrorCode   ierr;
1745   PetscReal        x,y,per;
1746   PetscViewer      v = (PetscViewer)monctx;
1747   static PetscReal prev; /* should be in the context */
1748   PetscDraw        draw;
1749   PetscFunctionBegin;
1750   if (!monctx) {
1751     MPI_Comm    comm;
1752 
1753     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1754     v      = PETSC_VIEWER_DRAW_(comm);
1755   }
1756   ierr   = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr);
1757   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1758   ierr   = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1759   ierr   = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr);
1760   x = (PetscReal) n;
1761   if (rnorm > 0.0) y = log10(rnorm); else y = -15.0;
1762   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1763   if (n < 20 || !(n % 5)) {
1764     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1765   }
1766 
1767   ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr);
1768   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1769   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1770   ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr);
1771   ierr =  SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr);
1772   x = (PetscReal) n;
1773   y = 100.0*per;
1774   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1775   if (n < 20 || !(n % 5)) {
1776     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1777   }
1778 
1779   ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr);
1780   if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1781   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1782   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr);
1783   x = (PetscReal) n;
1784   y = (prev - rnorm)/prev;
1785   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1786   if (n < 20 || !(n % 5)) {
1787     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1788   }
1789 
1790   ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr);
1791   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1792   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1793   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
1794   x = (PetscReal) n;
1795   y = (prev - rnorm)/(prev*per);
1796   if (n > 2) { /*skip initial crazy value */
1797     ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1798   }
1799   if (n < 20 || !(n % 5)) {
1800     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1801   }
1802   prev = rnorm;
1803   PetscFunctionReturn(0);
1804 }
1805 
1806 #undef __FUNCT__
1807 #define __FUNCT__ "SNESMonitorLGRangeCreate"
1808 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1809 {
1810   PetscErrorCode ierr;
1811 
1812   PetscFunctionBegin;
1813   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1814   PetscFunctionReturn(0);
1815 }
1816 
1817 #undef __FUNCT__
1818 #define __FUNCT__ "SNESMonitorLGRangeDestroy"
1819 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeDestroy(PetscDrawLG draw)
1820 {
1821   PetscErrorCode ierr;
1822 
1823   PetscFunctionBegin;
1824   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1825   PetscFunctionReturn(0);
1826 }
1827 
1828 /* ------------ Routines to set performance monitoring options ----------- */
1829 
1830 #undef __FUNCT__
1831 #define __FUNCT__ "SNESMonitorSet"
1832 /*@C
1833    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
1834    iteration of the nonlinear solver to display the iteration's
1835    progress.
1836 
1837    Collective on SNES
1838 
1839    Input Parameters:
1840 +  snes - the SNES context
1841 .  func - monitoring routine
1842 .  mctx - [optional] user-defined context for private data for the
1843           monitor routine (use PETSC_NULL if no context is desired)
1844 -  monitordestroy - [optional] routine that frees monitor context
1845           (may be PETSC_NULL)
1846 
1847    Calling sequence of func:
1848 $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1849 
1850 +    snes - the SNES context
1851 .    its - iteration number
1852 .    norm - 2-norm function value (may be estimated)
1853 -    mctx - [optional] monitoring context
1854 
1855    Options Database Keys:
1856 +    -snes_monitor        - sets SNESMonitorDefault()
1857 .    -snes_monitor_draw    - sets line graph monitor,
1858                             uses SNESMonitorLGCreate()
1859 _    -snes_monitor_cancel - cancels all monitors that have
1860                             been hardwired into a code by
1861                             calls to SNESMonitorSet(), but
1862                             does not cancel those set via
1863                             the options database.
1864 
1865    Notes:
1866    Several different monitoring routines may be set by calling
1867    SNESMonitorSet() multiple times; all will be called in the
1868    order in which they were set.
1869 
1870    Fortran notes: Only a single monitor function can be set for each SNES object
1871 
1872    Level: intermediate
1873 
1874 .keywords: SNES, nonlinear, set, monitor
1875 
1876 .seealso: SNESMonitorDefault(), SNESMonitorCancel()
1877 @*/
1878 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
1879 {
1880   PetscInt i;
1881 
1882   PetscFunctionBegin;
1883   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1884   if (snes->numbermonitors >= MAXSNESMONITORS) {
1885     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1886   }
1887   for (i=0; i<snes->numbermonitors;i++) {
1888     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0);
1889 
1890     /* check if both default monitors that share common ASCII viewer */
1891     if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) {
1892       if (mctx && snes->monitorcontext[i]) {
1893         PetscErrorCode          ierr;
1894         PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx;
1895         PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i];
1896         if (viewer1->viewer == viewer2->viewer) {
1897           ierr = (*monitordestroy)(mctx);CHKERRQ(ierr);
1898           PetscFunctionReturn(0);
1899         }
1900       }
1901     }
1902   }
1903   snes->monitor[snes->numbermonitors]           = monitor;
1904   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1905   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1906   PetscFunctionReturn(0);
1907 }
1908 
1909 #undef __FUNCT__
1910 #define __FUNCT__ "SNESMonitorCancel"
1911 /*@C
1912    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
1913 
1914    Collective on SNES
1915 
1916    Input Parameters:
1917 .  snes - the SNES context
1918 
1919    Options Database Key:
1920 .  -snes_monitor_cancel - cancels all monitors that have been hardwired
1921     into a code by calls to SNESMonitorSet(), but does not cancel those
1922     set via the options database
1923 
1924    Notes:
1925    There is no way to clear one specific monitor from a SNES object.
1926 
1927    Level: intermediate
1928 
1929 .keywords: SNES, nonlinear, set, monitor
1930 
1931 .seealso: SNESMonitorDefault(), SNESMonitorSet()
1932 @*/
1933 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes)
1934 {
1935   PetscErrorCode ierr;
1936   PetscInt       i;
1937 
1938   PetscFunctionBegin;
1939   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1940   for (i=0; i<snes->numbermonitors; i++) {
1941     if (snes->monitordestroy[i]) {
1942       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1943     }
1944   }
1945   snes->numbermonitors = 0;
1946   PetscFunctionReturn(0);
1947 }
1948 
1949 #undef __FUNCT__
1950 #define __FUNCT__ "SNESSetConvergenceTest"
1951 /*@C
1952    SNESSetConvergenceTest - Sets the function that is to be used
1953    to test for convergence of the nonlinear iterative solution.
1954 
1955    Collective on SNES
1956 
1957    Input Parameters:
1958 +  snes - the SNES context
1959 .  func - routine to test for convergence
1960 .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
1961 -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
1962 
1963    Calling sequence of func:
1964 $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1965 
1966 +    snes - the SNES context
1967 .    it - current iteration (0 is the first and is before any Newton step)
1968 .    cctx - [optional] convergence context
1969 .    reason - reason for convergence/divergence
1970 .    xnorm - 2-norm of current iterate
1971 .    gnorm - 2-norm of current step
1972 -    f - 2-norm of function
1973 
1974    Level: advanced
1975 
1976 .keywords: SNES, nonlinear, set, convergence, test
1977 
1978 .seealso: SNESDefaultConverged(), SNESSkipConverged()
1979 @*/
1980 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1981 {
1982   PetscErrorCode ierr;
1983 
1984   PetscFunctionBegin;
1985   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1986   if (!func) func = SNESSkipConverged;
1987   if (snes->ops->convergeddestroy) {
1988     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
1989   }
1990   snes->ops->converged        = func;
1991   snes->ops->convergeddestroy = destroy;
1992   snes->cnvP                  = cctx;
1993   PetscFunctionReturn(0);
1994 }
1995 
1996 #undef __FUNCT__
1997 #define __FUNCT__ "SNESGetConvergedReason"
1998 /*@
1999    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
2000 
2001    Not Collective
2002 
2003    Input Parameter:
2004 .  snes - the SNES context
2005 
2006    Output Parameter:
2007 .  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
2008             manual pages for the individual convergence tests for complete lists
2009 
2010    Level: intermediate
2011 
2012    Notes: Can only be called after the call the SNESSolve() is complete.
2013 
2014 .keywords: SNES, nonlinear, set, convergence, test
2015 
2016 .seealso: SNESSetConvergenceTest(), SNESConvergedReason
2017 @*/
2018 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
2019 {
2020   PetscFunctionBegin;
2021   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2022   PetscValidPointer(reason,2);
2023   *reason = snes->reason;
2024   PetscFunctionReturn(0);
2025 }
2026 
2027 #undef __FUNCT__
2028 #define __FUNCT__ "SNESSetConvergenceHistory"
2029 /*@
2030    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
2031 
2032    Collective on SNES
2033 
2034    Input Parameters:
2035 +  snes - iterative context obtained from SNESCreate()
2036 .  a   - array to hold history, this array will contain the function norms computed at each step
2037 .  its - integer array holds the number of linear iterations for each solve.
2038 .  na  - size of a and its
2039 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
2040            else it continues storing new values for new nonlinear solves after the old ones
2041 
2042    This routine is useful, e.g., when running a code for purposes
2043    of accurate performance monitoring, when no I/O should be done
2044    during the section of code that is being timed.
2045 
2046    Level: intermediate
2047 
2048 .keywords: SNES, set, convergence, history
2049 
2050 .seealso: SNESGetConvergenceHistory()
2051 
2052 @*/
2053 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset)
2054 {
2055   PetscFunctionBegin;
2056   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2057   if (na)  PetscValidScalarPointer(a,2);
2058   if (its) PetscValidIntPointer(its,3);
2059   snes->conv_hist       = a;
2060   snes->conv_hist_its   = its;
2061   snes->conv_hist_max   = na;
2062   snes->conv_hist_len   = 0;
2063   snes->conv_hist_reset = reset;
2064   PetscFunctionReturn(0);
2065 }
2066 
2067 #undef __FUNCT__
2068 #define __FUNCT__ "SNESGetConvergenceHistory"
2069 /*@C
2070    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
2071 
2072    Collective on SNES
2073 
2074    Input Parameter:
2075 .  snes - iterative context obtained from SNESCreate()
2076 
2077    Output Parameters:
2078 .  a   - array to hold history
2079 .  its - integer array holds the number of linear iterations (or
2080          negative if not converged) for each solve.
2081 -  na  - size of a and its
2082 
2083    Notes:
2084     The calling sequence for this routine in Fortran is
2085 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
2086 
2087    This routine is useful, e.g., when running a code for purposes
2088    of accurate performance monitoring, when no I/O should be done
2089    during the section of code that is being timed.
2090 
2091    Level: intermediate
2092 
2093 .keywords: SNES, get, convergence, history
2094 
2095 .seealso: SNESSetConvergencHistory()
2096 
2097 @*/
2098 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
2099 {
2100   PetscFunctionBegin;
2101   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2102   if (a)   *a   = snes->conv_hist;
2103   if (its) *its = snes->conv_hist_its;
2104   if (na)  *na  = snes->conv_hist_len;
2105   PetscFunctionReturn(0);
2106 }
2107 
2108 #undef __FUNCT__
2109 #define __FUNCT__ "SNESSetUpdate"
2110 /*@C
2111   SNESSetUpdate - Sets the general-purpose update function called
2112   at the beginning o every iteration of the nonlinear solve. Specifically
2113   it is called just before the Jacobian is "evaluated".
2114 
2115   Collective on SNES
2116 
2117   Input Parameters:
2118 . snes - The nonlinear solver context
2119 . func - The function
2120 
2121   Calling sequence of func:
2122 . func (SNES snes, PetscInt step);
2123 
2124 . step - The current step of the iteration
2125 
2126   Level: advanced
2127 
2128   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()
2129         This is not used by most users.
2130 
2131 .keywords: SNES, update
2132 
2133 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
2134 @*/
2135 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
2136 {
2137   PetscFunctionBegin;
2138   PetscValidHeaderSpecific(snes, SNES_CLASSID,1);
2139   snes->ops->update = func;
2140   PetscFunctionReturn(0);
2141 }
2142 
2143 #undef __FUNCT__
2144 #define __FUNCT__ "SNESDefaultUpdate"
2145 /*@
2146   SNESDefaultUpdate - The default update function which does nothing.
2147 
2148   Not collective
2149 
2150   Input Parameters:
2151 . snes - The nonlinear solver context
2152 . step - The current step of the iteration
2153 
2154   Level: intermediate
2155 
2156 .keywords: SNES, update
2157 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
2158 @*/
2159 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
2160 {
2161   PetscFunctionBegin;
2162   PetscFunctionReturn(0);
2163 }
2164 
2165 #undef __FUNCT__
2166 #define __FUNCT__ "SNESScaleStep_Private"
2167 /*
2168    SNESScaleStep_Private - Scales a step so that its length is less than the
2169    positive parameter delta.
2170 
2171     Input Parameters:
2172 +   snes - the SNES context
2173 .   y - approximate solution of linear system
2174 .   fnorm - 2-norm of current function
2175 -   delta - trust region size
2176 
2177     Output Parameters:
2178 +   gpnorm - predicted function norm at the new point, assuming local
2179     linearization.  The value is zero if the step lies within the trust
2180     region, and exceeds zero otherwise.
2181 -   ynorm - 2-norm of the step
2182 
2183     Note:
2184     For non-trust region methods such as SNESLS, the parameter delta
2185     is set to be the maximum allowable step size.
2186 
2187 .keywords: SNES, nonlinear, scale, step
2188 */
2189 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
2190 {
2191   PetscReal      nrm;
2192   PetscScalar    cnorm;
2193   PetscErrorCode ierr;
2194 
2195   PetscFunctionBegin;
2196   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2197   PetscValidHeaderSpecific(y,VEC_CLASSID,2);
2198   PetscCheckSameComm(snes,1,y,2);
2199 
2200   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
2201   if (nrm > *delta) {
2202      nrm = *delta/nrm;
2203      *gpnorm = (1.0 - nrm)*(*fnorm);
2204      cnorm = nrm;
2205      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
2206      *ynorm = *delta;
2207   } else {
2208      *gpnorm = 0.0;
2209      *ynorm = nrm;
2210   }
2211   PetscFunctionReturn(0);
2212 }
2213 
2214 #undef __FUNCT__
2215 #define __FUNCT__ "SNESSolve"
2216 /*@C
2217    SNESSolve - Solves a nonlinear system F(x) = b.
2218    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
2219 
2220    Collective on SNES
2221 
2222    Input Parameters:
2223 +  snes - the SNES context
2224 .  b - the constant part of the equation, or PETSC_NULL to use zero.
2225 -  x - the solution vector.
2226 
2227    Notes:
2228    The user should initialize the vector,x, with the initial guess
2229    for the nonlinear solve prior to calling SNESSolve.  In particular,
2230    to employ an initial guess of zero, the user should explicitly set
2231    this vector to zero by calling VecSet().
2232 
2233    Level: beginner
2234 
2235 .keywords: SNES, nonlinear, solve
2236 
2237 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian()
2238 @*/
2239 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
2240 {
2241   PetscErrorCode ierr;
2242   PetscTruth     flg;
2243   char           filename[PETSC_MAX_PATH_LEN];
2244   PetscViewer    viewer;
2245 
2246   PetscFunctionBegin;
2247   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2248   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2249   PetscCheckSameComm(snes,1,x,3);
2250   if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2);
2251   if (b) PetscCheckSameComm(snes,1,b,2);
2252 
2253   /* set solution vector */
2254   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
2255   if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); }
2256   snes->vec_sol = x;
2257   /* set afine vector if provided */
2258   if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
2259   if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); }
2260   snes->vec_rhs = b;
2261 
2262   if (!snes->vec_func && snes->vec_rhs) {
2263     ierr = VecDuplicate(b, &snes->vec_func);CHKERRQ(ierr);
2264   }
2265 
2266   ierr = SNESSetUp(snes);CHKERRQ(ierr);
2267 
2268   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
2269   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
2270 
2271   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2272   ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
2273   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2274   if (snes->domainerror){
2275     snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
2276     snes->domainerror = PETSC_FALSE;
2277   }
2278 
2279   if (!snes->reason) {
2280     SETERRQ(PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
2281   }
2282 
2283   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
2284   if (flg && !PetscPreLoadingOn) {
2285     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
2286     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
2287     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
2288   }
2289 
2290   flg  = PETSC_FALSE;
2291   ierr = PetscOptionsGetTruth(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,PETSC_NULL);CHKERRQ(ierr);
2292   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
2293   if (snes->printreason) {
2294     if (snes->reason > 0) {
2295       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2296     } else {
2297       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2298     }
2299   }
2300 
2301   PetscFunctionReturn(0);
2302 }
2303 
2304 /* --------- Internal routines for SNES Package --------- */
2305 
2306 #undef __FUNCT__
2307 #define __FUNCT__ "SNESSetType"
2308 /*@C
2309    SNESSetType - Sets the method for the nonlinear solver.
2310 
2311    Collective on SNES
2312 
2313    Input Parameters:
2314 +  snes - the SNES context
2315 -  type - a known method
2316 
2317    Options Database Key:
2318 .  -snes_type <type> - Sets the method; use -help for a list
2319    of available methods (for instance, ls or tr)
2320 
2321    Notes:
2322    See "petsc/include/petscsnes.h" for available methods (for instance)
2323 +    SNESLS - Newton's method with line search
2324      (systems of nonlinear equations)
2325 .    SNESTR - Newton's method with trust region
2326      (systems of nonlinear equations)
2327 
2328   Normally, it is best to use the SNESSetFromOptions() command and then
2329   set the SNES solver type from the options database rather than by using
2330   this routine.  Using the options database provides the user with
2331   maximum flexibility in evaluating the many nonlinear solvers.
2332   The SNESSetType() routine is provided for those situations where it
2333   is necessary to set the nonlinear solver independently of the command
2334   line or options database.  This might be the case, for example, when
2335   the choice of solver changes during the execution of the program,
2336   and the user's application is taking responsibility for choosing the
2337   appropriate method.
2338 
2339   Level: intermediate
2340 
2341 .keywords: SNES, set, type
2342 
2343 .seealso: SNESType, SNESCreate()
2344 
2345 @*/
2346 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type)
2347 {
2348   PetscErrorCode ierr,(*r)(SNES);
2349   PetscTruth     match;
2350 
2351   PetscFunctionBegin;
2352   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2353   PetscValidCharPointer(type,2);
2354 
2355   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
2356   if (match) PetscFunctionReturn(0);
2357 
2358   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
2359   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
2360   /* Destroy the previous private SNES context */
2361   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
2362   /* Reinitialize function pointers in SNESOps structure */
2363   snes->ops->setup          = 0;
2364   snes->ops->solve          = 0;
2365   snes->ops->view           = 0;
2366   snes->ops->setfromoptions = 0;
2367   snes->ops->destroy        = 0;
2368   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
2369   snes->setupcalled = PETSC_FALSE;
2370   ierr = (*r)(snes);CHKERRQ(ierr);
2371   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
2372   PetscFunctionReturn(0);
2373 }
2374 
2375 
2376 /* --------------------------------------------------------------------- */
2377 #undef __FUNCT__
2378 #define __FUNCT__ "SNESRegisterDestroy"
2379 /*@
2380    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2381    registered by SNESRegisterDynamic().
2382 
2383    Not Collective
2384 
2385    Level: advanced
2386 
2387 .keywords: SNES, nonlinear, register, destroy
2388 
2389 .seealso: SNESRegisterAll(), SNESRegisterAll()
2390 @*/
2391 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
2392 {
2393   PetscErrorCode ierr;
2394 
2395   PetscFunctionBegin;
2396   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
2397   SNESRegisterAllCalled = PETSC_FALSE;
2398   PetscFunctionReturn(0);
2399 }
2400 
2401 #undef __FUNCT__
2402 #define __FUNCT__ "SNESGetType"
2403 /*@C
2404    SNESGetType - Gets the SNES method type and name (as a string).
2405 
2406    Not Collective
2407 
2408    Input Parameter:
2409 .  snes - nonlinear solver context
2410 
2411    Output Parameter:
2412 .  type - SNES method (a character string)
2413 
2414    Level: intermediate
2415 
2416 .keywords: SNES, nonlinear, get, type, name
2417 @*/
2418 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type)
2419 {
2420   PetscFunctionBegin;
2421   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2422   PetscValidPointer(type,2);
2423   *type = ((PetscObject)snes)->type_name;
2424   PetscFunctionReturn(0);
2425 }
2426 
2427 #undef __FUNCT__
2428 #define __FUNCT__ "SNESGetSolution"
2429 /*@
2430    SNESGetSolution - Returns the vector where the approximate solution is
2431    stored.
2432 
2433    Not Collective, but Vec is parallel if SNES is parallel
2434 
2435    Input Parameter:
2436 .  snes - the SNES context
2437 
2438    Output Parameter:
2439 .  x - the solution
2440 
2441    Level: intermediate
2442 
2443 .keywords: SNES, nonlinear, get, solution
2444 
2445 .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
2446 @*/
2447 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
2448 {
2449   PetscFunctionBegin;
2450   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2451   PetscValidPointer(x,2);
2452   *x = snes->vec_sol;
2453   PetscFunctionReturn(0);
2454 }
2455 
2456 #undef __FUNCT__
2457 #define __FUNCT__ "SNESGetSolutionUpdate"
2458 /*@
2459    SNESGetSolutionUpdate - Returns the vector where the solution update is
2460    stored.
2461 
2462    Not Collective, but Vec is parallel if SNES is parallel
2463 
2464    Input Parameter:
2465 .  snes - the SNES context
2466 
2467    Output Parameter:
2468 .  x - the solution update
2469 
2470    Level: advanced
2471 
2472 .keywords: SNES, nonlinear, get, solution, update
2473 
2474 .seealso: SNESGetSolution(), SNESGetFunction()
2475 @*/
2476 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
2477 {
2478   PetscFunctionBegin;
2479   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2480   PetscValidPointer(x,2);
2481   *x = snes->vec_sol_update;
2482   PetscFunctionReturn(0);
2483 }
2484 
2485 #undef __FUNCT__
2486 #define __FUNCT__ "SNESGetFunction"
2487 /*@C
2488    SNESGetFunction - Returns the vector where the function is stored.
2489 
2490    Not Collective, but Vec is parallel if SNES is parallel
2491 
2492    Input Parameter:
2493 .  snes - the SNES context
2494 
2495    Output Parameter:
2496 +  r - the function (or PETSC_NULL)
2497 .  func - the function (or PETSC_NULL)
2498 -  ctx - the function context (or PETSC_NULL)
2499 
2500    Level: advanced
2501 
2502 .keywords: SNES, nonlinear, get, function
2503 
2504 .seealso: SNESSetFunction(), SNESGetSolution()
2505 @*/
2506 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
2507 {
2508   PetscFunctionBegin;
2509   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2510   if (r)    *r    = snes->vec_func;
2511   if (func) *func = snes->ops->computefunction;
2512   if (ctx)  *ctx  = snes->funP;
2513   PetscFunctionReturn(0);
2514 }
2515 
2516 #undef __FUNCT__
2517 #define __FUNCT__ "SNESSetOptionsPrefix"
2518 /*@C
2519    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2520    SNES options in the database.
2521 
2522    Collective on SNES
2523 
2524    Input Parameter:
2525 +  snes - the SNES context
2526 -  prefix - the prefix to prepend to all option names
2527 
2528    Notes:
2529    A hyphen (-) must NOT be given at the beginning of the prefix name.
2530    The first character of all runtime options is AUTOMATICALLY the hyphen.
2531 
2532    Level: advanced
2533 
2534 .keywords: SNES, set, options, prefix, database
2535 
2536 .seealso: SNESSetFromOptions()
2537 @*/
2538 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
2539 {
2540   PetscErrorCode ierr;
2541 
2542   PetscFunctionBegin;
2543   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2544   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2545   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2546   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2547   PetscFunctionReturn(0);
2548 }
2549 
2550 #undef __FUNCT__
2551 #define __FUNCT__ "SNESAppendOptionsPrefix"
2552 /*@C
2553    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2554    SNES options in the database.
2555 
2556    Collective on SNES
2557 
2558    Input Parameters:
2559 +  snes - the SNES context
2560 -  prefix - the prefix to prepend to all option names
2561 
2562    Notes:
2563    A hyphen (-) must NOT be given at the beginning of the prefix name.
2564    The first character of all runtime options is AUTOMATICALLY the hyphen.
2565 
2566    Level: advanced
2567 
2568 .keywords: SNES, append, options, prefix, database
2569 
2570 .seealso: SNESGetOptionsPrefix()
2571 @*/
2572 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
2573 {
2574   PetscErrorCode ierr;
2575 
2576   PetscFunctionBegin;
2577   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2578   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2579   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2580   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2581   PetscFunctionReturn(0);
2582 }
2583 
2584 #undef __FUNCT__
2585 #define __FUNCT__ "SNESGetOptionsPrefix"
2586 /*@C
2587    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2588    SNES options in the database.
2589 
2590    Not Collective
2591 
2592    Input Parameter:
2593 .  snes - the SNES context
2594 
2595    Output Parameter:
2596 .  prefix - pointer to the prefix string used
2597 
2598    Notes: On the fortran side, the user should pass in a string 'prefix' of
2599    sufficient length to hold the prefix.
2600 
2601    Level: advanced
2602 
2603 .keywords: SNES, get, options, prefix, database
2604 
2605 .seealso: SNESAppendOptionsPrefix()
2606 @*/
2607 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
2608 {
2609   PetscErrorCode ierr;
2610 
2611   PetscFunctionBegin;
2612   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2613   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2614   PetscFunctionReturn(0);
2615 }
2616 
2617 
2618 #undef __FUNCT__
2619 #define __FUNCT__ "SNESRegister"
2620 /*@C
2621   SNESRegister - See SNESRegisterDynamic()
2622 
2623   Level: advanced
2624 @*/
2625 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2626 {
2627   char           fullname[PETSC_MAX_PATH_LEN];
2628   PetscErrorCode ierr;
2629 
2630   PetscFunctionBegin;
2631   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2632   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2633   PetscFunctionReturn(0);
2634 }
2635 
2636 #undef __FUNCT__
2637 #define __FUNCT__ "SNESTestLocalMin"
2638 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2639 {
2640   PetscErrorCode ierr;
2641   PetscInt       N,i,j;
2642   Vec            u,uh,fh;
2643   PetscScalar    value;
2644   PetscReal      norm;
2645 
2646   PetscFunctionBegin;
2647   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2648   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2649   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2650 
2651   /* currently only works for sequential */
2652   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2653   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2654   for (i=0; i<N; i++) {
2655     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2656     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2657     for (j=-10; j<11; j++) {
2658       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2659       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2660       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2661       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2662       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2663       value = -value;
2664       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2665     }
2666   }
2667   ierr = VecDestroy(uh);CHKERRQ(ierr);
2668   ierr = VecDestroy(fh);CHKERRQ(ierr);
2669   PetscFunctionReturn(0);
2670 }
2671 
2672 #undef __FUNCT__
2673 #define __FUNCT__ "SNESKSPSetUseEW"
2674 /*@
2675    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
2676    computing relative tolerance for linear solvers within an inexact
2677    Newton method.
2678 
2679    Collective on SNES
2680 
2681    Input Parameters:
2682 +  snes - SNES context
2683 -  flag - PETSC_TRUE or PETSC_FALSE
2684 
2685     Options Database:
2686 +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
2687 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
2688 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
2689 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
2690 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
2691 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
2692 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
2693 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
2694 
2695    Notes:
2696    Currently, the default is to use a constant relative tolerance for
2697    the inner linear solvers.  Alternatively, one can use the
2698    Eisenstat-Walker method, where the relative convergence tolerance
2699    is reset at each Newton iteration according progress of the nonlinear
2700    solver.
2701 
2702    Level: advanced
2703 
2704    Reference:
2705    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2706    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2707 
2708 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2709 
2710 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2711 @*/
2712 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag)
2713 {
2714   PetscFunctionBegin;
2715   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2716   snes->ksp_ewconv = flag;
2717   PetscFunctionReturn(0);
2718 }
2719 
2720 #undef __FUNCT__
2721 #define __FUNCT__ "SNESKSPGetUseEW"
2722 /*@
2723    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
2724    for computing relative tolerance for linear solvers within an
2725    inexact Newton method.
2726 
2727    Not Collective
2728 
2729    Input Parameter:
2730 .  snes - SNES context
2731 
2732    Output Parameter:
2733 .  flag - PETSC_TRUE or PETSC_FALSE
2734 
2735    Notes:
2736    Currently, the default is to use a constant relative tolerance for
2737    the inner linear solvers.  Alternatively, one can use the
2738    Eisenstat-Walker method, where the relative convergence tolerance
2739    is reset at each Newton iteration according progress of the nonlinear
2740    solver.
2741 
2742    Level: advanced
2743 
2744    Reference:
2745    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2746    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2747 
2748 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2749 
2750 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2751 @*/
2752 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag)
2753 {
2754   PetscFunctionBegin;
2755   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2756   PetscValidPointer(flag,2);
2757   *flag = snes->ksp_ewconv;
2758   PetscFunctionReturn(0);
2759 }
2760 
2761 #undef __FUNCT__
2762 #define __FUNCT__ "SNESKSPSetParametersEW"
2763 /*@
2764    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
2765    convergence criteria for the linear solvers within an inexact
2766    Newton method.
2767 
2768    Collective on SNES
2769 
2770    Input Parameters:
2771 +    snes - SNES context
2772 .    version - version 1, 2 (default is 2) or 3
2773 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2774 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2775 .    gamma - multiplicative factor for version 2 rtol computation
2776              (0 <= gamma2 <= 1)
2777 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2778 .    alpha2 - power for safeguard
2779 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2780 
2781    Note:
2782    Version 3 was contributed by Luis Chacon, June 2006.
2783 
2784    Use PETSC_DEFAULT to retain the default for any of the parameters.
2785 
2786    Level: advanced
2787 
2788    Reference:
2789    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2790    inexact Newton method", Utah State University Math. Stat. Dept. Res.
2791    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
2792 
2793 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
2794 
2795 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
2796 @*/
2797 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
2798 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
2799 {
2800   SNESKSPEW *kctx;
2801   PetscFunctionBegin;
2802   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2803   kctx = (SNESKSPEW*)snes->kspconvctx;
2804   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2805 
2806   if (version != PETSC_DEFAULT)   kctx->version   = version;
2807   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
2808   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
2809   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
2810   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
2811   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
2812   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
2813 
2814   if (kctx->version < 1 || kctx->version > 3) {
2815     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
2816   }
2817   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
2818     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
2819   }
2820   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
2821     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
2822   }
2823   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
2824     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
2825   }
2826   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
2827     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
2828   }
2829   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
2830     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
2831   }
2832   PetscFunctionReturn(0);
2833 }
2834 
2835 #undef __FUNCT__
2836 #define __FUNCT__ "SNESKSPGetParametersEW"
2837 /*@
2838    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
2839    convergence criteria for the linear solvers within an inexact
2840    Newton method.
2841 
2842    Not Collective
2843 
2844    Input Parameters:
2845      snes - SNES context
2846 
2847    Output Parameters:
2848 +    version - version 1, 2 (default is 2) or 3
2849 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2850 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2851 .    gamma - multiplicative factor for version 2 rtol computation
2852              (0 <= gamma2 <= 1)
2853 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2854 .    alpha2 - power for safeguard
2855 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2856 
2857    Level: advanced
2858 
2859 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
2860 
2861 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
2862 @*/
2863 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
2864 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
2865 {
2866   SNESKSPEW *kctx;
2867   PetscFunctionBegin;
2868   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2869   kctx = (SNESKSPEW*)snes->kspconvctx;
2870   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2871   if(version)   *version   = kctx->version;
2872   if(rtol_0)    *rtol_0    = kctx->rtol_0;
2873   if(rtol_max)  *rtol_max  = kctx->rtol_max;
2874   if(gamma)     *gamma     = kctx->gamma;
2875   if(alpha)     *alpha     = kctx->alpha;
2876   if(alpha2)    *alpha2    = kctx->alpha2;
2877   if(threshold) *threshold = kctx->threshold;
2878   PetscFunctionReturn(0);
2879 }
2880 
2881 #undef __FUNCT__
2882 #define __FUNCT__ "SNESKSPEW_PreSolve"
2883 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
2884 {
2885   PetscErrorCode ierr;
2886   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2887   PetscReal      rtol=PETSC_DEFAULT,stol;
2888 
2889   PetscFunctionBegin;
2890   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2891   if (!snes->iter) { /* first time in, so use the original user rtol */
2892     rtol = kctx->rtol_0;
2893   } else {
2894     if (kctx->version == 1) {
2895       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
2896       if (rtol < 0.0) rtol = -rtol;
2897       stol = pow(kctx->rtol_last,kctx->alpha2);
2898       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2899     } else if (kctx->version == 2) {
2900       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2901       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
2902       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2903     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
2904       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2905       /* safeguard: avoid sharp decrease of rtol */
2906       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
2907       stol = PetscMax(rtol,stol);
2908       rtol = PetscMin(kctx->rtol_0,stol);
2909       /* safeguard: avoid oversolving */
2910       stol = kctx->gamma*(snes->ttol)/snes->norm;
2911       stol = PetscMax(rtol,stol);
2912       rtol = PetscMin(kctx->rtol_0,stol);
2913     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
2914   }
2915   /* safeguard: avoid rtol greater than one */
2916   rtol = PetscMin(rtol,kctx->rtol_max);
2917   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2918   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
2919   PetscFunctionReturn(0);
2920 }
2921 
2922 #undef __FUNCT__
2923 #define __FUNCT__ "SNESKSPEW_PostSolve"
2924 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
2925 {
2926   PetscErrorCode ierr;
2927   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2928   PCSide         pcside;
2929   Vec            lres;
2930 
2931   PetscFunctionBegin;
2932   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2933   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
2934   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
2935   if (kctx->version == 1) {
2936     ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr);
2937     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
2938       /* KSP residual is true linear residual */
2939       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
2940     } else {
2941       /* KSP residual is preconditioned residual */
2942       /* compute true linear residual norm */
2943       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
2944       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
2945       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
2946       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
2947       ierr = VecDestroy(lres);CHKERRQ(ierr);
2948     }
2949   }
2950   PetscFunctionReturn(0);
2951 }
2952 
2953 #undef __FUNCT__
2954 #define __FUNCT__ "SNES_KSPSolve"
2955 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
2956 {
2957   PetscErrorCode ierr;
2958 
2959   PetscFunctionBegin;
2960   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
2961   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
2962   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
2963   PetscFunctionReturn(0);
2964 }
2965 
2966 #undef __FUNCT__
2967 #define __FUNCT__ "SNESSetDM"
2968 /*@
2969    SNESSetDM - Sets the DM that may be used by some preconditioners
2970 
2971    Collective on SNES
2972 
2973    Input Parameters:
2974 +  snes - the preconditioner context
2975 -  dm - the dm
2976 
2977    Level: intermediate
2978 
2979 
2980 .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
2981 @*/
2982 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetDM(SNES snes,DM dm)
2983 {
2984   PetscErrorCode ierr;
2985   KSP            ksp;
2986 
2987   PetscFunctionBegin;
2988   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2989   if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);}
2990   snes->dm = dm;
2991   ierr = PetscObjectReference((PetscObject)snes->dm);CHKERRQ(ierr);
2992   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
2993   ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr);
2994   ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr);
2995   PetscFunctionReturn(0);
2996 }
2997 
2998 #undef __FUNCT__
2999 #define __FUNCT__ "SNESGetDM"
3000 /*@
3001    SNESGetDM - Gets the DM that may be used by some preconditioners
3002 
3003    Collective on SNES
3004 
3005    Input Parameter:
3006 . snes - the preconditioner context
3007 
3008    Output Parameter:
3009 .  dm - the dm
3010 
3011    Level: intermediate
3012 
3013 
3014 .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
3015 @*/
3016 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetDM(SNES snes,DM *dm)
3017 {
3018   PetscFunctionBegin;
3019   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3020   *dm = snes->dm;
3021   PetscFunctionReturn(0);
3022 }
3023