xref: /petsc/src/snes/interface/snes.c (revision e24b481b4a15e258c09dabebcc19da3246e723e1)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: snes.c,v 1.160 1998/10/24 02:25:31 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 
554   /* publish children of object */
555   ierr = PetscObjectPublish((PetscObject)v->sles);CHKERRQ(ierr);
556 
557 #else
558   PetscFunctionBegin;
559 #endif
560   PetscFunctionReturn(0);
561 }
562 
563 /* -----------------------------------------------------------*/
564 #undef __FUNC__
565 #define __FUNC__ "SNESCreate"
566 /*@C
567    SNESCreate - Creates a nonlinear solver context.
568 
569    Collective on MPI_Comm
570 
571    Input Parameters:
572 +  comm - MPI communicator
573 -  type - type of method, either
574    SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations)
575    or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization)
576 
577    Output Parameter:
578 .  outsnes - the new SNES context
579 
580    Options Database Keys:
581 +   -snes_mf - Activates default matrix-free Jacobian-vector products,
582                and no preconditioning matrix
583 .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
584                products, and a user-provided preconditioning matrix
585                as set by SNESSetJacobian()
586 -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
587 
588 .keywords: SNES, nonlinear, create, context
589 
590 .seealso: SNESSolve(), SNESDestroy()
591 @*/
592 int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
593 {
594   int                 ierr;
595   SNES                snes;
596   SNES_KSP_EW_ConvCtx *kctx;
597 
598   PetscFunctionBegin;
599   *outsnes = 0;
600   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
601     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
602   }
603   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,comm,SNESDestroy,SNESView);
604   PLogObjectCreate(snes);
605   snes->bops->publish     = SNESPublish_Petsc;
606   snes->max_its           = 50;
607   snes->max_funcs	  = 10000;
608   snes->norm		  = 0.0;
609   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
610     snes->rtol		  = 1.e-8;
611     snes->ttol            = 0.0;
612     snes->atol		  = 1.e-10;
613   } else {
614     snes->rtol		  = 1.e-8;
615     snes->ttol            = 0.0;
616     snes->atol		  = 1.e-50;
617   }
618   snes->xtol		  = 1.e-8;
619   snes->trunctol	  = 1.e-12; /* no longer used */
620   snes->nfuncs            = 0;
621   snes->nfailures         = 0;
622   snes->linear_its        = 0;
623   snes->numbermonitors    = 0;
624   snes->data              = 0;
625   snes->view              = 0;
626   snes->computeumfunction = 0;
627   snes->umfunP            = 0;
628   snes->fc                = 0;
629   snes->deltatol          = 1.e-12;
630   snes->fmin              = -1.e30;
631   snes->method_class      = type;
632   snes->set_method_called = 0;
633   snes->setupcalled      = 0;
634   snes->ksp_ewconv        = 0;
635   snes->xmonitor          = 0;
636   snes->vwork             = 0;
637   snes->nwork             = 0;
638   snes->conv_hist_size    = 0;
639   snes->conv_act_size     = 0;
640   snes->conv_hist         = 0;
641 
642   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
643   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
644   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
645   snes->kspconvctx  = (void*)kctx;
646   kctx->version     = 2;
647   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
648                              this was too large for some test cases */
649   kctx->rtol_last   = 0;
650   kctx->rtol_max    = .9;
651   kctx->gamma       = 1.0;
652   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
653   kctx->alpha       = kctx->alpha2;
654   kctx->threshold   = .1;
655   kctx->lresid_last = 0;
656   kctx->norm_last   = 0;
657 
658   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
659   PLogObjectParent(snes,snes->sles)
660 
661   *outsnes = snes;
662   PetscPublishAll(snes);
663   PetscFunctionReturn(0);
664 }
665 
666 /* --------------------------------------------------------------- */
667 #undef __FUNC__
668 #define __FUNC__ "SNESSetFunction"
669 /*@C
670    SNESSetFunction - Sets the function evaluation routine and function
671    vector for use by the SNES routines in solving systems of nonlinear
672    equations.
673 
674    Collective on SNES
675 
676    Input Parameters:
677 +  snes - the SNES context
678 .  func - function evaluation routine
679 .  r - vector to store function value
680 -  ctx - [optional] user-defined context for private data for the
681          function evaluation routine (may be PETSC_NULL)
682 
683    Calling sequence of func:
684 $    func (SNES snes,Vec x,Vec f,void *ctx);
685 
686 .  f - function vector
687 -  ctx - optional user-defined function context
688 
689    Notes:
690    The Newton-like methods typically solve linear systems of the form
691 $      f'(x) x = -f(x),
692    where f'(x) denotes the Jacobian matrix and f(x) is the function.
693 
694    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
695    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
696    SNESSetMinimizationFunction() and SNESSetGradient();
697 
698 .keywords: SNES, nonlinear, set, function
699 
700 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
701 @*/
702 int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
703 {
704   PetscFunctionBegin;
705   PetscValidHeaderSpecific(snes,SNES_COOKIE);
706   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
707     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
708   }
709   snes->computefunction     = func;
710   snes->vec_func            = snes->vec_func_always = r;
711   snes->funP                = ctx;
712   PetscFunctionReturn(0);
713 }
714 
715 #undef __FUNC__
716 #define __FUNC__ "SNESComputeFunction"
717 /*@
718    SNESComputeFunction - Computes the function that has been set with
719    SNESSetFunction().
720 
721    Collective on SNES
722 
723    Input Parameters:
724 +  snes - the SNES context
725 -  x - input vector
726 
727    Output Parameter:
728 .  y - function vector, as set by SNESSetFunction()
729 
730    Notes:
731    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
732    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
733    SNESComputeMinimizationFunction() and SNESComputeGradient();
734 
735 .keywords: SNES, nonlinear, compute, function
736 
737 .seealso: SNESSetFunction(), SNESGetFunction()
738 @*/
739 int SNESComputeFunction(SNES snes,Vec x, Vec y)
740 {
741   int    ierr;
742 
743   PetscFunctionBegin;
744   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
745     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
746   }
747   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
748   PetscStackPush("SNES user function");
749   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
750   PetscStackPop;
751   snes->nfuncs++;
752   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
753   PetscFunctionReturn(0);
754 }
755 
756 #undef __FUNC__
757 #define __FUNC__ "SNESSetMinimizationFunction"
758 /*@C
759    SNESSetMinimizationFunction - Sets the function evaluation routine for
760    unconstrained minimization.
761 
762    Collective on SNES
763 
764    Input Parameters:
765 +  snes - the SNES context
766 .  func - function evaluation routine
767 -  ctx - [optional] user-defined context for private data for the
768          function evaluation routine (may be PETSC_NULL)
769 
770    Calling sequence of func:
771 $     func (SNES snes,Vec x,double *f,void *ctx);
772 
773 +  x - input vector
774 .  f - function
775 -  ctx - [optional] user-defined function context
776 
777    Notes:
778    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
779    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
780    SNESSetFunction().
781 
782 .keywords: SNES, nonlinear, set, minimization, function
783 
784 .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
785            SNESSetHessian(), SNESSetGradient()
786 @*/
787 int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
788                       void *ctx)
789 {
790   PetscFunctionBegin;
791   PetscValidHeaderSpecific(snes,SNES_COOKIE);
792   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
793     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
794   }
795   snes->computeumfunction   = func;
796   snes->umfunP              = ctx;
797   PetscFunctionReturn(0);
798 }
799 
800 #undef __FUNC__
801 #define __FUNC__ "SNESComputeMinimizationFunction"
802 /*@
803    SNESComputeMinimizationFunction - Computes the function that has been
804    set with SNESSetMinimizationFunction().
805 
806    Collective on SNES
807 
808    Input Parameters:
809 +  snes - the SNES context
810 -  x - input vector
811 
812    Output Parameter:
813 .  y - function value
814 
815    Notes:
816    SNESComputeMinimizationFunction() is valid only for
817    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
818    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
819 
820 .keywords: SNES, nonlinear, compute, minimization, function
821 
822 .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
823           SNESComputeGradient(), SNESComputeHessian()
824 @*/
825 int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
826 {
827   int    ierr;
828 
829   PetscFunctionBegin;
830   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
831     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
832   }
833   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
834   PetscStackPush("SNES user minimzation function");
835   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
836   PetscStackPop;
837   snes->nfuncs++;
838   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
839   PetscFunctionReturn(0);
840 }
841 
842 #undef __FUNC__
843 #define __FUNC__ "SNESSetGradient"
844 /*@C
845    SNESSetGradient - Sets the gradient evaluation routine and gradient
846    vector for use by the SNES routines.
847 
848    Collective on SNES
849 
850    Input Parameters:
851 +  snes - the SNES context
852 .  func - function evaluation routine
853 .  ctx - optional user-defined context for private data for the
854          gradient evaluation routine (may be PETSC_NULL)
855 -  r - vector to store gradient value
856 
857    Calling sequence of func:
858 $     func (SNES, Vec x, Vec g, void *ctx);
859 
860 +  x - input vector
861 .  g - gradient vector
862 -  ctx - optional user-defined gradient context
863 
864    Notes:
865    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
866    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
867    SNESSetFunction().
868 
869 .keywords: SNES, nonlinear, set, function
870 
871 .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
872           SNESSetMinimizationFunction(),
873 @*/
874 int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
875 {
876   PetscFunctionBegin;
877   PetscValidHeaderSpecific(snes,SNES_COOKIE);
878   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
879     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
880   }
881   snes->computefunction     = func;
882   snes->vec_func            = snes->vec_func_always = r;
883   snes->funP                = ctx;
884   PetscFunctionReturn(0);
885 }
886 
887 #undef __FUNC__
888 #define __FUNC__ "SNESComputeGradient"
889 /*@
890    SNESComputeGradient - Computes the gradient that has been set with
891    SNESSetGradient().
892 
893    Collective on SNES
894 
895    Input Parameters:
896 +  snes - the SNES context
897 -  x - input vector
898 
899    Output Parameter:
900 .  y - gradient vector
901 
902    Notes:
903    SNESComputeGradient() is valid only for
904    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
905    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
906 
907 .keywords: SNES, nonlinear, compute, gradient
908 
909 .seealso:  SNESSetGradient(), SNESGetGradient(),
910            SNESComputeMinimizationFunction(), SNESComputeHessian()
911 @*/
912 int SNESComputeGradient(SNES snes,Vec x, Vec y)
913 {
914   int    ierr;
915 
916   PetscFunctionBegin;
917   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
918     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
919   }
920   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
921   PetscStackPush("SNES user gradient function");
922   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
923   PetscStackPop;
924   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
925   PetscFunctionReturn(0);
926 }
927 
928 #undef __FUNC__
929 #define __FUNC__ "SNESComputeJacobian"
930 /*@
931    SNESComputeJacobian - Computes the Jacobian matrix that has been
932    set with SNESSetJacobian().
933 
934    Collective on SNES and Mat
935 
936    Input Parameters:
937 +  snes - the SNES context
938 -  x - input vector
939 
940    Output Parameters:
941 +  A - Jacobian matrix
942 .  B - optional preconditioning matrix
943 -  flag - flag indicating matrix structure
944 
945    Notes:
946    Most users should not need to explicitly call this routine, as it
947    is used internally within the nonlinear solvers.
948 
949    See SLESSetOperators() for important information about setting the
950    flag parameter.
951 
952    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
953    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
954    methods is SNESComputeHessian().
955 
956 .keywords: SNES, compute, Jacobian, matrix
957 
958 .seealso:  SNESSetJacobian(), SLESSetOperators()
959 @*/
960 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
961 {
962   int    ierr;
963 
964   PetscFunctionBegin;
965   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
966     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
967   }
968   if (!snes->computejacobian) PetscFunctionReturn(0);
969   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
970   *flg = DIFFERENT_NONZERO_PATTERN;
971   PetscStackPush("SNES user Jacobian function");
972   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
973   PetscStackPop;
974   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
975   /* make sure user returned a correct Jacobian and preconditioner */
976   PetscValidHeaderSpecific(*A,MAT_COOKIE);
977   PetscValidHeaderSpecific(*B,MAT_COOKIE);
978   PetscFunctionReturn(0);
979 }
980 
981 #undef __FUNC__
982 #define __FUNC__ "SNESComputeHessian"
983 /*@
984    SNESComputeHessian - Computes the Hessian matrix that has been
985    set with SNESSetHessian().
986 
987    Collective on SNES and Mat
988 
989    Input Parameters:
990 +  snes - the SNES context
991 -  x - input vector
992 
993    Output Parameters:
994 +  A - Hessian matrix
995 .  B - optional preconditioning matrix
996 -  flag - flag indicating matrix structure
997 
998    Notes:
999    Most users should not need to explicitly call this routine, as it
1000    is used internally within the nonlinear solvers.
1001 
1002    See SLESSetOperators() for important information about setting the
1003    flag parameter.
1004 
1005    SNESComputeHessian() is valid only for
1006    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
1007    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
1008 
1009 .keywords: SNES, compute, Hessian, matrix
1010 
1011 .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
1012            SNESComputeMinimizationFunction()
1013 @*/
1014 int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
1015 {
1016   int    ierr;
1017 
1018   PetscFunctionBegin;
1019   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1020     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1021   }
1022   if (!snes->computejacobian) PetscFunctionReturn(0);
1023   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
1024   *flag = DIFFERENT_NONZERO_PATTERN;
1025   PetscStackPush("SNES user Hessian function");
1026   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
1027   PetscStackPop;
1028   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
1029   /* make sure user returned a correct Jacobian and preconditioner */
1030   PetscValidHeaderSpecific(*A,MAT_COOKIE);
1031   PetscValidHeaderSpecific(*B,MAT_COOKIE);
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 #undef __FUNC__
1036 #define __FUNC__ "SNESSetJacobian"
1037 /*@C
1038    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1039    location to store the matrix.
1040 
1041    Collective on SNES and Mat
1042 
1043    Input Parameters:
1044 +  snes - the SNES context
1045 .  A - Jacobian matrix
1046 .  B - preconditioner matrix (usually same as the Jacobian)
1047 .  func - Jacobian evaluation routine
1048 -  ctx - [optional] user-defined context for private data for the
1049          Jacobian evaluation routine (may be PETSC_NULL)
1050 
1051    Calling sequence of func:
1052 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1053 
1054 +  x - input vector
1055 .  A - Jacobian matrix
1056 .  B - preconditioner matrix, usually the same as A
1057 .  flag - flag indicating information about the preconditioner matrix
1058    structure (same as flag in SLESSetOperators())
1059 -  ctx - [optional] user-defined Jacobian context
1060 
1061    Notes:
1062    See SLESSetOperators() for important information about setting the flag
1063    output parameter in the routine func().  Be sure to read this information!
1064 
1065    The routine func() takes Mat * as the matrix arguments rather than Mat.
1066    This allows the Jacobian evaluation routine to replace A and/or B with a
1067    completely new new matrix structure (not just different matrix elements)
1068    when appropriate, for instance, if the nonzero structure is changing
1069    throughout the global iterations.
1070 
1071 .keywords: SNES, nonlinear, set, Jacobian, matrix
1072 
1073 .seealso: SLESSetOperators(), SNESSetFunction()
1074 @*/
1075 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
1076                     MatStructure*,void*),void *ctx)
1077 {
1078   PetscFunctionBegin;
1079   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1080   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1081     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1082   }
1083   snes->computejacobian = func;
1084   snes->jacP            = ctx;
1085   snes->jacobian        = A;
1086   snes->jacobian_pre    = B;
1087   PetscFunctionReturn(0);
1088 }
1089 
1090 #undef __FUNC__
1091 #define __FUNC__ "SNESGetJacobian"
1092 /*@
1093    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1094    provided context for evaluating the Jacobian.
1095 
1096    Not Collective, but Mat object will be parallel if SNES object is
1097 
1098    Input Parameter:
1099 .  snes - the nonlinear solver context
1100 
1101    Output Parameters:
1102 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1103 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1104 -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1105 
1106 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1107 @*/
1108 int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
1109 {
1110   PetscFunctionBegin;
1111   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1112     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1113   }
1114   if (A)   *A = snes->jacobian;
1115   if (B)   *B = snes->jacobian_pre;
1116   if (ctx) *ctx = snes->jacP;
1117   PetscFunctionReturn(0);
1118 }
1119 
1120 #undef __FUNC__
1121 #define __FUNC__ "SNESSetHessian"
1122 /*@C
1123    SNESSetHessian - Sets the function to compute Hessian as well as the
1124    location to store the matrix.
1125 
1126    Collective on SNES and Mat
1127 
1128    Input Parameters:
1129 +  snes - the SNES context
1130 .  A - Hessian matrix
1131 .  B - preconditioner matrix (usually same as the Hessian)
1132 .  func - Jacobian evaluation routine
1133 -  ctx - [optional] user-defined context for private data for the
1134          Hessian evaluation routine (may be PETSC_NULL)
1135 
1136    Calling sequence of func:
1137 $    func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1138 
1139 +  x - input vector
1140 .  A - Hessian matrix
1141 .  B - preconditioner matrix, usually the same as A
1142 .  flag - flag indicating information about the preconditioner matrix
1143    structure (same as flag in SLESSetOperators())
1144 -  ctx - [optional] user-defined Hessian context
1145 
1146    Notes:
1147    See SLESSetOperators() for important information about setting the flag
1148    output parameter in the routine func().  Be sure to read this information!
1149 
1150    The function func() takes Mat * as the matrix arguments rather than Mat.
1151    This allows the Hessian evaluation routine to replace A and/or B with a
1152    completely new new matrix structure (not just different matrix elements)
1153    when appropriate, for instance, if the nonzero structure is changing
1154    throughout the global iterations.
1155 
1156 .keywords: SNES, nonlinear, set, Hessian, matrix
1157 
1158 .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
1159 @*/
1160 int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
1161                     MatStructure*,void*),void *ctx)
1162 {
1163   PetscFunctionBegin;
1164   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1165   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1166     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1167   }
1168   snes->computejacobian = func;
1169   snes->jacP            = ctx;
1170   snes->jacobian        = A;
1171   snes->jacobian_pre    = B;
1172   PetscFunctionReturn(0);
1173 }
1174 
1175 #undef __FUNC__
1176 #define __FUNC__ "SNESGetHessian"
1177 /*@
1178    SNESGetHessian - Returns the Hessian matrix and optionally the user
1179    provided context for evaluating the Hessian.
1180 
1181    Not Collective, but Mat object is parallel if SNES object is parallel
1182 
1183    Input Parameter:
1184 .  snes - the nonlinear solver context
1185 
1186    Output Parameters:
1187 +  A - location to stash Hessian matrix (or PETSC_NULL)
1188 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1189 -  ctx - location to stash Hessian ctx (or PETSC_NULL)
1190 
1191 .seealso: SNESSetHessian(), SNESComputeHessian()
1192 
1193 .keywords: SNES, get, Hessian
1194 @*/
1195 int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
1196 {
1197   PetscFunctionBegin;
1198   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1199     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1200   }
1201   if (A)   *A = snes->jacobian;
1202   if (B)   *B = snes->jacobian_pre;
1203   if (ctx) *ctx = snes->jacP;
1204   PetscFunctionReturn(0);
1205 }
1206 
1207 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1208 
1209 #undef __FUNC__
1210 #define __FUNC__ "SNESSetUp"
1211 /*@
1212    SNESSetUp - Sets up the internal data structures for the later use
1213    of a nonlinear solver.
1214 
1215    Collective on SNES
1216 
1217    Input Parameters:
1218 +  snes - the SNES context
1219 -  x - the solution vector
1220 
1221    Notes:
1222    For basic use of the SNES solvers the user need not explicitly call
1223    SNESSetUp(), since these actions will automatically occur during
1224    the call to SNESSolve().  However, if one wishes to control this
1225    phase separately, SNESSetUp() should be called after SNESCreate()
1226    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1227 
1228 .keywords: SNES, nonlinear, setup
1229 
1230 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1231 @*/
1232 int SNESSetUp(SNES snes,Vec x)
1233 {
1234   int ierr, flg;
1235 
1236   PetscFunctionBegin;
1237   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1238   PetscValidHeaderSpecific(x,VEC_COOKIE);
1239   snes->vec_sol = snes->vec_sol_always = x;
1240 
1241   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1242   /*
1243       This version replaces the user provided Jacobian matrix with a
1244       matrix-free version but still employs the user-provided preconditioner matrix
1245   */
1246   if (flg) {
1247     Mat J;
1248     ierr = MatCreateSNESFDMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1249     PLogObjectParent(snes,J);
1250     snes->mfshell = J;
1251     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1252       snes->jacobian = J;
1253       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1254     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1255       snes->jacobian = J;
1256       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1257     } else {
1258       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option");
1259     }
1260     ierr = MatSNESFDMFSetFromOptions(J);CHKERRQ(ierr);
1261   }
1262   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1263   /*
1264       This version replaces both the user-provided Jacobian and the user-
1265       provided preconditioner matrix with the default matrix free version.
1266    */
1267   if (flg) {
1268     Mat J;
1269     ierr = MatCreateSNESFDMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1270     PLogObjectParent(snes,J);
1271     snes->mfshell = J;
1272     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1273       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
1274       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1275     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1276       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
1277       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1278     } else {
1279       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option");
1280     }
1281     ierr = MatSNESFDMFSetFromOptions(J);CHKERRQ(ierr);
1282   }
1283   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1284     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1285     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1286     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first");
1287     if (snes->vec_func == snes->vec_sol) {
1288       SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector");
1289     }
1290 
1291     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
1292     if (snes->ksp_ewconv && PetscStrcmp(snes->type_name,SNES_EQ_TR)) {
1293       SLES sles; KSP ksp;
1294       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1295       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1296       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1297              (void *)snes); CHKERRQ(ierr);
1298     }
1299   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1300     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1301     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1302     if (!snes->computeumfunction) {
1303       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first");
1304     }
1305     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first");
1306   } else {
1307     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
1308   }
1309   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1310   snes->setupcalled = 1;
1311   PetscFunctionReturn(0);
1312 }
1313 
1314 #undef __FUNC__
1315 #define __FUNC__ "SNESDestroy"
1316 /*@C
1317    SNESDestroy - Destroys the nonlinear solver context that was created
1318    with SNESCreate().
1319 
1320    Collective on SNES
1321 
1322    Input Parameter:
1323 .  snes - the SNES context
1324 
1325 .keywords: SNES, nonlinear, destroy
1326 
1327 .seealso: SNESCreate(), SNESSolve()
1328 @*/
1329 int SNESDestroy(SNES snes)
1330 {
1331   int ierr;
1332 
1333   PetscFunctionBegin;
1334   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1335   if (--snes->refct > 0) PetscFunctionReturn(0);
1336 
1337   if (snes->destroy) {ierr = (*(snes)->destroy)(snes); CHKERRQ(ierr);}
1338   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
1339   if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);}
1340   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
1341   if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);}
1342   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1343   PLogObjectDestroy((PetscObject)snes);
1344   PetscHeaderDestroy((PetscObject)snes);
1345   PetscFunctionReturn(0);
1346 }
1347 
1348 /* ----------- Routines to set solver parameters ---------- */
1349 
1350 #undef __FUNC__
1351 #define __FUNC__ "SNESSetTolerances"
1352 /*@
1353    SNESSetTolerances - Sets various parameters used in convergence tests.
1354 
1355    Collective on SNES
1356 
1357    Input Parameters:
1358 +  snes - the SNES context
1359 .  atol - absolute convergence tolerance
1360 .  rtol - relative convergence tolerance
1361 .  stol -  convergence tolerance in terms of the norm
1362            of the change in the solution between steps
1363 .  maxit - maximum number of iterations
1364 -  maxf - maximum number of function evaluations
1365 
1366    Options Database Keys:
1367 +    -snes_atol <atol> - Sets atol
1368 .    -snes_rtol <rtol> - Sets rtol
1369 .    -snes_stol <stol> - Sets stol
1370 .    -snes_max_it <maxit> - Sets maxit
1371 -    -snes_max_funcs <maxf> - Sets maxf
1372 
1373    Notes:
1374    The default maximum number of iterations is 50.
1375    The default maximum number of function evaluations is 1000.
1376 
1377 .keywords: SNES, nonlinear, set, convergence, tolerances
1378 
1379 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
1380 @*/
1381 int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
1382 {
1383   PetscFunctionBegin;
1384   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1385   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1386   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1387   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1388   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1389   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1390   PetscFunctionReturn(0);
1391 }
1392 
1393 #undef __FUNC__
1394 #define __FUNC__ "SNESGetTolerances"
1395 /*@
1396    SNESGetTolerances - Gets various parameters used in convergence tests.
1397 
1398    Not Collective
1399 
1400    Input Parameters:
1401 +  snes - the SNES context
1402 .  atol - absolute convergence tolerance
1403 .  rtol - relative convergence tolerance
1404 .  stol -  convergence tolerance in terms of the norm
1405            of the change in the solution between steps
1406 .  maxit - maximum number of iterations
1407 -  maxf - maximum number of function evaluations
1408 
1409    Notes:
1410    The user can specify PETSC_NULL for any parameter that is not needed.
1411 
1412 .keywords: SNES, nonlinear, get, convergence, tolerances
1413 
1414 .seealso: SNESSetTolerances()
1415 @*/
1416 int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
1417 {
1418   PetscFunctionBegin;
1419   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1420   if (atol)  *atol  = snes->atol;
1421   if (rtol)  *rtol  = snes->rtol;
1422   if (stol)  *stol  = snes->xtol;
1423   if (maxit) *maxit = snes->max_its;
1424   if (maxf)  *maxf  = snes->max_funcs;
1425   PetscFunctionReturn(0);
1426 }
1427 
1428 #undef __FUNC__
1429 #define __FUNC__ "SNESSetTrustRegionTolerance"
1430 /*@
1431    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1432 
1433    Collective on SNES
1434 
1435    Input Parameters:
1436 +  snes - the SNES context
1437 -  tol - tolerance
1438 
1439    Options Database Key:
1440 .  -snes_trtol <tol> - Sets tol
1441 
1442 .keywords: SNES, nonlinear, set, trust region, tolerance
1443 
1444 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
1445 @*/
1446 int SNESSetTrustRegionTolerance(SNES snes,double tol)
1447 {
1448   PetscFunctionBegin;
1449   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1450   snes->deltatol = tol;
1451   PetscFunctionReturn(0);
1452 }
1453 
1454 #undef __FUNC__
1455 #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
1456 /*@
1457    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
1458    for unconstrained minimization solvers.
1459 
1460    Collective on SNES
1461 
1462    Input Parameters:
1463 +  snes - the SNES context
1464 -  ftol - minimum function tolerance
1465 
1466    Options Database Key:
1467 .  -snes_fmin <ftol> - Sets ftol
1468 
1469    Note:
1470    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
1471    methods only.
1472 
1473 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
1474 
1475 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
1476 @*/
1477 int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
1478 {
1479   PetscFunctionBegin;
1480   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1481   snes->fmin = ftol;
1482   PetscFunctionReturn(0);
1483 }
1484 
1485 /* ------------ Routines to set performance monitoring options ----------- */
1486 
1487 #undef __FUNC__
1488 #define __FUNC__ "SNESSetMonitor"
1489 /*@C
1490    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
1491    iteration of the nonlinear solver to display the iteration's
1492    progress.
1493 
1494    Collective on SNES
1495 
1496    Input Parameters:
1497 +  snes - the SNES context
1498 .  func - monitoring routine
1499 -  mctx - [optional] user-defined context for private data for the
1500           monitor routine (may be PETSC_NULL)
1501 
1502    Calling sequence of func:
1503 $     int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
1504 
1505 +    snes - the SNES context
1506 .    its - iteration number
1507 .    mctx - [optional] monitoring context
1508 .    norm - 2-norm function value (may be estimated)
1509             for SNES_NONLINEAR_EQUATIONS methods
1510 -    norm - 2-norm gradient value (may be estimated)
1511             for SNES_UNCONSTRAINED_MINIMIZATION methods
1512 
1513    Options Database Keys:
1514 +    -snes_monitor        - sets SNESDefaultMonitor()
1515 .    -snes_xmonitor       - sets line graph monitor,
1516                             uses SNESLGMonitorCreate()
1517 _    -snes_cancelmonitors - cancels all monitors that have
1518                             been hardwired into a code by
1519                             calls to SNESSetMonitor(), but
1520                             does not cancel those set via
1521                             the options database.
1522 
1523    Notes:
1524    Several different monitoring routines may be set by calling
1525    SNESSetMonitor() multiple times; all will be called in the
1526    order in which they were set.
1527 
1528 .keywords: SNES, nonlinear, set, monitor
1529 
1530 .seealso: SNESDefaultMonitor(), SNESClearMonitor()
1531 @*/
1532 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
1533 {
1534   PetscFunctionBegin;
1535   if (snes->numbermonitors >= MAXSNESMONITORS) {
1536     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
1537   }
1538 
1539   snes->monitor[snes->numbermonitors]           = func;
1540   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1541   PetscFunctionReturn(0);
1542 }
1543 
1544 #undef __FUNC__
1545 #define __FUNC__ "SNESClearMonitor"
1546 /*@C
1547    SNESClearMonitor - Clears all the monitor functions for a SNES object.
1548 
1549    Collective on SNES
1550 
1551    Input Parameters:
1552 .  snes - the SNES context
1553 
1554    Options Database:
1555 .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1556     into a code by calls to SNESSetMonitor(), but does not cancel those
1557     set via the options database
1558 
1559    Notes:
1560    There is no way to clear one specific monitor from a SNES object.
1561 
1562 .keywords: SNES, nonlinear, set, monitor
1563 
1564 .seealso: SNESDefaultMonitor(), SNESSetMonitor()
1565 @*/
1566 int SNESClearMonitor( SNES snes )
1567 {
1568   PetscFunctionBegin;
1569   snes->numbermonitors = 0;
1570   PetscFunctionReturn(0);
1571 }
1572 
1573 #undef __FUNC__
1574 #define __FUNC__ "SNESSetConvergenceTest"
1575 /*@C
1576    SNESSetConvergenceTest - Sets the function that is to be used
1577    to test for convergence of the nonlinear iterative solution.
1578 
1579    Collective on SNES
1580 
1581    Input Parameters:
1582 +  snes - the SNES context
1583 .  func - routine to test for convergence
1584 -  cctx - [optional] context for private data for the convergence routine
1585           (may be PETSC_NULL)
1586 
1587    Calling sequence of func:
1588 $     int func (SNES snes,double xnorm,double gnorm,double f,void *cctx)
1589 
1590 +    snes - the SNES context
1591 .    cctx - [optional] convergence context
1592 .    xnorm - 2-norm of current iterate
1593 .    gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods)
1594 .    f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods)
1595 .    gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods)
1596 -    f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods)
1597 
1598 .keywords: SNES, nonlinear, set, convergence, test
1599 
1600 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1601           SNESConverged_UM_LS(), SNESConverged_UM_TR()
1602 @*/
1603 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
1604 {
1605   PetscFunctionBegin;
1606   (snes)->converged = func;
1607   (snes)->cnvP      = cctx;
1608   PetscFunctionReturn(0);
1609 }
1610 
1611 #undef __FUNC__
1612 #define __FUNC__ "SNESSetConvergenceHistory"
1613 /*@
1614    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1615 
1616    Collective on SNES
1617 
1618    Input Parameters:
1619 +  snes - iterative context obtained from SNESCreate()
1620 .  a   - array to hold history
1621 -  na  - size of a
1622 
1623    Notes:
1624    If set, this array will contain the function norms (for
1625    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1626    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1627    at each step.
1628 
1629    This routine is useful, e.g., when running a code for purposes
1630    of accurate performance monitoring, when no I/O should be done
1631    during the section of code that is being timed.
1632 
1633 .keywords: SNES, set, convergence, history
1634 @*/
1635 int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1636 {
1637   PetscFunctionBegin;
1638   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1639   if (na) PetscValidScalarPointer(a);
1640   snes->conv_hist      = a;
1641   snes->conv_hist_size = na;
1642   PetscFunctionReturn(0);
1643 }
1644 
1645 #undef __FUNC__
1646 #define __FUNC__ "SNESScaleStep_Private"
1647 /*
1648    SNESScaleStep_Private - Scales a step so that its length is less than the
1649    positive parameter delta.
1650 
1651     Input Parameters:
1652 +   snes - the SNES context
1653 .   y - approximate solution of linear system
1654 .   fnorm - 2-norm of current function
1655 -   delta - trust region size
1656 
1657     Output Parameters:
1658 +   gpnorm - predicted function norm at the new point, assuming local
1659     linearization.  The value is zero if the step lies within the trust
1660     region, and exceeds zero otherwise.
1661 -   ynorm - 2-norm of the step
1662 
1663     Note:
1664     For non-trust region methods such as SNES_EQ_LS, the parameter delta
1665     is set to be the maximum allowable step size.
1666 
1667 .keywords: SNES, nonlinear, scale, step
1668 */
1669 int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
1670                           double *gpnorm,double *ynorm)
1671 {
1672   double norm;
1673   Scalar cnorm;
1674   int    ierr;
1675 
1676   PetscFunctionBegin;
1677   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
1678   if (norm > *delta) {
1679      norm = *delta/norm;
1680      *gpnorm = (1.0 - norm)*(*fnorm);
1681      cnorm = norm;
1682      VecScale( &cnorm, y );
1683      *ynorm = *delta;
1684   } else {
1685      *gpnorm = 0.0;
1686      *ynorm = norm;
1687   }
1688   PetscFunctionReturn(0);
1689 }
1690 
1691 #undef __FUNC__
1692 #define __FUNC__ "SNESSolve"
1693 /*@
1694    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
1695    SNESCreate() and optional routines of the form SNESSetXXX().
1696 
1697    Collective on SNES
1698 
1699    Input Parameters:
1700 +  snes - the SNES context
1701 -  x - the solution vector
1702 
1703    Output Parameter:
1704 .  its - number of iterations until termination
1705 
1706    Notes:
1707    The user should initialize the vector, x, with the initial guess
1708    for the nonlinear solve prior to calling SNESSolve.  In particular,
1709    to employ an initial guess of zero, the user should explicitly set
1710    this vector to zero by calling VecSet().
1711 
1712 .keywords: SNES, nonlinear, solve
1713 
1714 .seealso: SNESCreate(), SNESDestroy()
1715 @*/
1716 int SNESSolve(SNES snes,Vec x,int *its)
1717 {
1718   int ierr, flg;
1719 
1720   PetscFunctionBegin;
1721   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1722   PetscValidIntPointer(its);
1723   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1724   else {snes->vec_sol = snes->vec_sol_always = x;}
1725   PLogEventBegin(SNES_Solve,snes,0,0,0);
1726   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
1727   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
1728   PLogEventEnd(SNES_Solve,snes,0,0,0);
1729   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
1730   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
1731   PetscFunctionReturn(0);
1732 }
1733 
1734 /* --------- Internal routines for SNES Package --------- */
1735 
1736 #undef __FUNC__
1737 #define __FUNC__ "SNESSetType"
1738 /*@C
1739    SNESSetType - Sets the method for the nonlinear solver.
1740 
1741    Collective on SNES
1742 
1743    Input Parameters:
1744 +  snes - the SNES context
1745 -  method - a known method
1746 
1747    Options Database Key:
1748 .  -snes_type <method> - Sets the method; use -help for a list
1749    of available methods (for instance, ls or tr)
1750 
1751    Notes:
1752    See "petsc/include/snes.h" for available methods (for instance)
1753 .    SNES_EQ_LS - Newton's method with line search
1754      (systems of nonlinear equations)
1755 .    SNES_EQ_TR - Newton's method with trust region
1756      (systems of nonlinear equations)
1757 .    SNES_UM_TR - Newton's method with trust region
1758      (unconstrained minimization)
1759 .    SNES_UM_LS - Newton's method with line search
1760      (unconstrained minimization)
1761 
1762   Normally, it is best to use the SNESSetFromOptions() command and then
1763   set the SNES solver type from the options database rather than by using
1764   this routine.  Using the options database provides the user with
1765   maximum flexibility in evaluating the many nonlinear solvers.
1766   The SNESSetType() routine is provided for those situations where it
1767   is necessary to set the nonlinear solver independently of the command
1768   line or options database.  This might be the case, for example, when
1769   the choice of solver changes during the execution of the program,
1770   and the user's application is taking responsibility for choosing the
1771   appropriate method.  In other words, this routine is for the advanced user.
1772 
1773 .keywords: SNES, set, method
1774 @*/
1775 int SNESSetType(SNES snes,SNESType method)
1776 {
1777   int ierr;
1778   int (*r)(SNES);
1779 
1780   PetscFunctionBegin;
1781   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1782 
1783   if (!PetscStrcmp(snes->type_name,method)) PetscFunctionReturn(0);
1784 
1785   if (snes->setupcalled) {
1786     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
1787     snes->data = 0;
1788   }
1789 
1790   /* Get the function pointers for the iterative method requested */
1791   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);}
1792 
1793   ierr =  FListFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr);
1794 
1795   if (!r) SETERRQ(1,1,"Unable to find requested SNES type");
1796 
1797   if (snes->data) PetscFree(snes->data);
1798   snes->data = 0;
1799   ierr = (*r)(snes); CHKERRQ(ierr);
1800 
1801   if (snes->type_name) PetscFree(snes->type_name);
1802   snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name);
1803   PetscStrcpy(snes->type_name,method);
1804   snes->set_method_called = 1;
1805 
1806   PetscFunctionReturn(0);
1807 }
1808 
1809 
1810 /* --------------------------------------------------------------------- */
1811 #undef __FUNC__
1812 #define __FUNC__ "SNESRegisterDestroy"
1813 /*@C
1814    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1815    registered by SNESRegister().
1816 
1817    Not Collective
1818 
1819 .keywords: SNES, nonlinear, register, destroy
1820 
1821 .seealso: SNESRegisterAll(), SNESRegisterAll()
1822 @*/
1823 int SNESRegisterDestroy(void)
1824 {
1825   int ierr;
1826 
1827   PetscFunctionBegin;
1828   if (SNESList) {
1829     ierr = FListDestroy( SNESList );CHKERRQ(ierr);
1830     SNESList = 0;
1831   }
1832   SNESRegisterAllCalled = 0;
1833   PetscFunctionReturn(0);
1834 }
1835 
1836 #undef __FUNC__
1837 #define __FUNC__ "SNESGetType"
1838 /*@C
1839    SNESGetType - Gets the SNES method type and name (as a string).
1840 
1841    Not Collective
1842 
1843    Input Parameter:
1844 .  snes - nonlinear solver context
1845 
1846    Output Parameter:
1847 .  method - SNES method (a charactor string)
1848 
1849 .keywords: SNES, nonlinear, get, method, name
1850 @*/
1851 int SNESGetType(SNES snes, SNESType *method)
1852 {
1853   PetscFunctionBegin;
1854   *method = snes->type_name;
1855   PetscFunctionReturn(0);
1856 }
1857 
1858 #undef __FUNC__
1859 #define __FUNC__ "SNESGetSolution"
1860 /*@C
1861    SNESGetSolution - Returns the vector where the approximate solution is
1862    stored.
1863 
1864    Not Collective, but Vec is parallel if SNES is parallel
1865 
1866    Input Parameter:
1867 .  snes - the SNES context
1868 
1869    Output Parameter:
1870 .  x - the solution
1871 
1872 .keywords: SNES, nonlinear, get, solution
1873 
1874 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
1875 @*/
1876 int SNESGetSolution(SNES snes,Vec *x)
1877 {
1878   PetscFunctionBegin;
1879   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1880   *x = snes->vec_sol_always;
1881   PetscFunctionReturn(0);
1882 }
1883 
1884 #undef __FUNC__
1885 #define __FUNC__ "SNESGetSolutionUpdate"
1886 /*@C
1887    SNESGetSolutionUpdate - Returns the vector where the solution update is
1888    stored.
1889 
1890    Not Collective, but Vec is parallel if SNES is parallel
1891 
1892    Input Parameter:
1893 .  snes - the SNES context
1894 
1895    Output Parameter:
1896 .  x - the solution update
1897 
1898 .keywords: SNES, nonlinear, get, solution, update
1899 
1900 .seealso: SNESGetSolution(), SNESGetFunction
1901 @*/
1902 int SNESGetSolutionUpdate(SNES snes,Vec *x)
1903 {
1904   PetscFunctionBegin;
1905   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1906   *x = snes->vec_sol_update_always;
1907   PetscFunctionReturn(0);
1908 }
1909 
1910 #undef __FUNC__
1911 #define __FUNC__ "SNESGetFunction"
1912 /*@C
1913    SNESGetFunction - Returns the vector where the function is stored.
1914 
1915    Not Collective, but Vec is parallel if SNES is parallel
1916 
1917    Input Parameter:
1918 .  snes - the SNES context
1919 
1920    Output Parameter:
1921 .  r - the function
1922 
1923    Notes:
1924    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
1925    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
1926    SNESGetMinimizationFunction() and SNESGetGradient();
1927 
1928 .keywords: SNES, nonlinear, get, function
1929 
1930 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
1931           SNESGetGradient()
1932 @*/
1933 int SNESGetFunction(SNES snes,Vec *r)
1934 {
1935   PetscFunctionBegin;
1936   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1937   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1938     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1939   }
1940   *r = snes->vec_func_always;
1941   PetscFunctionReturn(0);
1942 }
1943 
1944 #undef __FUNC__
1945 #define __FUNC__ "SNESGetGradient"
1946 /*@C
1947    SNESGetGradient - Returns the vector where the gradient is stored.
1948 
1949    Not Collective, but Vec is parallel if SNES is parallel
1950 
1951    Input Parameter:
1952 .  snes - the SNES context
1953 
1954    Output Parameter:
1955 .  r - the gradient
1956 
1957    Notes:
1958    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
1959    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
1960    SNESGetFunction().
1961 
1962 .keywords: SNES, nonlinear, get, gradient
1963 
1964 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
1965 @*/
1966 int SNESGetGradient(SNES snes,Vec *r)
1967 {
1968   PetscFunctionBegin;
1969   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1970   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1971     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1972   }
1973   *r = snes->vec_func_always;
1974   PetscFunctionReturn(0);
1975 }
1976 
1977 #undef __FUNC__
1978 #define __FUNC__ "SNESGetMinimizationFunction"
1979 /*@
1980    SNESGetMinimizationFunction - Returns the scalar function value for
1981    unconstrained minimization problems.
1982 
1983    Not Collective
1984 
1985    Input Parameter:
1986 .  snes - the SNES context
1987 
1988    Output Parameter:
1989 .  r - the function
1990 
1991    Notes:
1992    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
1993    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
1994    SNESGetFunction().
1995 
1996 .keywords: SNES, nonlinear, get, function
1997 
1998 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
1999 @*/
2000 int SNESGetMinimizationFunction(SNES snes,double *r)
2001 {
2002   PetscFunctionBegin;
2003   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2004   PetscValidScalarPointer(r);
2005   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
2006     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
2007   }
2008   *r = snes->fc;
2009   PetscFunctionReturn(0);
2010 }
2011 
2012 #undef __FUNC__
2013 #define __FUNC__ "SNESSetOptionsPrefix"
2014 /*@C
2015    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2016    SNES options in the database.
2017 
2018    Collective on SNES
2019 
2020    Input Parameter:
2021 +  snes - the SNES context
2022 -  prefix - the prefix to prepend to all option names
2023 
2024    Notes:
2025    A hyphen (-) must NOT be given at the beginning of the prefix name.
2026    The first character of all runtime options is AUTOMATICALLY the hyphen.
2027 
2028 .keywords: SNES, set, options, prefix, database
2029 
2030 .seealso: SNESSetFromOptions()
2031 @*/
2032 int SNESSetOptionsPrefix(SNES snes,char *prefix)
2033 {
2034   int ierr;
2035 
2036   PetscFunctionBegin;
2037   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2038   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
2039   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2040   PetscFunctionReturn(0);
2041 }
2042 
2043 #undef __FUNC__
2044 #define __FUNC__ "SNESAppendOptionsPrefix"
2045 /*@C
2046    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2047    SNES options in the database.
2048 
2049    Collective on SNES
2050 
2051    Input Parameters:
2052 +  snes - the SNES context
2053 -  prefix - the prefix to prepend to all option names
2054 
2055    Notes:
2056    A hyphen (-) must NOT be given at the beginning of the prefix name.
2057    The first character of all runtime options is AUTOMATICALLY the hyphen.
2058 
2059 .keywords: SNES, append, options, prefix, database
2060 
2061 .seealso: SNESGetOptionsPrefix()
2062 @*/
2063 int SNESAppendOptionsPrefix(SNES snes,char *prefix)
2064 {
2065   int ierr;
2066 
2067   PetscFunctionBegin;
2068   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2069   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
2070   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2071   PetscFunctionReturn(0);
2072 }
2073 
2074 #undef __FUNC__
2075 #define __FUNC__ "SNESGetOptionsPrefix"
2076 /*@
2077    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2078    SNES options in the database.
2079 
2080    Not Collective
2081 
2082    Input Parameter:
2083 .  snes - the SNES context
2084 
2085    Output Parameter:
2086 .  prefix - pointer to the prefix string used
2087 
2088 .keywords: SNES, get, options, prefix, database
2089 
2090 .seealso: SNESAppendOptionsPrefix()
2091 @*/
2092 int SNESGetOptionsPrefix(SNES snes,char **prefix)
2093 {
2094   int ierr;
2095 
2096   PetscFunctionBegin;
2097   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2098   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
2099   PetscFunctionReturn(0);
2100 }
2101 
2102 #undef __FUNC__
2103 #define __FUNC__ "SNESPrintHelp"
2104 /*@
2105    SNESPrintHelp - Prints all options for the SNES component.
2106 
2107    Collective on SNES
2108 
2109    Input Parameter:
2110 .  snes - the SNES context
2111 
2112    Options Database Keys:
2113 +  -help - Prints SNES options
2114 -  -h - Prints SNES options
2115 
2116 .keywords: SNES, nonlinear, help
2117 
2118 .seealso: SNESSetFromOptions()
2119 @*/
2120 int SNESPrintHelp(SNES snes)
2121 {
2122   char                p[64];
2123   SNES_KSP_EW_ConvCtx *kctx;
2124   int                 ierr;
2125 
2126   PetscFunctionBegin;
2127   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2128 
2129   PetscStrcpy(p,"-");
2130   if (snes->prefix) PetscStrcat(p, snes->prefix);
2131 
2132   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2133 
2134   (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");
2135   ierr = FListPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr);
2136   (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
2137   (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
2138   (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
2139   (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
2140   (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
2141   (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
2142   (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
2143   (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");
2144   (*PetscHelpPrintf)(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
2145   (*PetscHelpPrintf)(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
2146     residual norm at each iteration.\n",p);
2147   (*PetscHelpPrintf)(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
2148     residual norm for small residual norms. This is useful to conceal\n\
2149     meaningless digits that may be different on different machines.\n",p);
2150   (*PetscHelpPrintf)(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
2151   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
2152     (*PetscHelpPrintf)(snes->comm,
2153      " Options for solving systems of nonlinear equations only:\n");
2154     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
2155     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
2156     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
2157     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_ksp_monitor - if using matrix-free multiply then prints h at each KSP iteration\n",p);
2158     (*PetscHelpPrintf)(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
2159     (*PetscHelpPrintf)(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
2160     (*PetscHelpPrintf)(snes->comm,
2161      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
2162     (*PetscHelpPrintf)(snes->comm,
2163      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
2164     (*PetscHelpPrintf)(snes->comm,
2165      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
2166     (*PetscHelpPrintf)(snes->comm,
2167      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
2168     (*PetscHelpPrintf)(snes->comm,
2169      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
2170     (*PetscHelpPrintf)(snes->comm,
2171      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
2172     (*PetscHelpPrintf)(snes->comm,
2173      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
2174   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
2175     (*PetscHelpPrintf)(snes->comm,
2176      " Options for solving unconstrained minimization problems only:\n");
2177     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
2178     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
2179     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
2180   }
2181   (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
2182   (*PetscHelpPrintf)(snes->comm,"a particular method\n");
2183   if (snes->printhelp) {
2184     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2185   }
2186   PetscFunctionReturn(0);
2187 }
2188 
2189 /*MC
2190    SNESRegister - Adds a method to the nonlinear solver package.
2191 
2192    Synopsis:
2193    SNESRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
2194 
2195    Not collective
2196 
2197    Input Parameters:
2198 +  name_solver - name of a new user-defined solver
2199 .  path - path (either absolute or relative) the library containing this solver
2200 .  name_create - name of routine to create method context
2201 -  routine_create - routine to create method context
2202 
2203    Notes:
2204    SNESRegister() may be called multiple times to add several user-defined solvers.
2205 
2206    If dynamic libraries are used, then the fourth input argument (routine_create)
2207    is ignored.
2208 
2209    Sample usage:
2210 .vb
2211    SNESRegister("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2212                 "MySolverCreate",MySolverCreate);
2213 .ve
2214 
2215    Then, your solver can be chosen with the procedural interface via
2216 $     SNESSetType(snes,"my_solver")
2217    or at runtime via the option
2218 $     -snes_type my_solver
2219 
2220 .keywords: SNES, nonlinear, register
2221 
2222 .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2223 M*/
2224 
2225 #undef __FUNC__
2226 #define __FUNC__ "SNESRegister_Private"
2227 int SNESRegister_Private(char *sname,char *path,char *name,int (*function)(SNES))
2228 {
2229   char fullname[256];
2230   int  ierr;
2231 
2232   PetscFunctionBegin;
2233   PetscStrcpy(fullname,path); PetscStrcat(fullname,":");PetscStrcat(fullname,name);
2234   ierr = FListAdd_Private(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr);
2235   PetscFunctionReturn(0);
2236 }
2237