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