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