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