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