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