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