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