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