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