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