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