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