xref: /petsc/src/snes/interface/snes.c (revision bdad233f16b1da2fc2320f4673f853efec1db8b6)
1 /*$Id: snes.c,v 1.235 2001/08/21 21:03:49 bsmith Exp $*/
2 
3 #include "src/snes/snesimpl.h"      /*I "petscsnes.h"  I*/
4 
5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
6 PetscFList SNESList              = PETSC_NULL;
7 
8 /* Logging support */
9 int SNES_COOKIE;
10 int MATSNESMFCTX_COOKIE;
11 int SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval, SNES_MinimizationFunctionEval, SNES_GradientEval;
12 int SNES_HessianEval;
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "SNESGetProblemType"
16 /*@C
17    SNESGetProblemType -Indicates if SNES is solving a nonlinear system or a minimization
18 
19    Not Collective
20 
21    Input Parameter:
22 .  SNES - the SNES context
23 
24    Output Parameter:
25 .   type - SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations)
26    or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization)
27 
28    Level: intermediate
29 
30 .keywords: SNES, problem type
31 
32 .seealso: SNESCreate()
33 @*/
34 int SNESGetProblemType(SNES snes,SNESProblemType *type)
35 {
36   PetscFunctionBegin;
37   PetscValidHeaderSpecific(snes,SNES_COOKIE);
38   *type = snes->method_class;
39   PetscFunctionReturn(0);
40 }
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "SNESView"
44 /*@C
45    SNESView - Prints the SNES data structure.
46 
47    Collective on SNES
48 
49    Input Parameters:
50 +  SNES - the SNES context
51 -  viewer - visualization context
52 
53    Options Database Key:
54 .  -snes_view - Calls SNESView() at end of SNESSolve()
55 
56    Notes:
57    The available visualization contexts include
58 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
59 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
60          output where only the first processor opens
61          the file.  All other processors send their
62          data to the first processor to print.
63 
64    The user can open an alternative visualization context with
65    PetscViewerASCIIOpen() - output to a specified file.
66 
67    Level: beginner
68 
69 .keywords: SNES, view
70 
71 .seealso: PetscViewerASCIIOpen()
72 @*/
73 int SNESView(SNES snes,PetscViewer viewer)
74 {
75   SNES_KSP_EW_ConvCtx *kctx;
76   int                 ierr;
77   SLES                sles;
78   char                *type;
79   PetscTruth          isascii,isstring;
80 
81   PetscFunctionBegin;
82   PetscValidHeaderSpecific(snes,SNES_COOKIE);
83   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm);
84   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
85   PetscCheckSameComm(snes,viewer);
86 
87   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
88   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
89   if (isascii) {
90     if (snes->prefix) {
91       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr);
92     } else {
93       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
94     }
95     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
96     if (type) {
97       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
98     } else {
99       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
100     }
101     if (snes->view) {
102       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
103       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
104       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
105     }
106     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
107     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, solution=%g\n",
108                  snes->rtol,snes->atol,snes->xtol);CHKERRQ(ierr);
109     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr);
110     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr);
111     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
112       ierr = PetscViewerASCIIPrintf(viewer,"  min function tolerance=%g\n",snes->fmin);CHKERRQ(ierr);
113     }
114     if (snes->ksp_ewconv) {
115       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
116       if (kctx) {
117         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr);
118         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
119         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
120       }
121     }
122   } else if (isstring) {
123     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
124     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
125   }
126   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
127   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
128   ierr = SLESView(sles,viewer);CHKERRQ(ierr);
129   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
130   PetscFunctionReturn(0);
131 }
132 
133 /*
134   We retain a list of functions that also take SNES command
135   line options. These are called at the end SNESSetFromOptions()
136 */
137 #define MAXSETFROMOPTIONS 5
138 static int numberofsetfromoptions;
139 static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
140 
141 #undef __FUNCT__
142 #define __FUNCT__ "SNESAddOptionsChecker"
143 /*@
144   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
145 
146   Not Collective
147 
148   Input Parameter:
149 . snescheck - function that checks for options
150 
151   Level: developer
152 
153 .seealso: SNESSetFromOptions()
154 @*/
155 int SNESAddOptionsChecker(int (*snescheck)(SNES))
156 {
157   PetscFunctionBegin;
158   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
159     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %d allowed", MAXSETFROMOPTIONS);
160   }
161 
162   othersetfromoptions[numberofsetfromoptions++] = snescheck;
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "SNESSetFromOptions"
168 /*@
169    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
170 
171    Collective on SNES
172 
173    Input Parameter:
174 .  snes - the SNES context
175 
176    Options Database Keys:
177 +  -snes_type <type> - ls, tr, umls, umtr, test
178 .  -snes_stol - convergence tolerance in terms of the norm
179                 of the change in the solution between steps
180 .  -snes_atol <atol> - absolute tolerance of residual norm
181 .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
182 .  -snes_max_it <max_it> - maximum number of iterations
183 .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
184 .  -snes_max_fail <max_fail> - maximum number of failures
185 .  -snes_trtol <trtol> - trust region tolerance
186 .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
187                                solver; hence iterations will continue until max_it
188                                or some other criterion is reached. Saves expense
189                                of convergence test
190 .  -snes_monitor - prints residual norm at each iteration
191 .  -snes_vecmonitor - plots solution at each iteration
192 .  -snes_vecmonitor_update - plots update to solution at each iteration
193 .  -snes_xmonitor - plots residual norm at each iteration
194 .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
195 -  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
196 
197     Options Database for Eisenstat-Walker method:
198 +  -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence
199 .  -snes_ksp_eq_version ver - version of  Eisenstat-Walker method
200 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
201 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
202 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
203 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
204 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
205 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
206 
207    Notes:
208    To see all options, run your program with the -help option or consult
209    the users manual.
210 
211    Level: beginner
212 
213 .keywords: SNES, nonlinear, set, options, database
214 
215 .seealso: SNESSetOptionsPrefix()
216 @*/
217 int SNESSetFromOptions(SNES snes)
218 {
219   SLES                sles;
220   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
221   PetscTruth          flg;
222   int                 ierr, i;
223   char                *deft,type[256];
224 
225   PetscFunctionBegin;
226   PetscValidHeaderSpecific(snes,SNES_COOKIE);
227 
228   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
229     if (snes->type_name) {
230       deft = snes->type_name;
231     } else {
232       if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
233         deft = SNESEQLS;
234       } else {
235         deft = SNESUMTR;
236       }
237     }
238 
239     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
240     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
241     if (flg) {
242       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
243     } else if (!snes->type_name) {
244       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
245     }
246 
247     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
248     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr);
249 
250     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
251     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
252     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
253     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
254     ierr = PetscOptionsReal("-snes_fmin","Minimization function tolerance","SNESSetMinimizationFunctionTolerance",snes->fmin,&snes->fmin,0);CHKERRQ(ierr);
255 
256     ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr);
257 
258     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
259     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
260     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
261     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
262     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
263     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
264     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
265 
266     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
267     if (flg) {snes->converged = 0;}
268     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
269     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
270     ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
271     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
272     ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr);
273     if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);}
274     ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
275     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
276     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
277     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
278     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
279     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
280     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
281     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
282 
283     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
284     if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
285       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
286       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
287     } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
288       ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeHessian,snes->funP);CHKERRQ(ierr);
289       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
290     }
291 
292     for(i = 0; i < numberofsetfromoptions; i++) {
293       ierr = (*othersetfromoptions[i])(snes);                                                             CHKERRQ(ierr);
294     }
295 
296     if (snes->setfromoptions) {
297       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
298     }
299 
300   ierr = PetscOptionsEnd();CHKERRQ(ierr);
301 
302   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
303   ierr = SLESSetFromOptions(sles);CHKERRQ(ierr);
304 
305   PetscFunctionReturn(0);
306 }
307 
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "SNESSetApplicationContext"
311 /*@
312    SNESSetApplicationContext - Sets the optional user-defined context for
313    the nonlinear solvers.
314 
315    Collective on SNES
316 
317    Input Parameters:
318 +  snes - the SNES context
319 -  usrP - optional user context
320 
321    Level: intermediate
322 
323 .keywords: SNES, nonlinear, set, application, context
324 
325 .seealso: SNESGetApplicationContext()
326 @*/
327 int SNESSetApplicationContext(SNES snes,void *usrP)
328 {
329   PetscFunctionBegin;
330   PetscValidHeaderSpecific(snes,SNES_COOKIE);
331   snes->user		= usrP;
332   PetscFunctionReturn(0);
333 }
334 
335 #undef __FUNCT__
336 #define __FUNCT__ "SNESGetApplicationContext"
337 /*@C
338    SNESGetApplicationContext - Gets the user-defined context for the
339    nonlinear solvers.
340 
341    Not Collective
342 
343    Input Parameter:
344 .  snes - SNES context
345 
346    Output Parameter:
347 .  usrP - user context
348 
349    Level: intermediate
350 
351 .keywords: SNES, nonlinear, get, application, context
352 
353 .seealso: SNESSetApplicationContext()
354 @*/
355 int SNESGetApplicationContext(SNES snes,void **usrP)
356 {
357   PetscFunctionBegin;
358   PetscValidHeaderSpecific(snes,SNES_COOKIE);
359   *usrP = snes->user;
360   PetscFunctionReturn(0);
361 }
362 
363 #undef __FUNCT__
364 #define __FUNCT__ "SNESGetIterationNumber"
365 /*@
366    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
367    at this time.
368 
369    Not Collective
370 
371    Input Parameter:
372 .  snes - SNES context
373 
374    Output Parameter:
375 .  iter - iteration number
376 
377    Notes:
378    For example, during the computation of iteration 2 this would return 1.
379 
380    This is useful for using lagged Jacobians (where one does not recompute the
381    Jacobian at each SNES iteration). For example, the code
382 .vb
383       ierr = SNESGetIterationNumber(snes,&it);
384       if (!(it % 2)) {
385         [compute Jacobian here]
386       }
387 .ve
388    can be used in your ComputeJacobian() function to cause the Jacobian to be
389    recomputed every second SNES iteration.
390 
391    Level: intermediate
392 
393 .keywords: SNES, nonlinear, get, iteration, number
394 @*/
395 int SNESGetIterationNumber(SNES snes,int* iter)
396 {
397   PetscFunctionBegin;
398   PetscValidHeaderSpecific(snes,SNES_COOKIE);
399   PetscValidIntPointer(iter);
400   *iter = snes->iter;
401   PetscFunctionReturn(0);
402 }
403 
404 #undef __FUNCT__
405 #define __FUNCT__ "SNESGetFunctionNorm"
406 /*@
407    SNESGetFunctionNorm - Gets the norm of the current function that was set
408    with SNESSSetFunction().
409 
410    Collective on SNES
411 
412    Input Parameter:
413 .  snes - SNES context
414 
415    Output Parameter:
416 .  fnorm - 2-norm of function
417 
418    Note:
419    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
420    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
421    SNESGetGradientNorm().
422 
423    Level: intermediate
424 
425 .keywords: SNES, nonlinear, get, function, norm
426 
427 .seealso: SNESGetFunction()
428 @*/
429 int SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
430 {
431   PetscFunctionBegin;
432   PetscValidHeaderSpecific(snes,SNES_COOKIE);
433   PetscValidScalarPointer(fnorm);
434   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
435     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_NONLINEAR_EQUATIONS only");
436   }
437   *fnorm = snes->norm;
438   PetscFunctionReturn(0);
439 }
440 
441 #undef __FUNCT__
442 #define __FUNCT__ "SNESGetGradientNorm"
443 /*@
444    SNESGetGradientNorm - Gets the norm of the current gradient that was set
445    with SNESSSetGradient().
446 
447    Collective on SNES
448 
449    Input Parameter:
450 .  snes - SNES context
451 
452    Output Parameter:
453 .  fnorm - 2-norm of gradient
454 
455    Note:
456    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
457    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
458    is SNESGetFunctionNorm().
459 
460    Level: intermediate
461 
462 .keywords: SNES, nonlinear, get, gradient, norm
463 
464 .seelso: SNESSetGradient()
465 @*/
466 int SNESGetGradientNorm(SNES snes,PetscScalar *gnorm)
467 {
468   PetscFunctionBegin;
469   PetscValidHeaderSpecific(snes,SNES_COOKIE);
470   PetscValidScalarPointer(gnorm);
471   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
472     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_UNCONSTRAINED_MINIMIZATION only");
473   }
474   *gnorm = snes->norm;
475   PetscFunctionReturn(0);
476 }
477 
478 #undef __FUNCT__
479 #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
480 /*@
481    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
482    attempted by the nonlinear solver.
483 
484    Not Collective
485 
486    Input Parameter:
487 .  snes - SNES context
488 
489    Output Parameter:
490 .  nfails - number of unsuccessful steps attempted
491 
492    Notes:
493    This counter is reset to zero for each successive call to SNESSolve().
494 
495    Level: intermediate
496 
497 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
498 @*/
499 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
500 {
501   PetscFunctionBegin;
502   PetscValidHeaderSpecific(snes,SNES_COOKIE);
503   PetscValidIntPointer(nfails);
504   *nfails = snes->numFailures;
505   PetscFunctionReturn(0);
506 }
507 
508 #undef __FUNCT__
509 #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps"
510 /*@
511    SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps
512    attempted by the nonlinear solver before it gives up.
513 
514    Not Collective
515 
516    Input Parameters:
517 +  snes     - SNES context
518 -  maxFails - maximum of unsuccessful steps
519 
520    Level: intermediate
521 
522 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
523 @*/
524 int SNESSetMaximumUnsuccessfulSteps(SNES snes, int maxFails)
525 {
526   PetscFunctionBegin;
527   PetscValidHeaderSpecific(snes,SNES_COOKIE);
528   snes->maxFailures = maxFails;
529   PetscFunctionReturn(0);
530 }
531 
532 #undef __FUNCT__
533 #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps"
534 /*@
535    SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps
536    attempted by the nonlinear solver before it gives up.
537 
538    Not Collective
539 
540    Input Parameter:
541 .  snes     - SNES context
542 
543    Output Parameter:
544 .  maxFails - maximum of unsuccessful steps
545 
546    Level: intermediate
547 
548 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
549 @*/
550 int SNESGetMaximumUnsuccessfulSteps(SNES snes, int *maxFails)
551 {
552   PetscFunctionBegin;
553   PetscValidHeaderSpecific(snes,SNES_COOKIE);
554   PetscValidIntPointer(maxFails);
555   *maxFails = snes->maxFailures;
556   PetscFunctionReturn(0);
557 }
558 
559 #undef __FUNCT__
560 #define __FUNCT__ "SNESGetNumberLinearIterations"
561 /*@
562    SNESGetNumberLinearIterations - Gets the total number of linear iterations
563    used by the nonlinear solver.
564 
565    Not Collective
566 
567    Input Parameter:
568 .  snes - SNES context
569 
570    Output Parameter:
571 .  lits - number of linear iterations
572 
573    Notes:
574    This counter is reset to zero for each successive call to SNESSolve().
575 
576    Level: intermediate
577 
578 .keywords: SNES, nonlinear, get, number, linear, iterations
579 @*/
580 int SNESGetNumberLinearIterations(SNES snes,int* lits)
581 {
582   PetscFunctionBegin;
583   PetscValidHeaderSpecific(snes,SNES_COOKIE);
584   PetscValidIntPointer(lits);
585   *lits = snes->linear_its;
586   PetscFunctionReturn(0);
587 }
588 
589 #undef __FUNCT__
590 #define __FUNCT__ "SNESGetSLES"
591 /*@C
592    SNESGetSLES - Returns the SLES context for a SNES solver.
593 
594    Not Collective, but if SNES object is parallel, then SLES object is parallel
595 
596    Input Parameter:
597 .  snes - the SNES context
598 
599    Output Parameter:
600 .  sles - the SLES context
601 
602    Notes:
603    The user can then directly manipulate the SLES context to set various
604    options, etc.  Likewise, the user can then extract and manipulate the
605    KSP and PC contexts as well.
606 
607    Level: beginner
608 
609 .keywords: SNES, nonlinear, get, SLES, context
610 
611 .seealso: SLESGetPC(), SLESGetKSP()
612 @*/
613 int SNESGetSLES(SNES snes,SLES *sles)
614 {
615   PetscFunctionBegin;
616   PetscValidHeaderSpecific(snes,SNES_COOKIE);
617   *sles = snes->sles;
618   PetscFunctionReturn(0);
619 }
620 
621 #undef __FUNCT__
622 #define __FUNCT__ "SNESPublish_Petsc"
623 static int SNESPublish_Petsc(PetscObject obj)
624 {
625 #if defined(PETSC_HAVE_AMS)
626   SNES          v = (SNES) obj;
627   int          ierr;
628 #endif
629 
630   PetscFunctionBegin;
631 
632 #if defined(PETSC_HAVE_AMS)
633   /* if it is already published then return */
634   if (v->amem >=0) PetscFunctionReturn(0);
635 
636   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
637   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
638                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
639   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
640                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
641   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
642 #endif
643   PetscFunctionReturn(0);
644 }
645 
646 /* -----------------------------------------------------------*/
647 #undef __FUNCT__
648 #define __FUNCT__ "SNESCreate"
649 /*@C
650    SNESCreate - Creates a nonlinear solver context.
651 
652    Collective on MPI_Comm
653 
654    Input Parameters:
655 +  comm - MPI communicator
656 -  type - type of method, either
657    SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations)
658    or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization)
659 
660    Output Parameter:
661 .  outsnes - the new SNES context
662 
663    Options Database Keys:
664 +   -snes_mf - Activates default matrix-free Jacobian-vector products,
665                and no preconditioning matrix
666 .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
667                products, and a user-provided preconditioning matrix
668                as set by SNESSetJacobian()
669 -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
670 
671    Level: beginner
672 
673 .keywords: SNES, nonlinear, create, context
674 
675 .seealso: SNESSolve(), SNESDestroy(), SNESProblemType, SNES
676 @*/
677 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
678 {
679   int                 ierr;
680   SNES                snes;
681   SNES_KSP_EW_ConvCtx *kctx;
682 
683   PetscFunctionBegin;
684   PetscValidPointer(outsnes);
685   *outsnes = PETSC_NULL;
686 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
687   ierr = SNESInitializePackage(PETSC_NULL);                                                               CHKERRQ(ierr);
688 #endif
689 
690   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
691     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"incorrect method type");
692   }
693   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);
694   PetscLogObjectCreate(snes);
695   snes->bops->publish     = SNESPublish_Petsc;
696   snes->max_its           = 50;
697   snes->max_funcs	  = 10000;
698   snes->norm		  = 0.0;
699   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
700     snes->rtol		  = 1.e-8;
701     snes->ttol            = 0.0;
702     snes->atol		  = 1.e-10;
703   } else {
704     snes->rtol		  = 1.e-8;
705     snes->ttol            = 0.0;
706     snes->atol		  = 1.e-50;
707   }
708   snes->xtol		  = 1.e-8;
709   snes->nfuncs            = 0;
710   snes->numFailures       = 0;
711   snes->maxFailures       = 1;
712   snes->linear_its        = 0;
713   snes->numbermonitors    = 0;
714   snes->data              = 0;
715   snes->view              = 0;
716   snes->computeumfunction = 0;
717   snes->umfunP            = 0;
718   snes->fc                = 0;
719   snes->deltatol          = 1.e-12;
720   snes->fmin              = -1.e30;
721   snes->method_class      = type;
722   snes->set_method_called = 0;
723   snes->setupcalled       = 0;
724   snes->ksp_ewconv        = PETSC_FALSE;
725   snes->vwork             = 0;
726   snes->nwork             = 0;
727   snes->conv_hist_len     = 0;
728   snes->conv_hist_max     = 0;
729   snes->conv_hist         = PETSC_NULL;
730   snes->conv_hist_its     = PETSC_NULL;
731   snes->conv_hist_reset   = PETSC_TRUE;
732   snes->reason            = SNES_CONVERGED_ITERATING;
733 
734   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
735   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
736   PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
737   snes->kspconvctx  = (void*)kctx;
738   kctx->version     = 2;
739   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
740                              this was too large for some test cases */
741   kctx->rtol_last   = 0;
742   kctx->rtol_max    = .9;
743   kctx->gamma       = 1.0;
744   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
745   kctx->alpha       = kctx->alpha2;
746   kctx->threshold   = .1;
747   kctx->lresid_last = 0;
748   kctx->norm_last   = 0;
749 
750   ierr = SLESCreate(comm,&snes->sles);CHKERRQ(ierr);
751   PetscLogObjectParent(snes,snes->sles)
752 
753   *outsnes = snes;
754   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
755   PetscFunctionReturn(0);
756 }
757 
758 /* --------------------------------------------------------------- */
759 #undef __FUNCT__
760 #define __FUNCT__ "SNESSetFunction"
761 /*@C
762    SNESSetFunction - Sets the function evaluation routine and function
763    vector for use by the SNES routines in solving systems of nonlinear
764    equations.
765 
766    Collective on SNES
767 
768    Input Parameters:
769 +  snes - the SNES context
770 .  func - function evaluation routine
771 .  r - vector to store function value
772 -  ctx - [optional] user-defined context for private data for the
773          function evaluation routine (may be PETSC_NULL)
774 
775    Calling sequence of func:
776 $    func (SNES snes,Vec x,Vec f,void *ctx);
777 
778 .  f - function vector
779 -  ctx - optional user-defined function context
780 
781    Notes:
782    The Newton-like methods typically solve linear systems of the form
783 $      f'(x) x = -f(x),
784    where f'(x) denotes the Jacobian matrix and f(x) is the function.
785 
786    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
787    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
788    SNESSetMinimizationFunction() and SNESSetGradient();
789 
790    Level: beginner
791 
792 .keywords: SNES, nonlinear, set, function
793 
794 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
795 @*/
796 int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
797 {
798   PetscFunctionBegin;
799   PetscValidHeaderSpecific(snes,SNES_COOKIE);
800   PetscValidHeaderSpecific(r,VEC_COOKIE);
801   PetscCheckSameComm(snes,r);
802   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
803     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
804   }
805 
806   snes->computefunction     = func;
807   snes->vec_func            = snes->vec_func_always = r;
808   snes->funP                = ctx;
809   PetscFunctionReturn(0);
810 }
811 
812 #undef __FUNCT__
813 #define __FUNCT__ "SNESComputeFunction"
814 /*@
815    SNESComputeFunction - Calls the function that has been set with
816                          SNESSetFunction().
817 
818    Collective on SNES
819 
820    Input Parameters:
821 +  snes - the SNES context
822 -  x - input vector
823 
824    Output Parameter:
825 .  y - function vector, as set by SNESSetFunction()
826 
827    Notes:
828    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
829    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
830    SNESComputeMinimizationFunction() and SNESComputeGradient();
831 
832    SNESComputeFunction() is typically used within nonlinear solvers
833    implementations, so most users would not generally call this routine
834    themselves.
835 
836    Level: developer
837 
838 .keywords: SNES, nonlinear, compute, function
839 
840 .seealso: SNESSetFunction(), SNESGetFunction()
841 @*/
842 int SNESComputeFunction(SNES snes,Vec x,Vec y)
843 {
844   int    ierr;
845 
846   PetscFunctionBegin;
847   PetscValidHeaderSpecific(snes,SNES_COOKIE);
848   PetscValidHeaderSpecific(x,VEC_COOKIE);
849   PetscValidHeaderSpecific(y,VEC_COOKIE);
850   PetscCheckSameComm(snes,x);
851   PetscCheckSameComm(snes,y);
852   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
853     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
854   }
855 
856   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
857   PetscStackPush("SNES user function");
858   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
859   PetscStackPop;
860   snes->nfuncs++;
861   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
862   PetscFunctionReturn(0);
863 }
864 
865 #undef __FUNCT__
866 #define __FUNCT__ "SNESSetMinimizationFunction"
867 /*@C
868    SNESSetMinimizationFunction - Sets the function evaluation routine for
869    unconstrained minimization.
870 
871    Collective on SNES
872 
873    Input Parameters:
874 +  snes - the SNES context
875 .  func - function evaluation routine
876 -  ctx - [optional] user-defined context for private data for the
877          function evaluation routine (may be PETSC_NULL)
878 
879    Calling sequence of func:
880 $     func (SNES snes,Vec x,PetscReal *f,void *ctx);
881 
882 +  x - input vector
883 .  f - function
884 -  ctx - [optional] user-defined function context
885 
886    Level: beginner
887 
888    Notes:
889    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
890    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
891    SNESSetFunction().
892 
893 .keywords: SNES, nonlinear, set, minimization, function
894 
895 .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
896            SNESSetHessian(), SNESSetGradient()
897 @*/
898 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,PetscReal*,void*),void *ctx)
899 {
900   PetscFunctionBegin;
901   PetscValidHeaderSpecific(snes,SNES_COOKIE);
902   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
903     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
904   }
905   snes->computeumfunction   = func;
906   snes->umfunP              = ctx;
907   PetscFunctionReturn(0);
908 }
909 
910 #undef __FUNCT__
911 #define __FUNCT__ "SNESComputeMinimizationFunction"
912 /*@
913    SNESComputeMinimizationFunction - Computes the function that has been
914    set with SNESSetMinimizationFunction().
915 
916    Collective on SNES
917 
918    Input Parameters:
919 +  snes - the SNES context
920 -  x - input vector
921 
922    Output Parameter:
923 .  y - function value
924 
925    Notes:
926    SNESComputeMinimizationFunction() is valid only for
927    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
928    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
929 
930    SNESComputeMinimizationFunction() is typically used within minimization
931    implementations, so most users would not generally call this routine
932    themselves.
933 
934    Level: developer
935 
936 .keywords: SNES, nonlinear, compute, minimization, function
937 
938 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
939           SNESComputeGradient(), SNESComputeHessian()
940 @*/
941 int SNESComputeMinimizationFunction(SNES snes,Vec x,PetscReal *y)
942 {
943   int    ierr;
944 
945   PetscFunctionBegin;
946   PetscValidHeaderSpecific(snes,SNES_COOKIE);
947   PetscValidHeaderSpecific(x,VEC_COOKIE);
948   PetscCheckSameComm(snes,x);
949   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
950     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
951   }
952 
953   ierr = PetscLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr);
954   PetscStackPush("SNES user minimzation function");
955   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP);CHKERRQ(ierr);
956   PetscStackPop;
957   snes->nfuncs++;
958   ierr = PetscLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr);
959   PetscFunctionReturn(0);
960 }
961 
962 #undef __FUNCT__
963 #define __FUNCT__ "SNESSetGradient"
964 /*@C
965    SNESSetGradient - Sets the gradient evaluation routine and gradient
966    vector for use by the SNES routines.
967 
968    Collective on SNES
969 
970    Input Parameters:
971 +  snes - the SNES context
972 .  func - function evaluation routine
973 .  ctx - optional user-defined context for private data for the
974          gradient evaluation routine (may be PETSC_NULL)
975 -  r - vector to store gradient value
976 
977    Calling sequence of func:
978 $     func (SNES, Vec x, Vec g, void *ctx);
979 
980 +  x - input vector
981 .  g - gradient vector
982 -  ctx - optional user-defined gradient context
983 
984    Notes:
985    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
986    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
987    SNESSetFunction().
988 
989    Level: beginner
990 
991 .keywords: SNES, nonlinear, set, function
992 
993 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
994           SNESSetMinimizationFunction(),
995 @*/
996 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
997 {
998   PetscFunctionBegin;
999   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1000   PetscValidHeaderSpecific(r,VEC_COOKIE);
1001   PetscCheckSameComm(snes,r);
1002   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1003     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1004   }
1005   snes->computefunction     = func;
1006   snes->vec_func            = snes->vec_func_always = r;
1007   snes->funP                = ctx;
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 #undef __FUNCT__
1012 #define __FUNCT__ "SNESComputeGradient"
1013 /*@
1014    SNESComputeGradient - Computes the gradient that has been set with
1015    SNESSetGradient().
1016 
1017    Collective on SNES
1018 
1019    Input Parameters:
1020 +  snes - the SNES context
1021 -  x - input vector
1022 
1023    Output Parameter:
1024 .  y - gradient vector
1025 
1026    Notes:
1027    SNESComputeGradient() is valid only for
1028    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
1029    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
1030 
1031    SNESComputeGradient() is typically used within minimization
1032    implementations, so most users would not generally call this routine
1033    themselves.
1034 
1035    Level: developer
1036 
1037 .keywords: SNES, nonlinear, compute, gradient
1038 
1039 .seealso:  SNESSetGradient(), SNESGetGradient(),
1040            SNESComputeMinimizationFunction(), SNESComputeHessian()
1041 @*/
1042 int SNESComputeGradient(SNES snes,Vec x,Vec y)
1043 {
1044   int    ierr;
1045 
1046   PetscFunctionBegin;
1047   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1048   PetscValidHeaderSpecific(x,VEC_COOKIE);
1049   PetscValidHeaderSpecific(y,VEC_COOKIE);
1050   PetscCheckSameComm(snes,x);
1051   PetscCheckSameComm(snes,y);
1052   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1053     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1054   }
1055   ierr = PetscLogEventBegin(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr);
1056   PetscStackPush("SNES user gradient function");
1057   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
1058   PetscStackPop;
1059   ierr = PetscLogEventEnd(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr);
1060   PetscFunctionReturn(0);
1061 }
1062 
1063 #undef __FUNCT__
1064 #define __FUNCT__ "SNESComputeJacobian"
1065 /*@
1066    SNESComputeJacobian - Computes the Jacobian matrix that has been
1067    set with SNESSetJacobian().
1068 
1069    Collective on SNES and Mat
1070 
1071    Input Parameters:
1072 +  snes - the SNES context
1073 -  x - input vector
1074 
1075    Output Parameters:
1076 +  A - Jacobian matrix
1077 .  B - optional preconditioning matrix
1078 -  flag - flag indicating matrix structure
1079 
1080    Notes:
1081    Most users should not need to explicitly call this routine, as it
1082    is used internally within the nonlinear solvers.
1083 
1084    See SLESSetOperators() for important information about setting the
1085    flag parameter.
1086 
1087    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
1088    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
1089    methods is SNESComputeHessian().
1090 
1091    Level: developer
1092 
1093 .keywords: SNES, compute, Jacobian, matrix
1094 
1095 .seealso:  SNESSetJacobian(), SLESSetOperators()
1096 @*/
1097 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
1098 {
1099   int    ierr;
1100 
1101   PetscFunctionBegin;
1102   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1103   PetscValidHeaderSpecific(X,VEC_COOKIE);
1104   PetscCheckSameComm(snes,X);
1105   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1106     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1107   }
1108   if (!snes->computejacobian) PetscFunctionReturn(0);
1109   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1110   *flg = DIFFERENT_NONZERO_PATTERN;
1111   PetscStackPush("SNES user Jacobian function");
1112   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1113   PetscStackPop;
1114   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1115   /* make sure user returned a correct Jacobian and preconditioner */
1116   PetscValidHeaderSpecific(*A,MAT_COOKIE);
1117   PetscValidHeaderSpecific(*B,MAT_COOKIE);
1118   PetscFunctionReturn(0);
1119 }
1120 
1121 #undef __FUNCT__
1122 #define __FUNCT__ "SNESComputeHessian"
1123 /*@
1124    SNESComputeHessian - Computes the Hessian matrix that has been
1125    set with SNESSetHessian().
1126 
1127    Collective on SNES and Mat
1128 
1129    Input Parameters:
1130 +  snes - the SNES context
1131 -  x - input vector
1132 
1133    Output Parameters:
1134 +  A - Hessian matrix
1135 .  B - optional preconditioning matrix
1136 -  flag - flag indicating matrix structure
1137 
1138    Notes:
1139    Most users should not need to explicitly call this routine, as it
1140    is used internally within the nonlinear solvers.
1141 
1142    See SLESSetOperators() for important information about setting the
1143    flag parameter.
1144 
1145    SNESComputeHessian() is valid only for
1146    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
1147    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
1148 
1149    SNESComputeHessian() is typically used within minimization
1150    implementations, so most users would not generally call this routine
1151    themselves.
1152 
1153    Level: developer
1154 
1155 .keywords: SNES, compute, Hessian, matrix
1156 
1157 .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
1158            SNESComputeMinimizationFunction()
1159 @*/
1160 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
1161 {
1162   int    ierr;
1163 
1164   PetscFunctionBegin;
1165   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1166   PetscValidHeaderSpecific(x,VEC_COOKIE);
1167   PetscCheckSameComm(snes,x);
1168   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1169     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1170   }
1171   if (!snes->computejacobian) PetscFunctionReturn(0);
1172   ierr = PetscLogEventBegin(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr);
1173   *flag = DIFFERENT_NONZERO_PATTERN;
1174   PetscStackPush("SNES user Hessian function");
1175   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP);CHKERRQ(ierr);
1176   PetscStackPop;
1177   ierr = PetscLogEventEnd(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr);
1178   /* make sure user returned a correct Jacobian and preconditioner */
1179   PetscValidHeaderSpecific(*A,MAT_COOKIE);
1180   PetscValidHeaderSpecific(*B,MAT_COOKIE);
1181   PetscFunctionReturn(0);
1182 }
1183 
1184 #undef __FUNCT__
1185 #define __FUNCT__ "SNESSetJacobian"
1186 /*@C
1187    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1188    location to store the matrix.
1189 
1190    Collective on SNES and Mat
1191 
1192    Input Parameters:
1193 +  snes - the SNES context
1194 .  A - Jacobian matrix
1195 .  B - preconditioner matrix (usually same as the Jacobian)
1196 .  func - Jacobian evaluation routine
1197 -  ctx - [optional] user-defined context for private data for the
1198          Jacobian evaluation routine (may be PETSC_NULL)
1199 
1200    Calling sequence of func:
1201 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1202 
1203 +  x - input vector
1204 .  A - Jacobian matrix
1205 .  B - preconditioner matrix, usually the same as A
1206 .  flag - flag indicating information about the preconditioner matrix
1207    structure (same as flag in SLESSetOperators())
1208 -  ctx - [optional] user-defined Jacobian context
1209 
1210    Notes:
1211    See SLESSetOperators() for important information about setting the flag
1212    output parameter in the routine func().  Be sure to read this information!
1213 
1214    The routine func() takes Mat * as the matrix arguments rather than Mat.
1215    This allows the Jacobian evaluation routine to replace A and/or B with a
1216    completely new new matrix structure (not just different matrix elements)
1217    when appropriate, for instance, if the nonzero structure is changing
1218    throughout the global iterations.
1219 
1220    Level: beginner
1221 
1222 .keywords: SNES, nonlinear, set, Jacobian, matrix
1223 
1224 .seealso: SLESSetOperators(), SNESSetFunction()
1225 @*/
1226 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1227 {
1228   int ierr;
1229 
1230   PetscFunctionBegin;
1231   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1232   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE);
1233   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE);
1234   if (A) PetscCheckSameComm(snes,A);
1235   if (B) PetscCheckSameComm(snes,B);
1236   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1237     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1238   }
1239 
1240   if (func) snes->computejacobian = func;
1241   if (ctx)  snes->jacP            = ctx;
1242   if (A) {
1243     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1244     snes->jacobian = A;
1245     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1246   }
1247   if (B) {
1248     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1249     snes->jacobian_pre = B;
1250     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
1251   }
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 #undef __FUNCT__
1256 #define __FUNCT__ "SNESGetJacobian"
1257 /*@C
1258    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1259    provided context for evaluating the Jacobian.
1260 
1261    Not Collective, but Mat object will be parallel if SNES object is
1262 
1263    Input Parameter:
1264 .  snes - the nonlinear solver context
1265 
1266    Output Parameters:
1267 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1268 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1269 .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1270 -  func - location to put Jacobian function (or PETSC_NULL)
1271 
1272    Level: advanced
1273 
1274 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1275 @*/
1276 int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
1277 {
1278   PetscFunctionBegin;
1279   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1280   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1281     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1282   }
1283   if (A)    *A    = snes->jacobian;
1284   if (B)    *B    = snes->jacobian_pre;
1285   if (ctx)  *ctx  = snes->jacP;
1286   if (func) *func = snes->computejacobian;
1287   PetscFunctionReturn(0);
1288 }
1289 
1290 #undef __FUNCT__
1291 #define __FUNCT__ "SNESSetHessian"
1292 /*@C
1293    SNESSetHessian - Sets the function to compute Hessian as well as the
1294    location to store the matrix.
1295 
1296    Collective on SNES and Mat
1297 
1298    Input Parameters:
1299 +  snes - the SNES context
1300 .  A - Hessian matrix
1301 .  B - preconditioner matrix (usually same as the Hessian)
1302 .  func - Jacobian evaluation routine
1303 -  ctx - [optional] user-defined context for private data for the
1304          Hessian evaluation routine (may be PETSC_NULL)
1305 
1306    Calling sequence of func:
1307 $    func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1308 
1309 +  x - input vector
1310 .  A - Hessian matrix
1311 .  B - preconditioner matrix, usually the same as A
1312 .  flag - flag indicating information about the preconditioner matrix
1313    structure (same as flag in SLESSetOperators())
1314 -  ctx - [optional] user-defined Hessian context
1315 
1316    Notes:
1317    See SLESSetOperators() for important information about setting the flag
1318    output parameter in the routine func().  Be sure to read this information!
1319 
1320    The function func() takes Mat * as the matrix arguments rather than Mat.
1321    This allows the Hessian evaluation routine to replace A and/or B with a
1322    completely new new matrix structure (not just different matrix elements)
1323    when appropriate, for instance, if the nonzero structure is changing
1324    throughout the global iterations.
1325 
1326    Level: beginner
1327 
1328 .keywords: SNES, nonlinear, set, Hessian, matrix
1329 
1330 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
1331 @*/
1332 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1333 {
1334   int ierr;
1335 
1336   PetscFunctionBegin;
1337   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1338   PetscValidHeaderSpecific(A,MAT_COOKIE);
1339   PetscValidHeaderSpecific(B,MAT_COOKIE);
1340   PetscCheckSameComm(snes,A);
1341   PetscCheckSameComm(snes,B);
1342   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1343     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1344   }
1345   if (func) snes->computejacobian = func;
1346   if (ctx)  snes->jacP            = ctx;
1347   if (A) {
1348     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1349     snes->jacobian = A;
1350     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1351   }
1352   if (B) {
1353     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1354     snes->jacobian_pre = B;
1355     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
1356   }
1357   PetscFunctionReturn(0);
1358 }
1359 
1360 #undef __FUNCT__
1361 #define __FUNCT__ "SNESGetHessian"
1362 /*@
1363    SNESGetHessian - Returns the Hessian matrix and optionally the user
1364    provided context for evaluating the Hessian.
1365 
1366    Not Collective, but Mat object is parallel if SNES object is parallel
1367 
1368    Input Parameter:
1369 .  snes - the nonlinear solver context
1370 
1371    Output Parameters:
1372 +  A - location to stash Hessian matrix (or PETSC_NULL)
1373 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1374 -  ctx - location to stash Hessian ctx (or PETSC_NULL)
1375 
1376    Level: advanced
1377 
1378 .seealso: SNESSetHessian(), SNESComputeHessian()
1379 
1380 .keywords: SNES, get, Hessian
1381 @*/
1382 int SNESGetHessian(SNES snes,Mat *A,Mat *B,void **ctx)
1383 {
1384   PetscFunctionBegin;
1385   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1386   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1387     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1388   }
1389   if (A)   *A = snes->jacobian;
1390   if (B)   *B = snes->jacobian_pre;
1391   if (ctx) *ctx = snes->jacP;
1392   PetscFunctionReturn(0);
1393 }
1394 
1395 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1396 
1397 #undef __FUNCT__
1398 #define __FUNCT__ "SNESSetUp"
1399 /*@
1400    SNESSetUp - Sets up the internal data structures for the later use
1401    of a nonlinear solver.
1402 
1403    Collective on SNES
1404 
1405    Input Parameters:
1406 +  snes - the SNES context
1407 -  x - the solution vector
1408 
1409    Notes:
1410    For basic use of the SNES solvers the user need not explicitly call
1411    SNESSetUp(), since these actions will automatically occur during
1412    the call to SNESSolve().  However, if one wishes to control this
1413    phase separately, SNESSetUp() should be called after SNESCreate()
1414    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1415 
1416    Level: advanced
1417 
1418 .keywords: SNES, nonlinear, setup
1419 
1420 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1421 @*/
1422 int SNESSetUp(SNES snes,Vec x)
1423 {
1424   int        ierr;
1425   PetscTruth flg;
1426 
1427   PetscFunctionBegin;
1428   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1429   PetscValidHeaderSpecific(x,VEC_COOKIE);
1430   PetscCheckSameComm(snes,x);
1431   snes->vec_sol = snes->vec_sol_always = x;
1432 
1433   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1434   /*
1435       This version replaces the user provided Jacobian matrix with a
1436       matrix-free version but still employs the user-provided preconditioner matrix
1437   */
1438   if (flg) {
1439     Mat J;
1440     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1441     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1442     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n");
1443     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
1444     ierr = MatDestroy(J);CHKERRQ(ierr);
1445   }
1446   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1447   /*
1448       This version replaces both the user-provided Jacobian and the user-
1449       provided preconditioner matrix with the default matrix free version.
1450    */
1451   if (flg) {
1452     Mat  J;
1453     SLES sles;
1454     PC   pc;
1455 
1456     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1457     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1458     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n");
1459     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1460       ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
1461     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1462       ierr = SNESSetHessian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
1463     } else {
1464       SETERRQ(PETSC_ERR_SUP,"Method class doesn't support matrix-free option");
1465     }
1466     ierr = MatDestroy(J);CHKERRQ(ierr);
1467 
1468     /* force no preconditioner */
1469     ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1470     ierr = SLESGetPC(sles,&pc);CHKERRQ(ierr);
1471     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1472   }
1473 
1474   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1475     PetscTruth iseqtr;
1476 
1477     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1478     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1479     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1480     if (snes->vec_func == snes->vec_sol) {
1481       SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1482     }
1483 
1484     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
1485     ierr = PetscTypeCompare((PetscObject)snes,SNESEQTR,&iseqtr);CHKERRQ(ierr);
1486     if (snes->ksp_ewconv && !iseqtr) {
1487       SLES sles; KSP ksp;
1488       ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1489       ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
1490       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1491     }
1492   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1493     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first");
1494     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first");
1495     if (!snes->computeumfunction) {
1496       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetMinimizationFunction() first");
1497     }
1498     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetHessian()");
1499   } else {
1500     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class");
1501   }
1502   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
1503   snes->setupcalled = 1;
1504   PetscFunctionReturn(0);
1505 }
1506 
1507 #undef __FUNCT__
1508 #define __FUNCT__ "SNESDestroy"
1509 /*@C
1510    SNESDestroy - Destroys the nonlinear solver context that was created
1511    with SNESCreate().
1512 
1513    Collective on SNES
1514 
1515    Input Parameter:
1516 .  snes - the SNES context
1517 
1518    Level: beginner
1519 
1520 .keywords: SNES, nonlinear, destroy
1521 
1522 .seealso: SNESCreate(), SNESSolve()
1523 @*/
1524 int SNESDestroy(SNES snes)
1525 {
1526   int i,ierr;
1527 
1528   PetscFunctionBegin;
1529   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1530   if (--snes->refct > 0) PetscFunctionReturn(0);
1531 
1532   /* if memory was published with AMS then destroy it */
1533   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1534 
1535   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1536   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
1537   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1538   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1539   ierr = SLESDestroy(snes->sles);CHKERRQ(ierr);
1540   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1541   for (i=0; i<snes->numbermonitors; i++) {
1542     if (snes->monitordestroy[i]) {
1543       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1544     }
1545   }
1546   PetscLogObjectDestroy((PetscObject)snes);
1547   PetscHeaderDestroy((PetscObject)snes);
1548   PetscFunctionReturn(0);
1549 }
1550 
1551 /* ----------- Routines to set solver parameters ---------- */
1552 
1553 #undef __FUNCT__
1554 #define __FUNCT__ "SNESSetTolerances"
1555 /*@
1556    SNESSetTolerances - Sets various parameters used in convergence tests.
1557 
1558    Collective on SNES
1559 
1560    Input Parameters:
1561 +  snes - the SNES context
1562 .  atol - absolute convergence tolerance
1563 .  rtol - relative convergence tolerance
1564 .  stol -  convergence tolerance in terms of the norm
1565            of the change in the solution between steps
1566 .  maxit - maximum number of iterations
1567 -  maxf - maximum number of function evaluations
1568 
1569    Options Database Keys:
1570 +    -snes_atol <atol> - Sets atol
1571 .    -snes_rtol <rtol> - Sets rtol
1572 .    -snes_stol <stol> - Sets stol
1573 .    -snes_max_it <maxit> - Sets maxit
1574 -    -snes_max_funcs <maxf> - Sets maxf
1575 
1576    Notes:
1577    The default maximum number of iterations is 50.
1578    The default maximum number of function evaluations is 1000.
1579 
1580    Level: intermediate
1581 
1582 .keywords: SNES, nonlinear, set, convergence, tolerances
1583 
1584 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
1585 @*/
1586 int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
1587 {
1588   PetscFunctionBegin;
1589   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1590   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1591   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1592   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1593   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1594   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1595   PetscFunctionReturn(0);
1596 }
1597 
1598 #undef __FUNCT__
1599 #define __FUNCT__ "SNESGetTolerances"
1600 /*@
1601    SNESGetTolerances - Gets various parameters used in convergence tests.
1602 
1603    Not Collective
1604 
1605    Input Parameters:
1606 +  snes - the SNES context
1607 .  atol - absolute convergence tolerance
1608 .  rtol - relative convergence tolerance
1609 .  stol -  convergence tolerance in terms of the norm
1610            of the change in the solution between steps
1611 .  maxit - maximum number of iterations
1612 -  maxf - maximum number of function evaluations
1613 
1614    Notes:
1615    The user can specify PETSC_NULL for any parameter that is not needed.
1616 
1617    Level: intermediate
1618 
1619 .keywords: SNES, nonlinear, get, convergence, tolerances
1620 
1621 .seealso: SNESSetTolerances()
1622 @*/
1623 int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
1624 {
1625   PetscFunctionBegin;
1626   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1627   if (atol)  *atol  = snes->atol;
1628   if (rtol)  *rtol  = snes->rtol;
1629   if (stol)  *stol  = snes->xtol;
1630   if (maxit) *maxit = snes->max_its;
1631   if (maxf)  *maxf  = snes->max_funcs;
1632   PetscFunctionReturn(0);
1633 }
1634 
1635 #undef __FUNCT__
1636 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1637 /*@
1638    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1639 
1640    Collective on SNES
1641 
1642    Input Parameters:
1643 +  snes - the SNES context
1644 -  tol - tolerance
1645 
1646    Options Database Key:
1647 .  -snes_trtol <tol> - Sets tol
1648 
1649    Level: intermediate
1650 
1651 .keywords: SNES, nonlinear, set, trust region, tolerance
1652 
1653 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
1654 @*/
1655 int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1656 {
1657   PetscFunctionBegin;
1658   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1659   snes->deltatol = tol;
1660   PetscFunctionReturn(0);
1661 }
1662 
1663 #undef __FUNCT__
1664 #define __FUNCT__ "SNESSetMinimizationFunctionTolerance"
1665 /*@
1666    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
1667    for unconstrained minimization solvers.
1668 
1669    Collective on SNES
1670 
1671    Input Parameters:
1672 +  snes - the SNES context
1673 -  ftol - minimum function tolerance
1674 
1675    Options Database Key:
1676 .  -snes_fmin <ftol> - Sets ftol
1677 
1678    Note:
1679    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
1680    methods only.
1681 
1682    Level: intermediate
1683 
1684 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
1685 
1686 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
1687 @*/
1688 int SNESSetMinimizationFunctionTolerance(SNES snes,PetscReal ftol)
1689 {
1690   PetscFunctionBegin;
1691   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1692   snes->fmin = ftol;
1693   PetscFunctionReturn(0);
1694 }
1695 /*
1696    Duplicate the lg monitors for SNES from KSP; for some reason with
1697    dynamic libraries things don't work under Sun4 if we just use
1698    macros instead of functions
1699 */
1700 #undef __FUNCT__
1701 #define __FUNCT__ "SNESLGMonitor"
1702 int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1703 {
1704   int ierr;
1705 
1706   PetscFunctionBegin;
1707   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1708   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1709   PetscFunctionReturn(0);
1710 }
1711 
1712 #undef __FUNCT__
1713 #define __FUNCT__ "SNESLGMonitorCreate"
1714 int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw)
1715 {
1716   int ierr;
1717 
1718   PetscFunctionBegin;
1719   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1720   PetscFunctionReturn(0);
1721 }
1722 
1723 #undef __FUNCT__
1724 #define __FUNCT__ "SNESLGMonitorDestroy"
1725 int SNESLGMonitorDestroy(PetscDrawLG draw)
1726 {
1727   int ierr;
1728 
1729   PetscFunctionBegin;
1730   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1731   PetscFunctionReturn(0);
1732 }
1733 
1734 /* ------------ Routines to set performance monitoring options ----------- */
1735 
1736 #undef __FUNCT__
1737 #define __FUNCT__ "SNESSetMonitor"
1738 /*@C
1739    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
1740    iteration of the nonlinear solver to display the iteration's
1741    progress.
1742 
1743    Collective on SNES
1744 
1745    Input Parameters:
1746 +  snes - the SNES context
1747 .  func - monitoring routine
1748 .  mctx - [optional] user-defined context for private data for the
1749           monitor routine (use PETSC_NULL if no context is desitre)
1750 -  monitordestroy - [optional] routine that frees monitor context
1751           (may be PETSC_NULL)
1752 
1753    Calling sequence of func:
1754 $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1755 
1756 +    snes - the SNES context
1757 .    its - iteration number
1758 .    norm - 2-norm function value (may be estimated)
1759             for SNES_NONLINEAR_EQUATIONS methods
1760 .    norm - 2-norm gradient value (may be estimated)
1761             for SNES_UNCONSTRAINED_MINIMIZATION methods
1762 -    mctx - [optional] monitoring context
1763 
1764    Options Database Keys:
1765 +    -snes_monitor        - sets SNESDefaultMonitor()
1766 .    -snes_xmonitor       - sets line graph monitor,
1767                             uses SNESLGMonitorCreate()
1768 _    -snes_cancelmonitors - cancels all monitors that have
1769                             been hardwired into a code by
1770                             calls to SNESSetMonitor(), but
1771                             does not cancel those set via
1772                             the options database.
1773 
1774    Notes:
1775    Several different monitoring routines may be set by calling
1776    SNESSetMonitor() multiple times; all will be called in the
1777    order in which they were set.
1778 
1779    Level: intermediate
1780 
1781 .keywords: SNES, nonlinear, set, monitor
1782 
1783 .seealso: SNESDefaultMonitor(), SNESClearMonitor()
1784 @*/
1785 int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *))
1786 {
1787   PetscFunctionBegin;
1788   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1789   if (snes->numbermonitors >= MAXSNESMONITORS) {
1790     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1791   }
1792 
1793   snes->monitor[snes->numbermonitors]           = func;
1794   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1795   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1796   PetscFunctionReturn(0);
1797 }
1798 
1799 #undef __FUNCT__
1800 #define __FUNCT__ "SNESClearMonitor"
1801 /*@C
1802    SNESClearMonitor - Clears all the monitor functions for a SNES object.
1803 
1804    Collective on SNES
1805 
1806    Input Parameters:
1807 .  snes - the SNES context
1808 
1809    Options Database:
1810 .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1811     into a code by calls to SNESSetMonitor(), but does not cancel those
1812     set via the options database
1813 
1814    Notes:
1815    There is no way to clear one specific monitor from a SNES object.
1816 
1817    Level: intermediate
1818 
1819 .keywords: SNES, nonlinear, set, monitor
1820 
1821 .seealso: SNESDefaultMonitor(), SNESSetMonitor()
1822 @*/
1823 int SNESClearMonitor(SNES snes)
1824 {
1825   PetscFunctionBegin;
1826   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1827   snes->numbermonitors = 0;
1828   PetscFunctionReturn(0);
1829 }
1830 
1831 #undef __FUNCT__
1832 #define __FUNCT__ "SNESSetConvergenceTest"
1833 /*@C
1834    SNESSetConvergenceTest - Sets the function that is to be used
1835    to test for convergence of the nonlinear iterative solution.
1836 
1837    Collective on SNES
1838 
1839    Input Parameters:
1840 +  snes - the SNES context
1841 .  func - routine to test for convergence
1842 -  cctx - [optional] context for private data for the convergence routine
1843           (may be PETSC_NULL)
1844 
1845    Calling sequence of func:
1846 $     int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1847 
1848 +    snes - the SNES context
1849 .    cctx - [optional] convergence context
1850 .    reason - reason for convergence/divergence
1851 .    xnorm - 2-norm of current iterate
1852 .    gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods)
1853 .    f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods)
1854 .    gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods)
1855 -    f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods)
1856 
1857    Level: advanced
1858 
1859 .keywords: SNES, nonlinear, set, convergence, test
1860 
1861 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1862           SNESConverged_UM_LS(), SNESConverged_UM_TR()
1863 @*/
1864 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
1865 {
1866   PetscFunctionBegin;
1867   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1868   (snes)->converged = func;
1869   (snes)->cnvP      = cctx;
1870   PetscFunctionReturn(0);
1871 }
1872 
1873 #undef __FUNCT__
1874 #define __FUNCT__ "SNESGetConvergedReason"
1875 /*@C
1876    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1877 
1878    Not Collective
1879 
1880    Input Parameter:
1881 .  snes - the SNES context
1882 
1883    Output Parameter:
1884 .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1885             manual pages for the individual convergence tests for complete lists
1886 
1887    Level: intermediate
1888 
1889    Notes: Can only be called after the call the SNESSolve() is complete.
1890 
1891 .keywords: SNES, nonlinear, set, convergence, test
1892 
1893 .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1894           SNESConverged_UM_LS(), SNESConverged_UM_TR(), SNESConvergedReason
1895 @*/
1896 int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1897 {
1898   PetscFunctionBegin;
1899   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1900   *reason = snes->reason;
1901   PetscFunctionReturn(0);
1902 }
1903 
1904 #undef __FUNCT__
1905 #define __FUNCT__ "SNESSetConvergenceHistory"
1906 /*@
1907    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1908 
1909    Collective on SNES
1910 
1911    Input Parameters:
1912 +  snes - iterative context obtained from SNESCreate()
1913 .  a   - array to hold history
1914 .  its - integer array holds the number of linear iterations for each solve.
1915 .  na  - size of a and its
1916 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1917            else it continues storing new values for new nonlinear solves after the old ones
1918 
1919    Notes:
1920    If set, this array will contain the function norms (for
1921    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1922    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1923    at each step.
1924 
1925    This routine is useful, e.g., when running a code for purposes
1926    of accurate performance monitoring, when no I/O should be done
1927    during the section of code that is being timed.
1928 
1929    Level: intermediate
1930 
1931 .keywords: SNES, set, convergence, history
1932 
1933 .seealso: SNESGetConvergenceHistory()
1934 
1935 @*/
1936 int SNESSetConvergenceHistory(SNES snes,PetscReal *a,int *its,int na,PetscTruth reset)
1937 {
1938   PetscFunctionBegin;
1939   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1940   if (na) PetscValidScalarPointer(a);
1941   snes->conv_hist       = a;
1942   snes->conv_hist_its   = its;
1943   snes->conv_hist_max   = na;
1944   snes->conv_hist_reset = reset;
1945   PetscFunctionReturn(0);
1946 }
1947 
1948 #undef __FUNCT__
1949 #define __FUNCT__ "SNESGetConvergenceHistory"
1950 /*@C
1951    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1952 
1953    Collective on SNES
1954 
1955    Input Parameter:
1956 .  snes - iterative context obtained from SNESCreate()
1957 
1958    Output Parameters:
1959 .  a   - array to hold history
1960 .  its - integer array holds the number of linear iterations (or
1961          negative if not converged) for each solve.
1962 -  na  - size of a and its
1963 
1964    Notes:
1965     The calling sequence for this routine in Fortran is
1966 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1967 
1968    This routine is useful, e.g., when running a code for purposes
1969    of accurate performance monitoring, when no I/O should be done
1970    during the section of code that is being timed.
1971 
1972    Level: intermediate
1973 
1974 .keywords: SNES, get, convergence, history
1975 
1976 .seealso: SNESSetConvergencHistory()
1977 
1978 @*/
1979 int SNESGetConvergenceHistory(SNES snes,PetscReal **a,int **its,int *na)
1980 {
1981   PetscFunctionBegin;
1982   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1983   if (a)   *a   = snes->conv_hist;
1984   if (its) *its = snes->conv_hist_its;
1985   if (na) *na   = snes->conv_hist_len;
1986   PetscFunctionReturn(0);
1987 }
1988 
1989 #undef __FUNCT__
1990 #define __FUNCT__ "SNESSetRhsBC"
1991 /*@
1992   SNESSetRhsBC - Sets the function which applies boundary conditions
1993   to the Rhs of each system.
1994 
1995   Collective on SNES
1996 
1997   Input Parameters:
1998 . snes - The nonlinear solver context
1999 . func - The function
2000 
2001   Calling sequence of func:
2002 . func (SNES snes, Vec rhs, void *ctx);
2003 
2004 . rhs - The current rhs vector
2005 . ctx - The user-context
2006 
2007   Level: intermediate
2008 
2009 .keywords: SNES, Rhs, boundary conditions
2010 .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate()
2011 @*/
2012 int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *))
2013 {
2014   PetscFunctionBegin;
2015   PetscValidHeaderSpecific(snes, SNES_COOKIE);
2016   snes->applyrhsbc = func;
2017   PetscFunctionReturn(0);
2018 }
2019 
2020 #undef __FUNCT__
2021 #define __FUNCT__ "SNESDefaultRhsBC"
2022 /*@
2023   SNESDefaultRhsBC - The default boundary condition function which does nothing.
2024 
2025   Not collective
2026 
2027   Input Parameters:
2028 . snes - The nonlinear solver context
2029 . rhs  - The Rhs
2030 . ctx  - The user-context
2031 
2032   Level: beginner
2033 
2034 .keywords: SNES, Rhs, boundary conditions
2035 .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate()
2036 @*/
2037 int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx)
2038 {
2039   PetscFunctionBegin;
2040   PetscFunctionReturn(0);
2041 }
2042 
2043 #undef __FUNCT__
2044 #define __FUNCT__ "SNESSetSolutionBC"
2045 /*@
2046   SNESSetSolutionBC - Sets the function which applies boundary conditions
2047   to the solution of each system.
2048 
2049   Collective on SNES
2050 
2051   Input Parameters:
2052 . snes - The nonlinear solver context
2053 . func - The function
2054 
2055   Calling sequence of func:
2056 . func (TS ts, Vec rsol, void *ctx);
2057 
2058 . sol - The current solution vector
2059 . ctx - The user-context
2060 
2061   Level: intermediate
2062 
2063 .keywords: SNES, solution, boundary conditions
2064 .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate()
2065 @*/
2066 int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *))
2067 {
2068   PetscFunctionBegin;
2069   PetscValidHeaderSpecific(snes, SNES_COOKIE);
2070   snes->applysolbc = func;
2071   PetscFunctionReturn(0);
2072 }
2073 
2074 #undef __FUNCT__
2075 #define __FUNCT__ "SNESDefaultSolutionBC"
2076 /*@
2077   SNESDefaultSolutionBC - The default boundary condition function which does nothing.
2078 
2079   Not collective
2080 
2081   Input Parameters:
2082 . snes - The nonlinear solver context
2083 . sol  - The solution
2084 . ctx  - The user-context
2085 
2086   Level: beginner
2087 
2088 .keywords: SNES, solution, boundary conditions
2089 .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate()
2090 @*/
2091 int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx)
2092 {
2093   PetscFunctionBegin;
2094   PetscFunctionReturn(0);
2095 }
2096 
2097 #undef __FUNCT__
2098 #define __FUNCT__ "SNESSetUpdate"
2099 /*@
2100   SNESSetUpdate - Sets the general-purpose update function called
2101   at the beginning of every step of the iteration.
2102 
2103   Collective on SNES
2104 
2105   Input Parameters:
2106 . snes - The nonlinear solver context
2107 . func - The function
2108 
2109   Calling sequence of func:
2110 . func (TS ts, int step);
2111 
2112 . step - The current step of the iteration
2113 
2114   Level: intermediate
2115 
2116 .keywords: SNES, update
2117 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
2118 @*/
2119 int SNESSetUpdate(SNES snes, int (*func)(SNES, int))
2120 {
2121   PetscFunctionBegin;
2122   PetscValidHeaderSpecific(snes, SNES_COOKIE);
2123   snes->update = func;
2124   PetscFunctionReturn(0);
2125 }
2126 
2127 #undef __FUNCT__
2128 #define __FUNCT__ "SNESDefaultUpdate"
2129 /*@
2130   SNESDefaultUpdate - The default update function which does nothing.
2131 
2132   Not collective
2133 
2134   Input Parameters:
2135 . snes - The nonlinear solver context
2136 . step - The current step of the iteration
2137 
2138 .keywords: SNES, update
2139 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
2140 @*/
2141 int SNESDefaultUpdate(SNES snes, int step)
2142 {
2143   PetscFunctionBegin;
2144   PetscFunctionReturn(0);
2145 }
2146 
2147 #undef __FUNCT__
2148 #define __FUNCT__ "SNESScaleStep_Private"
2149 /*
2150    SNESScaleStep_Private - Scales a step so that its length is less than the
2151    positive parameter delta.
2152 
2153     Input Parameters:
2154 +   snes - the SNES context
2155 .   y - approximate solution of linear system
2156 .   fnorm - 2-norm of current function
2157 -   delta - trust region size
2158 
2159     Output Parameters:
2160 +   gpnorm - predicted function norm at the new point, assuming local
2161     linearization.  The value is zero if the step lies within the trust
2162     region, and exceeds zero otherwise.
2163 -   ynorm - 2-norm of the step
2164 
2165     Note:
2166     For non-trust region methods such as SNESEQLS, the parameter delta
2167     is set to be the maximum allowable step size.
2168 
2169 .keywords: SNES, nonlinear, scale, step
2170 */
2171 int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,
2172                           PetscReal *gpnorm,PetscReal *ynorm)
2173 {
2174   PetscReal norm;
2175   PetscScalar cnorm;
2176   int    ierr;
2177 
2178   PetscFunctionBegin;
2179   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2180   PetscValidHeaderSpecific(y,VEC_COOKIE);
2181   PetscCheckSameComm(snes,y);
2182 
2183   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
2184   if (norm > *delta) {
2185      norm = *delta/norm;
2186      *gpnorm = (1.0 - norm)*(*fnorm);
2187      cnorm = norm;
2188      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
2189      *ynorm = *delta;
2190   } else {
2191      *gpnorm = 0.0;
2192      *ynorm = norm;
2193   }
2194   PetscFunctionReturn(0);
2195 }
2196 
2197 #undef __FUNCT__
2198 #define __FUNCT__ "SNESSolve"
2199 /*@
2200    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
2201    SNESCreate() and optional routines of the form SNESSetXXX().
2202 
2203    Collective on SNES
2204 
2205    Input Parameters:
2206 +  snes - the SNES context
2207 -  x - the solution vector
2208 
2209    Output Parameter:
2210 .  its - number of iterations until termination
2211 
2212    Notes:
2213    The user should initialize the vector,x, with the initial guess
2214    for the nonlinear solve prior to calling SNESSolve.  In particular,
2215    to employ an initial guess of zero, the user should explicitly set
2216    this vector to zero by calling VecSet().
2217 
2218    Level: beginner
2219 
2220 .keywords: SNES, nonlinear, solve
2221 
2222 .seealso: SNESCreate(), SNESDestroy()
2223 @*/
2224 int SNESSolve(SNES snes,Vec x,int *its)
2225 {
2226   int        ierr;
2227   PetscTruth flg;
2228 
2229   PetscFunctionBegin;
2230   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2231   PetscValidHeaderSpecific(x,VEC_COOKIE);
2232   PetscCheckSameComm(snes,x);
2233   PetscValidIntPointer(its);
2234   if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
2235 
2236   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
2237   else {snes->vec_sol = snes->vec_sol_always = x;}
2238   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
2239   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2240   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
2241   ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr);
2242   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2243   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
2244   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); }
2245   PetscFunctionReturn(0);
2246 }
2247 
2248 /* --------- Internal routines for SNES Package --------- */
2249 
2250 #undef __FUNCT__
2251 #define __FUNCT__ "SNESSetType"
2252 /*@C
2253    SNESSetType - Sets the method for the nonlinear solver.
2254 
2255    Collective on SNES
2256 
2257    Input Parameters:
2258 +  snes - the SNES context
2259 -  type - a known method
2260 
2261    Options Database Key:
2262 .  -snes_type <type> - Sets the method; use -help for a list
2263    of available methods (for instance, ls or tr)
2264 
2265    Notes:
2266    See "petsc/include/petscsnes.h" for available methods (for instance)
2267 +    SNESEQLS - Newton's method with line search
2268      (systems of nonlinear equations)
2269 .    SNESEQTR - Newton's method with trust region
2270      (systems of nonlinear equations)
2271 .    SNESUMTR - Newton's method with trust region
2272      (unconstrained minimization)
2273 -    SNESUMLS - Newton's method with line search
2274      (unconstrained minimization)
2275 
2276   Normally, it is best to use the SNESSetFromOptions() command and then
2277   set the SNES solver type from the options database rather than by using
2278   this routine.  Using the options database provides the user with
2279   maximum flexibility in evaluating the many nonlinear solvers.
2280   The SNESSetType() routine is provided for those situations where it
2281   is necessary to set the nonlinear solver independently of the command
2282   line or options database.  This might be the case, for example, when
2283   the choice of solver changes during the execution of the program,
2284   and the user's application is taking responsibility for choosing the
2285   appropriate method.
2286 
2287   Level: intermediate
2288 
2289 .keywords: SNES, set, type
2290 
2291 .seealso: SNESType, SNESCreate()
2292 
2293 @*/
2294 int SNESSetType(SNES snes,SNESType type)
2295 {
2296   int        ierr,(*r)(SNES);
2297   PetscTruth match;
2298 
2299   PetscFunctionBegin;
2300   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2301   PetscValidCharPointer(type);
2302 
2303   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
2304   if (match) PetscFunctionReturn(0);
2305 
2306   if (snes->setupcalled) {
2307     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
2308     snes->data = 0;
2309   }
2310 
2311   /* Get the function pointers for the iterative method requested */
2312   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
2313 
2314   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
2315 
2316   if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type);
2317 
2318   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
2319   snes->data = 0;
2320   ierr = (*r)(snes);CHKERRQ(ierr);
2321 
2322   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
2323   snes->set_method_called = 1;
2324 
2325   PetscFunctionReturn(0);
2326 }
2327 
2328 
2329 /* --------------------------------------------------------------------- */
2330 #undef __FUNCT__
2331 #define __FUNCT__ "SNESRegisterDestroy"
2332 /*@C
2333    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2334    registered by SNESRegisterDynamic().
2335 
2336    Not Collective
2337 
2338    Level: advanced
2339 
2340 .keywords: SNES, nonlinear, register, destroy
2341 
2342 .seealso: SNESRegisterAll(), SNESRegisterAll()
2343 @*/
2344 int SNESRegisterDestroy(void)
2345 {
2346   int ierr;
2347 
2348   PetscFunctionBegin;
2349   if (SNESList) {
2350     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
2351     SNESList = 0;
2352   }
2353   SNESRegisterAllCalled = PETSC_FALSE;
2354   PetscFunctionReturn(0);
2355 }
2356 
2357 #undef __FUNCT__
2358 #define __FUNCT__ "SNESGetType"
2359 /*@C
2360    SNESGetType - Gets the SNES method type and name (as a string).
2361 
2362    Not Collective
2363 
2364    Input Parameter:
2365 .  snes - nonlinear solver context
2366 
2367    Output Parameter:
2368 .  type - SNES method (a character string)
2369 
2370    Level: intermediate
2371 
2372 .keywords: SNES, nonlinear, get, type, name
2373 @*/
2374 int SNESGetType(SNES snes,SNESType *type)
2375 {
2376   PetscFunctionBegin;
2377   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2378   *type = snes->type_name;
2379   PetscFunctionReturn(0);
2380 }
2381 
2382 #undef __FUNCT__
2383 #define __FUNCT__ "SNESGetSolution"
2384 /*@C
2385    SNESGetSolution - Returns the vector where the approximate solution is
2386    stored.
2387 
2388    Not Collective, but Vec is parallel if SNES is parallel
2389 
2390    Input Parameter:
2391 .  snes - the SNES context
2392 
2393    Output Parameter:
2394 .  x - the solution
2395 
2396    Level: advanced
2397 
2398 .keywords: SNES, nonlinear, get, solution
2399 
2400 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
2401 @*/
2402 int SNESGetSolution(SNES snes,Vec *x)
2403 {
2404   PetscFunctionBegin;
2405   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2406   *x = snes->vec_sol_always;
2407   PetscFunctionReturn(0);
2408 }
2409 
2410 #undef __FUNCT__
2411 #define __FUNCT__ "SNESGetSolutionUpdate"
2412 /*@C
2413    SNESGetSolutionUpdate - Returns the vector where the solution update is
2414    stored.
2415 
2416    Not Collective, but Vec is parallel if SNES is parallel
2417 
2418    Input Parameter:
2419 .  snes - the SNES context
2420 
2421    Output Parameter:
2422 .  x - the solution update
2423 
2424    Level: advanced
2425 
2426 .keywords: SNES, nonlinear, get, solution, update
2427 
2428 .seealso: SNESGetSolution(), SNESGetFunction
2429 @*/
2430 int SNESGetSolutionUpdate(SNES snes,Vec *x)
2431 {
2432   PetscFunctionBegin;
2433   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2434   *x = snes->vec_sol_update_always;
2435   PetscFunctionReturn(0);
2436 }
2437 
2438 #undef __FUNCT__
2439 #define __FUNCT__ "SNESGetFunction"
2440 /*@C
2441    SNESGetFunction - Returns the vector where the function is stored.
2442 
2443    Not Collective, but Vec is parallel if SNES is parallel
2444 
2445    Input Parameter:
2446 .  snes - the SNES context
2447 
2448    Output Parameter:
2449 +  r - the function (or PETSC_NULL)
2450 .  ctx - the function context (or PETSC_NULL)
2451 -  func - the function (or PETSC_NULL)
2452 
2453    Notes:
2454    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
2455    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
2456    SNESGetMinimizationFunction() and SNESGetGradient();
2457 
2458    Level: advanced
2459 
2460 .keywords: SNES, nonlinear, get, function
2461 
2462 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
2463           SNESGetGradient()
2464 
2465 @*/
2466 int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*))
2467 {
2468   PetscFunctionBegin;
2469   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2470   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
2471     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
2472   }
2473   if (r)    *r    = snes->vec_func_always;
2474   if (ctx)  *ctx  = snes->funP;
2475   if (func) *func = snes->computefunction;
2476   PetscFunctionReturn(0);
2477 }
2478 
2479 #undef __FUNCT__
2480 #define __FUNCT__ "SNESGetGradient"
2481 /*@C
2482    SNESGetGradient - Returns the vector where the gradient is stored.
2483 
2484    Not Collective, but Vec is parallel if SNES is parallel
2485 
2486    Input Parameter:
2487 .  snes - the SNES context
2488 
2489    Output Parameter:
2490 +  r - the gradient (or PETSC_NULL)
2491 -  ctx - the gradient context (or PETSC_NULL)
2492 
2493    Notes:
2494    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
2495    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
2496    SNESGetFunction().
2497 
2498    Level: advanced
2499 
2500 .keywords: SNES, nonlinear, get, gradient
2501 
2502 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction(),
2503           SNESSetGradient(), SNESSetFunction()
2504 
2505 @*/
2506 int SNESGetGradient(SNES snes,Vec *r,void **ctx)
2507 {
2508   PetscFunctionBegin;
2509   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2510   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
2511     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
2512   }
2513   if (r)   *r = snes->vec_func_always;
2514   if (ctx) *ctx = snes->funP;
2515   PetscFunctionReturn(0);
2516 }
2517 
2518 #undef __FUNCT__
2519 #define __FUNCT__ "SNESGetMinimizationFunction"
2520 /*@C
2521    SNESGetMinimizationFunction - Returns the scalar function value for
2522    unconstrained minimization problems.
2523 
2524    Not Collective
2525 
2526    Input Parameter:
2527 .  snes - the SNES context
2528 
2529    Output Parameter:
2530 +  r - the function (or PETSC_NULL)
2531 -  ctx - the function context (or PETSC_NULL)
2532 
2533    Notes:
2534    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
2535    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
2536    SNESGetFunction().
2537 
2538    Level: advanced
2539 
2540 .keywords: SNES, nonlinear, get, function
2541 
2542 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction(), SNESSetFunction()
2543 
2544 @*/
2545 int SNESGetMinimizationFunction(SNES snes,PetscReal *r,void **ctx)
2546 {
2547   PetscFunctionBegin;
2548   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2549   PetscValidScalarPointer(r);
2550   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
2551     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
2552   }
2553   if (r)   *r = snes->fc;
2554   if (ctx) *ctx = snes->umfunP;
2555   PetscFunctionReturn(0);
2556 }
2557 
2558 #undef __FUNCT__
2559 #define __FUNCT__ "SNESSetOptionsPrefix"
2560 /*@C
2561    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2562    SNES options in the database.
2563 
2564    Collective on SNES
2565 
2566    Input Parameter:
2567 +  snes - the SNES context
2568 -  prefix - the prefix to prepend to all option names
2569 
2570    Notes:
2571    A hyphen (-) must NOT be given at the beginning of the prefix name.
2572    The first character of all runtime options is AUTOMATICALLY the hyphen.
2573 
2574    Level: advanced
2575 
2576 .keywords: SNES, set, options, prefix, database
2577 
2578 .seealso: SNESSetFromOptions()
2579 @*/
2580 int SNESSetOptionsPrefix(SNES snes,char *prefix)
2581 {
2582   int ierr;
2583 
2584   PetscFunctionBegin;
2585   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2586   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2587   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2588   PetscFunctionReturn(0);
2589 }
2590 
2591 #undef __FUNCT__
2592 #define __FUNCT__ "SNESAppendOptionsPrefix"
2593 /*@C
2594    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2595    SNES options in the database.
2596 
2597    Collective on SNES
2598 
2599    Input Parameters:
2600 +  snes - the SNES context
2601 -  prefix - the prefix to prepend to all option names
2602 
2603    Notes:
2604    A hyphen (-) must NOT be given at the beginning of the prefix name.
2605    The first character of all runtime options is AUTOMATICALLY the hyphen.
2606 
2607    Level: advanced
2608 
2609 .keywords: SNES, append, options, prefix, database
2610 
2611 .seealso: SNESGetOptionsPrefix()
2612 @*/
2613 int SNESAppendOptionsPrefix(SNES snes,char *prefix)
2614 {
2615   int ierr;
2616 
2617   PetscFunctionBegin;
2618   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2619   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2620   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2621   PetscFunctionReturn(0);
2622 }
2623 
2624 #undef __FUNCT__
2625 #define __FUNCT__ "SNESGetOptionsPrefix"
2626 /*@C
2627    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2628    SNES options in the database.
2629 
2630    Not Collective
2631 
2632    Input Parameter:
2633 .  snes - the SNES context
2634 
2635    Output Parameter:
2636 .  prefix - pointer to the prefix string used
2637 
2638    Notes: On the fortran side, the user should pass in a string 'prifix' of
2639    sufficient length to hold the prefix.
2640 
2641    Level: advanced
2642 
2643 .keywords: SNES, get, options, prefix, database
2644 
2645 .seealso: SNESAppendOptionsPrefix()
2646 @*/
2647 int SNESGetOptionsPrefix(SNES snes,char **prefix)
2648 {
2649   int ierr;
2650 
2651   PetscFunctionBegin;
2652   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2653   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2654   PetscFunctionReturn(0);
2655 }
2656 
2657 /*MC
2658    SNESRegisterDynamic - Adds a method to the nonlinear solver package.
2659 
2660    Synopsis:
2661    int SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
2662 
2663    Not collective
2664 
2665    Input Parameters:
2666 +  name_solver - name of a new user-defined solver
2667 .  path - path (either absolute or relative) the library containing this solver
2668 .  name_create - name of routine to create method context
2669 -  routine_create - routine to create method context
2670 
2671    Notes:
2672    SNESRegisterDynamic() may be called multiple times to add several user-defined solvers.
2673 
2674    If dynamic libraries are used, then the fourth input argument (routine_create)
2675    is ignored.
2676 
2677    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
2678    and others of the form ${any_environmental_variable} occuring in pathname will be
2679    replaced with appropriate values.
2680 
2681    Sample usage:
2682 .vb
2683    SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2684                 "MySolverCreate",MySolverCreate);
2685 .ve
2686 
2687    Then, your solver can be chosen with the procedural interface via
2688 $     SNESSetType(snes,"my_solver")
2689    or at runtime via the option
2690 $     -snes_type my_solver
2691 
2692    Level: advanced
2693 
2694 .keywords: SNES, nonlinear, register
2695 
2696 .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2697 M*/
2698 
2699 #undef __FUNCT__
2700 #define __FUNCT__ "SNESRegister"
2701 int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES))
2702 {
2703   char fullname[256];
2704   int  ierr;
2705 
2706   PetscFunctionBegin;
2707   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2708   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2709   PetscFunctionReturn(0);
2710 }
2711