xref: /petsc/src/snes/interface/snes.c (revision 00ff320a6701aa587ed3bd591b944da12a89b695)
1 /*$Id: snes.c,v 1.215 2000/07/13 17:28:22 bsmith Exp bsmith $*/
2 
3 #include "src/snes/snesimpl.h"      /*I "petscsnes.h"  I*/
4 
5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
6 FList      SNESList = 0;
7 
8 #undef __FUNC__
9 #define __FUNC__ /*<a name=""></a>*/"SNESView"
10 /*@C
11    SNESView - Prints the SNES data structure.
12 
13    Collective on SNES
14 
15    Input Parameters:
16 +  SNES - the SNES context
17 -  viewer - visualization context
18 
19    Options Database Key:
20 .  -snes_view - Calls SNESView() at end of SNESSolve()
21 
22    Notes:
23    The available visualization contexts include
24 +     VIEWER_STDOUT_SELF - standard output (default)
25 -     VIEWER_STDOUT_WORLD - synchronized standard
26          output where only the first processor opens
27          the file.  All other processors send their
28          data to the first processor to print.
29 
30    The user can open an alternative visualization context with
31    ViewerASCIIOpen() - output to a specified file.
32 
33    Level: beginner
34 
35 .keywords: SNES, view
36 
37 .seealso: ViewerASCIIOpen()
38 @*/
39 int SNESView(SNES snes,Viewer viewer)
40 {
41   SNES_KSP_EW_ConvCtx *kctx;
42   int                 ierr;
43   SLES                sles;
44   char                *type;
45   PetscTruth          isascii,isstring;
46 
47   PetscFunctionBegin;
48   PetscValidHeaderSpecific(snes,SNES_COOKIE);
49   if (!viewer) viewer = VIEWER_STDOUT_(snes->comm);
50   PetscValidHeaderSpecific(viewer,VIEWER_COOKIE);
51   PetscCheckSameComm(snes,viewer);
52 
53   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
54   ierr = PetscTypeCompare((PetscObject)viewer,STRING_VIEWER,&isstring);CHKERRQ(ierr);
55   if (isascii) {
56     ierr = ViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
57     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
58     if (type) {
59       ierr = ViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
60     } else {
61       ierr = ViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
62     }
63     if (snes->view) {
64       ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
65       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
66       ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
67     }
68     ierr = ViewerASCIIPrintf(viewer,"  maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
69     ierr = ViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
70                  snes->rtol,snes->atol,snes->trunctol,snes->xtol);CHKERRQ(ierr);
71     ierr = ViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr);
72     ierr = ViewerASCIIPrintf(viewer,"  total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr);
73     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
74       ierr = ViewerASCIIPrintf(viewer,"  min function tolerance=%g\n",snes->fmin);CHKERRQ(ierr);
75     }
76     if (snes->ksp_ewconv) {
77       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
78       if (kctx) {
79         ierr = ViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr);
80         ierr = ViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
81         ierr = ViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
82       }
83     }
84   } else if (isstring) {
85     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
86     ierr = ViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
87   }
88   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
89   ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
90   ierr = SLESView(sles,viewer);CHKERRQ(ierr);
91   ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
92   PetscFunctionReturn(0);
93 }
94 
95 /*
96        We retain a list of functions that also take SNES command
97     line options. These are called at the end SNESSetFromOptions()
98 */
99 #define MAXSETFROMOPTIONS 5
100 static int numberofsetfromoptions;
101 static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
102 
103 #undef __FUNC__
104 #define __FUNC__ /*<a name=""></a>*/"SNESAddOptionsChecker"
105 /*@
106     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
107 
108     Not Collective
109 
110     Input Parameter:
111 .   snescheck - function that checks for options
112 
113     Level: developer
114 
115 .seealso: SNESSetFromOptions()
116 @*/
117 int SNESAddOptionsChecker(int (*snescheck)(SNES))
118 {
119   PetscFunctionBegin;
120   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
121     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed");
122   }
123 
124   othersetfromoptions[numberofsetfromoptions++] = snescheck;
125   PetscFunctionReturn(0);
126 }
127 
128 #undef __FUNC__
129 #define __FUNC__ /*<a name=""></a>*/"SNESSetTypeFromOptions"
130 /*@
131    SNESSetTypeFromOptions - Sets the SNES solver type from the options database,
132         or sets a default if none is give.
133 
134    Collective on SNES
135 
136    Input Parameter:
137 .  snes - the SNES context
138 
139    Options Database Keys:
140 .  -snes_type <type> - ls, tr, umls, umtr, test
141 
142    Level: beginner
143 
144 .keywords: SNES, nonlinear, set, options, database
145 
146 .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetFromOptions()
147 @*/
148 int SNESSetTypeFromOptions(SNES snes)
149 {
150   char       type[256];
151   int        ierr;
152   PetscTruth flg;
153 
154   PetscFunctionBegin;
155   PetscValidHeaderSpecific(snes,SNES_COOKIE);
156   if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp()");
157   ierr = OptionsGetString(snes->prefix,"-snes_type",type,256,&flg);CHKERRQ(ierr);
158   if (flg) {
159     ierr = SNESSetType(snes,(SNESType) type);CHKERRQ(ierr);
160   }
161   /*
162       If SNES type has not yet been set, set it now
163   */
164   if (!snes->type_name) {
165     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
166       ierr = SNESSetType(snes,SNESEQLS);CHKERRQ(ierr);
167     } else {
168       ierr = SNESSetType(snes,SNESUMTR);CHKERRQ(ierr);
169     }
170   }
171   PetscFunctionReturn(0);
172 }
173 
174 #undef __FUNC__
175 #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions"
176 /*@
177    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
178 
179    Collective on SNES
180 
181    Input Parameter:
182 .  snes - the SNES context
183 
184    Options Database Keys:
185 +  -snes_type <type> - ls, tr, umls, umtr, test
186 .  -snes_stol - convergence tolerance in terms of the norm
187                 of the change in the solution between steps
188 .  -snes_atol <atol> - absolute tolerance of residual norm
189 .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
190 .  -snes_max_it <max_it> - maximum number of iterations
191 .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
192 .  -snes_trtol <trtol> - trust region tolerance
193 .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
194                                solver; hence iterations will continue until max_it
195                                or some other criterion is reached. Saves expense
196                                of convergence test
197 .  -snes_monitor - prints residual norm at each iteration
198 .  -snes_vecmonitor - plots solution at each iteration
199 .  -snes_vecmonitor_update - plots update to solution at each iteration
200 .  -snes_xmonitor - plots residual norm at each iteration
201 .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
202 -  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
203 
204     Options Database for Eisenstat-Walker method:
205 +  -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence
206 .  -snes_ksp_eq_version ver - version of  Eisenstat-Walker method
207 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
208 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
209 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
210 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
211 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
212 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
213 
214    Notes:
215    To see all options, run your program with the -help option or consult
216    the users manual.
217 
218    Level: beginner
219 
220 .keywords: SNES, nonlinear, set, options, database
221 
222 .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetTypeFromOptions()
223 @*/
224 int SNESSetFromOptions(SNES snes)
225 {
226   PetscReal  tmp;
227   SLES       sles;
228   PetscTruth flg;
229   int        ierr,i,loc[4],nmax = 4;
230   int        version   = PETSC_DEFAULT;
231   PetscReal  rtol_0    = PETSC_DEFAULT;
232   PetscReal  rtol_max  = PETSC_DEFAULT;
233   PetscReal  gamma2    = PETSC_DEFAULT;
234   PetscReal  alpha     = PETSC_DEFAULT;
235   PetscReal  alpha2    = PETSC_DEFAULT;
236   PetscReal  threshold = PETSC_DEFAULT;
237 
238   PetscFunctionBegin;
239   PetscValidHeaderSpecific(snes,SNES_COOKIE);
240   ierr = SNESSetTypeFromOptions(snes);CHKERRQ(ierr);
241 
242   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
243   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp,&flg);CHKERRQ(ierr);
244   if (flg) {
245     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
246   }
247   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp,&flg);CHKERRQ(ierr);
248   if (flg) {
249     ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
250   }
251   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp,&flg);CHKERRQ(ierr);
252   if (flg) {
253     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
254   }
255   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
256   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
257   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp,&flg);CHKERRQ(ierr);
258   if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp);CHKERRQ(ierr); }
259   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp,&flg);CHKERRQ(ierr);
260   if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);CHKERRQ(ierr);}
261   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv",&flg);CHKERRQ(ierr);
262   if (flg) { snes->ksp_ewconv = 1; }
263   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,PETSC_NULL);CHKERRQ(ierr);
264   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,PETSC_NULL);CHKERRQ(ierr);
265   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,PETSC_NULL);CHKERRQ(ierr);
266   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,PETSC_NULL);CHKERRQ(ierr);
267   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,PETSC_NULL);CHKERRQ(ierr);
268   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,PETSC_NULL);CHKERRQ(ierr);
269   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,PETSC_NULL);CHKERRQ(ierr);
270 
271   ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg);CHKERRQ(ierr);
272   if (flg) {snes->converged = 0;}
273 
274   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
275                                   alpha2,threshold);CHKERRQ(ierr);
276   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg);CHKERRQ(ierr);
277   if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
278   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg);CHKERRQ(ierr);
279   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
280   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg);CHKERRQ(ierr);
281   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
282   ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor",&flg);CHKERRQ(ierr);
283   if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
284   ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor_update",&flg);CHKERRQ(ierr);
285   if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
286   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
287   if (flg) {
288     ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
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,MatSNESMFKSPMonitor,PETSC_NULL,0);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__ /*<a name=""></a>*/"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__ /*<a name=""></a>*/"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__ /*<a name=""></a>*/"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__ /*<a name=""></a>*/"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__ /*<a name=""></a>*/"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__ /*<a name=""></a>*/"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__ /*<a name=""></a>*/"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__ /*<a name=""></a>*/"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__ /*<a name=""></a>*/"SNESPublish_Petsc"
590 static int SNESPublish_Petsc(PetscObject obj)
591 {
592 #if defined(PETSC_HAVE_AMS)
593   SNES          v = (SNES) obj;
594   int          ierr;
595 #endif
596 
597   PetscFunctionBegin;
598 
599 #if defined(PETSC_HAVE_AMS)
600   /* if it is already published then return */
601   if (v->amem >=0) PetscFunctionReturn(0);
602 
603   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
604   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
605                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
606   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
607                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
608   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
609 #endif
610   PetscFunctionReturn(0);
611 }
612 
613 /* -----------------------------------------------------------*/
614 #undef __FUNC__
615 #define __FUNC__ /*<a name=""></a>*/"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->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   snes->reason            = SNES_CONVERGED_ITERATING;
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   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
717   PetscFunctionReturn(0);
718 }
719 
720 /* --------------------------------------------------------------- */
721 #undef __FUNC__
722 #define __FUNC__ /*<a name=""></a>*/"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   PetscValidHeaderSpecific(r,VEC_COOKIE);
763   PetscCheckSameComm(snes,r);
764   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
765     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
766   }
767 
768   snes->computefunction     = func;
769   snes->vec_func            = snes->vec_func_always = r;
770   snes->funP                = ctx;
771   PetscFunctionReturn(0);
772 }
773 
774 #undef __FUNC__
775 #define __FUNC__ /*<a name=""></a>*/"SNESComputeFunction"
776 /*@
777    SNESComputeFunction - Calls the function that has been set with
778                          SNESSetFunction().
779 
780    Collective on SNES
781 
782    Input Parameters:
783 +  snes - the SNES context
784 -  x - input vector
785 
786    Output Parameter:
787 .  y - function vector, as set by SNESSetFunction()
788 
789    Notes:
790    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
791    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
792    SNESComputeMinimizationFunction() and SNESComputeGradient();
793 
794    SNESComputeFunction() is typically used within nonlinear solvers
795    implementations, so most users would not generally call this routine
796    themselves.
797 
798    Level: developer
799 
800 .keywords: SNES, nonlinear, compute, function
801 
802 .seealso: SNESSetFunction(), SNESGetFunction()
803 @*/
804 int SNESComputeFunction(SNES snes,Vec x,Vec y)
805 {
806   int    ierr;
807 
808   PetscFunctionBegin;
809   PetscValidHeaderSpecific(snes,SNES_COOKIE);
810   PetscValidHeaderSpecific(x,VEC_COOKIE);
811   PetscValidHeaderSpecific(y,VEC_COOKIE);
812   PetscCheckSameComm(snes,x);
813   PetscCheckSameComm(snes,y);
814   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
815     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
816   }
817 
818   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
819   PetscStackPush("SNES user function");
820   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
821   PetscStackPop;
822   snes->nfuncs++;
823   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
824   PetscFunctionReturn(0);
825 }
826 
827 #undef __FUNC__
828 #define __FUNC__ /*<a name=""></a>*/"SNESSetMinimizationFunction"
829 /*@C
830    SNESSetMinimizationFunction - Sets the function evaluation routine for
831    unconstrained minimization.
832 
833    Collective on SNES
834 
835    Input Parameters:
836 +  snes - the SNES context
837 .  func - function evaluation routine
838 -  ctx - [optional] user-defined context for private data for the
839          function evaluation routine (may be PETSC_NULL)
840 
841    Calling sequence of func:
842 $     func (SNES snes,Vec x,PetscReal *f,void *ctx);
843 
844 +  x - input vector
845 .  f - function
846 -  ctx - [optional] user-defined function context
847 
848    Level: beginner
849 
850    Notes:
851    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
852    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
853    SNESSetFunction().
854 
855 .keywords: SNES, nonlinear, set, minimization, function
856 
857 .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
858            SNESSetHessian(), SNESSetGradient()
859 @*/
860 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,PetscReal*,void*),void *ctx)
861 {
862   PetscFunctionBegin;
863   PetscValidHeaderSpecific(snes,SNES_COOKIE);
864   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
865     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
866   }
867   snes->computeumfunction   = func;
868   snes->umfunP              = ctx;
869   PetscFunctionReturn(0);
870 }
871 
872 #undef __FUNC__
873 #define __FUNC__ /*<a name=""></a>*/"SNESComputeMinimizationFunction"
874 /*@
875    SNESComputeMinimizationFunction - Computes the function that has been
876    set with SNESSetMinimizationFunction().
877 
878    Collective on SNES
879 
880    Input Parameters:
881 +  snes - the SNES context
882 -  x - input vector
883 
884    Output Parameter:
885 .  y - function value
886 
887    Notes:
888    SNESComputeMinimizationFunction() is valid only for
889    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
890    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
891 
892    SNESComputeMinimizationFunction() is typically used within minimization
893    implementations, so most users would not generally call this routine
894    themselves.
895 
896    Level: developer
897 
898 .keywords: SNES, nonlinear, compute, minimization, function
899 
900 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
901           SNESComputeGradient(), SNESComputeHessian()
902 @*/
903 int SNESComputeMinimizationFunction(SNES snes,Vec x,PetscReal *y)
904 {
905   int    ierr;
906 
907   PetscFunctionBegin;
908   PetscValidHeaderSpecific(snes,SNES_COOKIE);
909   PetscValidHeaderSpecific(x,VEC_COOKIE);
910   PetscCheckSameComm(snes,x);
911   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
912     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
913   }
914 
915   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
916   PetscStackPush("SNES user minimzation function");
917   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP);CHKERRQ(ierr);
918   PetscStackPop;
919   snes->nfuncs++;
920   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
921   PetscFunctionReturn(0);
922 }
923 
924 #undef __FUNC__
925 #define __FUNC__ /*<a name=""></a>*/"SNESSetGradient"
926 /*@C
927    SNESSetGradient - Sets the gradient evaluation routine and gradient
928    vector for use by the SNES routines.
929 
930    Collective on SNES
931 
932    Input Parameters:
933 +  snes - the SNES context
934 .  func - function evaluation routine
935 .  ctx - optional user-defined context for private data for the
936          gradient evaluation routine (may be PETSC_NULL)
937 -  r - vector to store gradient value
938 
939    Calling sequence of func:
940 $     func (SNES, Vec x, Vec g, void *ctx);
941 
942 +  x - input vector
943 .  g - gradient vector
944 -  ctx - optional user-defined gradient context
945 
946    Notes:
947    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
948    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
949    SNESSetFunction().
950 
951    Level: beginner
952 
953 .keywords: SNES, nonlinear, set, function
954 
955 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
956           SNESSetMinimizationFunction(),
957 @*/
958 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
959 {
960   PetscFunctionBegin;
961   PetscValidHeaderSpecific(snes,SNES_COOKIE);
962   PetscValidHeaderSpecific(r,VEC_COOKIE);
963   PetscCheckSameComm(snes,r);
964   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
965     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
966   }
967   snes->computefunction     = func;
968   snes->vec_func            = snes->vec_func_always = r;
969   snes->funP                = ctx;
970   PetscFunctionReturn(0);
971 }
972 
973 #undef __FUNC__
974 #define __FUNC__ /*<a name=""></a>*/"SNESComputeGradient"
975 /*@
976    SNESComputeGradient - Computes the gradient that has been set with
977    SNESSetGradient().
978 
979    Collective on SNES
980 
981    Input Parameters:
982 +  snes - the SNES context
983 -  x - input vector
984 
985    Output Parameter:
986 .  y - gradient vector
987 
988    Notes:
989    SNESComputeGradient() is valid only for
990    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
991    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
992 
993    SNESComputeGradient() is typically used within minimization
994    implementations, so most users would not generally call this routine
995    themselves.
996 
997    Level: developer
998 
999 .keywords: SNES, nonlinear, compute, gradient
1000 
1001 .seealso:  SNESSetGradient(), SNESGetGradient(),
1002            SNESComputeMinimizationFunction(), SNESComputeHessian()
1003 @*/
1004 int SNESComputeGradient(SNES snes,Vec x,Vec y)
1005 {
1006   int    ierr;
1007 
1008   PetscFunctionBegin;
1009   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1010   PetscValidHeaderSpecific(x,VEC_COOKIE);
1011   PetscValidHeaderSpecific(y,VEC_COOKIE);
1012   PetscCheckSameComm(snes,x);
1013   PetscCheckSameComm(snes,y);
1014   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1015     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1016   }
1017   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
1018   PetscStackPush("SNES user gradient function");
1019   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
1020   PetscStackPop;
1021   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
1022   PetscFunctionReturn(0);
1023 }
1024 
1025 #undef __FUNC__
1026 #define __FUNC__ /*<a name=""></a>*/"SNESComputeJacobian"
1027 /*@
1028    SNESComputeJacobian - Computes the Jacobian matrix that has been
1029    set with SNESSetJacobian().
1030 
1031    Collective on SNES and Mat
1032 
1033    Input Parameters:
1034 +  snes - the SNES context
1035 -  x - input vector
1036 
1037    Output Parameters:
1038 +  A - Jacobian matrix
1039 .  B - optional preconditioning matrix
1040 -  flag - flag indicating matrix structure
1041 
1042    Notes:
1043    Most users should not need to explicitly call this routine, as it
1044    is used internally within the nonlinear solvers.
1045 
1046    See SLESSetOperators() for important information about setting the
1047    flag parameter.
1048 
1049    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
1050    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
1051    methods is SNESComputeHessian().
1052 
1053    Level: developer
1054 
1055 .keywords: SNES, compute, Jacobian, matrix
1056 
1057 .seealso:  SNESSetJacobian(), SLESSetOperators()
1058 @*/
1059 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
1060 {
1061   int    ierr;
1062 
1063   PetscFunctionBegin;
1064   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1065   PetscValidHeaderSpecific(X,VEC_COOKIE);
1066   PetscCheckSameComm(snes,X);
1067   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1068     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1069   }
1070   if (!snes->computejacobian) PetscFunctionReturn(0);
1071   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
1072   *flg = DIFFERENT_NONZERO_PATTERN;
1073   PetscStackPush("SNES user Jacobian function");
1074   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1075   PetscStackPop;
1076   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
1077   /* make sure user returned a correct Jacobian and preconditioner */
1078   PetscValidHeaderSpecific(*A,MAT_COOKIE);
1079   PetscValidHeaderSpecific(*B,MAT_COOKIE);
1080   PetscFunctionReturn(0);
1081 }
1082 
1083 #undef __FUNC__
1084 #define __FUNC__ /*<a name=""></a>*/"SNESComputeHessian"
1085 /*@
1086    SNESComputeHessian - Computes the Hessian matrix that has been
1087    set with SNESSetHessian().
1088 
1089    Collective on SNES and Mat
1090 
1091    Input Parameters:
1092 +  snes - the SNES context
1093 -  x - input vector
1094 
1095    Output Parameters:
1096 +  A - Hessian matrix
1097 .  B - optional preconditioning matrix
1098 -  flag - flag indicating matrix structure
1099 
1100    Notes:
1101    Most users should not need to explicitly call this routine, as it
1102    is used internally within the nonlinear solvers.
1103 
1104    See SLESSetOperators() for important information about setting the
1105    flag parameter.
1106 
1107    SNESComputeHessian() is valid only for
1108    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
1109    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
1110 
1111    SNESComputeHessian() is typically used within minimization
1112    implementations, so most users would not generally call this routine
1113    themselves.
1114 
1115    Level: developer
1116 
1117 .keywords: SNES, compute, Hessian, matrix
1118 
1119 .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
1120            SNESComputeMinimizationFunction()
1121 @*/
1122 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
1123 {
1124   int    ierr;
1125 
1126   PetscFunctionBegin;
1127   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1128   PetscValidHeaderSpecific(x,VEC_COOKIE);
1129   PetscCheckSameComm(snes,x);
1130   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1131     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1132   }
1133   if (!snes->computejacobian) PetscFunctionReturn(0);
1134   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
1135   *flag = DIFFERENT_NONZERO_PATTERN;
1136   PetscStackPush("SNES user Hessian function");
1137   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP);CHKERRQ(ierr);
1138   PetscStackPop;
1139   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
1140   /* make sure user returned a correct Jacobian and preconditioner */
1141   PetscValidHeaderSpecific(*A,MAT_COOKIE);
1142   PetscValidHeaderSpecific(*B,MAT_COOKIE);
1143   PetscFunctionReturn(0);
1144 }
1145 
1146 #undef __FUNC__
1147 #define __FUNC__ /*<a name=""></a>*/"SNESSetJacobian"
1148 /*@C
1149    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1150    location to store the matrix.
1151 
1152    Collective on SNES and Mat
1153 
1154    Input Parameters:
1155 +  snes - the SNES context
1156 .  A - Jacobian matrix
1157 .  B - preconditioner matrix (usually same as the Jacobian)
1158 .  func - Jacobian evaluation routine
1159 -  ctx - [optional] user-defined context for private data for the
1160          Jacobian evaluation routine (may be PETSC_NULL)
1161 
1162    Calling sequence of func:
1163 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1164 
1165 +  x - input vector
1166 .  A - Jacobian matrix
1167 .  B - preconditioner matrix, usually the same as A
1168 .  flag - flag indicating information about the preconditioner matrix
1169    structure (same as flag in SLESSetOperators())
1170 -  ctx - [optional] user-defined Jacobian context
1171 
1172    Notes:
1173    See SLESSetOperators() for important information about setting the flag
1174    output parameter in the routine func().  Be sure to read this information!
1175 
1176    The routine func() takes Mat * as the matrix arguments rather than Mat.
1177    This allows the Jacobian evaluation routine to replace A and/or B with a
1178    completely new new matrix structure (not just different matrix elements)
1179    when appropriate, for instance, if the nonzero structure is changing
1180    throughout the global iterations.
1181 
1182    Level: beginner
1183 
1184 .keywords: SNES, nonlinear, set, Jacobian, matrix
1185 
1186 .seealso: SLESSetOperators(), SNESSetFunction()
1187 @*/
1188 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1189 {
1190   PetscFunctionBegin;
1191   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1192   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE);
1193   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE);
1194   if (A) PetscCheckSameComm(snes,A);
1195   if (B) PetscCheckSameComm(snes,B);
1196   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1197     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1198   }
1199 
1200   snes->computejacobian = func;
1201   snes->jacP            = ctx;
1202   snes->jacobian        = A;
1203   snes->jacobian_pre    = B;
1204   PetscFunctionReturn(0);
1205 }
1206 
1207 #undef __FUNC__
1208 #define __FUNC__ /*<a name=""></a>*/"SNESGetJacobian"
1209 /*@C
1210    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1211    provided context for evaluating the Jacobian.
1212 
1213    Not Collective, but Mat object will be parallel if SNES object is
1214 
1215    Input Parameter:
1216 .  snes - the nonlinear solver context
1217 
1218    Output Parameters:
1219 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1220 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1221 .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1222 -  func - location to put Jacobian function (or PETSC_NULL)
1223 
1224    Level: advanced
1225 
1226 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1227 @*/
1228 int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
1229 {
1230   PetscFunctionBegin;
1231   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1232   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1233     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1234   }
1235   if (A)    *A    = snes->jacobian;
1236   if (B)    *B    = snes->jacobian_pre;
1237   if (ctx)  *ctx  = snes->jacP;
1238   if (func) *func = snes->computejacobian;
1239   PetscFunctionReturn(0);
1240 }
1241 
1242 #undef __FUNC__
1243 #define __FUNC__ /*<a name=""></a>*/"SNESSetHessian"
1244 /*@C
1245    SNESSetHessian - Sets the function to compute Hessian as well as the
1246    location to store the matrix.
1247 
1248    Collective on SNES and Mat
1249 
1250    Input Parameters:
1251 +  snes - the SNES context
1252 .  A - Hessian matrix
1253 .  B - preconditioner matrix (usually same as the Hessian)
1254 .  func - Jacobian evaluation routine
1255 -  ctx - [optional] user-defined context for private data for the
1256          Hessian evaluation routine (may be PETSC_NULL)
1257 
1258    Calling sequence of func:
1259 $    func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1260 
1261 +  x - input vector
1262 .  A - Hessian matrix
1263 .  B - preconditioner matrix, usually the same as A
1264 .  flag - flag indicating information about the preconditioner matrix
1265    structure (same as flag in SLESSetOperators())
1266 -  ctx - [optional] user-defined Hessian context
1267 
1268    Notes:
1269    See SLESSetOperators() for important information about setting the flag
1270    output parameter in the routine func().  Be sure to read this information!
1271 
1272    The function func() takes Mat * as the matrix arguments rather than Mat.
1273    This allows the Hessian evaluation routine to replace A and/or B with a
1274    completely new new matrix structure (not just different matrix elements)
1275    when appropriate, for instance, if the nonzero structure is changing
1276    throughout the global iterations.
1277 
1278    Level: beginner
1279 
1280 .keywords: SNES, nonlinear, set, Hessian, matrix
1281 
1282 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
1283 @*/
1284 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1285 {
1286   PetscFunctionBegin;
1287   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1288   PetscValidHeaderSpecific(A,MAT_COOKIE);
1289   PetscValidHeaderSpecific(B,MAT_COOKIE);
1290   PetscCheckSameComm(snes,A);
1291   PetscCheckSameComm(snes,B);
1292   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1293     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1294   }
1295   snes->computejacobian = func;
1296   snes->jacP            = ctx;
1297   snes->jacobian        = A;
1298   snes->jacobian_pre    = B;
1299   PetscFunctionReturn(0);
1300 }
1301 
1302 #undef __FUNC__
1303 #define __FUNC__ /*<a name=""></a>*/"SNESGetHessian"
1304 /*@
1305    SNESGetHessian - Returns the Hessian matrix and optionally the user
1306    provided context for evaluating the Hessian.
1307 
1308    Not Collective, but Mat object is parallel if SNES object is parallel
1309 
1310    Input Parameter:
1311 .  snes - the nonlinear solver context
1312 
1313    Output Parameters:
1314 +  A - location to stash Hessian matrix (or PETSC_NULL)
1315 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1316 -  ctx - location to stash Hessian ctx (or PETSC_NULL)
1317 
1318    Level: advanced
1319 
1320 .seealso: SNESSetHessian(), SNESComputeHessian()
1321 
1322 .keywords: SNES, get, Hessian
1323 @*/
1324 int SNESGetHessian(SNES snes,Mat *A,Mat *B,void **ctx)
1325 {
1326   PetscFunctionBegin;
1327   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1328   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1329     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1330   }
1331   if (A)   *A = snes->jacobian;
1332   if (B)   *B = snes->jacobian_pre;
1333   if (ctx) *ctx = snes->jacP;
1334   PetscFunctionReturn(0);
1335 }
1336 
1337 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1338 
1339 #undef __FUNC__
1340 #define __FUNC__ /*<a name=""></a>*/"SNESSetUp"
1341 /*@
1342    SNESSetUp - Sets up the internal data structures for the later use
1343    of a nonlinear solver.
1344 
1345    Collective on SNES
1346 
1347    Input Parameters:
1348 +  snes - the SNES context
1349 -  x - the solution vector
1350 
1351    Notes:
1352    For basic use of the SNES solvers the user need not explicitly call
1353    SNESSetUp(), since these actions will automatically occur during
1354    the call to SNESSolve().  However, if one wishes to control this
1355    phase separately, SNESSetUp() should be called after SNESCreate()
1356    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1357 
1358    Level: advanced
1359 
1360 .keywords: SNES, nonlinear, setup
1361 
1362 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1363 @*/
1364 int SNESSetUp(SNES snes,Vec x)
1365 {
1366   int        ierr;
1367   PetscTruth flg;
1368 
1369   PetscFunctionBegin;
1370   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1371   PetscValidHeaderSpecific(x,VEC_COOKIE);
1372   PetscCheckSameComm(snes,x);
1373   snes->vec_sol = snes->vec_sol_always = x;
1374 
1375   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1376   /*
1377       This version replaces the user provided Jacobian matrix with a
1378       matrix-free version but still employs the user-provided preconditioner matrix
1379   */
1380   if (flg) {
1381     Mat J;
1382     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1383     PLogObjectParent(snes,J);
1384     snes->mfshell = J;
1385     snes->jacobian = J;
1386     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1387       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1388     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1389       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1390     } else {
1391       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option");
1392     }
1393     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1394   }
1395   ierr = OptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1396   /*
1397       This version replaces both the user-provided Jacobian and the user-
1398       provided preconditioner matrix with the default matrix free version.
1399    */
1400   if (flg) {
1401     Mat J;
1402     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1403     PLogObjectParent(snes,J);
1404     snes->mfshell = J;
1405     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1406       ierr = SNESSetJacobian(snes,J,J,0,snes->funP);CHKERRQ(ierr);
1407       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1408     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1409       ierr = SNESSetHessian(snes,J,J,0,snes->funP);CHKERRQ(ierr);
1410       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1411     } else {
1412       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option");
1413     }
1414     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1415   }
1416   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1417     PetscTruth iseqtr;
1418 
1419     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1420     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1421     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1422     if (snes->vec_func == snes->vec_sol) {
1423       SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector");
1424     }
1425 
1426     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
1427     ierr = PetscTypeCompare((PetscObject)snes,SNESEQTR,&iseqtr);CHKERRQ(ierr);
1428     if (snes->ksp_ewconv && !iseqtr) {
1429       SLES sles; KSP ksp;
1430       ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1431       ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
1432       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,(void *)snes);CHKERRQ(ierr);
1433     }
1434   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1435     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1436     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1437     if (!snes->computeumfunction) {
1438       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first");
1439     }
1440     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian()");
1441   } else {
1442     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
1443   }
1444   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
1445   snes->setupcalled = 1;
1446   PetscFunctionReturn(0);
1447 }
1448 
1449 #undef __FUNC__
1450 #define __FUNC__ /*<a name=""></a>*/"SNESDestroy"
1451 /*@C
1452    SNESDestroy - Destroys the nonlinear solver context that was created
1453    with SNESCreate().
1454 
1455    Collective on SNES
1456 
1457    Input Parameter:
1458 .  snes - the SNES context
1459 
1460    Level: beginner
1461 
1462 .keywords: SNES, nonlinear, destroy
1463 
1464 .seealso: SNESCreate(), SNESSolve()
1465 @*/
1466 int SNESDestroy(SNES snes)
1467 {
1468   int i,ierr;
1469 
1470   PetscFunctionBegin;
1471   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1472   if (--snes->refct > 0) PetscFunctionReturn(0);
1473 
1474   /* if memory was published with AMS then destroy it */
1475   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1476 
1477   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1478   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
1479   if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);}
1480   ierr = SLESDestroy(snes->sles);CHKERRQ(ierr);
1481   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1482   for (i=0; i<snes->numbermonitors; i++) {
1483     if (snes->monitordestroy[i]) {
1484       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1485     }
1486   }
1487   PLogObjectDestroy((PetscObject)snes);
1488   PetscHeaderDestroy((PetscObject)snes);
1489   PetscFunctionReturn(0);
1490 }
1491 
1492 /* ----------- Routines to set solver parameters ---------- */
1493 
1494 #undef __FUNC__
1495 #define __FUNC__ /*<a name=""></a>*/"SNESSetTolerances"
1496 /*@
1497    SNESSetTolerances - Sets various parameters used in convergence tests.
1498 
1499    Collective on SNES
1500 
1501    Input Parameters:
1502 +  snes - the SNES context
1503 .  atol - absolute convergence tolerance
1504 .  rtol - relative convergence tolerance
1505 .  stol -  convergence tolerance in terms of the norm
1506            of the change in the solution between steps
1507 .  maxit - maximum number of iterations
1508 -  maxf - maximum number of function evaluations
1509 
1510    Options Database Keys:
1511 +    -snes_atol <atol> - Sets atol
1512 .    -snes_rtol <rtol> - Sets rtol
1513 .    -snes_stol <stol> - Sets stol
1514 .    -snes_max_it <maxit> - Sets maxit
1515 -    -snes_max_funcs <maxf> - Sets maxf
1516 
1517    Notes:
1518    The default maximum number of iterations is 50.
1519    The default maximum number of function evaluations is 1000.
1520 
1521    Level: intermediate
1522 
1523 .keywords: SNES, nonlinear, set, convergence, tolerances
1524 
1525 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
1526 @*/
1527 int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
1528 {
1529   PetscFunctionBegin;
1530   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1531   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1532   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1533   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1534   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1535   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1536   PetscFunctionReturn(0);
1537 }
1538 
1539 #undef __FUNC__
1540 #define __FUNC__ /*<a name=""></a>*/"SNESGetTolerances"
1541 /*@
1542    SNESGetTolerances - Gets various parameters used in convergence tests.
1543 
1544    Not Collective
1545 
1546    Input Parameters:
1547 +  snes - the SNES context
1548 .  atol - absolute convergence tolerance
1549 .  rtol - relative convergence tolerance
1550 .  stol -  convergence tolerance in terms of the norm
1551            of the change in the solution between steps
1552 .  maxit - maximum number of iterations
1553 -  maxf - maximum number of function evaluations
1554 
1555    Notes:
1556    The user can specify PETSC_NULL for any parameter that is not needed.
1557 
1558    Level: intermediate
1559 
1560 .keywords: SNES, nonlinear, get, convergence, tolerances
1561 
1562 .seealso: SNESSetTolerances()
1563 @*/
1564 int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
1565 {
1566   PetscFunctionBegin;
1567   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1568   if (atol)  *atol  = snes->atol;
1569   if (rtol)  *rtol  = snes->rtol;
1570   if (stol)  *stol  = snes->xtol;
1571   if (maxit) *maxit = snes->max_its;
1572   if (maxf)  *maxf  = snes->max_funcs;
1573   PetscFunctionReturn(0);
1574 }
1575 
1576 #undef __FUNC__
1577 #define __FUNC__ /*<a name=""></a>*/"SNESSetTrustRegionTolerance"
1578 /*@
1579    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1580 
1581    Collective on SNES
1582 
1583    Input Parameters:
1584 +  snes - the SNES context
1585 -  tol - tolerance
1586 
1587    Options Database Key:
1588 .  -snes_trtol <tol> - Sets tol
1589 
1590    Level: intermediate
1591 
1592 .keywords: SNES, nonlinear, set, trust region, tolerance
1593 
1594 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
1595 @*/
1596 int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1597 {
1598   PetscFunctionBegin;
1599   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1600   snes->deltatol = tol;
1601   PetscFunctionReturn(0);
1602 }
1603 
1604 #undef __FUNC__
1605 #define __FUNC__ /*<a name=""></a>*/"SNESSetMinimizationFunctionTolerance"
1606 /*@
1607    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
1608    for unconstrained minimization solvers.
1609 
1610    Collective on SNES
1611 
1612    Input Parameters:
1613 +  snes - the SNES context
1614 -  ftol - minimum function tolerance
1615 
1616    Options Database Key:
1617 .  -snes_fmin <ftol> - Sets ftol
1618 
1619    Note:
1620    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
1621    methods only.
1622 
1623    Level: intermediate
1624 
1625 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
1626 
1627 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
1628 @*/
1629 int SNESSetMinimizationFunctionTolerance(SNES snes,PetscReal ftol)
1630 {
1631   PetscFunctionBegin;
1632   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1633   snes->fmin = ftol;
1634   PetscFunctionReturn(0);
1635 }
1636 /*
1637    Duplicate the lg monitors for SNES from KSP; for some reason with
1638    dynamic libraries things don't work under Sun4 if we just use
1639    macros instead of functions
1640 */
1641 #undef __FUNC__
1642 #define __FUNC__ /*<a name=""></a>*/"SNESLGMonitor"
1643 int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1644 {
1645   int ierr;
1646 
1647   PetscFunctionBegin;
1648   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1649   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1650   PetscFunctionReturn(0);
1651 }
1652 
1653 #undef __FUNC__
1654 #define __FUNC__ /*<a name=""></a>*/"SNESLGMonitorCreate"
1655 int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,DrawLG *draw)
1656 {
1657   int ierr;
1658 
1659   PetscFunctionBegin;
1660   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1661   PetscFunctionReturn(0);
1662 }
1663 
1664 #undef __FUNC__
1665 #define __FUNC__ /*<a name=""></a>*/"SNESLGMonitorDestroy"
1666 int SNESLGMonitorDestroy(DrawLG draw)
1667 {
1668   int ierr;
1669 
1670   PetscFunctionBegin;
1671   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1672   PetscFunctionReturn(0);
1673 }
1674 
1675 /* ------------ Routines to set performance monitoring options ----------- */
1676 
1677 #undef __FUNC__
1678 #define __FUNC__ /*<a name=""></a>*/"SNESSetMonitor"
1679 /*@C
1680    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
1681    iteration of the nonlinear solver to display the iteration's
1682    progress.
1683 
1684    Collective on SNES
1685 
1686    Input Parameters:
1687 +  snes - the SNES context
1688 .  func - monitoring routine
1689 .  mctx - [optional] user-defined context for private data for the
1690           monitor routine (use PETSC_NULL if no context is desitre)
1691 -  monitordestroy - [optional] routine that frees monitor context
1692           (may be PETSC_NULL)
1693 
1694    Calling sequence of func:
1695 $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1696 
1697 +    snes - the SNES context
1698 .    its - iteration number
1699 .    norm - 2-norm function value (may be estimated)
1700             for SNES_NONLINEAR_EQUATIONS methods
1701 .    norm - 2-norm gradient value (may be estimated)
1702             for SNES_UNCONSTRAINED_MINIMIZATION methods
1703 -    mctx - [optional] monitoring context
1704 
1705    Options Database Keys:
1706 +    -snes_monitor        - sets SNESDefaultMonitor()
1707 .    -snes_xmonitor       - sets line graph monitor,
1708                             uses SNESLGMonitorCreate()
1709 _    -snes_cancelmonitors - cancels all monitors that have
1710                             been hardwired into a code by
1711                             calls to SNESSetMonitor(), but
1712                             does not cancel those set via
1713                             the options database.
1714 
1715    Notes:
1716    Several different monitoring routines may be set by calling
1717    SNESSetMonitor() multiple times; all will be called in the
1718    order in which they were set.
1719 
1720    Level: intermediate
1721 
1722 .keywords: SNES, nonlinear, set, monitor
1723 
1724 .seealso: SNESDefaultMonitor(), SNESClearMonitor()
1725 @*/
1726 int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *))
1727 {
1728   PetscFunctionBegin;
1729   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1730   if (snes->numbermonitors >= MAXSNESMONITORS) {
1731     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
1732   }
1733 
1734   snes->monitor[snes->numbermonitors]           = func;
1735   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1736   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1737   PetscFunctionReturn(0);
1738 }
1739 
1740 #undef __FUNC__
1741 #define __FUNC__ /*<a name=""></a>*/"SNESClearMonitor"
1742 /*@C
1743    SNESClearMonitor - Clears all the monitor functions for a SNES object.
1744 
1745    Collective on SNES
1746 
1747    Input Parameters:
1748 .  snes - the SNES context
1749 
1750    Options Database:
1751 .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1752     into a code by calls to SNESSetMonitor(), but does not cancel those
1753     set via the options database
1754 
1755    Notes:
1756    There is no way to clear one specific monitor from a SNES object.
1757 
1758    Level: intermediate
1759 
1760 .keywords: SNES, nonlinear, set, monitor
1761 
1762 .seealso: SNESDefaultMonitor(), SNESSetMonitor()
1763 @*/
1764 int SNESClearMonitor(SNES snes)
1765 {
1766   PetscFunctionBegin;
1767   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1768   snes->numbermonitors = 0;
1769   PetscFunctionReturn(0);
1770 }
1771 
1772 #undef __FUNC__
1773 #define __FUNC__ /*<a name=""></a>*/"SNESSetConvergenceTest"
1774 /*@C
1775    SNESSetConvergenceTest - Sets the function that is to be used
1776    to test for convergence of the nonlinear iterative solution.
1777 
1778    Collective on SNES
1779 
1780    Input Parameters:
1781 +  snes - the SNES context
1782 .  func - routine to test for convergence
1783 -  cctx - [optional] context for private data for the convergence routine
1784           (may be PETSC_NULL)
1785 
1786    Calling sequence of func:
1787 $     int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1788 
1789 +    snes - the SNES context
1790 .    cctx - [optional] convergence context
1791 .    reason - reason for convergence/divergence
1792 .    xnorm - 2-norm of current iterate
1793 .    gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods)
1794 .    f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods)
1795 .    gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods)
1796 -    f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods)
1797 
1798    Level: advanced
1799 
1800 .keywords: SNES, nonlinear, set, convergence, test
1801 
1802 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1803           SNESConverged_UM_LS(), SNESConverged_UM_TR()
1804 @*/
1805 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
1806 {
1807   PetscFunctionBegin;
1808   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1809   (snes)->converged = func;
1810   (snes)->cnvP      = cctx;
1811   PetscFunctionReturn(0);
1812 }
1813 
1814 #undef __FUNC__
1815 #define __FUNC__ /*<a name=""></a>*/"SNESGetConvergedReason"
1816 /*@C
1817    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1818 
1819    Not Collective
1820 
1821    Input Parameter:
1822 .  snes - the SNES context
1823 
1824    Output Parameter:
1825 .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1826             manual pages for the individual convergence tests for complete lists
1827 
1828    Level: intermediate
1829 
1830    Notes: Can only be called after the call the SNESSolve() is complete.
1831 
1832 .keywords: SNES, nonlinear, set, convergence, test
1833 
1834 .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1835           SNESConverged_UM_LS(), SNESConverged_UM_TR()
1836 @*/
1837 int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1838 {
1839   PetscFunctionBegin;
1840   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1841   *reason = snes->reason;
1842   PetscFunctionReturn(0);
1843 }
1844 
1845 #undef __FUNC__
1846 #define __FUNC__ /*<a name=""></a>*/"SNESSetConvergenceHistory"
1847 /*@
1848    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1849 
1850    Collective on SNES
1851 
1852    Input Parameters:
1853 +  snes - iterative context obtained from SNESCreate()
1854 .  a   - array to hold history
1855 .  its - integer array holds the number of linear iterations (or
1856          negative if not converged) for each solve.
1857 .  na  - size of a and its
1858 -  reset - PETSC_TRUTH indicates each new nonlinear solve resets the history counter to zero,
1859            else it continues storing new values for new nonlinear solves after the old ones
1860 
1861    Notes:
1862    If set, this array will contain the function norms (for
1863    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1864    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1865    at each step.
1866 
1867    This routine is useful, e.g., when running a code for purposes
1868    of accurate performance monitoring, when no I/O should be done
1869    during the section of code that is being timed.
1870 
1871    Level: intermediate
1872 
1873 .keywords: SNES, set, convergence, history
1874 
1875 .seealso: SNESGetConvergenceHistory()
1876 
1877 @*/
1878 int SNESSetConvergenceHistory(SNES snes,PetscReal *a,int *its,int na,PetscTruth reset)
1879 {
1880   PetscFunctionBegin;
1881   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1882   if (na) PetscValidScalarPointer(a);
1883   snes->conv_hist       = a;
1884   snes->conv_hist_its   = its;
1885   snes->conv_hist_max   = na;
1886   snes->conv_hist_reset = reset;
1887   PetscFunctionReturn(0);
1888 }
1889 
1890 #undef __FUNC__
1891 #define __FUNC__ /*<a name=""></a>*/"SNESGetConvergenceHistory"
1892 /*@C
1893    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1894 
1895    Collective on SNES
1896 
1897    Input Parameter:
1898 .  snes - iterative context obtained from SNESCreate()
1899 
1900    Output Parameters:
1901 .  a   - array to hold history
1902 .  its - integer array holds the number of linear iterations (or
1903          negative if not converged) for each solve.
1904 -  na  - size of a and its
1905 
1906    Notes:
1907     The calling sequence for this routine in Fortran is
1908 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1909 
1910    This routine is useful, e.g., when running a code for purposes
1911    of accurate performance monitoring, when no I/O should be done
1912    during the section of code that is being timed.
1913 
1914    Level: intermediate
1915 
1916 .keywords: SNES, get, convergence, history
1917 
1918 .seealso: SNESSetConvergencHistory()
1919 
1920 @*/
1921 int SNESGetConvergenceHistory(SNES snes,PetscReal **a,int **its,int *na)
1922 {
1923   PetscFunctionBegin;
1924   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1925   if (a)   *a   = snes->conv_hist;
1926   if (its) *its = snes->conv_hist_its;
1927   if (na) *na   = snes->conv_hist_len;
1928   PetscFunctionReturn(0);
1929 }
1930 
1931 #undef __FUNC__
1932 #define __FUNC__ /*<a name=""></a>*/"SNESScaleStep_Private"
1933 /*
1934    SNESScaleStep_Private - Scales a step so that its length is less than the
1935    positive parameter delta.
1936 
1937     Input Parameters:
1938 +   snes - the SNES context
1939 .   y - approximate solution of linear system
1940 .   fnorm - 2-norm of current function
1941 -   delta - trust region size
1942 
1943     Output Parameters:
1944 +   gpnorm - predicted function norm at the new point, assuming local
1945     linearization.  The value is zero if the step lies within the trust
1946     region, and exceeds zero otherwise.
1947 -   ynorm - 2-norm of the step
1948 
1949     Note:
1950     For non-trust region methods such as SNESEQLS, the parameter delta
1951     is set to be the maximum allowable step size.
1952 
1953 .keywords: SNES, nonlinear, scale, step
1954 */
1955 int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,
1956                           PetscReal *gpnorm,PetscReal *ynorm)
1957 {
1958   PetscReal norm;
1959   Scalar cnorm;
1960   int    ierr;
1961 
1962   PetscFunctionBegin;
1963   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1964   PetscValidHeaderSpecific(y,VEC_COOKIE);
1965   PetscCheckSameComm(snes,y);
1966 
1967   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
1968   if (norm > *delta) {
1969      norm = *delta/norm;
1970      *gpnorm = (1.0 - norm)*(*fnorm);
1971      cnorm = norm;
1972      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
1973      *ynorm = *delta;
1974   } else {
1975      *gpnorm = 0.0;
1976      *ynorm = norm;
1977   }
1978   PetscFunctionReturn(0);
1979 }
1980 
1981 #undef __FUNC__
1982 #define __FUNC__ /*<a name=""></a>*/"SNESSolve"
1983 /*@
1984    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
1985    SNESCreate() and optional routines of the form SNESSetXXX().
1986 
1987    Collective on SNES
1988 
1989    Input Parameters:
1990 +  snes - the SNES context
1991 -  x - the solution vector
1992 
1993    Output Parameter:
1994 .  its - number of iterations until termination
1995 
1996    Notes:
1997    The user should initialize the vector,x, with the initial guess
1998    for the nonlinear solve prior to calling SNESSolve.  In particular,
1999    to employ an initial guess of zero, the user should explicitly set
2000    this vector to zero by calling VecSet().
2001 
2002    Level: beginner
2003 
2004 .keywords: SNES, nonlinear, solve
2005 
2006 .seealso: SNESCreate(), SNESDestroy()
2007 @*/
2008 int SNESSolve(SNES snes,Vec x,int *its)
2009 {
2010   int        ierr;
2011   PetscTruth flg;
2012 
2013   PetscFunctionBegin;
2014   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2015   PetscValidHeaderSpecific(x,VEC_COOKIE);
2016   PetscCheckSameComm(snes,x);
2017   PetscValidIntPointer(its);
2018   if (!snes->solve) SETERRQ(1,1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
2019 
2020   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
2021   else {snes->vec_sol = snes->vec_sol_always = x;}
2022   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
2023   PLogEventBegin(SNES_Solve,snes,0,0,0);
2024   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
2025   ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr);
2026   PLogEventEnd(SNES_Solve,snes,0,0,0);
2027   ierr = OptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
2028   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD);CHKERRQ(ierr); }
2029   PetscFunctionReturn(0);
2030 }
2031 
2032 /* --------- Internal routines for SNES Package --------- */
2033 
2034 #undef __FUNC__
2035 #define __FUNC__ /*<a name=""></a>*/"SNESSetType"
2036 /*@C
2037    SNESSetType - Sets the method for the nonlinear solver.
2038 
2039    Collective on SNES
2040 
2041    Input Parameters:
2042 +  snes - the SNES context
2043 -  type - a known method
2044 
2045    Options Database Key:
2046 .  -snes_type <type> - Sets the method; use -help for a list
2047    of available methods (for instance, ls or tr)
2048 
2049    Notes:
2050    See "petsc/include/petscsnes.h" for available methods (for instance)
2051 +    SNESEQLS - Newton's method with line search
2052      (systems of nonlinear equations)
2053 .    SNESEQTR - Newton's method with trust region
2054      (systems of nonlinear equations)
2055 .    SNESUMTR - Newton's method with trust region
2056      (unconstrained minimization)
2057 -    SNESUMLS - Newton's method with line search
2058      (unconstrained minimization)
2059 
2060   Normally, it is best to use the SNESSetFromOptions() command and then
2061   set the SNES solver type from the options database rather than by using
2062   this routine.  Using the options database provides the user with
2063   maximum flexibility in evaluating the many nonlinear solvers.
2064   The SNESSetType() routine is provided for those situations where it
2065   is necessary to set the nonlinear solver independently of the command
2066   line or options database.  This might be the case, for example, when
2067   the choice of solver changes during the execution of the program,
2068   and the user's application is taking responsibility for choosing the
2069   appropriate method.  In other words, this routine is not for beginners.
2070 
2071   Level: intermediate
2072 
2073 .keywords: SNES, set, type
2074 @*/
2075 int SNESSetType(SNES snes,SNESType type)
2076 {
2077   int        ierr,(*r)(SNES);
2078   PetscTruth match;
2079 
2080   PetscFunctionBegin;
2081   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2082   PetscValidCharPointer(type);
2083 
2084   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
2085   if (match) PetscFunctionReturn(0);
2086 
2087   if (snes->setupcalled) {
2088     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
2089     snes->data = 0;
2090   }
2091 
2092   /* Get the function pointers for the iterative method requested */
2093   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
2094 
2095   ierr =  FListFind(snes->comm,SNESList,type,(int (**)(void *)) &r);CHKERRQ(ierr);
2096 
2097   if (!r) SETERRQ1(1,1,"Unable to find requested SNES type %s",type);
2098 
2099   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
2100   snes->data = 0;
2101   ierr = (*r)(snes);CHKERRQ(ierr);
2102 
2103   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
2104   snes->set_method_called = 1;
2105 
2106   PetscFunctionReturn(0);
2107 }
2108 
2109 
2110 /* --------------------------------------------------------------------- */
2111 #undef __FUNC__
2112 #define __FUNC__ /*<a name=""></a>*/"SNESRegisterDestroy"
2113 /*@C
2114    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2115    registered by SNESRegisterDynamic().
2116 
2117    Not Collective
2118 
2119    Level: advanced
2120 
2121 .keywords: SNES, nonlinear, register, destroy
2122 
2123 .seealso: SNESRegisterAll(), SNESRegisterAll()
2124 @*/
2125 int SNESRegisterDestroy(void)
2126 {
2127   int ierr;
2128 
2129   PetscFunctionBegin;
2130   if (SNESList) {
2131     ierr = FListDestroy(SNESList);CHKERRQ(ierr);
2132     SNESList = 0;
2133   }
2134   SNESRegisterAllCalled = PETSC_FALSE;
2135   PetscFunctionReturn(0);
2136 }
2137 
2138 #undef __FUNC__
2139 #define __FUNC__ /*<a name=""></a>*/"SNESGetType"
2140 /*@C
2141    SNESGetType - Gets the SNES method type and name (as a string).
2142 
2143    Not Collective
2144 
2145    Input Parameter:
2146 .  snes - nonlinear solver context
2147 
2148    Output Parameter:
2149 .  type - SNES method (a charactor string)
2150 
2151    Level: intermediate
2152 
2153 .keywords: SNES, nonlinear, get, type, name
2154 @*/
2155 int SNESGetType(SNES snes,SNESType *type)
2156 {
2157   PetscFunctionBegin;
2158   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2159   *type = snes->type_name;
2160   PetscFunctionReturn(0);
2161 }
2162 
2163 #undef __FUNC__
2164 #define __FUNC__ /*<a name=""></a>*/"SNESGetSolution"
2165 /*@C
2166    SNESGetSolution - Returns the vector where the approximate solution is
2167    stored.
2168 
2169    Not Collective, but Vec is parallel if SNES is parallel
2170 
2171    Input Parameter:
2172 .  snes - the SNES context
2173 
2174    Output Parameter:
2175 .  x - the solution
2176 
2177    Level: advanced
2178 
2179 .keywords: SNES, nonlinear, get, solution
2180 
2181 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
2182 @*/
2183 int SNESGetSolution(SNES snes,Vec *x)
2184 {
2185   PetscFunctionBegin;
2186   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2187   *x = snes->vec_sol_always;
2188   PetscFunctionReturn(0);
2189 }
2190 
2191 #undef __FUNC__
2192 #define __FUNC__ /*<a name=""></a>*/"SNESGetSolutionUpdate"
2193 /*@C
2194    SNESGetSolutionUpdate - Returns the vector where the solution update is
2195    stored.
2196 
2197    Not Collective, but Vec is parallel if SNES is parallel
2198 
2199    Input Parameter:
2200 .  snes - the SNES context
2201 
2202    Output Parameter:
2203 .  x - the solution update
2204 
2205    Level: advanced
2206 
2207 .keywords: SNES, nonlinear, get, solution, update
2208 
2209 .seealso: SNESGetSolution(), SNESGetFunction
2210 @*/
2211 int SNESGetSolutionUpdate(SNES snes,Vec *x)
2212 {
2213   PetscFunctionBegin;
2214   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2215   *x = snes->vec_sol_update_always;
2216   PetscFunctionReturn(0);
2217 }
2218 
2219 #undef __FUNC__
2220 #define __FUNC__ /*<a name=""></a>*/"SNESGetFunction"
2221 /*@C
2222    SNESGetFunction - Returns the vector where the function is stored.
2223 
2224    Not Collective, but Vec is parallel if SNES is parallel
2225 
2226    Input Parameter:
2227 .  snes - the SNES context
2228 
2229    Output Parameter:
2230 +  r - the function (or PETSC_NULL)
2231 .  ctx - the function context (or PETSC_NULL)
2232 -  func - the function (or PETSC_NULL)
2233 
2234    Notes:
2235    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
2236    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
2237    SNESGetMinimizationFunction() and SNESGetGradient();
2238 
2239    Level: advanced
2240 
2241 .keywords: SNES, nonlinear, get, function
2242 
2243 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
2244           SNESGetGradient()
2245 
2246 @*/
2247 int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*))
2248 {
2249   PetscFunctionBegin;
2250   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2251   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
2252     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
2253   }
2254   if (r)    *r    = snes->vec_func_always;
2255   if (ctx)  *ctx  = snes->funP;
2256   if (func) *func = snes->computefunction;
2257   PetscFunctionReturn(0);
2258 }
2259 
2260 #undef __FUNC__
2261 #define __FUNC__ /*<a name=""></a>*/"SNESGetGradient"
2262 /*@C
2263    SNESGetGradient - Returns the vector where the gradient is stored.
2264 
2265    Not Collective, but Vec is parallel if SNES is parallel
2266 
2267    Input Parameter:
2268 .  snes - the SNES context
2269 
2270    Output Parameter:
2271 +  r - the gradient (or PETSC_NULL)
2272 -  ctx - the gradient context (or PETSC_NULL)
2273 
2274    Notes:
2275    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
2276    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
2277    SNESGetFunction().
2278 
2279    Level: advanced
2280 
2281 .keywords: SNES, nonlinear, get, gradient
2282 
2283 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction(),
2284           SNESSetGradient(), SNESSetFunction()
2285 
2286 @*/
2287 int SNESGetGradient(SNES snes,Vec *r,void **ctx)
2288 {
2289   PetscFunctionBegin;
2290   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2291   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
2292     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
2293   }
2294   if (r)   *r = snes->vec_func_always;
2295   if (ctx) *ctx = snes->funP;
2296   PetscFunctionReturn(0);
2297 }
2298 
2299 #undef __FUNC__
2300 #define __FUNC__ /*<a name=""></a>*/"SNESGetMinimizationFunction"
2301 /*@C
2302    SNESGetMinimizationFunction - Returns the scalar function value for
2303    unconstrained minimization problems.
2304 
2305    Not Collective
2306 
2307    Input Parameter:
2308 .  snes - the SNES context
2309 
2310    Output Parameter:
2311 +  r - the function (or PETSC_NULL)
2312 -  ctx - the function context (or PETSC_NULL)
2313 
2314    Notes:
2315    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
2316    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
2317    SNESGetFunction().
2318 
2319    Level: advanced
2320 
2321 .keywords: SNES, nonlinear, get, function
2322 
2323 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction(), SNESSetFunction()
2324 
2325 @*/
2326 int SNESGetMinimizationFunction(SNES snes,PetscReal *r,void **ctx)
2327 {
2328   PetscFunctionBegin;
2329   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2330   PetscValidScalarPointer(r);
2331   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
2332     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
2333   }
2334   if (r)   *r = snes->fc;
2335   if (ctx) *ctx = snes->umfunP;
2336   PetscFunctionReturn(0);
2337 }
2338 
2339 #undef __FUNC__
2340 #define __FUNC__ /*<a name=""></a>*/"SNESSetOptionsPrefix"
2341 /*@C
2342    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2343    SNES options in the database.
2344 
2345    Collective on SNES
2346 
2347    Input Parameter:
2348 +  snes - the SNES context
2349 -  prefix - the prefix to prepend to all option names
2350 
2351    Notes:
2352    A hyphen (-) must NOT be given at the beginning of the prefix name.
2353    The first character of all runtime options is AUTOMATICALLY the hyphen.
2354 
2355    Level: advanced
2356 
2357 .keywords: SNES, set, options, prefix, database
2358 
2359 .seealso: SNESSetFromOptions()
2360 @*/
2361 int SNESSetOptionsPrefix(SNES snes,char *prefix)
2362 {
2363   int ierr;
2364 
2365   PetscFunctionBegin;
2366   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2367   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2368   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2369   PetscFunctionReturn(0);
2370 }
2371 
2372 #undef __FUNC__
2373 #define __FUNC__ /*<a name=""></a>*/"SNESAppendOptionsPrefix"
2374 /*@C
2375    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2376    SNES options in the database.
2377 
2378    Collective on SNES
2379 
2380    Input Parameters:
2381 +  snes - the SNES context
2382 -  prefix - the prefix to prepend to all option names
2383 
2384    Notes:
2385    A hyphen (-) must NOT be given at the beginning of the prefix name.
2386    The first character of all runtime options is AUTOMATICALLY the hyphen.
2387 
2388    Level: advanced
2389 
2390 .keywords: SNES, append, options, prefix, database
2391 
2392 .seealso: SNESGetOptionsPrefix()
2393 @*/
2394 int SNESAppendOptionsPrefix(SNES snes,char *prefix)
2395 {
2396   int ierr;
2397 
2398   PetscFunctionBegin;
2399   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2400   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2401   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2402   PetscFunctionReturn(0);
2403 }
2404 
2405 #undef __FUNC__
2406 #define __FUNC__ /*<a name=""></a>*/"SNESGetOptionsPrefix"
2407 /*@C
2408    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2409    SNES options in the database.
2410 
2411    Not Collective
2412 
2413    Input Parameter:
2414 .  snes - the SNES context
2415 
2416    Output Parameter:
2417 .  prefix - pointer to the prefix string used
2418 
2419    Notes: On the fortran side, the user should pass in a string 'prifix' of
2420    sufficient length to hold the prefix.
2421 
2422    Level: advanced
2423 
2424 .keywords: SNES, get, options, prefix, database
2425 
2426 .seealso: SNESAppendOptionsPrefix()
2427 @*/
2428 int SNESGetOptionsPrefix(SNES snes,char **prefix)
2429 {
2430   int ierr;
2431 
2432   PetscFunctionBegin;
2433   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2434   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2435   PetscFunctionReturn(0);
2436 }
2437 
2438 #undef __FUNC__
2439 #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp"
2440 /*@
2441    SNESPrintHelp - Prints all options for the SNES component.
2442 
2443    Collective on SNES
2444 
2445    Input Parameter:
2446 .  snes - the SNES context
2447 
2448    Options Database Keys:
2449 +  -help - Prints SNES options
2450 -  -h - Prints SNES options
2451 
2452    Level: beginner
2453 
2454 .keywords: SNES, nonlinear, help
2455 
2456 .seealso: SNESSetFromOptions()
2457 @*/
2458 int SNESPrintHelp(SNES snes)
2459 {
2460   char                p[64];
2461   SNES_KSP_EW_ConvCtx *kctx;
2462   int                 ierr;
2463 
2464   PetscFunctionBegin;
2465   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2466 
2467   ierr = PetscStrcpy(p,"-");CHKERRQ(ierr);
2468   if (snes->prefix) PetscStrcat(p,snes->prefix);
2469 
2470   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2471 
2472   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
2473   ierr = (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");CHKERRQ(ierr);
2474   ierr = FListPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr);
2475   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);CHKERRQ(ierr);
2476   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);CHKERRQ(ierr);
2477   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);CHKERRQ(ierr);
2478   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);CHKERRQ(ierr);
2479   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);CHKERRQ(ierr);
2480   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);CHKERRQ(ierr);
2481   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);CHKERRQ(ierr);
2482   ierr = (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");CHKERRQ(ierr);
2483   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);CHKERRQ(ierr);
2484   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
2485     residual norm at each iteration.\n",p);CHKERRQ(ierr);
2486   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
2487     residual norm for small residual norms. This is useful to conceal\n\
2488     meaningless digits that may be different on different machines.\n",p);CHKERRQ(ierr);
2489   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);CHKERRQ(ierr);
2490   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_vecmonitor: plots solution at each iteration \n",p);CHKERRQ(ierr);
2491   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_vecmonitor_update: plots update to solution at each iteration \n",p);CHKERRQ(ierr);
2492   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
2493     ierr = (*PetscHelpPrintf)(snes->comm,
2494      " Options for solving systems of nonlinear equations only:\n");CHKERRQ(ierr);
2495     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);CHKERRQ(ierr);
2496     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);CHKERRQ(ierr);
2497     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);CHKERRQ(ierr);
2498     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_ksp_monitor - if using matrix-free multiply then prints h at each KSP iteration\n",p);CHKERRQ(ierr);
2499     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);CHKERRQ(ierr);
2500     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);CHKERRQ(ierr);
2501     ierr = (*PetscHelpPrintf)(snes->comm,
2502      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);CHKERRQ(ierr);
2503     ierr = (*PetscHelpPrintf)(snes->comm,
2504      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);CHKERRQ(ierr);
2505     ierr = (*PetscHelpPrintf)(snes->comm,
2506      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);CHKERRQ(ierr);
2507     ierr = (*PetscHelpPrintf)(snes->comm,
2508      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);CHKERRQ(ierr);
2509     ierr = (*PetscHelpPrintf)(snes->comm,
2510      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);CHKERRQ(ierr);
2511     ierr = (*PetscHelpPrintf)(snes->comm,
2512      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);CHKERRQ(ierr);
2513     ierr = (*PetscHelpPrintf)(snes->comm,
2514      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);CHKERRQ(ierr);
2515   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
2516     ierr = (*PetscHelpPrintf)(snes->comm," Options for solving unconstrained minimization problems only:\n");CHKERRQ(ierr);
2517     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);CHKERRQ(ierr);
2518     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);CHKERRQ(ierr);
2519     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);CHKERRQ(ierr);
2520   }
2521   ierr = (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <type> for help on ",p);CHKERRQ(ierr);
2522   ierr = (*PetscHelpPrintf)(snes->comm,"a particular method\n");CHKERRQ(ierr);
2523   if (snes->printhelp) {
2524     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2525   }
2526   PetscFunctionReturn(0);
2527 }
2528 
2529 /*MC
2530    SNESRegisterDynamic - Adds a method to the nonlinear solver package.
2531 
2532    Synopsis:
2533    SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
2534 
2535    Not collective
2536 
2537    Input Parameters:
2538 +  name_solver - name of a new user-defined solver
2539 .  path - path (either absolute or relative) the library containing this solver
2540 .  name_create - name of routine to create method context
2541 -  routine_create - routine to create method context
2542 
2543    Notes:
2544    SNESRegisterDynamic() may be called multiple times to add several user-defined solvers.
2545 
2546    If dynamic libraries are used, then the fourth input argument (routine_create)
2547    is ignored.
2548 
2549    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LDIR}, ${BOPT},
2550    and others of the form ${any_environmental_variable} occuring in pathname will be
2551    replaced with appropriate values.
2552 
2553    Sample usage:
2554 .vb
2555    SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2556                 "MySolverCreate",MySolverCreate);
2557 .ve
2558 
2559    Then, your solver can be chosen with the procedural interface via
2560 $     SNESSetType(snes,"my_solver")
2561    or at runtime via the option
2562 $     -snes_type my_solver
2563 
2564    Level: advanced
2565 
2566 .keywords: SNES, nonlinear, register
2567 
2568 .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2569 M*/
2570 
2571 #undef __FUNC__
2572 #define __FUNC__ /*<a name=""></a>*/"SNESRegister"
2573 int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES))
2574 {
2575   char fullname[256];
2576   int  ierr;
2577 
2578   PetscFunctionBegin;
2579   ierr = FListConcat(path,name,fullname);CHKERRQ(ierr);
2580   ierr = FListAdd(&SNESList,sname,fullname,(int (*)(void*))function);CHKERRQ(ierr);
2581   PetscFunctionReturn(0);
2582 }
2583