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