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