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