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