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