xref: /petsc/src/snes/interface/snes.c (revision 639f9d9dbbc54d6ac4e42e98283c540b41bb2cee)
1 #ifndef lint
2 static char vcid[] = "$Id: snes.c,v 1.94 1996/10/15 21:58:38 balay Exp bsmith $";
3 #endif
4 
5 #include "draw.h"          /*I "draw.h"  I*/
6 #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
7 #include "src/sys/nreg.h"
8 #include "pinclude/pviewer.h"
9 #include <math.h>
10 
11 extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*);
12 extern int SNESPrintTypes_Private(MPI_Comm,char*,char*);
13 
14 /*@
15    SNESView - Prints the SNES data structure.
16 
17    Input Parameters:
18 .  SNES - the SNES context
19 .  viewer - visualization context
20 
21    Options Database Key:
22 $  -snes_view : calls SNESView() at end of SNESSolve()
23 
24    Notes:
25    The available visualization contexts include
26 $     VIEWER_STDOUT_SELF - standard output (default)
27 $     VIEWER_STDOUT_WORLD - synchronized standard
28 $       output where only the first processor opens
29 $       the file.  All other processors send their
30 $       data to the first processor to print.
31 
32    The user can open alternative vistualization contexts with
33 $    ViewerFileOpenASCII() - output to a specified file
34 
35 .keywords: SNES, view
36 
37 .seealso: ViewerFileOpenASCII()
38 @*/
39 int SNESView(SNES snes,Viewer viewer)
40 {
41   SNES_KSP_EW_ConvCtx *kctx;
42   FILE                *fd;
43   int                 ierr;
44   SLES                sles;
45   char                *method;
46   ViewerType          vtype;
47 
48   PetscValidHeaderSpecific(snes,SNES_COOKIE);
49   if (viewer) {PetscValidHeader(viewer);}
50   else { viewer = VIEWER_STDOUT_SELF; }
51 
52   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
53   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
54     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
55     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
56     SNESGetType(snes,PETSC_NULL,&method);
57     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
58     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
59     PetscFPrintf(snes->comm,fd,
60       "  maximum iterations=%d, maximum function evaluations=%d\n",
61       snes->max_its,snes->max_funcs);
62     PetscFPrintf(snes->comm,fd,
63     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
64       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
65     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
66       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
67     if (snes->ksp_ewconv) {
68       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
69       if (kctx) {
70         PetscFPrintf(snes->comm,fd,
71      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
72         kctx->version);
73         PetscFPrintf(snes->comm,fd,
74           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
75           kctx->rtol_max,kctx->threshold);
76         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
77           kctx->gamma,kctx->alpha,kctx->alpha2);
78       }
79     }
80   } else if (vtype == STRING_VIEWER) {
81     SNESGetType(snes,PETSC_NULL,&method);
82     ViewerStringSPrintf(viewer," %-3.3s",method);
83   }
84   SNESGetSLES(snes,&sles);
85   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
86   return 0;
87 }
88 
89 /*
90        We retain a list of functions that also take SNES command
91     line options. These are called at the end SNESSetFromOptions()
92 */
93 #define MAXSETFROMOPTIONS 5
94 static int numberofsetfromoptions;
95 static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
96 
97 /*@
98     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
99 
100   Input Parameter:
101 .   snescheck - function that checks for options
102 
103 .seealso: SNESSetFromOptions()
104 @*/
105 int SNESAddOptionsChecker(int (*snescheck)(SNES) )
106 {
107   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
108     SETERRQ(1,"SNESAddOptionsChecker:Too many options checkers, only 5 allowed");
109   }
110 
111   othersetfromoptions[numberofsetfromoptions++] = snescheck;
112   return 0;
113 }
114 
115 /*@
116    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
117 
118    Input Parameter:
119 .  snes - the SNES context
120 
121 .keywords: SNES, nonlinear, set, options, database
122 
123 .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
124 @*/
125 int SNESSetFromOptions(SNES snes)
126 {
127   SNESType method;
128   double   tmp;
129   SLES     sles;
130   int      ierr, flg,i;
131   int      version   = PETSC_DEFAULT;
132   double   rtol_0    = PETSC_DEFAULT;
133   double   rtol_max  = PETSC_DEFAULT;
134   double   gamma2    = PETSC_DEFAULT;
135   double   alpha     = PETSC_DEFAULT;
136   double   alpha2    = PETSC_DEFAULT;
137   double   threshold = PETSC_DEFAULT;
138 
139   PetscValidHeaderSpecific(snes,SNES_COOKIE);
140   if (snes->setup_called) SETERRQ(1,"SNESSetFromOptions:Must call prior to SNESSetUp");
141   ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr);
142   if (flg) {
143     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
144   }
145   else if (!snes->set_method_called) {
146     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
147       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
148     }
149     else {
150       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
151     }
152   }
153 
154   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
155   if (flg) { SNESPrintHelp(snes); }
156   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
157   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);}
158   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
159   if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);}
160   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
161   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); }
162   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
163   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
164   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
165   if (flg) { SNESSetTrustRegionTolerance(snes,tmp); }
166   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
167   if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); }
168   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
169   if (flg) { snes->ksp_ewconv = 1; }
170   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
171   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
172   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
173   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
174   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
175   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
176   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
177   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
178                                   alpha2,threshold); CHKERRQ(ierr);
179   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
180   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
181   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
182   if (flg) {
183    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
184   }
185   ierr = OptionsHasName(snes->prefix,"-snes_xmonitor",&flg); CHKERRQ(ierr);
186   if (flg) {
187     int    rank = 0;
188     DrawLG lg;
189     MPI_Initialized(&rank);
190     if (rank) MPI_Comm_rank(snes->comm,&rank);
191     if (!rank) {
192       ierr = SNESLGMonitorCreate(0,0,0,0,300,300,&lg); CHKERRQ(ierr);
193       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
194       snes->xmonitor = lg;
195       PLogObjectParent(snes,lg);
196     }
197   }
198   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
199   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
200     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
201                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
202     PLogInfo(snes,
203       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
204   }
205   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
206     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
207                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
208     PLogInfo(snes,
209       "SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
210   }
211 
212   for ( i=0; i<numberofsetfromoptions; i++ ) {
213     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
214   }
215 
216   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
217   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
218   if (!snes->setfromoptions) return 0;
219   return (*snes->setfromoptions)(snes);
220 }
221 
222 /*@
223    SNESPrintHelp - Prints all options for the SNES component.
224 
225    Input Parameter:
226 .  snes - the SNES context
227 
228    Options Database Keys:
229 $  -help, -h
230 
231 .keywords: SNES, nonlinear, help
232 
233 .seealso: SNESSetFromOptions()
234 @*/
235 int SNESPrintHelp(SNES snes)
236 {
237   char                p[64];
238   SNES_KSP_EW_ConvCtx *kctx;
239 
240   PetscValidHeaderSpecific(snes,SNES_COOKIE);
241 
242   PetscStrcpy(p,"-");
243   if (snes->prefix) PetscStrcat(p, snes->prefix);
244 
245   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
246 
247   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
248   SNESPrintTypes_Private(snes->comm,p,"snes_type");
249   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
250   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
251   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
252   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
253   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
254   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
255   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
256   PetscPrintf(snes->comm," SNES Monitoring Options: Choose 1 of the following\n");
257   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor\n",p);
258   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
259   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
260     PetscPrintf(snes->comm,
261      " Options for solving systems of nonlinear equations only:\n");
262     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
263     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
264     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
265     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p);
266     PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p);
267     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
268     PetscPrintf(snes->comm,
269      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
270     PetscPrintf(snes->comm,
271      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
272     PetscPrintf(snes->comm,
273      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
274     PetscPrintf(snes->comm,
275      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
276     PetscPrintf(snes->comm,
277      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
278     PetscPrintf(snes->comm,
279      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
280     PetscPrintf(snes->comm,
281      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
282   }
283   else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
284     PetscPrintf(snes->comm,
285      " Options for solving unconstrained minimization problems only:\n");
286     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
287     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
288     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
289     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p);
290      PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p);
291   }
292   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
293   PetscPrintf(snes->comm,"a particular method\n");
294   if (snes->printhelp) (*snes->printhelp)(snes,p);
295   return 0;
296 }
297 /*@
298    SNESSetApplicationContext - Sets the optional user-defined context for
299    the nonlinear solvers.
300 
301    Input Parameters:
302 .  snes - the SNES context
303 .  usrP - optional user context
304 
305 .keywords: SNES, nonlinear, set, application, context
306 
307 .seealso: SNESGetApplicationContext()
308 @*/
309 int SNESSetApplicationContext(SNES snes,void *usrP)
310 {
311   PetscValidHeaderSpecific(snes,SNES_COOKIE);
312   snes->user		= usrP;
313   return 0;
314 }
315 
316 /*@C
317    SNESGetApplicationContext - Gets the user-defined context for the
318    nonlinear solvers.
319 
320    Input Parameter:
321 .  snes - SNES context
322 
323    Output Parameter:
324 .  usrP - user context
325 
326 .keywords: SNES, nonlinear, get, application, context
327 
328 .seealso: SNESSetApplicationContext()
329 @*/
330 int SNESGetApplicationContext( SNES snes,  void **usrP )
331 {
332   PetscValidHeaderSpecific(snes,SNES_COOKIE);
333   *usrP = snes->user;
334   return 0;
335 }
336 
337 /*@
338    SNESGetIterationNumber - Gets the current iteration number of the
339    nonlinear solver.
340 
341    Input Parameter:
342 .  snes - SNES context
343 
344    Output Parameter:
345 .  iter - iteration number
346 
347 .keywords: SNES, nonlinear, get, iteration, number
348 @*/
349 int SNESGetIterationNumber(SNES snes,int* iter)
350 {
351   PetscValidHeaderSpecific(snes,SNES_COOKIE);
352   PetscValidIntPointer(iter);
353   *iter = snes->iter;
354   return 0;
355 }
356 
357 /*@
358    SNESGetFunctionNorm - Gets the norm of the current function that was set
359    with SNESSSetFunction().
360 
361    Input Parameter:
362 .  snes - SNES context
363 
364    Output Parameter:
365 .  fnorm - 2-norm of function
366 
367    Note:
368    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
369    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
370    SNESGetGradientNorm().
371 
372 .keywords: SNES, nonlinear, get, function, norm
373 
374 .seealso: SNESSetFunction()
375 @*/
376 int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
377 {
378   PetscValidHeaderSpecific(snes,SNES_COOKIE);
379   PetscValidScalarPointer(fnorm);
380   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
381     SETERRQ(1,"SNESGetFunctionNorm:For SNES_NONLINEAR_EQUATIONS only");
382   }
383   *fnorm = snes->norm;
384   return 0;
385 }
386 
387 /*@
388    SNESGetGradientNorm - Gets the norm of the current gradient that was set
389    with SNESSSetGradient().
390 
391    Input Parameter:
392 .  snes - SNES context
393 
394    Output Parameter:
395 .  fnorm - 2-norm of gradient
396 
397    Note:
398    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
399    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
400    is SNESGetFunctionNorm().
401 
402 .keywords: SNES, nonlinear, get, gradient, norm
403 
404 .seelso: SNESSetGradient()
405 @*/
406 int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
407 {
408   PetscValidHeaderSpecific(snes,SNES_COOKIE);
409   PetscValidScalarPointer(gnorm);
410   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
411     SETERRQ(1,"SNESGetGradientNorm:For SNES_UNCONSTRAINED_MINIMIZATION only");
412   }
413   *gnorm = snes->norm;
414   return 0;
415 }
416 
417 /*@
418    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
419    attempted by the nonlinear solver.
420 
421    Input Parameter:
422 .  snes - SNES context
423 
424    Output Parameter:
425 .  nfails - number of unsuccessful steps attempted
426 
427 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
428 @*/
429 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
430 {
431   PetscValidHeaderSpecific(snes,SNES_COOKIE);
432   PetscValidIntPointer(nfails);
433   *nfails = snes->nfailures;
434   return 0;
435 }
436 /*@C
437    SNESGetSLES - Returns the SLES context for a SNES solver.
438 
439    Input Parameter:
440 .  snes - the SNES context
441 
442    Output Parameter:
443 .  sles - the SLES context
444 
445    Notes:
446    The user can then directly manipulate the SLES context to set various
447    options, etc.  Likewise, the user can then extract and manipulate the
448    KSP and PC contexts as well.
449 
450 .keywords: SNES, nonlinear, get, SLES, context
451 
452 .seealso: SLESGetPC(), SLESGetKSP()
453 @*/
454 int SNESGetSLES(SNES snes,SLES *sles)
455 {
456   PetscValidHeaderSpecific(snes,SNES_COOKIE);
457   *sles = snes->sles;
458   return 0;
459 }
460 /* -----------------------------------------------------------*/
461 /*@C
462    SNESCreate - Creates a nonlinear solver context.
463 
464    Input Parameter:
465 .  comm - MPI communicator
466 .  type - type of method, one of
467 $    SNES_NONLINEAR_EQUATIONS
468 $      (for systems of nonlinear equations)
469 $    SNES_UNCONSTRAINED_MINIMIZATION
470 $      (for unconstrained minimization)
471 
472    Output Parameter:
473 .  outsnes - the new SNES context
474 
475    Options Database Key:
476 $   -snes_mf - use default matrix-free Jacobian-vector products,
477 $              and no preconditioning matrix
478 $   -snes_mf_operator - use default matrix-free Jacobian-vector
479 $             products, and a user-provided preconditioning matrix
480 $             as set by SNESSetJacobian()
481 $   -snes_fd - use (slow!) finite differences to compute Jacobian
482 
483 .keywords: SNES, nonlinear, create, context
484 
485 .seealso: SNESSolve(), SNESDestroy()
486 @*/
487 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
488 {
489   int                 ierr;
490   SNES                snes;
491   SNES_KSP_EW_ConvCtx *kctx;
492 
493   *outsnes = 0;
494   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS)
495     SETERRQ(1,"SNESCreate:incorrect method type");
496   PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm);
497   PLogObjectCreate(snes);
498   snes->max_its           = 50;
499   snes->max_funcs	  = 1000;
500   snes->norm		  = 0.0;
501   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
502     snes->rtol		  = 1.e-8;
503     snes->ttol            = 0.0;
504     snes->atol		  = 1.e-10;
505   }
506   else {
507     snes->rtol		  = 1.e-8;
508     snes->ttol            = 0.0;
509     snes->atol		  = 1.e-50;
510   }
511   snes->xtol		  = 1.e-8;
512   snes->trunctol	  = 1.e-12; /* no longer used */
513   snes->nfuncs            = 0;
514   snes->nfailures         = 0;
515   snes->numbermonitors    = 0;
516   snes->data              = 0;
517   snes->view              = 0;
518   snes->computeumfunction = 0;
519   snes->umfunP            = 0;
520   snes->fc                = 0;
521   snes->deltatol          = 1.e-12;
522   snes->fmin              = -1.e30;
523   snes->method_class      = type;
524   snes->set_method_called = 0;
525   snes->setup_called      = 0;
526   snes->ksp_ewconv        = 0;
527   snes->type              = -1;
528   snes->xmonitor          = 0;
529   snes->vwork             = 0;
530   snes->nwork             = 0;
531 
532   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
533   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
534   snes->kspconvctx  = (void*)kctx;
535   kctx->version     = 2;
536   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
537                              this was too large for some test cases */
538   kctx->rtol_last   = 0;
539   kctx->rtol_max    = .9;
540   kctx->gamma       = 1.0;
541   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
542   kctx->alpha       = kctx->alpha2;
543   kctx->threshold   = .1;
544   kctx->lresid_last = 0;
545   kctx->norm_last   = 0;
546 
547   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
548   PLogObjectParent(snes,snes->sles)
549 
550   *outsnes = snes;
551   return 0;
552 }
553 
554 /* --------------------------------------------------------------- */
555 /*@C
556    SNESSetFunction - Sets the function evaluation routine and function
557    vector for use by the SNES routines in solving systems of nonlinear
558    equations.
559 
560    Input Parameters:
561 .  snes - the SNES context
562 .  func - function evaluation routine
563 .  r - vector to store function value
564 .  ctx - [optional] user-defined context for private data for the
565          function evaluation routine (may be PETSC_NULL)
566 
567    Calling sequence of func:
568 .  func (SNES snes,Vec x,Vec f,void *ctx);
569 
570 .  x - input vector
571 .  f - function vector
572 .  ctx - optional user-defined function context
573 
574    Notes:
575    The Newton-like methods typically solve linear systems of the form
576 $      f'(x) x = -f(x),
577 $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
578 
579    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
580    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
581    SNESSetMinimizationFunction() and SNESSetGradient();
582 
583 .keywords: SNES, nonlinear, set, function
584 
585 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
586 @*/
587 int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
588 {
589   PetscValidHeaderSpecific(snes,SNES_COOKIE);
590   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
591     "SNESSetFunction:For SNES_NONLINEAR_EQUATIONS only");
592   snes->computefunction     = func;
593   snes->vec_func            = snes->vec_func_always = r;
594   snes->funP                = ctx;
595   return 0;
596 }
597 
598 /*@
599    SNESComputeFunction - Computes the function that has been set with
600    SNESSetFunction().
601 
602    Input Parameters:
603 .  snes - the SNES context
604 .  x - input vector
605 
606    Output Parameter:
607 .  y - function vector, as set by SNESSetFunction()
608 
609    Notes:
610    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
611    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
612    SNESComputeMinimizationFunction() and SNESComputeGradient();
613 
614 .keywords: SNES, nonlinear, compute, function
615 
616 .seealso: SNESSetFunction(), SNESGetFunction()
617 @*/
618 int SNESComputeFunction(SNES snes,Vec x, Vec y)
619 {
620   int    ierr;
621 
622   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
623     SETERRQ(1,"SNESComputeFunction: For SNES_NONLINEAR_EQUATIONS only");
624   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
625   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
626   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
627   return 0;
628 }
629 
630 /*@C
631    SNESSetMinimizationFunction - Sets the function evaluation routine for
632    unconstrained minimization.
633 
634    Input Parameters:
635 .  snes - the SNES context
636 .  func - function evaluation routine
637 .  ctx - [optional] user-defined context for private data for the
638          function evaluation routine (may be PETSC_NULL)
639 
640    Calling sequence of func:
641 .  func (SNES snes,Vec x,double *f,void *ctx);
642 
643 .  x - input vector
644 .  f - function
645 .  ctx - [optional] user-defined function context
646 
647    Notes:
648    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
649    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
650    SNESSetFunction().
651 
652 .keywords: SNES, nonlinear, set, minimization, function
653 
654 .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
655            SNESSetHessian(), SNESSetGradient()
656 @*/
657 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
658                       void *ctx)
659 {
660   PetscValidHeaderSpecific(snes,SNES_COOKIE);
661   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
662     "SNESSetMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION");
663   snes->computeumfunction   = func;
664   snes->umfunP              = ctx;
665   return 0;
666 }
667 
668 /*@
669    SNESComputeMinimizationFunction - Computes the function that has been
670    set with SNESSetMinimizationFunction().
671 
672    Input Parameters:
673 .  snes - the SNES context
674 .  x - input vector
675 
676    Output Parameter:
677 .  y - function value
678 
679    Notes:
680    SNESComputeMinimizationFunction() is valid only for
681    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
682    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
683 
684 .keywords: SNES, nonlinear, compute, minimization, function
685 
686 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
687           SNESComputeGradient(), SNESComputeHessian()
688 @*/
689 int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
690 {
691   int    ierr;
692   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
693     "SNESComputeMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION");
694   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
695   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
696   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
697   return 0;
698 }
699 
700 /*@C
701    SNESSetGradient - Sets the gradient evaluation routine and gradient
702    vector for use by the SNES routines.
703 
704    Input Parameters:
705 .  snes - the SNES context
706 .  func - function evaluation routine
707 .  ctx - optional user-defined context for private data for the
708          gradient evaluation routine (may be PETSC_NULL)
709 .  r - vector to store gradient value
710 
711    Calling sequence of func:
712 .  func (SNES, Vec x, Vec g, void *ctx);
713 
714 .  x - input vector
715 .  g - gradient vector
716 .  ctx - optional user-defined gradient context
717 
718    Notes:
719    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
720    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
721    SNESSetFunction().
722 
723 .keywords: SNES, nonlinear, set, function
724 
725 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
726           SNESSetMinimizationFunction(),
727 @*/
728 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
729 {
730   PetscValidHeaderSpecific(snes,SNES_COOKIE);
731   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
732     "SNESSetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
733   snes->computefunction     = func;
734   snes->vec_func            = snes->vec_func_always = r;
735   snes->funP                = ctx;
736   return 0;
737 }
738 
739 /*@
740    SNESComputeGradient - Computes the gradient that has been set with
741    SNESSetGradient().
742 
743    Input Parameters:
744 .  snes - the SNES context
745 .  x - input vector
746 
747    Output Parameter:
748 .  y - gradient vector
749 
750    Notes:
751    SNESComputeGradient() is valid only for
752    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
753    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
754 
755 .keywords: SNES, nonlinear, compute, gradient
756 
757 .seealso:  SNESSetGradient(), SNESGetGradient(),
758            SNESComputeMinimizationFunction(), SNESComputeHessian()
759 @*/
760 int SNESComputeGradient(SNES snes,Vec x, Vec y)
761 {
762   int    ierr;
763   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
764     SETERRQ(1,"SNESComputeGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
765   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
766   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
767   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
768   return 0;
769 }
770 
771 /*@
772    SNESComputeJacobian - Computes the Jacobian matrix that has been
773    set with SNESSetJacobian().
774 
775    Input Parameters:
776 .  snes - the SNES context
777 .  x - input vector
778 
779    Output Parameters:
780 .  A - Jacobian matrix
781 .  B - optional preconditioning matrix
782 .  flag - flag indicating matrix structure
783 
784    Notes:
785    Most users should not need to explicitly call this routine, as it
786    is used internally within the nonlinear solvers.
787 
788    See SLESSetOperators() for important information about setting the
789    flag parameter.
790 
791    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
792    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
793    methods is SNESComputeJacobian().
794 
795 .keywords: SNES, compute, Jacobian, matrix
796 
797 .seealso:  SNESSetJacobian(), SLESSetOperators()
798 @*/
799 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
800 {
801   int    ierr;
802   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
803     SETERRQ(1,"SNESComputeJacobian: For SNES_NONLINEAR_EQUATIONS only");
804   if (!snes->computejacobian) return 0;
805   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
806   *flg = DIFFERENT_NONZERO_PATTERN;
807   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
808   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
809   /* make sure user returned a correct Jacobian and preconditioner */
810   PetscValidHeaderSpecific(*A,MAT_COOKIE);
811   PetscValidHeaderSpecific(*B,MAT_COOKIE);
812   return 0;
813 }
814 
815 /*@
816    SNESComputeHessian - Computes the Hessian matrix that has been
817    set with SNESSetHessian().
818 
819    Input Parameters:
820 .  snes - the SNES context
821 .  x - input vector
822 
823    Output Parameters:
824 .  A - Hessian matrix
825 .  B - optional preconditioning matrix
826 .  flag - flag indicating matrix structure
827 
828    Notes:
829    Most users should not need to explicitly call this routine, as it
830    is used internally within the nonlinear solvers.
831 
832    See SLESSetOperators() for important information about setting the
833    flag parameter.
834 
835    SNESComputeHessian() is valid only for
836    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
837    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
838 
839 .keywords: SNES, compute, Hessian, matrix
840 
841 .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
842            SNESComputeMinimizationFunction()
843 @*/
844 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
845 {
846   int    ierr;
847   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
848     SETERRQ(1,"SNESComputeHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
849   if (!snes->computejacobian) return 0;
850   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
851   *flag = DIFFERENT_NONZERO_PATTERN;
852   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
853   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
854   /* make sure user returned a correct Jacobian and preconditioner */
855   PetscValidHeaderSpecific(*A,MAT_COOKIE);
856   PetscValidHeaderSpecific(*B,MAT_COOKIE);
857   return 0;
858 }
859 
860 /*@C
861    SNESSetJacobian - Sets the function to compute Jacobian as well as the
862    location to store the matrix.
863 
864    Input Parameters:
865 .  snes - the SNES context
866 .  A - Jacobian matrix
867 .  B - preconditioner matrix (usually same as the Jacobian)
868 .  func - Jacobian evaluation routine
869 .  ctx - [optional] user-defined context for private data for the
870          Jacobian evaluation routine (may be PETSC_NULL)
871 
872    Calling sequence of func:
873 .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
874 
875 .  x - input vector
876 .  A - Jacobian matrix
877 .  B - preconditioner matrix, usually the same as A
878 .  flag - flag indicating information about the preconditioner matrix
879    structure (same as flag in SLESSetOperators())
880 .  ctx - [optional] user-defined Jacobian context
881 
882    Notes:
883    See SLESSetOperators() for important information about setting the flag
884    output parameter in the routine func().  Be sure to read this information!
885 
886    The routine func() takes Mat * as the matrix arguments rather than Mat.
887    This allows the Jacobian evaluation routine to replace A and/or B with a
888    completely new new matrix structure (not just different matrix elements)
889    when appropriate, for instance, if the nonzero structure is changing
890    throughout the global iterations.
891 
892 .keywords: SNES, nonlinear, set, Jacobian, matrix
893 
894 .seealso: SLESSetOperators(), SNESSetFunction()
895 @*/
896 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
897                     MatStructure*,void*),void *ctx)
898 {
899   PetscValidHeaderSpecific(snes,SNES_COOKIE);
900   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
901     SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only");
902   snes->computejacobian = func;
903   snes->jacP            = ctx;
904   snes->jacobian        = A;
905   snes->jacobian_pre    = B;
906   return 0;
907 }
908 
909 /*@
910    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
911    provided context for evaluating the Jacobian.
912 
913    Input Parameter:
914 .  snes - the nonlinear solver context
915 
916    Output Parameters:
917 .  A - location to stash Jacobian matrix (or PETSC_NULL)
918 .  B - location to stash preconditioner matrix (or PETSC_NULL)
919 .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
920 
921 .seealso: SNESSetJacobian(), SNESComputeJacobian()
922 @*/
923 int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
924 {
925   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
926     SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only");
927   if (A)   *A = snes->jacobian;
928   if (B)   *B = snes->jacobian_pre;
929   if (ctx) *ctx = snes->jacP;
930   return 0;
931 }
932 
933 /*@C
934    SNESSetHessian - Sets the function to compute Hessian as well as the
935    location to store the matrix.
936 
937    Input Parameters:
938 .  snes - the SNES context
939 .  A - Hessian matrix
940 .  B - preconditioner matrix (usually same as the Hessian)
941 .  func - Jacobian evaluation routine
942 .  ctx - [optional] user-defined context for private data for the
943          Hessian evaluation routine (may be PETSC_NULL)
944 
945    Calling sequence of func:
946 .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
947 
948 .  x - input vector
949 .  A - Hessian matrix
950 .  B - preconditioner matrix, usually the same as A
951 .  flag - flag indicating information about the preconditioner matrix
952    structure (same as flag in SLESSetOperators())
953 .  ctx - [optional] user-defined Hessian context
954 
955    Notes:
956    See SLESSetOperators() for important information about setting the flag
957    output parameter in the routine func().  Be sure to read this information!
958 
959    The function func() takes Mat * as the matrix arguments rather than Mat.
960    This allows the Hessian evaluation routine to replace A and/or B with a
961    completely new new matrix structure (not just different matrix elements)
962    when appropriate, for instance, if the nonzero structure is changing
963    throughout the global iterations.
964 
965 .keywords: SNES, nonlinear, set, Hessian, matrix
966 
967 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
968 @*/
969 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
970                     MatStructure*,void*),void *ctx)
971 {
972   PetscValidHeaderSpecific(snes,SNES_COOKIE);
973   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
974     SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
975   snes->computejacobian = func;
976   snes->jacP            = ctx;
977   snes->jacobian        = A;
978   snes->jacobian_pre    = B;
979   return 0;
980 }
981 
982 /*@
983    SNESGetHessian - Returns the Hessian matrix and optionally the user
984    provided context for evaluating the Hessian.
985 
986    Input Parameter:
987 .  snes - the nonlinear solver context
988 
989    Output Parameters:
990 .  A - location to stash Hessian matrix (or PETSC_NULL)
991 .  B - location to stash preconditioner matrix (or PETSC_NULL)
992 .  ctx - location to stash Hessian ctx (or PETSC_NULL)
993 
994 .seealso: SNESSetHessian(), SNESComputeHessian()
995 @*/
996 int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
997 {
998   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
999     SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
1000   if (A)   *A = snes->jacobian;
1001   if (B)   *B = snes->jacobian_pre;
1002   if (ctx) *ctx = snes->jacP;
1003   return 0;
1004 }
1005 
1006 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1007 
1008 /*@
1009    SNESSetUp - Sets up the internal data structures for the later use
1010    of a nonlinear solver.
1011 
1012    Input Parameter:
1013 .  snes - the SNES context
1014 .  x - the solution vector
1015 
1016    Notes:
1017    For basic use of the SNES solvers the user need not explicitly call
1018    SNESSetUp(), since these actions will automatically occur during
1019    the call to SNESSolve().  However, if one wishes to control this
1020    phase separately, SNESSetUp() should be called after SNESCreate()
1021    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1022 
1023 .keywords: SNES, nonlinear, setup
1024 
1025 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1026 @*/
1027 int SNESSetUp(SNES snes,Vec x)
1028 {
1029   int ierr, flg;
1030   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1031   PetscValidHeaderSpecific(x,VEC_COOKIE);
1032   snes->vec_sol = snes->vec_sol_always = x;
1033 
1034   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1035   /*
1036       This version replaces the user provided Jacobian matrix with a
1037       matrix-free version but still employs the user-provided preconditioner matrix
1038   */
1039   if (flg) {
1040     Mat J;
1041     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1042     PLogObjectParent(snes,J);
1043     snes->mfshell = J;
1044     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1045       snes->jacobian = J;
1046       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1047     }
1048     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1049       snes->jacobian = J;
1050       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1051     } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free operator option");
1052   }
1053   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1054   /*
1055       This version replaces both the user-provided Jacobian and the user-
1056       provided preconditioner matrix with the default matrix free version.
1057    */
1058   if (flg) {
1059     Mat J;
1060     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1061     PLogObjectParent(snes,J);
1062     snes->mfshell = J;
1063     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1064       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
1065       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1066     }
1067     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1068       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
1069       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1070     } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free option");
1071   }
1072   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1073     if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first");
1074     if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first");
1075     if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetJacobian() first");
1076     if (snes->vec_func == snes->vec_sol) SETERRQ(1,"SNESSetUp:Solution vector cannot be function vector");
1077 
1078     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
1079     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1080       SLES sles; KSP ksp;
1081       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1082       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1083       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1084              (void *)snes); CHKERRQ(ierr);
1085     }
1086   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1087     if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first");
1088     if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first");
1089     if (!snes->computeumfunction)
1090       SETERRQ(1,"SNESSetUp:Must call SNESSetMinimizationFunction() first");
1091     if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetHessian() first");
1092   } else SETERRQ(1,"SNESSetUp:Unknown method class");
1093   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1094   snes->setup_called = 1;
1095   return 0;
1096 }
1097 
1098 /*@C
1099    SNESDestroy - Destroys the nonlinear solver context that was created
1100    with SNESCreate().
1101 
1102    Input Parameter:
1103 .  snes - the SNES context
1104 
1105 .keywords: SNES, nonlinear, destroy
1106 
1107 .seealso: SNESCreate(), SNESSolve()
1108 @*/
1109 int SNESDestroy(SNES snes)
1110 {
1111   int ierr;
1112   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1113   ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);
1114   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
1115   if (snes->mfshell) MatDestroy(snes->mfshell);
1116   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
1117   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
1118   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
1119   PLogObjectDestroy((PetscObject)snes);
1120   PetscHeaderDestroy((PetscObject)snes);
1121   return 0;
1122 }
1123 
1124 /* ----------- Routines to set solver parameters ---------- */
1125 
1126 /*@
1127    SNESSetTolerances - Sets various parameters used in convergence tests.
1128 
1129    Input Parameters:
1130 .  snes - the SNES context
1131 .  atol - absolute convergence tolerance
1132 .  rtol - relative convergence tolerance
1133 .  stol -  convergence tolerance in terms of the norm
1134            of the change in the solution between steps
1135 .  maxit - maximum number of iterations
1136 .  maxf - maximum number of function evaluations
1137 
1138    Options Database Keys:
1139 $    -snes_atol <atol>
1140 $    -snes_rtol <rtol>
1141 $    -snes_stol <stol>
1142 $    -snes_max_it <maxit>
1143 $    -snes_max_funcs <maxf>
1144 
1145    Notes:
1146    The default maximum number of iterations is 50.
1147    The default maximum number of function evaluations is 1000.
1148 
1149 .keywords: SNES, nonlinear, set, convergence, tolerances
1150 
1151 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
1152 @*/
1153 int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
1154 {
1155   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1156   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1157   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1158   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1159   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1160   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1161   return 0;
1162 }
1163 
1164 /*@
1165    SNESGetTolerances - Gets various parameters used in convergence tests.
1166 
1167    Input Parameters:
1168 .  snes - the SNES context
1169 .  atol - absolute convergence tolerance
1170 .  rtol - relative convergence tolerance
1171 .  stol -  convergence tolerance in terms of the norm
1172            of the change in the solution between steps
1173 .  maxit - maximum number of iterations
1174 .  maxf - maximum number of function evaluations
1175 
1176    Notes:
1177    The user can specify PETSC_NULL for any parameter that is not needed.
1178 
1179 .keywords: SNES, nonlinear, get, convergence, tolerances
1180 
1181 .seealso: SNESSetTolerances()
1182 @*/
1183 int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
1184 {
1185   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1186   if (atol)  *atol  = snes->atol;
1187   if (rtol)  *rtol  = snes->rtol;
1188   if (stol)  *stol  = snes->xtol;
1189   if (maxit) *maxit = snes->max_its;
1190   if (maxf)  *maxf  = snes->max_funcs;
1191   return 0;
1192 }
1193 
1194 /*@
1195    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1196 
1197    Input Parameters:
1198 .  snes - the SNES context
1199 .  tol - tolerance
1200 
1201    Options Database Key:
1202 $    -snes_trtol <tol>
1203 
1204 .keywords: SNES, nonlinear, set, trust region, tolerance
1205 
1206 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
1207 @*/
1208 int SNESSetTrustRegionTolerance(SNES snes,double tol)
1209 {
1210   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1211   snes->deltatol = tol;
1212   return 0;
1213 }
1214 
1215 /*@
1216    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
1217    for unconstrained minimization solvers.
1218 
1219    Input Parameters:
1220 .  snes - the SNES context
1221 .  ftol - minimum function tolerance
1222 
1223    Options Database Key:
1224 $    -snes_fmin <ftol>
1225 
1226    Note:
1227    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
1228    methods only.
1229 
1230 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
1231 
1232 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
1233 @*/
1234 int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
1235 {
1236   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1237   snes->fmin = ftol;
1238   return 0;
1239 }
1240 
1241 /* ------------ Routines to set performance monitoring options ----------- */
1242 
1243 /*@C
1244    SNESSetMonitor - Sets the function that is to be used at every
1245    iteration of the nonlinear solver to display the iteration's
1246    progress.
1247 
1248    Input Parameters:
1249 .  snes - the SNES context
1250 .  func - monitoring routine
1251 .  mctx - [optional] user-defined context for private data for the
1252           monitor routine (may be PETSC_NULL)
1253 
1254    Calling sequence of func:
1255    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
1256 
1257 $    snes - the SNES context
1258 $    its - iteration number
1259 $    mctx - [optional] monitoring context
1260 $
1261 $ SNES_NONLINEAR_EQUATIONS methods:
1262 $    norm - 2-norm function value (may be estimated)
1263 $
1264 $ SNES_UNCONSTRAINED_MINIMIZATION methods:
1265 $    norm - 2-norm gradient value (may be estimated)
1266 
1267    Notes:
1268    Several different monitor routines may be set and all will be called.
1269 
1270 .keywords: SNES, nonlinear, set, monitor
1271 
1272 .seealso: SNESDefaultMonitor()
1273 @*/
1274 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
1275 {
1276   if (!func) {
1277     snes->numbermonitors = 0;
1278     return 0;
1279   }
1280   if (snes->numbermonitors >= MAXSNESMONITORS) {
1281     SETERRQ(1,"SNESSetMonitor:Too many monitors set");
1282   }
1283 
1284   snes->monitor[snes->numbermonitors]           = func;
1285   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1286   return 0;
1287 }
1288 
1289 /*@C
1290    SNESSetConvergenceTest - Sets the function that is to be used
1291    to test for convergence of the nonlinear iterative solution.
1292 
1293    Input Parameters:
1294 .  snes - the SNES context
1295 .  func - routine to test for convergence
1296 .  cctx - [optional] context for private data for the convergence routine
1297           (may be PETSC_NULL)
1298 
1299    Calling sequence of func:
1300    int func (SNES snes,double xnorm,double gnorm,
1301              double f,void *cctx)
1302 
1303 $    snes - the SNES context
1304 $    cctx - [optional] convergence context
1305 $    xnorm - 2-norm of current iterate
1306 $
1307 $ SNES_NONLINEAR_EQUATIONS methods:
1308 $    gnorm - 2-norm of current step
1309 $    f - 2-norm of function
1310 $
1311 $ SNES_UNCONSTRAINED_MINIMIZATION methods:
1312 $    gnorm - 2-norm of current gradient
1313 $    f - function value
1314 
1315 .keywords: SNES, nonlinear, set, convergence, test
1316 
1317 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1318           SNESConverged_UM_LS(), SNESConverged_UM_TR()
1319 @*/
1320 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
1321 {
1322   (snes)->converged = func;
1323   (snes)->cnvP      = cctx;
1324   return 0;
1325 }
1326 
1327 /*
1328    SNESScaleStep_Private - Scales a step so that its length is less than the
1329    positive parameter delta.
1330 
1331     Input Parameters:
1332 .   snes - the SNES context
1333 .   y - approximate solution of linear system
1334 .   fnorm - 2-norm of current function
1335 .   delta - trust region size
1336 
1337     Output Parameters:
1338 .   gpnorm - predicted function norm at the new point, assuming local
1339     linearization.  The value is zero if the step lies within the trust
1340     region, and exceeds zero otherwise.
1341 .   ynorm - 2-norm of the step
1342 
1343     Note:
1344     For non-trust region methods such as SNES_EQ_LS, the parameter delta
1345     is set to be the maximum allowable step size.
1346 
1347 .keywords: SNES, nonlinear, scale, step
1348 */
1349 int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
1350                           double *gpnorm,double *ynorm)
1351 {
1352   double norm;
1353   Scalar cnorm;
1354   VecNorm(y,NORM_2, &norm );
1355   if (norm > *delta) {
1356      norm = *delta/norm;
1357      *gpnorm = (1.0 - norm)*(*fnorm);
1358      cnorm = norm;
1359      VecScale( &cnorm, y );
1360      *ynorm = *delta;
1361   } else {
1362      *gpnorm = 0.0;
1363      *ynorm = norm;
1364   }
1365   return 0;
1366 }
1367 
1368 /*@
1369    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
1370    SNESCreate() and optional routines of the form SNESSetXXX().
1371 
1372    Input Parameter:
1373 .  snes - the SNES context
1374 .  x - the solution vector
1375 
1376    Output Parameter:
1377    its - number of iterations until termination
1378 
1379    Note:
1380    The user should initialize the vector, x, with the initial guess
1381    for the nonlinear solve prior to calling SNESSolve.  In particular,
1382    to employ an initial guess of zero, the user should explicitly set
1383    this vector to zero by calling VecSet().
1384 
1385 .keywords: SNES, nonlinear, solve
1386 
1387 .seealso: SNESCreate(), SNESDestroy()
1388 @*/
1389 int SNESSolve(SNES snes,Vec x,int *its)
1390 {
1391   int ierr, flg;
1392 
1393   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1394   PetscValidIntPointer(its);
1395   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1396   else {snes->vec_sol = snes->vec_sol_always = x;}
1397   PLogEventBegin(SNES_Solve,snes,0,0,0);
1398   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
1399   PLogEventEnd(SNES_Solve,snes,0,0,0);
1400   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
1401   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
1402   return 0;
1403 }
1404 
1405 /* --------- Internal routines for SNES Package --------- */
1406 static NRList *__SNESList = 0;
1407 
1408 /*@
1409    SNESSetType - Sets the method for the nonlinear solver.
1410 
1411    Input Parameters:
1412 .  snes - the SNES context
1413 .  method - a known method
1414 
1415    Notes:
1416    See "petsc/include/snes.h" for available methods (for instance)
1417 $  Systems of nonlinear equations:
1418 $    SNES_EQ_LS - Newton's method with line search
1419 $    SNES_EQ_TR - Newton's method with trust region
1420 $  Unconstrained minimization:
1421 $    SNES_UM_TR - Newton's method with trust region
1422 $    SNES_UM_LS - Newton's method with line search
1423 
1424   Options Database Command:
1425 $ -snes_type  <method>
1426 $    Use -help for a list of available methods
1427 $    (for instance, ls or tr)
1428 
1429 .keywords: SNES, set, method
1430 @*/
1431 int SNESSetType(SNES snes,SNESType method)
1432 {
1433   int (*r)(SNES);
1434   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1435   if (snes->type == (int) method) return 0;
1436 
1437   /* Get the function pointers for the iterative method requested */
1438   if (!__SNESList) {SNESRegisterAll();}
1439   if (!__SNESList) {SETERRQ(1,"SNESSetType:Could not get methods");}
1440   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1441   if (!r) {SETERRQ(1,"SNESSetType:Unknown method");}
1442   if (snes->data) PetscFree(snes->data);
1443   snes->set_method_called = 1;
1444   return (*r)(snes);
1445 }
1446 
1447 /* --------------------------------------------------------------------- */
1448 /*@C
1449    SNESRegister - Adds the method to the nonlinear solver package, given
1450    a function pointer and a nonlinear solver name of the type SNESType.
1451 
1452    Input Parameters:
1453 .  name - for instance SNES_EQ_LS, SNES_EQ_TR, ...
1454 .  sname - corresponding string for name
1455 .  create - routine to create method context
1456 
1457 .keywords: SNES, nonlinear, register
1458 
1459 .seealso: SNESRegisterAll(), SNESRegisterDestroy()
1460 @*/
1461 int SNESRegister(int name, char *sname, int (*create)(SNES))
1462 {
1463   int ierr;
1464   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
1465   NRRegister( __SNESList, name, sname, (int (*)(void*))create );
1466   return 0;
1467 }
1468 /* --------------------------------------------------------------------- */
1469 /*@C
1470    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1471    registered by SNESRegister().
1472 
1473 .keywords: SNES, nonlinear, register, destroy
1474 
1475 .seealso: SNESRegisterAll(), SNESRegisterAll()
1476 @*/
1477 int SNESRegisterDestroy()
1478 {
1479   if (__SNESList) {
1480     NRDestroy( __SNESList );
1481     __SNESList = 0;
1482   }
1483   return 0;
1484 }
1485 
1486 /*
1487    SNESGetTypeFromOptions_Private - Sets the selected method from the
1488    options database.
1489 
1490    Input Parameter:
1491 .  ctx - the SNES context
1492 
1493    Output Parameter:
1494 .  method -  solver method
1495 
1496    Returns:
1497    Returns 1 if the method is found; 0 otherwise.
1498 
1499    Options Database Key:
1500 $  -snes_type  method
1501 */
1502 int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
1503 {
1504   int  ierr;
1505   char sbuf[50];
1506 
1507   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1508   if (*flg) {
1509     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
1510     *method = (SNESType)NRFindID( __SNESList, sbuf );
1511   }
1512   return 0;
1513 }
1514 
1515 /*@C
1516    SNESGetType - Gets the SNES method type and name (as a string).
1517 
1518    Input Parameter:
1519 .  snes - nonlinear solver context
1520 
1521    Output Parameter:
1522 .  method - SNES method (or use PETSC_NULL)
1523 .  name - name of SNES method (or use PETSC_NULL)
1524 
1525 .keywords: SNES, nonlinear, get, method, name
1526 @*/
1527 int SNESGetType(SNES snes, SNESType *method,char **name)
1528 {
1529   int ierr;
1530   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
1531   if (method) *method = (SNESType) snes->type;
1532   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
1533   return 0;
1534 }
1535 
1536 #include <stdio.h>
1537 /*
1538    SNESPrintTypes_Private - Prints the SNES methods available from the
1539    options database.
1540 
1541    Input Parameters:
1542 .  comm   - communicator (usually MPI_COMM_WORLD)
1543 .  prefix - prefix (usually "-")
1544 .  name   - the options database name (by default "snes_type")
1545 */
1546 int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
1547 {
1548   FuncList *entry;
1549   if (!__SNESList) {SNESRegisterAll();}
1550   entry = __SNESList->head;
1551   PetscPrintf(comm," %s%s (one of)",prefix,name);
1552   while (entry) {
1553     PetscPrintf(comm," %s",entry->name);
1554     entry = entry->next;
1555   }
1556   PetscPrintf(comm,"\n");
1557   return 0;
1558 }
1559 
1560 /*@C
1561    SNESGetSolution - Returns the vector where the approximate solution is
1562    stored.
1563 
1564    Input Parameter:
1565 .  snes - the SNES context
1566 
1567    Output Parameter:
1568 .  x - the solution
1569 
1570 .keywords: SNES, nonlinear, get, solution
1571 
1572 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
1573 @*/
1574 int SNESGetSolution(SNES snes,Vec *x)
1575 {
1576   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1577   *x = snes->vec_sol_always;
1578   return 0;
1579 }
1580 
1581 /*@C
1582    SNESGetSolutionUpdate - Returns the vector where the solution update is
1583    stored.
1584 
1585    Input Parameter:
1586 .  snes - the SNES context
1587 
1588    Output Parameter:
1589 .  x - the solution update
1590 
1591    Notes:
1592    This vector is implementation dependent.
1593 
1594 .keywords: SNES, nonlinear, get, solution, update
1595 
1596 .seealso: SNESGetSolution(), SNESGetFunction
1597 @*/
1598 int SNESGetSolutionUpdate(SNES snes,Vec *x)
1599 {
1600   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1601   *x = snes->vec_sol_update_always;
1602   return 0;
1603 }
1604 
1605 /*@C
1606    SNESGetFunction - Returns the vector where the function is stored.
1607 
1608    Input Parameter:
1609 .  snes - the SNES context
1610 
1611    Output Parameter:
1612 .  r - the function
1613 
1614    Notes:
1615    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
1616    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
1617    SNESGetMinimizationFunction() and SNESGetGradient();
1618 
1619 .keywords: SNES, nonlinear, get, function
1620 
1621 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
1622           SNESGetGradient()
1623 @*/
1624 int SNESGetFunction(SNES snes,Vec *r)
1625 {
1626   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1627   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
1628     "SNESGetFunction:For SNES_NONLINEAR_EQUATIONS only");
1629   *r = snes->vec_func_always;
1630   return 0;
1631 }
1632 
1633 /*@C
1634    SNESGetGradient - Returns the vector where the gradient is stored.
1635 
1636    Input Parameter:
1637 .  snes - the SNES context
1638 
1639    Output Parameter:
1640 .  r - the gradient
1641 
1642    Notes:
1643    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
1644    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
1645    SNESGetFunction().
1646 
1647 .keywords: SNES, nonlinear, get, gradient
1648 
1649 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
1650 @*/
1651 int SNESGetGradient(SNES snes,Vec *r)
1652 {
1653   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1654   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
1655     "SNESGetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
1656   *r = snes->vec_func_always;
1657   return 0;
1658 }
1659 
1660 /*@
1661    SNESGetMinimizationFunction - Returns the scalar function value for
1662    unconstrained minimization problems.
1663 
1664    Input Parameter:
1665 .  snes - the SNES context
1666 
1667    Output Parameter:
1668 .  r - the function
1669 
1670    Notes:
1671    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
1672    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
1673    SNESGetFunction().
1674 
1675 .keywords: SNES, nonlinear, get, function
1676 
1677 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
1678 @*/
1679 int SNESGetMinimizationFunction(SNES snes,double *r)
1680 {
1681   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1682   PetscValidScalarPointer(r);
1683   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
1684     "SNESGetMinimizationFunction:For SNES_UNCONSTRAINED_MINIMIZATION only");
1685   *r = snes->fc;
1686   return 0;
1687 }
1688 
1689 /*@C
1690    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1691    SNES options in the database. You must NOT include the - at the beginning of
1692    the prefix name.
1693 
1694    Input Parameter:
1695 .  snes - the SNES context
1696 .  prefix - the prefix to prepend to all option names
1697 
1698 .keywords: SNES, set, options, prefix, database
1699 
1700 .seealso: SNESSetFromOptions()
1701 @*/
1702 int SNESSetOptionsPrefix(SNES snes,char *prefix)
1703 {
1704   int ierr;
1705 
1706   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1707   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
1708   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
1709   return 0;
1710 }
1711 
1712 /*@C
1713    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1714    SNES options in the database. You must NOT include the - at the beginning of
1715    the prefix name.
1716 
1717    Input Parameter:
1718 .  snes - the SNES context
1719 .  prefix - the prefix to prepend to all option names
1720 
1721 .keywords: SNES, append, options, prefix, database
1722 
1723 .seealso: SNESGetOptionsPrefix()
1724 @*/
1725 int SNESAppendOptionsPrefix(SNES snes,char *prefix)
1726 {
1727   int ierr;
1728 
1729   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1730   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
1731   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
1732   return 0;
1733 }
1734 
1735 /*@
1736    SNESGetOptionsPrefix - Sets the prefix used for searching for all
1737    SNES options in the database.
1738 
1739    Input Parameter:
1740 .  snes - the SNES context
1741 
1742    Output Parameter:
1743 .  prefix - pointer to the prefix string used
1744 
1745 .keywords: SNES, get, options, prefix, database
1746 
1747 .seealso: SNESAppendOptionsPrefix()
1748 @*/
1749 int SNESGetOptionsPrefix(SNES snes,char **prefix)
1750 {
1751   int ierr;
1752 
1753   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1754   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
1755   return 0;
1756 }
1757 
1758 
1759 
1760 
1761 
1762