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