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