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