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