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