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