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