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