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