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