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