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