xref: /petsc/src/snes/interface/snes.c (revision 5a5d4f665b3034485fbcdfbe42345f740a4bd170)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: snes.c,v 1.132 1997/10/28 14:24:42 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(1,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(1,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) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
610   snes->computefunction     = func;
611   snes->vec_func            = snes->vec_func_always = r;
612   snes->funP                = ctx;
613   PetscFunctionReturn(0);
614 }
615 
616 #undef __FUNC__
617 #define __FUNC__ "SNESComputeFunction"
618 /*@
619    SNESComputeFunction - Computes the function that has been set with
620    SNESSetFunction().
621 
622    Input Parameters:
623 .  snes - the SNES context
624 .  x - input vector
625 
626    Output Parameter:
627 .  y - function vector, as set by SNESSetFunction()
628 
629    Notes:
630    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
631    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
632    SNESComputeMinimizationFunction() and SNESComputeGradient();
633 
634 .keywords: SNES, nonlinear, compute, function
635 
636 .seealso: SNESSetFunction(), SNESGetFunction()
637 @*/
638 int SNESComputeFunction(SNES snes,Vec x, Vec y)
639 {
640   int    ierr;
641 
642   PetscFunctionBegin;
643   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
644     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
645   }
646   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
647   PetscStackPush("SNES user function");
648   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
649   PetscStackPop;
650   snes->nfuncs++;
651   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
652   PetscFunctionReturn(0);
653 }
654 
655 #undef __FUNC__
656 #define __FUNC__ "SNESSetMinimizationFunction"
657 /*@C
658    SNESSetMinimizationFunction - Sets the function evaluation routine for
659    unconstrained minimization.
660 
661    Input Parameters:
662 .  snes - the SNES context
663 .  func - function evaluation routine
664 .  ctx - [optional] user-defined context for private data for the
665          function evaluation routine (may be PETSC_NULL)
666 
667    Calling sequence of func:
668 .  func (SNES snes,Vec x,double *f,void *ctx);
669 
670 .  x - input vector
671 .  f - function
672 .  ctx - [optional] user-defined function context
673 
674    Notes:
675    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
676    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
677    SNESSetFunction().
678 
679 .keywords: SNES, nonlinear, set, minimization, function
680 
681 .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
682            SNESSetHessian(), SNESSetGradient()
683 @*/
684 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
685                       void *ctx)
686 {
687   PetscFunctionBegin;
688   PetscValidHeaderSpecific(snes,SNES_COOKIE);
689   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
690   snes->computeumfunction   = func;
691   snes->umfunP              = ctx;
692   PetscFunctionReturn(0);
693 }
694 
695 #undef __FUNC__
696 #define __FUNC__ "SNESComputeMinimizationFunction"
697 /*@
698    SNESComputeMinimizationFunction - Computes the function that has been
699    set with SNESSetMinimizationFunction().
700 
701    Input Parameters:
702 .  snes - the SNES context
703 .  x - input vector
704 
705    Output Parameter:
706 .  y - function value
707 
708    Notes:
709    SNESComputeMinimizationFunction() is valid only for
710    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
711    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
712 
713 .keywords: SNES, nonlinear, compute, minimization, function
714 
715 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
716           SNESComputeGradient(), SNESComputeHessian()
717 @*/
718 int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
719 {
720   int    ierr;
721 
722   PetscFunctionBegin;
723   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
724   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
725   PetscStackPush("SNES user minimzation function");
726   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
727   PetscStackPop;
728   snes->nfuncs++;
729   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
730   PetscFunctionReturn(0);
731 }
732 
733 #undef __FUNC__
734 #define __FUNC__ "SNESSetGradient"
735 /*@C
736    SNESSetGradient - Sets the gradient evaluation routine and gradient
737    vector for use by the SNES routines.
738 
739    Input Parameters:
740 .  snes - the SNES context
741 .  func - function evaluation routine
742 .  ctx - optional user-defined context for private data for the
743          gradient evaluation routine (may be PETSC_NULL)
744 .  r - vector to store gradient value
745 
746    Calling sequence of func:
747 .  func (SNES, Vec x, Vec g, void *ctx);
748 
749 .  x - input vector
750 .  g - gradient vector
751 .  ctx - optional user-defined gradient context
752 
753    Notes:
754    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
755    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
756    SNESSetFunction().
757 
758 .keywords: SNES, nonlinear, set, function
759 
760 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
761           SNESSetMinimizationFunction(),
762 @*/
763 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
764 {
765   PetscFunctionBegin;
766   PetscValidHeaderSpecific(snes,SNES_COOKIE);
767   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
768   snes->computefunction     = func;
769   snes->vec_func            = snes->vec_func_always = r;
770   snes->funP                = ctx;
771   PetscFunctionReturn(0);
772 }
773 
774 #undef __FUNC__
775 #define __FUNC__ "SNESComputeGradient"
776 /*@
777    SNESComputeGradient - Computes the gradient that has been set with
778    SNESSetGradient().
779 
780    Input Parameters:
781 .  snes - the SNES context
782 .  x - input vector
783 
784    Output Parameter:
785 .  y - gradient vector
786 
787    Notes:
788    SNESComputeGradient() is valid only for
789    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
790    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
791 
792 .keywords: SNES, nonlinear, compute, gradient
793 
794 .seealso:  SNESSetGradient(), SNESGetGradient(),
795            SNESComputeMinimizationFunction(), SNESComputeHessian()
796 @*/
797 int SNESComputeGradient(SNES snes,Vec x, Vec y)
798 {
799   int    ierr;
800 
801   PetscFunctionBegin;
802   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
803     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
804   }
805   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
806   PetscStackPush("SNES user gradient function");
807   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
808   PetscStackPop;
809   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
810   PetscFunctionReturn(0);
811 }
812 
813 #undef __FUNC__
814 #define __FUNC__ "SNESComputeJacobian"
815 /*@
816    SNESComputeJacobian - Computes the Jacobian matrix that has been
817    set with SNESSetJacobian().
818 
819    Input Parameters:
820 .  snes - the SNES context
821 .  x - input vector
822 
823    Output Parameters:
824 .  A - Jacobian matrix
825 .  B - optional preconditioning matrix
826 .  flag - flag indicating matrix structure
827 
828    Notes:
829    Most users should not need to explicitly call this routine, as it
830    is used internally within the nonlinear solvers.
831 
832    See SLESSetOperators() for important information about setting the
833    flag parameter.
834 
835    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
836    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
837    methods is SNESComputeHessian().
838 
839 .keywords: SNES, compute, Jacobian, matrix
840 
841 .seealso:  SNESSetJacobian(), SLESSetOperators()
842 @*/
843 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
844 {
845   int    ierr;
846 
847   PetscFunctionBegin;
848   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
849     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
850   }
851   if (!snes->computejacobian) PetscFunctionReturn(0);
852   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
853   *flg = DIFFERENT_NONZERO_PATTERN;
854   PetscStackPush("SNES user Jacobian function");
855   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
856   PetscStackPop;
857   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
858   /* make sure user returned a correct Jacobian and preconditioner */
859   PetscValidHeaderSpecific(*A,MAT_COOKIE);
860   PetscValidHeaderSpecific(*B,MAT_COOKIE);
861   PetscFunctionReturn(0);
862 }
863 
864 #undef __FUNC__
865 #define __FUNC__ "SNESComputeHessian"
866 /*@
867    SNESComputeHessian - Computes the Hessian matrix that has been
868    set with SNESSetHessian().
869 
870    Input Parameters:
871 .  snes - the SNES context
872 .  x - input vector
873 
874    Output Parameters:
875 .  A - Hessian matrix
876 .  B - optional preconditioning matrix
877 .  flag - flag indicating matrix structure
878 
879    Notes:
880    Most users should not need to explicitly call this routine, as it
881    is used internally within the nonlinear solvers.
882 
883    See SLESSetOperators() for important information about setting the
884    flag parameter.
885 
886    SNESComputeHessian() is valid only for
887    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
888    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
889 
890 .keywords: SNES, compute, Hessian, matrix
891 
892 .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
893            SNESComputeMinimizationFunction()
894 @*/
895 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
896 {
897   int    ierr;
898 
899   PetscFunctionBegin;
900   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
901     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
902   }
903   if (!snes->computejacobian) PetscFunctionReturn(0);
904   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
905   *flag = DIFFERENT_NONZERO_PATTERN;
906   PetscStackPush("SNES user Hessian function");
907   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
908   PetscStackPop;
909   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
910   /* make sure user returned a correct Jacobian and preconditioner */
911   PetscValidHeaderSpecific(*A,MAT_COOKIE);
912   PetscValidHeaderSpecific(*B,MAT_COOKIE);
913   PetscFunctionReturn(0);
914 }
915 
916 #undef __FUNC__
917 #define __FUNC__ "SNESSetJacobian"
918 /*@C
919    SNESSetJacobian - Sets the function to compute Jacobian as well as the
920    location to store the matrix.
921 
922    Input Parameters:
923 .  snes - the SNES context
924 .  A - Jacobian matrix
925 .  B - preconditioner matrix (usually same as the Jacobian)
926 .  func - Jacobian evaluation routine
927 .  ctx - [optional] user-defined context for private data for the
928          Jacobian evaluation routine (may be PETSC_NULL)
929 
930    Calling sequence of func:
931 .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
932 
933 .  x - input vector
934 .  A - Jacobian matrix
935 .  B - preconditioner matrix, usually the same as A
936 .  flag - flag indicating information about the preconditioner matrix
937    structure (same as flag in SLESSetOperators())
938 .  ctx - [optional] user-defined Jacobian context
939 
940    Notes:
941    See SLESSetOperators() for important information about setting the flag
942    output parameter in the routine func().  Be sure to read this information!
943 
944    The routine func() takes Mat * as the matrix arguments rather than Mat.
945    This allows the Jacobian evaluation routine to replace A and/or B with a
946    completely new new matrix structure (not just different matrix elements)
947    when appropriate, for instance, if the nonzero structure is changing
948    throughout the global iterations.
949 
950 .keywords: SNES, nonlinear, set, Jacobian, matrix
951 
952 .seealso: SLESSetOperators(), SNESSetFunction()
953 @*/
954 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
955                     MatStructure*,void*),void *ctx)
956 {
957   PetscFunctionBegin;
958   PetscValidHeaderSpecific(snes,SNES_COOKIE);
959   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
960   snes->computejacobian = func;
961   snes->jacP            = ctx;
962   snes->jacobian        = A;
963   snes->jacobian_pre    = B;
964   PetscFunctionReturn(0);
965 }
966 
967 #undef __FUNC__
968 #define __FUNC__ "SNESGetJacobian"
969 /*@
970    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
971    provided context for evaluating the Jacobian.
972 
973    Input Parameter:
974 .  snes - the nonlinear solver context
975 
976    Output Parameters:
977 .  A - location to stash Jacobian matrix (or PETSC_NULL)
978 .  B - location to stash preconditioner matrix (or PETSC_NULL)
979 .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
980 
981 .seealso: SNESSetJacobian(), SNESComputeJacobian()
982 @*/
983 int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
984 {
985   PetscFunctionBegin;
986   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
987     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
988   }
989   if (A)   *A = snes->jacobian;
990   if (B)   *B = snes->jacobian_pre;
991   if (ctx) *ctx = snes->jacP;
992   PetscFunctionReturn(0);
993 }
994 
995 #undef __FUNC__
996 #define __FUNC__ "SNESSetHessian"
997 /*@C
998    SNESSetHessian - Sets the function to compute Hessian as well as the
999    location to store the matrix.
1000 
1001    Input Parameters:
1002 .  snes - the SNES context
1003 .  A - Hessian matrix
1004 .  B - preconditioner matrix (usually same as the Hessian)
1005 .  func - Jacobian evaluation routine
1006 .  ctx - [optional] user-defined context for private data for the
1007          Hessian evaluation routine (may be PETSC_NULL)
1008 
1009    Calling sequence of func:
1010 .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1011 
1012 .  x - input vector
1013 .  A - Hessian matrix
1014 .  B - preconditioner matrix, usually the same as A
1015 .  flag - flag indicating information about the preconditioner matrix
1016    structure (same as flag in SLESSetOperators())
1017 .  ctx - [optional] user-defined Hessian context
1018 
1019    Notes:
1020    See SLESSetOperators() for important information about setting the flag
1021    output parameter in the routine func().  Be sure to read this information!
1022 
1023    The function func() takes Mat * as the matrix arguments rather than Mat.
1024    This allows the Hessian evaluation routine to replace A and/or B with a
1025    completely new new matrix structure (not just different matrix elements)
1026    when appropriate, for instance, if the nonzero structure is changing
1027    throughout the global iterations.
1028 
1029 .keywords: SNES, nonlinear, set, Hessian, matrix
1030 
1031 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
1032 @*/
1033 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
1034                     MatStructure*,void*),void *ctx)
1035 {
1036   PetscFunctionBegin;
1037   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1038   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1039     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1040   }
1041   snes->computejacobian = func;
1042   snes->jacP            = ctx;
1043   snes->jacobian        = A;
1044   snes->jacobian_pre    = B;
1045   PetscFunctionReturn(0);
1046 }
1047 
1048 #undef __FUNC__
1049 #define __FUNC__ "SNESGetHessian"
1050 /*@
1051    SNESGetHessian - Returns the Hessian matrix and optionally the user
1052    provided context for evaluating the Hessian.
1053 
1054    Input Parameter:
1055 .  snes - the nonlinear solver context
1056 
1057    Output Parameters:
1058 .  A - location to stash Hessian matrix (or PETSC_NULL)
1059 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1060 .  ctx - location to stash Hessian ctx (or PETSC_NULL)
1061 
1062 .seealso: SNESSetHessian(), SNESComputeHessian()
1063 @*/
1064 int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
1065 {
1066   PetscFunctionBegin;
1067   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1068     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1069   }
1070   if (A)   *A = snes->jacobian;
1071   if (B)   *B = snes->jacobian_pre;
1072   if (ctx) *ctx = snes->jacP;
1073   PetscFunctionReturn(0);
1074 }
1075 
1076 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1077 
1078 #undef __FUNC__
1079 #define __FUNC__ "SNESSetUp"
1080 /*@
1081    SNESSetUp - Sets up the internal data structures for the later use
1082    of a nonlinear solver.
1083 
1084    Input Parameter:
1085 .  snes - the SNES context
1086 .  x - the solution vector
1087 
1088    Notes:
1089    For basic use of the SNES solvers the user need not explicitly call
1090    SNESSetUp(), since these actions will automatically occur during
1091    the call to SNESSolve().  However, if one wishes to control this
1092    phase separately, SNESSetUp() should be called after SNESCreate()
1093    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1094 
1095 .keywords: SNES, nonlinear, setup
1096 
1097 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1098 @*/
1099 int SNESSetUp(SNES snes,Vec x)
1100 {
1101   int ierr, flg;
1102 
1103   PetscFunctionBegin;
1104   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1105   PetscValidHeaderSpecific(x,VEC_COOKIE);
1106   snes->vec_sol = snes->vec_sol_always = x;
1107 
1108   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1109   /*
1110       This version replaces the user provided Jacobian matrix with a
1111       matrix-free version but still employs the user-provided preconditioner matrix
1112   */
1113   if (flg) {
1114     Mat J;
1115     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1116     PLogObjectParent(snes,J);
1117     snes->mfshell = J;
1118     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1119       snes->jacobian = J;
1120       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1121     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1122       snes->jacobian = J;
1123       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1124     } else {
1125       SETERRQ(1,0,"Method class doesn't support matrix-free operator option");
1126     }
1127   }
1128   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1129   /*
1130       This version replaces both the user-provided Jacobian and the user-
1131       provided preconditioner matrix with the default matrix free version.
1132    */
1133   if (flg) {
1134     Mat J;
1135     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1136     PLogObjectParent(snes,J);
1137     snes->mfshell = J;
1138     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1139       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
1140       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1141     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1142       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
1143       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1144     } else {
1145       SETERRQ(1,0,"Method class doesn't support matrix-free option");
1146     }
1147   }
1148   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1149     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first");
1150     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first");
1151     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first");
1152     if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector");
1153 
1154     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
1155     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1156       SLES sles; KSP ksp;
1157       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1158       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1159       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1160              (void *)snes); CHKERRQ(ierr);
1161     }
1162   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1163     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first");
1164     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first");
1165     if (!snes->computeumfunction) SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first");
1166     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first");
1167   } else {
1168     SETERRQ(1,0,"Unknown method class");
1169   }
1170   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1171   snes->setup_called = 1;
1172   PetscFunctionReturn(0);
1173 }
1174 
1175 #undef __FUNC__
1176 #define __FUNC__ "SNESDestroy"
1177 /*@C
1178    SNESDestroy - Destroys the nonlinear solver context that was created
1179    with SNESCreate().
1180 
1181    Input Parameter:
1182 .  snes - the SNES context
1183 
1184 .keywords: SNES, nonlinear, destroy
1185 
1186 .seealso: SNESCreate(), SNESSolve()
1187 @*/
1188 int SNESDestroy(SNES snes)
1189 {
1190   int ierr;
1191 
1192   PetscFunctionBegin;
1193   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1194   if (--snes->refct > 0) PetscFunctionReturn(0);
1195 
1196   if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);}
1197   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
1198   if (snes->mfshell) MatDestroy(snes->mfshell);
1199   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
1200   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
1201   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
1202   PLogObjectDestroy((PetscObject)snes);
1203   PetscHeaderDestroy((PetscObject)snes);
1204   PetscFunctionReturn(0);
1205 }
1206 
1207 /* ----------- Routines to set solver parameters ---------- */
1208 
1209 #undef __FUNC__
1210 #define __FUNC__ "SNESSetTolerances"
1211 /*@
1212    SNESSetTolerances - Sets various parameters used in convergence tests.
1213 
1214    Input Parameters:
1215 .  snes - the SNES context
1216 .  atol - absolute convergence tolerance
1217 .  rtol - relative convergence tolerance
1218 .  stol -  convergence tolerance in terms of the norm
1219            of the change in the solution between steps
1220 .  maxit - maximum number of iterations
1221 .  maxf - maximum number of function evaluations
1222 
1223    Options Database Keys:
1224 $    -snes_atol <atol>
1225 $    -snes_rtol <rtol>
1226 $    -snes_stol <stol>
1227 $    -snes_max_it <maxit>
1228 $    -snes_max_funcs <maxf>
1229 
1230    Notes:
1231    The default maximum number of iterations is 50.
1232    The default maximum number of function evaluations is 1000.
1233 
1234 .keywords: SNES, nonlinear, set, convergence, tolerances
1235 
1236 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
1237 @*/
1238 int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
1239 {
1240   PetscFunctionBegin;
1241   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1242   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1243   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1244   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1245   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1246   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1247   PetscFunctionReturn(0);
1248 }
1249 
1250 #undef __FUNC__
1251 #define __FUNC__ "SNESGetTolerances"
1252 /*@
1253    SNESGetTolerances - Gets various parameters used in convergence tests.
1254 
1255    Input Parameters:
1256 .  snes - the SNES context
1257 .  atol - absolute convergence tolerance
1258 .  rtol - relative convergence tolerance
1259 .  stol -  convergence tolerance in terms of the norm
1260            of the change in the solution between steps
1261 .  maxit - maximum number of iterations
1262 .  maxf - maximum number of function evaluations
1263 
1264    Notes:
1265    The user can specify PETSC_NULL for any parameter that is not needed.
1266 
1267 .keywords: SNES, nonlinear, get, convergence, tolerances
1268 
1269 .seealso: SNESSetTolerances()
1270 @*/
1271 int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
1272 {
1273   PetscFunctionBegin;
1274   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1275   if (atol)  *atol  = snes->atol;
1276   if (rtol)  *rtol  = snes->rtol;
1277   if (stol)  *stol  = snes->xtol;
1278   if (maxit) *maxit = snes->max_its;
1279   if (maxf)  *maxf  = snes->max_funcs;
1280   PetscFunctionReturn(0);
1281 }
1282 
1283 #undef __FUNC__
1284 #define __FUNC__ "SNESSetTrustRegionTolerance"
1285 /*@
1286    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1287 
1288    Input Parameters:
1289 .  snes - the SNES context
1290 .  tol - tolerance
1291 
1292    Options Database Key:
1293 $    -snes_trtol <tol>
1294 
1295 .keywords: SNES, nonlinear, set, trust region, tolerance
1296 
1297 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
1298 @*/
1299 int SNESSetTrustRegionTolerance(SNES snes,double tol)
1300 {
1301   PetscFunctionBegin;
1302   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1303   snes->deltatol = tol;
1304   PetscFunctionReturn(0);
1305 }
1306 
1307 #undef __FUNC__
1308 #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
1309 /*@
1310    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
1311    for unconstrained minimization solvers.
1312 
1313    Input Parameters:
1314 .  snes - the SNES context
1315 .  ftol - minimum function tolerance
1316 
1317    Options Database Key:
1318 $    -snes_fmin <ftol>
1319 
1320    Note:
1321    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
1322    methods only.
1323 
1324 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
1325 
1326 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
1327 @*/
1328 int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
1329 {
1330   PetscFunctionBegin;
1331   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1332   snes->fmin = ftol;
1333   PetscFunctionReturn(0);
1334 }
1335 
1336 /* ------------ Routines to set performance monitoring options ----------- */
1337 
1338 #undef __FUNC__
1339 #define __FUNC__ "SNESSetMonitor"
1340 /*@C
1341    SNESSetMonitor - Sets the function that is to be used at every
1342    iteration of the nonlinear solver to display the iteration's
1343    progress.
1344 
1345    Input Parameters:
1346 .  snes - the SNES context
1347 .  func - monitoring routine
1348 .  mctx - [optional] user-defined context for private data for the
1349           monitor routine (may be PETSC_NULL)
1350 
1351    Calling sequence of func:
1352    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
1353 
1354 $    snes - the SNES context
1355 $    its - iteration number
1356 $    mctx - [optional] monitoring context
1357 $
1358 $ SNES_NONLINEAR_EQUATIONS methods:
1359 $    norm - 2-norm function value (may be estimated)
1360 $
1361 $ SNES_UNCONSTRAINED_MINIMIZATION methods:
1362 $    norm - 2-norm gradient value (may be estimated)
1363 
1364    Options Database Keys:
1365 $    -snes_monitor        : sets SNESDefaultMonitor()
1366 $    -snes_xmonitor       : sets line graph monitor,
1367 $                           uses SNESLGMonitorCreate()
1368 $    -snes_cancelmonitors : cancels all monitors that have
1369 $                           been hardwired into a code by
1370 $                           calls to SNESSetMonitor(), but
1371 $                           does not cancel those set via
1372 $                           the options database.
1373 
1374 
1375    Notes:
1376    Several different monitoring routines may be set by calling
1377    SNESSetMonitor() multiple times; all will be called in the
1378    order in which they were set.
1379 
1380 .keywords: SNES, nonlinear, set, monitor
1381 
1382 .seealso: SNESDefaultMonitor()
1383 @*/
1384 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
1385 {
1386   PetscFunctionBegin;
1387   if (!func) {
1388     snes->numbermonitors = 0;
1389     PetscFunctionReturn(0);
1390   }
1391   if (snes->numbermonitors >= MAXSNESMONITORS) {
1392     SETERRQ(1,0,"Too many monitors set");
1393   }
1394 
1395   snes->monitor[snes->numbermonitors]           = func;
1396   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1397   PetscFunctionReturn(0);
1398 }
1399 
1400 #undef __FUNC__
1401 #define __FUNC__ "SNESSetConvergenceTest"
1402 /*@C
1403    SNESSetConvergenceTest - Sets the function that is to be used
1404    to test for convergence of the nonlinear iterative solution.
1405 
1406    Input Parameters:
1407 .  snes - the SNES context
1408 .  func - routine to test for convergence
1409 .  cctx - [optional] context for private data for the convergence routine
1410           (may be PETSC_NULL)
1411 
1412    Calling sequence of func:
1413    int func (SNES snes,double xnorm,double gnorm,
1414              double f,void *cctx)
1415 
1416 $    snes - the SNES context
1417 $    cctx - [optional] convergence context
1418 $    xnorm - 2-norm of current iterate
1419 $
1420 $ SNES_NONLINEAR_EQUATIONS methods:
1421 $    gnorm - 2-norm of current step
1422 $    f - 2-norm of function
1423 $
1424 $ SNES_UNCONSTRAINED_MINIMIZATION methods:
1425 $    gnorm - 2-norm of current gradient
1426 $    f - function value
1427 
1428 .keywords: SNES, nonlinear, set, convergence, test
1429 
1430 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1431           SNESConverged_UM_LS(), SNESConverged_UM_TR()
1432 @*/
1433 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
1434 {
1435   PetscFunctionBegin;
1436   (snes)->converged = func;
1437   (snes)->cnvP      = cctx;
1438   PetscFunctionReturn(0);
1439 }
1440 
1441 #undef __FUNC__
1442 #define __FUNC__ "SNESSetConvergenceHistory"
1443 /*@
1444    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1445 
1446    Input Parameters:
1447 .  snes - iterative context obtained from SNESCreate()
1448 .  a   - array to hold history
1449 .  na  - size of a
1450 
1451    Notes:
1452    If set, this array will contain the function norms (for
1453    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1454    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1455    at each step.
1456 
1457    This routine is useful, e.g., when running a code for purposes
1458    of accurate performance monitoring, when no I/O should be done
1459    during the section of code that is being timed.
1460 
1461 .keywords: SNES, set, convergence, history
1462 @*/
1463 int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1464 {
1465   PetscFunctionBegin;
1466   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1467   if (na) PetscValidScalarPointer(a);
1468   snes->conv_hist      = a;
1469   snes->conv_hist_size = na;
1470   PetscFunctionReturn(0);
1471 }
1472 
1473 #undef __FUNC__
1474 #define __FUNC__ "SNESScaleStep_Private"
1475 /*
1476    SNESScaleStep_Private - Scales a step so that its length is less than the
1477    positive parameter delta.
1478 
1479     Input Parameters:
1480 .   snes - the SNES context
1481 .   y - approximate solution of linear system
1482 .   fnorm - 2-norm of current function
1483 .   delta - trust region size
1484 
1485     Output Parameters:
1486 .   gpnorm - predicted function norm at the new point, assuming local
1487     linearization.  The value is zero if the step lies within the trust
1488     region, and exceeds zero otherwise.
1489 .   ynorm - 2-norm of the step
1490 
1491     Note:
1492     For non-trust region methods such as SNES_EQ_LS, the parameter delta
1493     is set to be the maximum allowable step size.
1494 
1495 .keywords: SNES, nonlinear, scale, step
1496 */
1497 int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
1498                           double *gpnorm,double *ynorm)
1499 {
1500   double norm;
1501   Scalar cnorm;
1502   int    ierr;
1503 
1504   PetscFunctionBegin;
1505   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
1506   if (norm > *delta) {
1507      norm = *delta/norm;
1508      *gpnorm = (1.0 - norm)*(*fnorm);
1509      cnorm = norm;
1510      VecScale( &cnorm, y );
1511      *ynorm = *delta;
1512   } else {
1513      *gpnorm = 0.0;
1514      *ynorm = norm;
1515   }
1516   PetscFunctionReturn(0);
1517 }
1518 
1519 #undef __FUNC__
1520 #define __FUNC__ "SNESSolve"
1521 /*@
1522    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
1523    SNESCreate() and optional routines of the form SNESSetXXX().
1524 
1525    Input Parameter:
1526 .  snes - the SNES context
1527 .  x - the solution vector
1528 
1529    Output Parameter:
1530    its - number of iterations until termination
1531 
1532    Note:
1533    The user should initialize the vector, x, with the initial guess
1534    for the nonlinear solve prior to calling SNESSolve.  In particular,
1535    to employ an initial guess of zero, the user should explicitly set
1536    this vector to zero by calling VecSet().
1537 
1538 .keywords: SNES, nonlinear, solve
1539 
1540 .seealso: SNESCreate(), SNESDestroy()
1541 @*/
1542 int SNESSolve(SNES snes,Vec x,int *its)
1543 {
1544   int ierr, flg;
1545 
1546   PetscFunctionBegin;
1547   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1548   PetscValidIntPointer(its);
1549   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1550   else {snes->vec_sol = snes->vec_sol_always = x;}
1551   PLogEventBegin(SNES_Solve,snes,0,0,0);
1552   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
1553   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
1554   PLogEventEnd(SNES_Solve,snes,0,0,0);
1555   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
1556   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
1557   PetscFunctionReturn(0);
1558 }
1559 
1560 /* --------- Internal routines for SNES Package --------- */
1561 
1562 #undef __FUNC__
1563 #define __FUNC__ "SNESSetType"
1564 /*@
1565    SNESSetType - Sets the method for the nonlinear solver.
1566 
1567    Input Parameters:
1568 .  snes - the SNES context
1569 .  method - a known method
1570 
1571   Options Database Command:
1572 $ -snes_type  <method>
1573 $    Use -help for a list of available methods
1574 $    (for instance, ls or tr)
1575 
1576    Notes:
1577    See "petsc/include/snes.h" for available methods (for instance)
1578 $  Systems of nonlinear equations:
1579 $    SNES_EQ_LS - Newton's method with line search
1580 $    SNES_EQ_TR - Newton's method with trust region
1581 $  Unconstrained minimization:
1582 $    SNES_UM_TR - Newton's method with trust region
1583 $    SNES_UM_LS - Newton's method with line search
1584 
1585   Normally, it is best to use the SNESSetFromOptions() command and then
1586   set the SNES solver type from the options database rather than by using
1587   this routine.  Using the options database provides the user with
1588   maximum flexibility in evaluating the many nonlinear solvers.
1589   The SNESSetType() routine is provided for those situations where it
1590   is necessary to set the nonlinear solver independently of the command
1591   line or options database.  This might be the case, for example, when
1592   the choice of solver changes during the execution of the program,
1593   and the user's application is taking responsibility for choosing the
1594   appropriate method.  In other words, this routine is for the advanced user.
1595 
1596 .keywords: SNES, set, method
1597 @*/
1598 int SNESSetType(SNES snes,SNESType method)
1599 {
1600   int ierr;
1601   int (*r)(SNES);
1602 
1603   PetscFunctionBegin;
1604   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1605   if (snes->type == (int) method) PetscFunctionReturn(0);
1606 
1607   /* Get the function pointers for the iterative method requested */
1608   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
1609   if (!__SNESList) {SETERRQ(1,0,"Could not get methods");}
1610   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1611   if (!r) {SETERRQ(1,0,"Unknown method");}
1612   if (snes->data) PetscFree(snes->data);
1613   snes->set_method_called = 1;
1614   ierr = (*r)(snes);CHKERRQ(ierr);
1615   PetscFunctionReturn(0);
1616 }
1617 
1618 /* --------------------------------------------------------------------- */
1619 #undef __FUNC__
1620 #define __FUNC__ "SNESRegister"
1621 /*@C
1622    SNESRegister - Adds the method to the nonlinear solver package, given
1623    a function pointer and a nonlinear solver name of the type SNESType.
1624 
1625    Input Parameters:
1626 .  name - either a predefined name such as SNES_EQ_LS, or SNES_NEW
1627           to indicate a new user-defined solver
1628 .  sname - corresponding string for name
1629 .  create - routine to create method context
1630 
1631    Output Parameter:
1632 .  oname - type associated with this new solver
1633 
1634    Notes:
1635    Multiple user-defined nonlinear solvers can be added by calling
1636    SNESRegister() with the input parameter "name" set to be SNES_NEW;
1637    each call will return a unique solver type in the output
1638    parameter "oname".
1639 
1640 .keywords: SNES, nonlinear, register
1641 
1642 .seealso: SNESRegisterAll(), SNESRegisterDestroy()
1643 @*/
1644 int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES))
1645 {
1646   int        ierr;
1647   static int numberregistered = 0;
1648 
1649   PetscFunctionBegin;
1650   if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++);
1651 
1652   if (oname) *oname = name;
1653   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
1654   NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create );
1655   PetscFunctionReturn(0);
1656 }
1657 
1658 /* --------------------------------------------------------------------- */
1659 #undef __FUNC__
1660 #define __FUNC__ "SNESRegisterDestroy"
1661 /*@C
1662    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1663    registered by SNESRegister().
1664 
1665 .keywords: SNES, nonlinear, register, destroy
1666 
1667 .seealso: SNESRegisterAll(), SNESRegisterAll()
1668 @*/
1669 int SNESRegisterDestroy()
1670 {
1671   PetscFunctionBegin;
1672   if (__SNESList) {
1673     NRDestroy( __SNESList );
1674     __SNESList = 0;
1675   }
1676   SNESRegisterAllCalled = 0;
1677   PetscFunctionReturn(0);
1678 }
1679 
1680 #undef __FUNC__
1681 #define __FUNC__ "SNESGetType"
1682 /*@C
1683    SNESGetType - Gets the SNES method type and name (as a string).
1684 
1685    Input Parameter:
1686 .  snes - nonlinear solver context
1687 
1688    Output Parameter:
1689 .  method - SNES method (or use PETSC_NULL)
1690 .  name - name of SNES method (or use PETSC_NULL)
1691 
1692 .keywords: SNES, nonlinear, get, method, name
1693 @*/
1694 int SNESGetType(SNES snes, SNESType *method,char **name)
1695 {
1696   int ierr;
1697 
1698   PetscFunctionBegin;
1699   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
1700   if (method) *method = (SNESType) snes->type;
1701   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
1702   PetscFunctionReturn(0);
1703 }
1704 
1705 #undef __FUNC__
1706 #define __FUNC__ "SNESGetSolution"
1707 /*@C
1708    SNESGetSolution - Returns the vector where the approximate solution is
1709    stored.
1710 
1711    Input Parameter:
1712 .  snes - the SNES context
1713 
1714    Output Parameter:
1715 .  x - the solution
1716 
1717 .keywords: SNES, nonlinear, get, solution
1718 
1719 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
1720 @*/
1721 int SNESGetSolution(SNES snes,Vec *x)
1722 {
1723   PetscFunctionBegin;
1724   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1725   *x = snes->vec_sol_always;
1726   PetscFunctionReturn(0);
1727 }
1728 
1729 #undef __FUNC__
1730 #define __FUNC__ "SNESGetSolutionUpdate"
1731 /*@C
1732    SNESGetSolutionUpdate - Returns the vector where the solution update is
1733    stored.
1734 
1735    Input Parameter:
1736 .  snes - the SNES context
1737 
1738    Output Parameter:
1739 .  x - the solution update
1740 
1741    Notes:
1742    This vector is implementation dependent.
1743 
1744 .keywords: SNES, nonlinear, get, solution, update
1745 
1746 .seealso: SNESGetSolution(), SNESGetFunction
1747 @*/
1748 int SNESGetSolutionUpdate(SNES snes,Vec *x)
1749 {
1750   PetscFunctionBegin;
1751   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1752   *x = snes->vec_sol_update_always;
1753   PetscFunctionReturn(0);
1754 }
1755 
1756 #undef __FUNC__
1757 #define __FUNC__ "SNESGetFunction"
1758 /*@C
1759    SNESGetFunction - Returns the vector where the function is stored.
1760 
1761    Input Parameter:
1762 .  snes - the SNES context
1763 
1764    Output Parameter:
1765 .  r - the function
1766 
1767    Notes:
1768    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
1769    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
1770    SNESGetMinimizationFunction() and SNESGetGradient();
1771 
1772 .keywords: SNES, nonlinear, get, function
1773 
1774 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
1775           SNESGetGradient()
1776 @*/
1777 int SNESGetFunction(SNES snes,Vec *r)
1778 {
1779   PetscFunctionBegin;
1780   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1781   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
1782   *r = snes->vec_func_always;
1783   PetscFunctionReturn(0);
1784 }
1785 
1786 #undef __FUNC__
1787 #define __FUNC__ "SNESGetGradient"
1788 /*@C
1789    SNESGetGradient - Returns the vector where the gradient is stored.
1790 
1791    Input Parameter:
1792 .  snes - the SNES context
1793 
1794    Output Parameter:
1795 .  r - the gradient
1796 
1797    Notes:
1798    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
1799    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
1800    SNESGetFunction().
1801 
1802 .keywords: SNES, nonlinear, get, gradient
1803 
1804 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
1805 @*/
1806 int SNESGetGradient(SNES snes,Vec *r)
1807 {
1808   PetscFunctionBegin;
1809   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1810   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1811     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1812   }
1813   *r = snes->vec_func_always;
1814   PetscFunctionReturn(0);
1815 }
1816 
1817 #undef __FUNC__
1818 #define __FUNC__ "SNESGetMinimizationFunction"
1819 /*@
1820    SNESGetMinimizationFunction - Returns the scalar function value for
1821    unconstrained minimization problems.
1822 
1823    Input Parameter:
1824 .  snes - the SNES context
1825 
1826    Output Parameter:
1827 .  r - the function
1828 
1829    Notes:
1830    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
1831    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
1832    SNESGetFunction().
1833 
1834 .keywords: SNES, nonlinear, get, function
1835 
1836 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
1837 @*/
1838 int SNESGetMinimizationFunction(SNES snes,double *r)
1839 {
1840   PetscFunctionBegin;
1841   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1842   PetscValidScalarPointer(r);
1843   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1844     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1845   }
1846   *r = snes->fc;
1847   PetscFunctionReturn(0);
1848 }
1849 
1850 #undef __FUNC__
1851 #define __FUNC__ "SNESSetOptionsPrefix"
1852 /*@C
1853    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1854    SNES options in the database.
1855 
1856    Input Parameter:
1857 .  snes - the SNES context
1858 .  prefix - the prefix to prepend to all option names
1859 
1860    Notes:
1861    A hyphen (-) must NOT be given at the beginning of the prefix name.
1862    The first character of all runtime options is AUTOMATICALLY the
1863    hyphen.
1864 
1865 .keywords: SNES, set, options, prefix, database
1866 
1867 .seealso: SNESSetFromOptions()
1868 @*/
1869 int SNESSetOptionsPrefix(SNES snes,char *prefix)
1870 {
1871   int ierr;
1872 
1873   PetscFunctionBegin;
1874   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1875   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
1876   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
1877   PetscFunctionReturn(0);
1878 }
1879 
1880 #undef __FUNC__
1881 #define __FUNC__ "SNESAppendOptionsPrefix"
1882 /*@C
1883    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1884    SNES options in the database.
1885 
1886    Input Parameter:
1887 .  snes - the SNES context
1888 .  prefix - the prefix to prepend to all option names
1889 
1890    Notes:
1891    A hyphen (-) must NOT be given at the beginning of the prefix name.
1892    The first character of all runtime options is AUTOMATICALLY the
1893    hyphen.
1894 
1895 .keywords: SNES, append, options, prefix, database
1896 
1897 .seealso: SNESGetOptionsPrefix()
1898 @*/
1899 int SNESAppendOptionsPrefix(SNES snes,char *prefix)
1900 {
1901   int ierr;
1902 
1903   PetscFunctionBegin;
1904   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1905   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
1906   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
1907   PetscFunctionReturn(0);
1908 }
1909 
1910 #undef __FUNC__
1911 #define __FUNC__ "SNESGetOptionsPrefix"
1912 /*@
1913    SNESGetOptionsPrefix - Sets the prefix used for searching for all
1914    SNES options in the database.
1915 
1916    Input Parameter:
1917 .  snes - the SNES context
1918 
1919    Output Parameter:
1920 .  prefix - pointer to the prefix string used
1921 
1922 .keywords: SNES, get, options, prefix, database
1923 
1924 .seealso: SNESAppendOptionsPrefix()
1925 @*/
1926 int SNESGetOptionsPrefix(SNES snes,char **prefix)
1927 {
1928   int ierr;
1929 
1930   PetscFunctionBegin;
1931   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1932   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
1933   PetscFunctionReturn(0);
1934 }
1935 
1936 #undef __FUNC__
1937 #define __FUNC__ "SNESPrintHelp"
1938 /*@
1939    SNESPrintHelp - Prints all options for the SNES component.
1940 
1941    Input Parameter:
1942 .  snes - the SNES context
1943 
1944    Options Database Keys:
1945 $  -help, -h
1946 
1947 .keywords: SNES, nonlinear, help
1948 
1949 .seealso: SNESSetFromOptions()
1950 @*/
1951 int SNESPrintHelp(SNES snes)
1952 {
1953   char                p[64];
1954   SNES_KSP_EW_ConvCtx *kctx;
1955   int                 ierr;
1956 
1957   PetscFunctionBegin;
1958   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1959 
1960   PetscStrcpy(p,"-");
1961   if (snes->prefix) PetscStrcat(p, snes->prefix);
1962 
1963   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
1964 
1965   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
1966   ierr = NRPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",__SNESList);CHKERRQ(ierr);
1967   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
1968   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
1969   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
1970   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
1971   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
1972   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
1973   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
1974   PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n");
1975   PetscPrintf(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
1976   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
1977     residual norm at each iteration.\n",p);
1978   PetscPrintf(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
1979     residual norm for small residual norms. This is useful to conceal\n\
1980     meaningless digits that may be different on different machines.\n",p);
1981   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
1982   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
1983     PetscPrintf(snes->comm,
1984      " Options for solving systems of nonlinear equations only:\n");
1985     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
1986     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
1987     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
1988     PetscPrintf(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
1989     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
1990     PetscPrintf(snes->comm,
1991      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
1992     PetscPrintf(snes->comm,
1993      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
1994     PetscPrintf(snes->comm,
1995      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
1996     PetscPrintf(snes->comm,
1997      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
1998     PetscPrintf(snes->comm,
1999      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
2000     PetscPrintf(snes->comm,
2001      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
2002     PetscPrintf(snes->comm,
2003      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
2004   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
2005     PetscPrintf(snes->comm,
2006      " Options for solving unconstrained minimization problems only:\n");
2007     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
2008     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
2009     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
2010   }
2011   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
2012   PetscPrintf(snes->comm,"a particular method\n");
2013   if (snes->printhelp) {
2014     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2015   }
2016   PetscFunctionReturn(0);
2017 }
2018 
2019 
2020 
2021 
2022 
2023