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