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