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