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