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