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