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