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