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