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