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