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