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