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