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