xref: /petsc/src/snes/interface/snes.c (revision 09573ac72a50d3e7ecd55a2b7f0ef28450cd0a8b)
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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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 PETSCSNES_DLLEXPORT 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);
1162     PetscStackPop;
1163     if (PetscExceptionValue(ierr)) {
1164       PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr);
1165     }
1166     CHKERRQ(ierr);
1167   } else if (snes->vec_rhs) {
1168     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
1169   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
1170   if (snes->vec_rhs) {
1171     ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr);
1172   }
1173   snes->nfuncs++;
1174   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1175   PetscFunctionReturn(0);
1176 }
1177 
1178 #undef __FUNCT__
1179 #define __FUNCT__ "SNESComputeJacobian"
1180 /*@
1181    SNESComputeJacobian - Computes the Jacobian matrix that has been
1182    set with SNESSetJacobian().
1183 
1184    Collective on SNES and Mat
1185 
1186    Input Parameters:
1187 +  snes - the SNES context
1188 -  x - input vector
1189 
1190    Output Parameters:
1191 +  A - Jacobian matrix
1192 .  B - optional preconditioning matrix
1193 -  flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER)
1194 
1195   Options Database Keys:
1196 +    -snes_lag_preconditioner <lag>
1197 -    -snes_lag_jacobian <lag>
1198 
1199    Notes:
1200    Most users should not need to explicitly call this routine, as it
1201    is used internally within the nonlinear solvers.
1202 
1203    See KSPSetOperators() for important information about setting the
1204    flag parameter.
1205 
1206    Level: developer
1207 
1208 .keywords: SNES, compute, Jacobian, matrix
1209 
1210 .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
1211 @*/
1212 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
1213 {
1214   PetscErrorCode ierr;
1215   PetscBool      flag;
1216 
1217   PetscFunctionBegin;
1218   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1219   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1220   PetscValidPointer(flg,5);
1221   PetscCheckSameComm(snes,1,X,2);
1222   if (!snes->ops->computejacobian) PetscFunctionReturn(0);
1223 
1224   /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
1225 
1226   if (snes->lagjacobian == -2) {
1227     snes->lagjacobian = -1;
1228     ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr);
1229   } else if (snes->lagjacobian == -1) {
1230     *flg = SAME_PRECONDITIONER;
1231     ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr);
1232     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1233     if (flag) {
1234       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1235       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1236     }
1237     PetscFunctionReturn(0);
1238   } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) {
1239     *flg = SAME_PRECONDITIONER;
1240     ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr);
1241     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1242     if (flag) {
1243       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1244       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1245     }
1246     PetscFunctionReturn(0);
1247   }
1248 
1249   *flg = DIFFERENT_NONZERO_PATTERN;
1250   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1251   PetscStackPush("SNES user Jacobian function");
1252   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1253   PetscStackPop;
1254   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1255 
1256   if (snes->lagpreconditioner == -2) {
1257     ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr);
1258     snes->lagpreconditioner = -1;
1259   } else if (snes->lagpreconditioner == -1) {
1260     *flg = SAME_PRECONDITIONER;
1261     ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr);
1262   } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) {
1263     *flg = SAME_PRECONDITIONER;
1264     ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr);
1265   }
1266 
1267   /* make sure user returned a correct Jacobian and preconditioner */
1268   /* PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
1269     PetscValidHeaderSpecific(*B,MAT_CLASSID,4);   */
1270   PetscFunctionReturn(0);
1271 }
1272 
1273 #undef __FUNCT__
1274 #define __FUNCT__ "SNESSetJacobian"
1275 /*@C
1276    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1277    location to store the matrix.
1278 
1279    Logically Collective on SNES and Mat
1280 
1281    Input Parameters:
1282 +  snes - the SNES context
1283 .  A - Jacobian matrix
1284 .  B - preconditioner matrix (usually same as the Jacobian)
1285 .  func - Jacobian evaluation routine
1286 -  ctx - [optional] user-defined context for private data for the
1287          Jacobian evaluation routine (may be PETSC_NULL)
1288 
1289    Calling sequence of func:
1290 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1291 
1292 +  x - input vector
1293 .  A - Jacobian matrix
1294 .  B - preconditioner matrix, usually the same as A
1295 .  flag - flag indicating information about the preconditioner matrix
1296    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1297 -  ctx - [optional] user-defined Jacobian context
1298 
1299    Notes:
1300    See KSPSetOperators() for important information about setting the flag
1301    output parameter in the routine func().  Be sure to read this information!
1302 
1303    The routine func() takes Mat * as the matrix arguments rather than Mat.
1304    This allows the Jacobian evaluation routine to replace A and/or B with a
1305    completely new new matrix structure (not just different matrix elements)
1306    when appropriate, for instance, if the nonzero structure is changing
1307    throughout the global iterations.
1308 
1309    If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on
1310    each matrix.
1311 
1312    Level: beginner
1313 
1314 .keywords: SNES, nonlinear, set, Jacobian, matrix
1315 
1316 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
1317 @*/
1318 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1319 {
1320   PetscErrorCode ierr;
1321 
1322   PetscFunctionBegin;
1323   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1324   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
1325   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
1326   if (A) PetscCheckSameComm(snes,1,A,2);
1327   if (B) PetscCheckSameComm(snes,1,B,3);
1328   if (func) snes->ops->computejacobian = func;
1329   if (ctx)  snes->jacP                 = ctx;
1330   if (A) {
1331     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1332     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1333     snes->jacobian = A;
1334   }
1335   if (B) {
1336     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
1337     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1338     snes->jacobian_pre = B;
1339   }
1340   PetscFunctionReturn(0);
1341 }
1342 
1343 #undef __FUNCT__
1344 #define __FUNCT__ "SNESGetJacobian"
1345 /*@C
1346    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1347    provided context for evaluating the Jacobian.
1348 
1349    Not Collective, but Mat object will be parallel if SNES object is
1350 
1351    Input Parameter:
1352 .  snes - the nonlinear solver context
1353 
1354    Output Parameters:
1355 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1356 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1357 .  func - location to put Jacobian function (or PETSC_NULL)
1358 -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1359 
1360    Level: advanced
1361 
1362 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1363 @*/
1364 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1365 {
1366   PetscFunctionBegin;
1367   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1368   if (A)    *A    = snes->jacobian;
1369   if (B)    *B    = snes->jacobian_pre;
1370   if (func) *func = snes->ops->computejacobian;
1371   if (ctx)  *ctx  = snes->jacP;
1372   PetscFunctionReturn(0);
1373 }
1374 
1375 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1376 
1377 #undef __FUNCT__
1378 #define __FUNCT__ "SNESSetUp"
1379 /*@
1380    SNESSetUp - Sets up the internal data structures for the later use
1381    of a nonlinear solver.
1382 
1383    Collective on SNES
1384 
1385    Input Parameters:
1386 .  snes - the SNES context
1387 
1388    Notes:
1389    For basic use of the SNES solvers the user need not explicitly call
1390    SNESSetUp(), since these actions will automatically occur during
1391    the call to SNESSolve().  However, if one wishes to control this
1392    phase separately, SNESSetUp() should be called after SNESCreate()
1393    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1394 
1395    Level: advanced
1396 
1397 .keywords: SNES, nonlinear, setup
1398 
1399 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1400 @*/
1401 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
1402 {
1403   PetscErrorCode ierr;
1404 
1405   PetscFunctionBegin;
1406   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1407   if (snes->setupcalled) PetscFunctionReturn(0);
1408 
1409   if (!((PetscObject)snes)->type_name) {
1410     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
1411   }
1412 
1413   if (!snes->vec_func && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1414   if (!snes->ops->computefunction && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1415   if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1416   if (!snes->ops->computejacobian && snes->dm) {
1417     Mat           J;
1418     ISColoring    coloring;
1419     MatFDColoring fd;
1420 
1421     ierr = DMGetMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr);
1422     ierr = DMGetColoring(snes->dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr);
1423     ierr = MatFDColoringCreate(J,coloring,&fd);CHKERRQ(ierr);
1424     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
1425     ierr = SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fd);CHKERRQ(ierr);
1426     ierr = ISColoringDestroy(coloring);CHKERRQ(ierr);
1427   }
1428   if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);}
1429 
1430   if (snes->ops->setup) {
1431     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
1432   }
1433   snes->setupcalled = PETSC_TRUE;
1434   PetscFunctionReturn(0);
1435 }
1436 
1437 #undef __FUNCT__
1438 #define __FUNCT__ "SNESDestroy"
1439 /*@
1440    SNESDestroy - Destroys the nonlinear solver context that was created
1441    with SNESCreate().
1442 
1443    Collective on SNES
1444 
1445    Input Parameter:
1446 .  snes - the SNES context
1447 
1448    Level: beginner
1449 
1450 .keywords: SNES, nonlinear, destroy
1451 
1452 .seealso: SNESCreate(), SNESSolve()
1453 @*/
1454 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
1455 {
1456   PetscErrorCode ierr;
1457 
1458   PetscFunctionBegin;
1459   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1460   if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0);
1461 
1462   /* if memory was published with AMS then destroy it */
1463   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1464 
1465   if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);}
1466   if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);}
1467 
1468   if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);}
1469   if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);}
1470   if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);}
1471   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1472   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1473   if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);}
1474   ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);
1475   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1476   ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
1477   if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);}
1478   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
1479   PetscFunctionReturn(0);
1480 }
1481 
1482 /* ----------- Routines to set solver parameters ---------- */
1483 
1484 #undef __FUNCT__
1485 #define __FUNCT__ "SNESSetLagPreconditioner"
1486 /*@
1487    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
1488 
1489    Logically Collective on SNES
1490 
1491    Input Parameters:
1492 +  snes - the SNES context
1493 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1494          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1495 
1496    Options Database Keys:
1497 .    -snes_lag_preconditioner <lag>
1498 
1499    Notes:
1500    The default is 1
1501    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1502    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
1503 
1504    Level: intermediate
1505 
1506 .keywords: SNES, nonlinear, set, convergence, tolerances
1507 
1508 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
1509 
1510 @*/
1511 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag)
1512 {
1513   PetscFunctionBegin;
1514   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1515   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1516   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1517   PetscValidLogicalCollectiveInt(snes,lag,2);
1518   snes->lagpreconditioner = lag;
1519   PetscFunctionReturn(0);
1520 }
1521 
1522 #undef __FUNCT__
1523 #define __FUNCT__ "SNESGetLagPreconditioner"
1524 /*@
1525    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
1526 
1527    Not Collective
1528 
1529    Input Parameter:
1530 .  snes - the SNES context
1531 
1532    Output Parameter:
1533 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1534          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1535 
1536    Options Database Keys:
1537 .    -snes_lag_preconditioner <lag>
1538 
1539    Notes:
1540    The default is 1
1541    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1542 
1543    Level: intermediate
1544 
1545 .keywords: SNES, nonlinear, set, convergence, tolerances
1546 
1547 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
1548 
1549 @*/
1550 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
1551 {
1552   PetscFunctionBegin;
1553   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1554   *lag = snes->lagpreconditioner;
1555   PetscFunctionReturn(0);
1556 }
1557 
1558 #undef __FUNCT__
1559 #define __FUNCT__ "SNESSetLagJacobian"
1560 /*@
1561    SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
1562      often the preconditioner is rebuilt.
1563 
1564    Logically Collective on SNES
1565 
1566    Input Parameters:
1567 +  snes - the SNES context
1568 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1569          the Jacobian is built etc. -2 means rebuild at next chance but then never again
1570 
1571    Options Database Keys:
1572 .    -snes_lag_jacobian <lag>
1573 
1574    Notes:
1575    The default is 1
1576    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1577    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
1578    at the next Newton step but never again (unless it is reset to another value)
1579 
1580    Level: intermediate
1581 
1582 .keywords: SNES, nonlinear, set, convergence, tolerances
1583 
1584 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
1585 
1586 @*/
1587 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagJacobian(SNES snes,PetscInt lag)
1588 {
1589   PetscFunctionBegin;
1590   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1591   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1592   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1593   PetscValidLogicalCollectiveInt(snes,lag,2);
1594   snes->lagjacobian = lag;
1595   PetscFunctionReturn(0);
1596 }
1597 
1598 #undef __FUNCT__
1599 #define __FUNCT__ "SNESGetLagJacobian"
1600 /*@
1601    SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
1602 
1603    Not Collective
1604 
1605    Input Parameter:
1606 .  snes - the SNES context
1607 
1608    Output Parameter:
1609 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1610          the Jacobian is built etc.
1611 
1612    Options Database Keys:
1613 .    -snes_lag_jacobian <lag>
1614 
1615    Notes:
1616    The default is 1
1617    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1618 
1619    Level: intermediate
1620 
1621 .keywords: SNES, nonlinear, set, convergence, tolerances
1622 
1623 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
1624 
1625 @*/
1626 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagJacobian(SNES snes,PetscInt *lag)
1627 {
1628   PetscFunctionBegin;
1629   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1630   *lag = snes->lagjacobian;
1631   PetscFunctionReturn(0);
1632 }
1633 
1634 #undef __FUNCT__
1635 #define __FUNCT__ "SNESSetTolerances"
1636 /*@
1637    SNESSetTolerances - Sets various parameters used in convergence tests.
1638 
1639    Logically Collective on SNES
1640 
1641    Input Parameters:
1642 +  snes - the SNES context
1643 .  abstol - absolute convergence tolerance
1644 .  rtol - relative convergence tolerance
1645 .  stol -  convergence tolerance in terms of the norm
1646            of the change in the solution between steps
1647 .  maxit - maximum number of iterations
1648 -  maxf - maximum number of function evaluations
1649 
1650    Options Database Keys:
1651 +    -snes_atol <abstol> - Sets abstol
1652 .    -snes_rtol <rtol> - Sets rtol
1653 .    -snes_stol <stol> - Sets stol
1654 .    -snes_max_it <maxit> - Sets maxit
1655 -    -snes_max_funcs <maxf> - Sets maxf
1656 
1657    Notes:
1658    The default maximum number of iterations is 50.
1659    The default maximum number of function evaluations is 1000.
1660 
1661    Level: intermediate
1662 
1663 .keywords: SNES, nonlinear, set, convergence, tolerances
1664 
1665 .seealso: SNESSetTrustRegionTolerance()
1666 @*/
1667 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
1668 {
1669   PetscFunctionBegin;
1670   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1671   PetscValidLogicalCollectiveReal(snes,abstol,2);
1672   PetscValidLogicalCollectiveReal(snes,rtol,3);
1673   PetscValidLogicalCollectiveReal(snes,stol,4);
1674   PetscValidLogicalCollectiveInt(snes,maxit,5);
1675   PetscValidLogicalCollectiveInt(snes,maxf,6);
1676 
1677   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1678   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1679   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1680   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1681   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1682   PetscFunctionReturn(0);
1683 }
1684 
1685 #undef __FUNCT__
1686 #define __FUNCT__ "SNESGetTolerances"
1687 /*@
1688    SNESGetTolerances - Gets various parameters used in convergence tests.
1689 
1690    Not Collective
1691 
1692    Input Parameters:
1693 +  snes - the SNES context
1694 .  atol - absolute convergence tolerance
1695 .  rtol - relative convergence tolerance
1696 .  stol -  convergence tolerance in terms of the norm
1697            of the change in the solution between steps
1698 .  maxit - maximum number of iterations
1699 -  maxf - maximum number of function evaluations
1700 
1701    Notes:
1702    The user can specify PETSC_NULL for any parameter that is not needed.
1703 
1704    Level: intermediate
1705 
1706 .keywords: SNES, nonlinear, get, convergence, tolerances
1707 
1708 .seealso: SNESSetTolerances()
1709 @*/
1710 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
1711 {
1712   PetscFunctionBegin;
1713   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1714   if (atol)  *atol  = snes->abstol;
1715   if (rtol)  *rtol  = snes->rtol;
1716   if (stol)  *stol  = snes->xtol;
1717   if (maxit) *maxit = snes->max_its;
1718   if (maxf)  *maxf  = snes->max_funcs;
1719   PetscFunctionReturn(0);
1720 }
1721 
1722 #undef __FUNCT__
1723 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1724 /*@
1725    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1726 
1727    Logically Collective on SNES
1728 
1729    Input Parameters:
1730 +  snes - the SNES context
1731 -  tol - tolerance
1732 
1733    Options Database Key:
1734 .  -snes_trtol <tol> - Sets tol
1735 
1736    Level: intermediate
1737 
1738 .keywords: SNES, nonlinear, set, trust region, tolerance
1739 
1740 .seealso: SNESSetTolerances()
1741 @*/
1742 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1743 {
1744   PetscFunctionBegin;
1745   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1746   PetscValidLogicalCollectiveReal(snes,tol,2);
1747   snes->deltatol = tol;
1748   PetscFunctionReturn(0);
1749 }
1750 
1751 /*
1752    Duplicate the lg monitors for SNES from KSP; for some reason with
1753    dynamic libraries things don't work under Sun4 if we just use
1754    macros instead of functions
1755 */
1756 #undef __FUNCT__
1757 #define __FUNCT__ "SNESMonitorLG"
1758 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1759 {
1760   PetscErrorCode ierr;
1761 
1762   PetscFunctionBegin;
1763   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1764   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1765   PetscFunctionReturn(0);
1766 }
1767 
1768 #undef __FUNCT__
1769 #define __FUNCT__ "SNESMonitorLGCreate"
1770 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1771 {
1772   PetscErrorCode ierr;
1773 
1774   PetscFunctionBegin;
1775   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1776   PetscFunctionReturn(0);
1777 }
1778 
1779 #undef __FUNCT__
1780 #define __FUNCT__ "SNESMonitorLGDestroy"
1781 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw)
1782 {
1783   PetscErrorCode ierr;
1784 
1785   PetscFunctionBegin;
1786   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1787   PetscFunctionReturn(0);
1788 }
1789 
1790 extern PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
1791 #undef __FUNCT__
1792 #define __FUNCT__ "SNESMonitorLGRange"
1793 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
1794 {
1795   PetscDrawLG      lg;
1796   PetscErrorCode   ierr;
1797   PetscReal        x,y,per;
1798   PetscViewer      v = (PetscViewer)monctx;
1799   static PetscReal prev; /* should be in the context */
1800   PetscDraw        draw;
1801   PetscFunctionBegin;
1802   if (!monctx) {
1803     MPI_Comm    comm;
1804 
1805     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1806     v      = PETSC_VIEWER_DRAW_(comm);
1807   }
1808   ierr   = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr);
1809   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1810   ierr   = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1811   ierr   = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr);
1812   x = (PetscReal) n;
1813   if (rnorm > 0.0) y = log10(rnorm); else y = -15.0;
1814   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1815   if (n < 20 || !(n % 5)) {
1816     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1817   }
1818 
1819   ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr);
1820   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1821   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1822   ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr);
1823   ierr =  SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr);
1824   x = (PetscReal) n;
1825   y = 100.0*per;
1826   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1827   if (n < 20 || !(n % 5)) {
1828     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1829   }
1830 
1831   ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr);
1832   if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1833   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1834   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr);
1835   x = (PetscReal) n;
1836   y = (prev - rnorm)/prev;
1837   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1838   if (n < 20 || !(n % 5)) {
1839     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1840   }
1841 
1842   ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr);
1843   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1844   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1845   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
1846   x = (PetscReal) n;
1847   y = (prev - rnorm)/(prev*per);
1848   if (n > 2) { /*skip initial crazy value */
1849     ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1850   }
1851   if (n < 20 || !(n % 5)) {
1852     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1853   }
1854   prev = rnorm;
1855   PetscFunctionReturn(0);
1856 }
1857 
1858 #undef __FUNCT__
1859 #define __FUNCT__ "SNESMonitorLGRangeCreate"
1860 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1861 {
1862   PetscErrorCode ierr;
1863 
1864   PetscFunctionBegin;
1865   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1866   PetscFunctionReturn(0);
1867 }
1868 
1869 #undef __FUNCT__
1870 #define __FUNCT__ "SNESMonitorLGRangeDestroy"
1871 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGRangeDestroy(PetscDrawLG draw)
1872 {
1873   PetscErrorCode ierr;
1874 
1875   PetscFunctionBegin;
1876   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1877   PetscFunctionReturn(0);
1878 }
1879 
1880 /* ------------ Routines to set performance monitoring options ----------- */
1881 
1882 #undef __FUNCT__
1883 #define __FUNCT__ "SNESMonitorSet"
1884 /*@C
1885    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
1886    iteration of the nonlinear solver to display the iteration's
1887    progress.
1888 
1889    Logically Collective on SNES
1890 
1891    Input Parameters:
1892 +  snes - the SNES context
1893 .  func - monitoring routine
1894 .  mctx - [optional] user-defined context for private data for the
1895           monitor routine (use PETSC_NULL if no context is desired)
1896 -  monitordestroy - [optional] routine that frees monitor context
1897           (may be PETSC_NULL)
1898 
1899    Calling sequence of func:
1900 $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1901 
1902 +    snes - the SNES context
1903 .    its - iteration number
1904 .    norm - 2-norm function value (may be estimated)
1905 -    mctx - [optional] monitoring context
1906 
1907    Options Database Keys:
1908 +    -snes_monitor        - sets SNESMonitorDefault()
1909 .    -snes_monitor_draw    - sets line graph monitor,
1910                             uses SNESMonitorLGCreate()
1911 _    -snes_monitor_cancel - cancels all monitors that have
1912                             been hardwired into a code by
1913                             calls to SNESMonitorSet(), but
1914                             does not cancel those set via
1915                             the options database.
1916 
1917    Notes:
1918    Several different monitoring routines may be set by calling
1919    SNESMonitorSet() multiple times; all will be called in the
1920    order in which they were set.
1921 
1922    Fortran notes: Only a single monitor function can be set for each SNES object
1923 
1924    Level: intermediate
1925 
1926 .keywords: SNES, nonlinear, set, monitor
1927 
1928 .seealso: SNESMonitorDefault(), SNESMonitorCancel()
1929 @*/
1930 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
1931 {
1932   PetscInt i;
1933 
1934   PetscFunctionBegin;
1935   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1936   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1937   for (i=0; i<snes->numbermonitors;i++) {
1938     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0);
1939 
1940     /* check if both default monitors that share common ASCII viewer */
1941     if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) {
1942       if (mctx && snes->monitorcontext[i]) {
1943         PetscErrorCode          ierr;
1944         PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx;
1945         PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i];
1946         if (viewer1->viewer == viewer2->viewer) {
1947           ierr = (*monitordestroy)(mctx);CHKERRQ(ierr);
1948           PetscFunctionReturn(0);
1949         }
1950       }
1951     }
1952   }
1953   snes->monitor[snes->numbermonitors]           = monitor;
1954   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1955   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1956   PetscFunctionReturn(0);
1957 }
1958 
1959 #undef __FUNCT__
1960 #define __FUNCT__ "SNESMonitorCancel"
1961 /*@C
1962    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
1963 
1964    Logically Collective on SNES
1965 
1966    Input Parameters:
1967 .  snes - the SNES context
1968 
1969    Options Database Key:
1970 .  -snes_monitor_cancel - cancels all monitors that have been hardwired
1971     into a code by calls to SNESMonitorSet(), but does not cancel those
1972     set via the options database
1973 
1974    Notes:
1975    There is no way to clear one specific monitor from a SNES object.
1976 
1977    Level: intermediate
1978 
1979 .keywords: SNES, nonlinear, set, monitor
1980 
1981 .seealso: SNESMonitorDefault(), SNESMonitorSet()
1982 @*/
1983 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes)
1984 {
1985   PetscErrorCode ierr;
1986   PetscInt       i;
1987 
1988   PetscFunctionBegin;
1989   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1990   for (i=0; i<snes->numbermonitors; i++) {
1991     if (snes->monitordestroy[i]) {
1992       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1993     }
1994   }
1995   snes->numbermonitors = 0;
1996   PetscFunctionReturn(0);
1997 }
1998 
1999 #undef __FUNCT__
2000 #define __FUNCT__ "SNESSetConvergenceTest"
2001 /*@C
2002    SNESSetConvergenceTest - Sets the function that is to be used
2003    to test for convergence of the nonlinear iterative solution.
2004 
2005    Logically Collective on SNES
2006 
2007    Input Parameters:
2008 +  snes - the SNES context
2009 .  func - routine to test for convergence
2010 .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
2011 -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
2012 
2013    Calling sequence of func:
2014 $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
2015 
2016 +    snes - the SNES context
2017 .    it - current iteration (0 is the first and is before any Newton step)
2018 .    cctx - [optional] convergence context
2019 .    reason - reason for convergence/divergence
2020 .    xnorm - 2-norm of current iterate
2021 .    gnorm - 2-norm of current step
2022 -    f - 2-norm of function
2023 
2024    Level: advanced
2025 
2026 .keywords: SNES, nonlinear, set, convergence, test
2027 
2028 .seealso: SNESDefaultConverged(), SNESSkipConverged()
2029 @*/
2030 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
2031 {
2032   PetscErrorCode ierr;
2033 
2034   PetscFunctionBegin;
2035   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2036   if (!func) func = SNESSkipConverged;
2037   if (snes->ops->convergeddestroy) {
2038     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
2039   }
2040   snes->ops->converged        = func;
2041   snes->ops->convergeddestroy = destroy;
2042   snes->cnvP                  = cctx;
2043   PetscFunctionReturn(0);
2044 }
2045 
2046 #undef __FUNCT__
2047 #define __FUNCT__ "SNESGetConvergedReason"
2048 /*@
2049    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
2050 
2051    Not Collective
2052 
2053    Input Parameter:
2054 .  snes - the SNES context
2055 
2056    Output Parameter:
2057 .  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
2058             manual pages for the individual convergence tests for complete lists
2059 
2060    Level: intermediate
2061 
2062    Notes: Can only be called after the call the SNESSolve() is complete.
2063 
2064 .keywords: SNES, nonlinear, set, convergence, test
2065 
2066 .seealso: SNESSetConvergenceTest(), SNESConvergedReason
2067 @*/
2068 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
2069 {
2070   PetscFunctionBegin;
2071   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2072   PetscValidPointer(reason,2);
2073   *reason = snes->reason;
2074   PetscFunctionReturn(0);
2075 }
2076 
2077 #undef __FUNCT__
2078 #define __FUNCT__ "SNESSetConvergenceHistory"
2079 /*@
2080    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
2081 
2082    Logically Collective on SNES
2083 
2084    Input Parameters:
2085 +  snes - iterative context obtained from SNESCreate()
2086 .  a   - array to hold history, this array will contain the function norms computed at each step
2087 .  its - integer array holds the number of linear iterations for each solve.
2088 .  na  - size of a and its
2089 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
2090            else it continues storing new values for new nonlinear solves after the old ones
2091 
2092    This routine is useful, e.g., when running a code for purposes
2093    of accurate performance monitoring, when no I/O should be done
2094    during the section of code that is being timed.
2095 
2096    Level: intermediate
2097 
2098 .keywords: SNES, set, convergence, history
2099 
2100 .seealso: SNESGetConvergenceHistory()
2101 
2102 @*/
2103 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool  reset)
2104 {
2105   PetscFunctionBegin;
2106   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2107   if (na)  PetscValidScalarPointer(a,2);
2108   if (its) PetscValidIntPointer(its,3);
2109   snes->conv_hist       = a;
2110   snes->conv_hist_its   = its;
2111   snes->conv_hist_max   = na;
2112   snes->conv_hist_len   = 0;
2113   snes->conv_hist_reset = reset;
2114   PetscFunctionReturn(0);
2115 }
2116 
2117 #undef __FUNCT__
2118 #define __FUNCT__ "SNESGetConvergenceHistory"
2119 /*@C
2120    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
2121 
2122    Not Collective
2123 
2124    Input Parameter:
2125 .  snes - iterative context obtained from SNESCreate()
2126 
2127    Output Parameters:
2128 .  a   - array to hold history
2129 .  its - integer array holds the number of linear iterations (or
2130          negative if not converged) for each solve.
2131 -  na  - size of a and its
2132 
2133    Notes:
2134     The calling sequence for this routine in Fortran is
2135 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
2136 
2137    This routine is useful, e.g., when running a code for purposes
2138    of accurate performance monitoring, when no I/O should be done
2139    during the section of code that is being timed.
2140 
2141    Level: intermediate
2142 
2143 .keywords: SNES, get, convergence, history
2144 
2145 .seealso: SNESSetConvergencHistory()
2146 
2147 @*/
2148 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
2149 {
2150   PetscFunctionBegin;
2151   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2152   if (a)   *a   = snes->conv_hist;
2153   if (its) *its = snes->conv_hist_its;
2154   if (na)  *na  = snes->conv_hist_len;
2155   PetscFunctionReturn(0);
2156 }
2157 
2158 #undef __FUNCT__
2159 #define __FUNCT__ "SNESSetUpdate"
2160 /*@C
2161   SNESSetUpdate - Sets the general-purpose update function called
2162   at the beginning of every iteration of the nonlinear solve. Specifically
2163   it is called just before the Jacobian is "evaluated".
2164 
2165   Logically Collective on SNES
2166 
2167   Input Parameters:
2168 . snes - The nonlinear solver context
2169 . func - The function
2170 
2171   Calling sequence of func:
2172 . func (SNES snes, PetscInt step);
2173 
2174 . step - The current step of the iteration
2175 
2176   Level: advanced
2177 
2178   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()
2179         This is not used by most users.
2180 
2181 .keywords: SNES, update
2182 
2183 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
2184 @*/
2185 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
2186 {
2187   PetscFunctionBegin;
2188   PetscValidHeaderSpecific(snes, SNES_CLASSID,1);
2189   snes->ops->update = func;
2190   PetscFunctionReturn(0);
2191 }
2192 
2193 #undef __FUNCT__
2194 #define __FUNCT__ "SNESDefaultUpdate"
2195 /*@
2196   SNESDefaultUpdate - The default update function which does nothing.
2197 
2198   Not collective
2199 
2200   Input Parameters:
2201 . snes - The nonlinear solver context
2202 . step - The current step of the iteration
2203 
2204   Level: intermediate
2205 
2206 .keywords: SNES, update
2207 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
2208 @*/
2209 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
2210 {
2211   PetscFunctionBegin;
2212   PetscFunctionReturn(0);
2213 }
2214 
2215 #undef __FUNCT__
2216 #define __FUNCT__ "SNESScaleStep_Private"
2217 /*
2218    SNESScaleStep_Private - Scales a step so that its length is less than the
2219    positive parameter delta.
2220 
2221     Input Parameters:
2222 +   snes - the SNES context
2223 .   y - approximate solution of linear system
2224 .   fnorm - 2-norm of current function
2225 -   delta - trust region size
2226 
2227     Output Parameters:
2228 +   gpnorm - predicted function norm at the new point, assuming local
2229     linearization.  The value is zero if the step lies within the trust
2230     region, and exceeds zero otherwise.
2231 -   ynorm - 2-norm of the step
2232 
2233     Note:
2234     For non-trust region methods such as SNESLS, the parameter delta
2235     is set to be the maximum allowable step size.
2236 
2237 .keywords: SNES, nonlinear, scale, step
2238 */
2239 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
2240 {
2241   PetscReal      nrm;
2242   PetscScalar    cnorm;
2243   PetscErrorCode ierr;
2244 
2245   PetscFunctionBegin;
2246   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2247   PetscValidHeaderSpecific(y,VEC_CLASSID,2);
2248   PetscCheckSameComm(snes,1,y,2);
2249 
2250   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
2251   if (nrm > *delta) {
2252      nrm = *delta/nrm;
2253      *gpnorm = (1.0 - nrm)*(*fnorm);
2254      cnorm = nrm;
2255      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
2256      *ynorm = *delta;
2257   } else {
2258      *gpnorm = 0.0;
2259      *ynorm = nrm;
2260   }
2261   PetscFunctionReturn(0);
2262 }
2263 
2264 #undef __FUNCT__
2265 #define __FUNCT__ "SNESSolve"
2266 /*@C
2267    SNESSolve - Solves a nonlinear system F(x) = b.
2268    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
2269 
2270    Collective on SNES
2271 
2272    Input Parameters:
2273 +  snes - the SNES context
2274 .  b - the constant part of the equation, or PETSC_NULL to use zero.
2275 -  x - the solution vector.
2276 
2277    Notes:
2278    The user should initialize the vector,x, with the initial guess
2279    for the nonlinear solve prior to calling SNESSolve.  In particular,
2280    to employ an initial guess of zero, the user should explicitly set
2281    this vector to zero by calling VecSet().
2282 
2283    Level: beginner
2284 
2285 .keywords: SNES, nonlinear, solve
2286 
2287 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian()
2288 @*/
2289 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
2290 {
2291   PetscErrorCode ierr;
2292   PetscBool      flg;
2293   char           filename[PETSC_MAX_PATH_LEN];
2294   PetscViewer    viewer;
2295 
2296   PetscFunctionBegin;
2297   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2298   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2299   PetscCheckSameComm(snes,1,x,3);
2300   if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2);
2301   if (b) PetscCheckSameComm(snes,1,b,2);
2302 
2303   /* set solution vector */
2304   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
2305   if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); }
2306   snes->vec_sol = x;
2307   /* set afine vector if provided */
2308   if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
2309   if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); }
2310   snes->vec_rhs = b;
2311 
2312   if (!snes->vec_func && snes->vec_rhs) {
2313     ierr = VecDuplicate(b, &snes->vec_func);CHKERRQ(ierr);
2314   }
2315 
2316   ierr = SNESSetUp(snes);CHKERRQ(ierr);
2317 
2318   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
2319   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
2320 
2321   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2322   ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
2323   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2324   if (snes->domainerror){
2325     snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
2326     snes->domainerror = PETSC_FALSE;
2327   }
2328   if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
2329 
2330   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
2331   if (flg && !PetscPreLoadingOn) {
2332     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
2333     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
2334     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
2335   }
2336 
2337   flg  = PETSC_FALSE;
2338   ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,PETSC_NULL);CHKERRQ(ierr);
2339   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
2340   if (snes->printreason) {
2341     if (snes->reason > 0) {
2342       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2343     } else {
2344       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2345     }
2346   }
2347 
2348   if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
2349   PetscFunctionReturn(0);
2350 }
2351 
2352 /* --------- Internal routines for SNES Package --------- */
2353 
2354 #undef __FUNCT__
2355 #define __FUNCT__ "SNESSetType"
2356 /*@C
2357    SNESSetType - Sets the method for the nonlinear solver.
2358 
2359    Collective on SNES
2360 
2361    Input Parameters:
2362 +  snes - the SNES context
2363 -  type - a known method
2364 
2365    Options Database Key:
2366 .  -snes_type <type> - Sets the method; use -help for a list
2367    of available methods (for instance, ls or tr)
2368 
2369    Notes:
2370    See "petsc/include/petscsnes.h" for available methods (for instance)
2371 +    SNESLS - Newton's method with line search
2372      (systems of nonlinear equations)
2373 .    SNESTR - Newton's method with trust region
2374      (systems of nonlinear equations)
2375 
2376   Normally, it is best to use the SNESSetFromOptions() command and then
2377   set the SNES solver type from the options database rather than by using
2378   this routine.  Using the options database provides the user with
2379   maximum flexibility in evaluating the many nonlinear solvers.
2380   The SNESSetType() routine is provided for those situations where it
2381   is necessary to set the nonlinear solver independently of the command
2382   line or options database.  This might be the case, for example, when
2383   the choice of solver changes during the execution of the program,
2384   and the user's application is taking responsibility for choosing the
2385   appropriate method.
2386 
2387   Level: intermediate
2388 
2389 .keywords: SNES, set, type
2390 
2391 .seealso: SNESType, SNESCreate()
2392 
2393 @*/
2394 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type)
2395 {
2396   PetscErrorCode ierr,(*r)(SNES);
2397   PetscBool      match;
2398 
2399   PetscFunctionBegin;
2400   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2401   PetscValidCharPointer(type,2);
2402 
2403   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
2404   if (match) PetscFunctionReturn(0);
2405 
2406   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
2407   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
2408   /* Destroy the previous private SNES context */
2409   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
2410   /* Reinitialize function pointers in SNESOps structure */
2411   snes->ops->setup          = 0;
2412   snes->ops->solve          = 0;
2413   snes->ops->view           = 0;
2414   snes->ops->setfromoptions = 0;
2415   snes->ops->destroy        = 0;
2416   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
2417   snes->setupcalled = PETSC_FALSE;
2418   ierr = (*r)(snes);CHKERRQ(ierr);
2419   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
2420 #if defined(PETSC_HAVE_AMS)
2421   if (PetscAMSPublishAll) {
2422     ierr = PetscObjectAMSPublish((PetscObject)snes);CHKERRQ(ierr);
2423   }
2424 #endif
2425   PetscFunctionReturn(0);
2426 }
2427 
2428 
2429 /* --------------------------------------------------------------------- */
2430 #undef __FUNCT__
2431 #define __FUNCT__ "SNESRegisterDestroy"
2432 /*@
2433    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2434    registered by SNESRegisterDynamic().
2435 
2436    Not Collective
2437 
2438    Level: advanced
2439 
2440 .keywords: SNES, nonlinear, register, destroy
2441 
2442 .seealso: SNESRegisterAll(), SNESRegisterAll()
2443 @*/
2444 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
2445 {
2446   PetscErrorCode ierr;
2447 
2448   PetscFunctionBegin;
2449   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
2450   SNESRegisterAllCalled = PETSC_FALSE;
2451   PetscFunctionReturn(0);
2452 }
2453 
2454 #undef __FUNCT__
2455 #define __FUNCT__ "SNESGetType"
2456 /*@C
2457    SNESGetType - Gets the SNES method type and name (as a string).
2458 
2459    Not Collective
2460 
2461    Input Parameter:
2462 .  snes - nonlinear solver context
2463 
2464    Output Parameter:
2465 .  type - SNES method (a character string)
2466 
2467    Level: intermediate
2468 
2469 .keywords: SNES, nonlinear, get, type, name
2470 @*/
2471 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type)
2472 {
2473   PetscFunctionBegin;
2474   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2475   PetscValidPointer(type,2);
2476   *type = ((PetscObject)snes)->type_name;
2477   PetscFunctionReturn(0);
2478 }
2479 
2480 #undef __FUNCT__
2481 #define __FUNCT__ "SNESGetSolution"
2482 /*@
2483    SNESGetSolution - Returns the vector where the approximate solution is
2484    stored.
2485 
2486    Not Collective, but Vec is parallel if SNES is parallel
2487 
2488    Input Parameter:
2489 .  snes - the SNES context
2490 
2491    Output Parameter:
2492 .  x - the solution
2493 
2494    Level: intermediate
2495 
2496 .keywords: SNES, nonlinear, get, solution
2497 
2498 .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
2499 @*/
2500 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
2501 {
2502   PetscFunctionBegin;
2503   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2504   PetscValidPointer(x,2);
2505   *x = snes->vec_sol;
2506   PetscFunctionReturn(0);
2507 }
2508 
2509 #undef __FUNCT__
2510 #define __FUNCT__ "SNESGetSolutionUpdate"
2511 /*@
2512    SNESGetSolutionUpdate - Returns the vector where the solution update is
2513    stored.
2514 
2515    Not Collective, but Vec is parallel if SNES is parallel
2516 
2517    Input Parameter:
2518 .  snes - the SNES context
2519 
2520    Output Parameter:
2521 .  x - the solution update
2522 
2523    Level: advanced
2524 
2525 .keywords: SNES, nonlinear, get, solution, update
2526 
2527 .seealso: SNESGetSolution(), SNESGetFunction()
2528 @*/
2529 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
2530 {
2531   PetscFunctionBegin;
2532   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2533   PetscValidPointer(x,2);
2534   *x = snes->vec_sol_update;
2535   PetscFunctionReturn(0);
2536 }
2537 
2538 #undef __FUNCT__
2539 #define __FUNCT__ "SNESGetFunction"
2540 /*@C
2541    SNESGetFunction - Returns the vector where the function is stored.
2542 
2543    Not Collective, but Vec is parallel if SNES is parallel
2544 
2545    Input Parameter:
2546 .  snes - the SNES context
2547 
2548    Output Parameter:
2549 +  r - the function (or PETSC_NULL)
2550 .  func - the function (or PETSC_NULL)
2551 -  ctx - the function context (or PETSC_NULL)
2552 
2553    Level: advanced
2554 
2555 .keywords: SNES, nonlinear, get, function
2556 
2557 .seealso: SNESSetFunction(), SNESGetSolution()
2558 @*/
2559 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
2560 {
2561   PetscFunctionBegin;
2562   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2563   if (r)    *r    = snes->vec_func;
2564   if (func) *func = snes->ops->computefunction;
2565   if (ctx)  *ctx  = snes->funP;
2566   PetscFunctionReturn(0);
2567 }
2568 
2569 #undef __FUNCT__
2570 #define __FUNCT__ "SNESSetOptionsPrefix"
2571 /*@C
2572    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2573    SNES options in the database.
2574 
2575    Logically Collective on SNES
2576 
2577    Input Parameter:
2578 +  snes - the SNES context
2579 -  prefix - the prefix to prepend to all option names
2580 
2581    Notes:
2582    A hyphen (-) must NOT be given at the beginning of the prefix name.
2583    The first character of all runtime options is AUTOMATICALLY the hyphen.
2584 
2585    Level: advanced
2586 
2587 .keywords: SNES, set, options, prefix, database
2588 
2589 .seealso: SNESSetFromOptions()
2590 @*/
2591 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
2592 {
2593   PetscErrorCode ierr;
2594 
2595   PetscFunctionBegin;
2596   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2597   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2598   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2599   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2600   PetscFunctionReturn(0);
2601 }
2602 
2603 #undef __FUNCT__
2604 #define __FUNCT__ "SNESAppendOptionsPrefix"
2605 /*@C
2606    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2607    SNES options in the database.
2608 
2609    Logically Collective on SNES
2610 
2611    Input Parameters:
2612 +  snes - the SNES context
2613 -  prefix - the prefix to prepend to all option names
2614 
2615    Notes:
2616    A hyphen (-) must NOT be given at the beginning of the prefix name.
2617    The first character of all runtime options is AUTOMATICALLY the hyphen.
2618 
2619    Level: advanced
2620 
2621 .keywords: SNES, append, options, prefix, database
2622 
2623 .seealso: SNESGetOptionsPrefix()
2624 @*/
2625 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
2626 {
2627   PetscErrorCode ierr;
2628 
2629   PetscFunctionBegin;
2630   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2631   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2632   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
2633   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2634   PetscFunctionReturn(0);
2635 }
2636 
2637 #undef __FUNCT__
2638 #define __FUNCT__ "SNESGetOptionsPrefix"
2639 /*@C
2640    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2641    SNES options in the database.
2642 
2643    Not Collective
2644 
2645    Input Parameter:
2646 .  snes - the SNES context
2647 
2648    Output Parameter:
2649 .  prefix - pointer to the prefix string used
2650 
2651    Notes: On the fortran side, the user should pass in a string 'prefix' of
2652    sufficient length to hold the prefix.
2653 
2654    Level: advanced
2655 
2656 .keywords: SNES, get, options, prefix, database
2657 
2658 .seealso: SNESAppendOptionsPrefix()
2659 @*/
2660 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
2661 {
2662   PetscErrorCode ierr;
2663 
2664   PetscFunctionBegin;
2665   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2666   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2667   PetscFunctionReturn(0);
2668 }
2669 
2670 
2671 #undef __FUNCT__
2672 #define __FUNCT__ "SNESRegister"
2673 /*@C
2674   SNESRegister - See SNESRegisterDynamic()
2675 
2676   Level: advanced
2677 @*/
2678 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2679 {
2680   char           fullname[PETSC_MAX_PATH_LEN];
2681   PetscErrorCode ierr;
2682 
2683   PetscFunctionBegin;
2684   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2685   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2686   PetscFunctionReturn(0);
2687 }
2688 
2689 #undef __FUNCT__
2690 #define __FUNCT__ "SNESTestLocalMin"
2691 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2692 {
2693   PetscErrorCode ierr;
2694   PetscInt       N,i,j;
2695   Vec            u,uh,fh;
2696   PetscScalar    value;
2697   PetscReal      norm;
2698 
2699   PetscFunctionBegin;
2700   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2701   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2702   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2703 
2704   /* currently only works for sequential */
2705   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2706   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2707   for (i=0; i<N; i++) {
2708     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2709     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2710     for (j=-10; j<11; j++) {
2711       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2712       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2713       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2714       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2715       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2716       value = -value;
2717       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2718     }
2719   }
2720   ierr = VecDestroy(uh);CHKERRQ(ierr);
2721   ierr = VecDestroy(fh);CHKERRQ(ierr);
2722   PetscFunctionReturn(0);
2723 }
2724 
2725 #undef __FUNCT__
2726 #define __FUNCT__ "SNESKSPSetUseEW"
2727 /*@
2728    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
2729    computing relative tolerance for linear solvers within an inexact
2730    Newton method.
2731 
2732    Logically Collective on SNES
2733 
2734    Input Parameters:
2735 +  snes - SNES context
2736 -  flag - PETSC_TRUE or PETSC_FALSE
2737 
2738     Options Database:
2739 +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
2740 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
2741 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
2742 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
2743 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
2744 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
2745 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
2746 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
2747 
2748    Notes:
2749    Currently, the default is to use a constant relative tolerance for
2750    the inner linear solvers.  Alternatively, one can use the
2751    Eisenstat-Walker method, where the relative convergence tolerance
2752    is reset at each Newton iteration according progress of the nonlinear
2753    solver.
2754 
2755    Level: advanced
2756 
2757    Reference:
2758    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2759    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2760 
2761 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2762 
2763 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2764 @*/
2765 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscBool  flag)
2766 {
2767   PetscFunctionBegin;
2768   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2769   PetscValidLogicalCollectiveBool(snes,flag,2);
2770   snes->ksp_ewconv = flag;
2771   PetscFunctionReturn(0);
2772 }
2773 
2774 #undef __FUNCT__
2775 #define __FUNCT__ "SNESKSPGetUseEW"
2776 /*@
2777    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
2778    for computing relative tolerance for linear solvers within an
2779    inexact Newton method.
2780 
2781    Not Collective
2782 
2783    Input Parameter:
2784 .  snes - SNES context
2785 
2786    Output Parameter:
2787 .  flag - PETSC_TRUE or PETSC_FALSE
2788 
2789    Notes:
2790    Currently, the default is to use a constant relative tolerance for
2791    the inner linear solvers.  Alternatively, one can use the
2792    Eisenstat-Walker method, where the relative convergence tolerance
2793    is reset at each Newton iteration according progress of the nonlinear
2794    solver.
2795 
2796    Level: advanced
2797 
2798    Reference:
2799    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2800    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2801 
2802 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2803 
2804 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2805 @*/
2806 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
2807 {
2808   PetscFunctionBegin;
2809   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2810   PetscValidPointer(flag,2);
2811   *flag = snes->ksp_ewconv;
2812   PetscFunctionReturn(0);
2813 }
2814 
2815 #undef __FUNCT__
2816 #define __FUNCT__ "SNESKSPSetParametersEW"
2817 /*@
2818    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
2819    convergence criteria for the linear solvers within an inexact
2820    Newton method.
2821 
2822    Logically Collective on SNES
2823 
2824    Input Parameters:
2825 +    snes - SNES context
2826 .    version - version 1, 2 (default is 2) or 3
2827 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2828 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2829 .    gamma - multiplicative factor for version 2 rtol computation
2830              (0 <= gamma2 <= 1)
2831 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2832 .    alpha2 - power for safeguard
2833 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2834 
2835    Note:
2836    Version 3 was contributed by Luis Chacon, June 2006.
2837 
2838    Use PETSC_DEFAULT to retain the default for any of the parameters.
2839 
2840    Level: advanced
2841 
2842    Reference:
2843    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2844    inexact Newton method", Utah State University Math. Stat. Dept. Res.
2845    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
2846 
2847 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
2848 
2849 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
2850 @*/
2851 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
2852 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
2853 {
2854   SNESKSPEW *kctx;
2855   PetscFunctionBegin;
2856   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2857   kctx = (SNESKSPEW*)snes->kspconvctx;
2858   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2859   PetscValidLogicalCollectiveInt(snes,version,2);
2860   PetscValidLogicalCollectiveReal(snes,rtol_0,3);
2861   PetscValidLogicalCollectiveReal(snes,rtol_max,4);
2862   PetscValidLogicalCollectiveReal(snes,gamma,5);
2863   PetscValidLogicalCollectiveReal(snes,alpha,6);
2864   PetscValidLogicalCollectiveReal(snes,alpha2,7);
2865   PetscValidLogicalCollectiveReal(snes,threshold,8);
2866 
2867   if (version != PETSC_DEFAULT)   kctx->version   = version;
2868   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
2869   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
2870   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
2871   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
2872   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
2873   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
2874 
2875   if (kctx->version < 1 || kctx->version > 3) {
2876     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
2877   }
2878   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
2879     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
2880   }
2881   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
2882     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
2883   }
2884   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
2885     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
2886   }
2887   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
2888     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
2889   }
2890   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
2891     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
2892   }
2893   PetscFunctionReturn(0);
2894 }
2895 
2896 #undef __FUNCT__
2897 #define __FUNCT__ "SNESKSPGetParametersEW"
2898 /*@
2899    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
2900    convergence criteria for the linear solvers within an inexact
2901    Newton method.
2902 
2903    Not Collective
2904 
2905    Input Parameters:
2906      snes - SNES context
2907 
2908    Output Parameters:
2909 +    version - version 1, 2 (default is 2) or 3
2910 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2911 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2912 .    gamma - multiplicative factor for version 2 rtol computation
2913              (0 <= gamma2 <= 1)
2914 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2915 .    alpha2 - power for safeguard
2916 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2917 
2918    Level: advanced
2919 
2920 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
2921 
2922 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
2923 @*/
2924 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
2925 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
2926 {
2927   SNESKSPEW *kctx;
2928   PetscFunctionBegin;
2929   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2930   kctx = (SNESKSPEW*)snes->kspconvctx;
2931   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2932   if(version)   *version   = kctx->version;
2933   if(rtol_0)    *rtol_0    = kctx->rtol_0;
2934   if(rtol_max)  *rtol_max  = kctx->rtol_max;
2935   if(gamma)     *gamma     = kctx->gamma;
2936   if(alpha)     *alpha     = kctx->alpha;
2937   if(alpha2)    *alpha2    = kctx->alpha2;
2938   if(threshold) *threshold = kctx->threshold;
2939   PetscFunctionReturn(0);
2940 }
2941 
2942 #undef __FUNCT__
2943 #define __FUNCT__ "SNESKSPEW_PreSolve"
2944 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
2945 {
2946   PetscErrorCode ierr;
2947   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2948   PetscReal      rtol=PETSC_DEFAULT,stol;
2949 
2950   PetscFunctionBegin;
2951   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2952   if (!snes->iter) { /* first time in, so use the original user rtol */
2953     rtol = kctx->rtol_0;
2954   } else {
2955     if (kctx->version == 1) {
2956       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
2957       if (rtol < 0.0) rtol = -rtol;
2958       stol = pow(kctx->rtol_last,kctx->alpha2);
2959       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2960     } else if (kctx->version == 2) {
2961       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2962       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
2963       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2964     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
2965       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2966       /* safeguard: avoid sharp decrease of rtol */
2967       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
2968       stol = PetscMax(rtol,stol);
2969       rtol = PetscMin(kctx->rtol_0,stol);
2970       /* safeguard: avoid oversolving */
2971       stol = kctx->gamma*(snes->ttol)/snes->norm;
2972       stol = PetscMax(rtol,stol);
2973       rtol = PetscMin(kctx->rtol_0,stol);
2974     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
2975   }
2976   /* safeguard: avoid rtol greater than one */
2977   rtol = PetscMin(rtol,kctx->rtol_max);
2978   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2979   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
2980   PetscFunctionReturn(0);
2981 }
2982 
2983 #undef __FUNCT__
2984 #define __FUNCT__ "SNESKSPEW_PostSolve"
2985 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
2986 {
2987   PetscErrorCode ierr;
2988   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2989   PCSide         pcside;
2990   Vec            lres;
2991 
2992   PetscFunctionBegin;
2993   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2994   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
2995   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
2996   if (kctx->version == 1) {
2997     ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr);
2998     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
2999       /* KSP residual is true linear residual */
3000       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
3001     } else {
3002       /* KSP residual is preconditioned residual */
3003       /* compute true linear residual norm */
3004       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
3005       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
3006       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
3007       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
3008       ierr = VecDestroy(lres);CHKERRQ(ierr);
3009     }
3010   }
3011   PetscFunctionReturn(0);
3012 }
3013 
3014 #undef __FUNCT__
3015 #define __FUNCT__ "SNES_KSPSolve"
3016 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
3017 {
3018   PetscErrorCode ierr;
3019 
3020   PetscFunctionBegin;
3021   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
3022   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
3023   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
3024   PetscFunctionReturn(0);
3025 }
3026 
3027 #undef __FUNCT__
3028 #define __FUNCT__ "SNESSetDM"
3029 /*@
3030    SNESSetDM - Sets the DM that may be used by some preconditioners
3031 
3032    Logically Collective on SNES
3033 
3034    Input Parameters:
3035 +  snes - the preconditioner context
3036 -  dm - the dm
3037 
3038    Level: intermediate
3039 
3040 
3041 .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
3042 @*/
3043 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetDM(SNES snes,DM dm)
3044 {
3045   PetscErrorCode ierr;
3046   KSP            ksp;
3047 
3048   PetscFunctionBegin;
3049   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3050   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
3051   if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);}
3052   snes->dm = dm;
3053   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
3054   ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr);
3055   ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr);
3056   PetscFunctionReturn(0);
3057 }
3058 
3059 #undef __FUNCT__
3060 #define __FUNCT__ "SNESGetDM"
3061 /*@
3062    SNESGetDM - Gets the DM that may be used by some preconditioners
3063 
3064    Not Collective but DM obtained is parallel on SNES
3065 
3066    Input Parameter:
3067 . snes - the preconditioner context
3068 
3069    Output Parameter:
3070 .  dm - the dm
3071 
3072    Level: intermediate
3073 
3074 
3075 .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
3076 @*/
3077 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetDM(SNES snes,DM *dm)
3078 {
3079   PetscFunctionBegin;
3080   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3081   *dm = snes->dm;
3082   PetscFunctionReturn(0);
3083 }
3084 
3085 #if defined(PETSC_HAVE_MATLAB_ENGINE)
3086 #include "mex.h"
3087 
3088 typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
3089 
3090 #undef __FUNCT__
3091 #define __FUNCT__ "SNESComputeFunction_Matlab"
3092 /*
3093    SNESComputeFunction_Matlab - Calls the function that has been set with
3094                          SNESSetFunctionMatlab().
3095 
3096    Collective on SNES
3097 
3098    Input Parameters:
3099 +  snes - the SNES context
3100 -  x - input vector
3101 
3102    Output Parameter:
3103 .  y - function vector, as set by SNESSetFunction()
3104 
3105    Notes:
3106    SNESComputeFunction() is typically used within nonlinear solvers
3107    implementations, so most users would not generally call this routine
3108    themselves.
3109 
3110    Level: developer
3111 
3112 .keywords: SNES, nonlinear, compute, function
3113 
3114 .seealso: SNESSetFunction(), SNESGetFunction()
3115 */
3116 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
3117 {
3118   PetscErrorCode    ierr;
3119   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
3120   int               nlhs = 1,nrhs = 5;
3121   mxArray	    *plhs[1],*prhs[5];
3122   long long int     lx = 0,ly = 0,ls = 0;
3123 
3124   PetscFunctionBegin;
3125   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3126   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
3127   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
3128   PetscCheckSameComm(snes,1,x,2);
3129   PetscCheckSameComm(snes,1,y,3);
3130 
3131   /* call Matlab function in ctx with arguments x and y */
3132 
3133   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
3134   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3135   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
3136   prhs[0] =  mxCreateDoubleScalar((double)ls);
3137   prhs[1] =  mxCreateDoubleScalar((double)lx);
3138   prhs[2] =  mxCreateDoubleScalar((double)ly);
3139   prhs[3] =  mxCreateString(sctx->funcname);
3140   prhs[4] =  sctx->ctx;
3141   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");CHKERRQ(ierr);
3142   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3143   mxDestroyArray(prhs[0]);
3144   mxDestroyArray(prhs[1]);
3145   mxDestroyArray(prhs[2]);
3146   mxDestroyArray(prhs[3]);
3147   mxDestroyArray(plhs[0]);
3148   PetscFunctionReturn(0);
3149 }
3150 
3151 
3152 #undef __FUNCT__
3153 #define __FUNCT__ "SNESSetFunctionMatlab"
3154 /*
3155    SNESSetFunctionMatlab - Sets the function evaluation routine and function
3156    vector for use by the SNES routines in solving systems of nonlinear
3157    equations from Matlab. Here the function is a string containing the name of a Matlab function
3158 
3159    Logically Collective on SNES
3160 
3161    Input Parameters:
3162 +  snes - the SNES context
3163 .  r - vector to store function value
3164 -  func - function evaluation routine
3165 
3166    Calling sequence of func:
3167 $    func (SNES snes,Vec x,Vec f,void *ctx);
3168 
3169 
3170    Notes:
3171    The Newton-like methods typically solve linear systems of the form
3172 $      f'(x) x = -f(x),
3173    where f'(x) denotes the Jacobian matrix and f(x) is the function.
3174 
3175    Level: beginner
3176 
3177 .keywords: SNES, nonlinear, set, function
3178 
3179 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
3180 */
3181 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunctionMatlab(SNES snes,Vec r,const char *func,mxArray *ctx)
3182 {
3183   PetscErrorCode    ierr;
3184   SNESMatlabContext *sctx;
3185 
3186   PetscFunctionBegin;
3187   /* currently sctx is memory bleed */
3188   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
3189   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3190   /*
3191      This should work, but it doesn't
3192   sctx->ctx = ctx;
3193   mexMakeArrayPersistent(sctx->ctx);
3194   */
3195   sctx->ctx = mxDuplicateArray(ctx);
3196   ierr = SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);CHKERRQ(ierr);
3197   PetscFunctionReturn(0);
3198 }
3199 
3200 #undef __FUNCT__
3201 #define __FUNCT__ "SNESComputeJacobian_Matlab"
3202 /*
3203    SNESComputeJacobian_Matlab - Calls the function that has been set with
3204                          SNESSetJacobianMatlab().
3205 
3206    Collective on SNES
3207 
3208    Input Parameters:
3209 +  snes - the SNES context
3210 .  x - input vector
3211 .  A, B - the matrices
3212 -  ctx - user context
3213 
3214    Output Parameter:
3215 .  flag - structure of the matrix
3216 
3217    Level: developer
3218 
3219 .keywords: SNES, nonlinear, compute, function
3220 
3221 .seealso: SNESSetFunction(), SNESGetFunction()
3222 @*/
3223 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag, void *ctx)
3224 {
3225   PetscErrorCode    ierr;
3226   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
3227   int               nlhs = 2,nrhs = 6;
3228   mxArray	    *plhs[2],*prhs[6];
3229   long long int     lx = 0,lA = 0,ls = 0, lB = 0;
3230 
3231   PetscFunctionBegin;
3232   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3233   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
3234 
3235   /* call Matlab function in ctx with arguments x and y */
3236 
3237   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
3238   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3239   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
3240   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
3241   prhs[0] =  mxCreateDoubleScalar((double)ls);
3242   prhs[1] =  mxCreateDoubleScalar((double)lx);
3243   prhs[2] =  mxCreateDoubleScalar((double)lA);
3244   prhs[3] =  mxCreateDoubleScalar((double)lB);
3245   prhs[4] =  mxCreateString(sctx->funcname);
3246   prhs[5] =  sctx->ctx;
3247   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");CHKERRQ(ierr);
3248   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3249   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
3250   mxDestroyArray(prhs[0]);
3251   mxDestroyArray(prhs[1]);
3252   mxDestroyArray(prhs[2]);
3253   mxDestroyArray(prhs[3]);
3254   mxDestroyArray(prhs[4]);
3255   mxDestroyArray(plhs[0]);
3256   mxDestroyArray(plhs[1]);
3257   PetscFunctionReturn(0);
3258 }
3259 
3260 
3261 #undef __FUNCT__
3262 #define __FUNCT__ "SNESSetJacobianMatlab"
3263 /*
3264    SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
3265    vector for use by the SNES routines in solving systems of nonlinear
3266    equations from Matlab. Here the function is a string containing the name of a Matlab function
3267 
3268    Logically Collective on SNES
3269 
3270    Input Parameters:
3271 +  snes - the SNES context
3272 .  A,B - Jacobian matrices
3273 .  func - function evaluation routine
3274 -  ctx - user context
3275 
3276    Calling sequence of func:
3277 $    flag = func (SNES snes,Vec x,Mat A,Mat B,void *ctx);
3278 
3279 
3280    Level: developer
3281 
3282 .keywords: SNES, nonlinear, set, function
3283 
3284 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
3285 */
3286 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *func,mxArray *ctx)
3287 {
3288   PetscErrorCode    ierr;
3289   SNESMatlabContext *sctx;
3290 
3291   PetscFunctionBegin;
3292   /* currently sctx is memory bleed */
3293   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
3294   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3295   /*
3296      This should work, but it doesn't
3297   sctx->ctx = ctx;
3298   mexMakeArrayPersistent(sctx->ctx);
3299   */
3300   sctx->ctx = mxDuplicateArray(ctx);
3301   ierr = SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
3302   PetscFunctionReturn(0);
3303 }
3304 
3305 #undef __FUNCT__
3306 #define __FUNCT__ "SNESMonitor_Matlab"
3307 /*
3308    SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().
3309 
3310    Collective on SNES
3311 
3312 .seealso: SNESSetFunction(), SNESGetFunction()
3313 @*/
3314 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
3315 {
3316   PetscErrorCode  ierr;
3317   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
3318   int             nlhs = 1,nrhs = 6;
3319   mxArray	  *plhs[1],*prhs[6];
3320   long long int   lx = 0,ls = 0;
3321   Vec             x=snes->vec_sol;
3322 
3323   PetscFunctionBegin;
3324   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3325 
3326   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
3327   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3328   prhs[0] =  mxCreateDoubleScalar((double)ls);
3329   prhs[1] =  mxCreateDoubleScalar((double)it);
3330   prhs[2] =  mxCreateDoubleScalar((double)fnorm);
3331   prhs[3] =  mxCreateDoubleScalar((double)lx);
3332   prhs[4] =  mxCreateString(sctx->funcname);
3333   prhs[5] =  sctx->ctx;
3334   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");CHKERRQ(ierr);
3335   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3336   mxDestroyArray(prhs[0]);
3337   mxDestroyArray(prhs[1]);
3338   mxDestroyArray(prhs[2]);
3339   mxDestroyArray(prhs[3]);
3340   mxDestroyArray(prhs[4]);
3341   mxDestroyArray(plhs[0]);
3342   PetscFunctionReturn(0);
3343 }
3344 
3345 
3346 #undef __FUNCT__
3347 #define __FUNCT__ "SNESMonitorSetMatlab"
3348 /*
3349    SNESMonitorSetMatlab - Sets the monitor function from Matlab
3350 
3351    Level: developer
3352 
3353 .keywords: SNES, nonlinear, set, function
3354 
3355 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
3356 */
3357 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSetMatlab(SNES snes,const char *func,mxArray *ctx)
3358 {
3359   PetscErrorCode    ierr;
3360   SNESMatlabContext *sctx;
3361 
3362   PetscFunctionBegin;
3363   /* currently sctx is memory bleed */
3364   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
3365   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3366   /*
3367      This should work, but it doesn't
3368   sctx->ctx = ctx;
3369   mexMakeArrayPersistent(sctx->ctx);
3370   */
3371   sctx->ctx = mxDuplicateArray(ctx);
3372   ierr = SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
3373   PetscFunctionReturn(0);
3374 }
3375 
3376 #endif
3377