xref: /petsc/src/snes/interface/snes.c (revision ecf371e4f11dc92d760a8a34d4f20e65e993c904)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*ecf371e4SBarry Smith static char vcid[] = "$Id: snes.c,v 1.139 1998/03/12 23:22:43 bsmith Exp bsmith $";
39b94acceSBarry Smith #endif
49b94acceSBarry Smith 
570f55243SBarry Smith #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
6f5eb4b81SSatish Balay #include "src/sys/nreg.h"
79b94acceSBarry Smith #include "pinclude/pviewer.h"
89b94acceSBarry Smith #include <math.h>
99b94acceSBarry Smith 
1084cb2905SBarry Smith int SNESRegisterAllCalled = 0;
1182bf6240SBarry Smith DLList SNESList = 0;
1282bf6240SBarry Smith 
1384cb2905SBarry Smith 
145615d1e5SSatish Balay #undef __FUNC__
15d4bb536fSBarry Smith #define __FUNC__ "SNESView"
169b94acceSBarry Smith /*@
179b94acceSBarry Smith    SNESView - Prints the SNES data structure.
189b94acceSBarry Smith 
199b94acceSBarry Smith    Input Parameters:
209b94acceSBarry Smith .  SNES - the SNES context
219b94acceSBarry Smith .  viewer - visualization context
229b94acceSBarry Smith 
239b94acceSBarry Smith    Options Database Key:
249b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
259b94acceSBarry Smith 
269b94acceSBarry Smith    Notes:
279b94acceSBarry Smith    The available visualization contexts include
286d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
296d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
309b94acceSBarry Smith $       output where only the first processor opens
319b94acceSBarry Smith $       the file.  All other processors send their
329b94acceSBarry Smith $       data to the first processor to print.
339b94acceSBarry Smith 
349b94acceSBarry Smith    The user can open alternative vistualization contexts with
35dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
369b94acceSBarry Smith 
379b94acceSBarry Smith .keywords: SNES, view
389b94acceSBarry Smith 
39dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
409b94acceSBarry Smith @*/
419b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
429b94acceSBarry Smith {
439b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
449b94acceSBarry Smith   FILE                *fd;
459b94acceSBarry Smith   int                 ierr;
469b94acceSBarry Smith   SLES                sles;
479b94acceSBarry Smith   char                *method;
4819bcc07fSBarry Smith   ViewerType          vtype;
499b94acceSBarry Smith 
503a40ed3dSBarry Smith   PetscFunctionBegin;
5174679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5274679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5374679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5474679c65SBarry Smith 
5519bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5619bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5790ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5877c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
5982bf6240SBarry Smith     SNESGetType(snes,&method);
6077c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
619b94acceSBarry Smith     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
6277c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
639b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
649b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6577c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
669b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
679b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
687a00f4a9SLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
697a00f4a9SLois Curfman McInnes     "  total number of linear solver iterations=%d\n",snes->linear_its);
70ae3c334cSLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
7168632166SLois Curfman McInnes      "  total number of function evaluations=%d\n",snes->nfuncs);
729b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
7377c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
749b94acceSBarry Smith     if (snes->ksp_ewconv) {
759b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
769b94acceSBarry Smith       if (kctx) {
7777c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
789b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
799b94acceSBarry Smith         kctx->version);
8077c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
819b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
829b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
8377c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
849b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
859b94acceSBarry Smith       }
869b94acceSBarry Smith     }
87c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
8882bf6240SBarry Smith     SNESGetType(snes,&method);
89c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
9019bcc07fSBarry Smith   }
919b94acceSBarry Smith   SNESGetSLES(snes,&sles);
929b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
933a40ed3dSBarry Smith   PetscFunctionReturn(0);
949b94acceSBarry Smith }
959b94acceSBarry Smith 
96639f9d9dSBarry Smith /*
97639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
98639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
99639f9d9dSBarry Smith */
100639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
101639f9d9dSBarry Smith static int numberofsetfromoptions;
102639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
103639f9d9dSBarry Smith 
1045615d1e5SSatish Balay #undef __FUNC__
105d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker"
106639f9d9dSBarry Smith /*@
107639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
108639f9d9dSBarry Smith 
109639f9d9dSBarry Smith     Input Parameter:
110639f9d9dSBarry Smith .   snescheck - function that checks for options
111639f9d9dSBarry Smith 
112639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
113639f9d9dSBarry Smith @*/
114639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
115639f9d9dSBarry Smith {
1163a40ed3dSBarry Smith   PetscFunctionBegin;
117639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
118a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed");
119639f9d9dSBarry Smith   }
120639f9d9dSBarry Smith 
121639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
1223a40ed3dSBarry Smith   PetscFunctionReturn(0);
123639f9d9dSBarry Smith }
124639f9d9dSBarry Smith 
1255615d1e5SSatish Balay #undef __FUNC__
1265615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1279b94acceSBarry Smith /*@
128682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1299b94acceSBarry Smith 
1309b94acceSBarry Smith    Input Parameter:
1319b94acceSBarry Smith .  snes - the SNES context
1329b94acceSBarry Smith 
13311ca99fdSLois Curfman McInnes    Notes:
13411ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
13511ca99fdSLois Curfman McInnes    the users manual.
13683e2fdc7SBarry Smith 
1379b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1389b94acceSBarry Smith 
139a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1409b94acceSBarry Smith @*/
1419b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1429b94acceSBarry Smith {
14382bf6240SBarry Smith   char     method[256];
1449b94acceSBarry Smith   double   tmp;
1459b94acceSBarry Smith   SLES     sles;
14681f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1473c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1489b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1499b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1509b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1519b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1529b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1539b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1549b94acceSBarry Smith 
1553a40ed3dSBarry Smith   PetscFunctionBegin;
15681f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
15781f57ec7SBarry Smith 
15877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15982bf6240SBarry Smith   if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp");
160ca161407SBarry Smith 
16182bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
16282bf6240SBarry Smith   ierr = OptionsGetString(snes->prefix,"-snes_type",method,256,&flg);
163052efed2SBarry Smith   if (flg) {
16482bf6240SBarry Smith     ierr = SNESSetType(snes,(SNESType) method); CHKERRQ(ierr);
1655334005bSBarry Smith   }
1663c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
167d64ed03dSBarry Smith   if (flg) {
168d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
169d64ed03dSBarry Smith   }
1703c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
171d64ed03dSBarry Smith   if (flg) {
172d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
173d64ed03dSBarry Smith   }
1743c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
175d64ed03dSBarry Smith   if (flg) {
176d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
177d64ed03dSBarry Smith   }
1783c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1793c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
180d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
181d64ed03dSBarry Smith   if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); }
182d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
183d64ed03dSBarry Smith   if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);  CHKERRQ(ierr);}
1843c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1853c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
186b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
187b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
188b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
189b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
190b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
191b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
192b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
19393c39befSBarry Smith 
19493c39befSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr);
19593c39befSBarry Smith   if (flg) {snes->converged = 0;}
19693c39befSBarry Smith 
1979b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
1989b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
1995bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
2005bbad29bSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);}
2013c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
202639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
2033c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
204d31a3109SLois Curfman McInnes   if (flg) {
205639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
206d31a3109SLois Curfman McInnes   }
20781f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
2083c7409f5SSatish Balay   if (flg) {
20917699dbbSLois Curfman McInnes     int    rank = 0;
210d7e8b826SBarry Smith     DrawLG lg;
21117699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
21217699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
21317699dbbSLois Curfman McInnes     if (!rank) {
21481f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2159b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
216f8c826e1SBarry Smith       snes->xmonitor = lg;
2179b94acceSBarry Smith       PLogObjectParent(snes,lg);
2189b94acceSBarry Smith     }
2199b94acceSBarry Smith   }
2203c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2213c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2229b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2239b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
224981c4779SBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
225981c4779SBarry Smith   } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
22631615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
22731615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
228d64ed03dSBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2299b94acceSBarry Smith   }
230639f9d9dSBarry Smith 
231639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
232639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
233639f9d9dSBarry Smith   }
234639f9d9dSBarry Smith 
2359b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2369b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
23782bf6240SBarry Smith   /*
23882bf6240SBarry Smith         Since the private setfromoptions requires the type to all ready have
23982bf6240SBarry Smith       been set we make sure a type is set by this time
24082bf6240SBarry Smith   */
24182bf6240SBarry Smith   if (!snes->type_name) {
24282bf6240SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
24382bf6240SBarry Smith       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
24482bf6240SBarry Smith     } else {
24582bf6240SBarry Smith       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
24682bf6240SBarry Smith     }
24782bf6240SBarry Smith   }
24882bf6240SBarry Smith 
24982bf6240SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
25082bf6240SBarry Smith   if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);}
25182bf6240SBarry Smith 
2523a40ed3dSBarry Smith   if (snes->setfromoptions) {
2533a40ed3dSBarry Smith     ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
2543a40ed3dSBarry Smith   }
2553a40ed3dSBarry Smith   PetscFunctionReturn(0);
2569b94acceSBarry Smith }
2579b94acceSBarry Smith 
258a847f771SSatish Balay 
2595615d1e5SSatish Balay #undef __FUNC__
260d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext"
2619b94acceSBarry Smith /*@
2629b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2639b94acceSBarry Smith    the nonlinear solvers.
2649b94acceSBarry Smith 
2659b94acceSBarry Smith    Input Parameters:
2669b94acceSBarry Smith .  snes - the SNES context
2679b94acceSBarry Smith .  usrP - optional user context
2689b94acceSBarry Smith 
2699b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2709b94acceSBarry Smith 
2719b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2729b94acceSBarry Smith @*/
2739b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
2749b94acceSBarry Smith {
2753a40ed3dSBarry Smith   PetscFunctionBegin;
27677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2779b94acceSBarry Smith   snes->user		= usrP;
2783a40ed3dSBarry Smith   PetscFunctionReturn(0);
2799b94acceSBarry Smith }
28074679c65SBarry Smith 
2815615d1e5SSatish Balay #undef __FUNC__
282d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext"
2839b94acceSBarry Smith /*@C
2849b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
2859b94acceSBarry Smith    nonlinear solvers.
2869b94acceSBarry Smith 
2879b94acceSBarry Smith    Input Parameter:
2889b94acceSBarry Smith .  snes - SNES context
2899b94acceSBarry Smith 
2909b94acceSBarry Smith    Output Parameter:
2919b94acceSBarry Smith .  usrP - user context
2929b94acceSBarry Smith 
2939b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
2949b94acceSBarry Smith 
2959b94acceSBarry Smith .seealso: SNESSetApplicationContext()
2969b94acceSBarry Smith @*/
2979b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
2989b94acceSBarry Smith {
2993a40ed3dSBarry Smith   PetscFunctionBegin;
30077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3019b94acceSBarry Smith   *usrP = snes->user;
3023a40ed3dSBarry Smith   PetscFunctionReturn(0);
3039b94acceSBarry Smith }
30474679c65SBarry Smith 
3055615d1e5SSatish Balay #undef __FUNC__
306d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber"
3079b94acceSBarry Smith /*@
3089b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3099b94acceSBarry Smith    nonlinear solver.
3109b94acceSBarry Smith 
3119b94acceSBarry Smith    Input Parameter:
3129b94acceSBarry Smith .  snes - SNES context
3139b94acceSBarry Smith 
3149b94acceSBarry Smith    Output Parameter:
3159b94acceSBarry Smith .  iter - iteration number
3169b94acceSBarry Smith 
3179b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3189b94acceSBarry Smith @*/
3199b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3209b94acceSBarry Smith {
3213a40ed3dSBarry Smith   PetscFunctionBegin;
32277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
32374679c65SBarry Smith   PetscValidIntPointer(iter);
3249b94acceSBarry Smith   *iter = snes->iter;
3253a40ed3dSBarry Smith   PetscFunctionReturn(0);
3269b94acceSBarry Smith }
32774679c65SBarry Smith 
3285615d1e5SSatish Balay #undef __FUNC__
3295615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
3309b94acceSBarry Smith /*@
3319b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3329b94acceSBarry Smith    with SNESSSetFunction().
3339b94acceSBarry Smith 
3349b94acceSBarry Smith    Input Parameter:
3359b94acceSBarry Smith .  snes - SNES context
3369b94acceSBarry Smith 
3379b94acceSBarry Smith    Output Parameter:
3389b94acceSBarry Smith .  fnorm - 2-norm of function
3399b94acceSBarry Smith 
3409b94acceSBarry Smith    Note:
3419b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
342a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
343a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3449b94acceSBarry Smith 
3459b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
346a86d99e1SLois Curfman McInnes 
347a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3489b94acceSBarry Smith @*/
3499b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3509b94acceSBarry Smith {
3513a40ed3dSBarry Smith   PetscFunctionBegin;
35277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
35374679c65SBarry Smith   PetscValidScalarPointer(fnorm);
35474679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
355d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only");
35674679c65SBarry Smith   }
3579b94acceSBarry Smith   *fnorm = snes->norm;
3583a40ed3dSBarry Smith   PetscFunctionReturn(0);
3599b94acceSBarry Smith }
36074679c65SBarry Smith 
3615615d1e5SSatish Balay #undef __FUNC__
3625615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
3639b94acceSBarry Smith /*@
3649b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
3659b94acceSBarry Smith    with SNESSSetGradient().
3669b94acceSBarry Smith 
3679b94acceSBarry Smith    Input Parameter:
3689b94acceSBarry Smith .  snes - SNES context
3699b94acceSBarry Smith 
3709b94acceSBarry Smith    Output Parameter:
3719b94acceSBarry Smith .  fnorm - 2-norm of gradient
3729b94acceSBarry Smith 
3739b94acceSBarry Smith    Note:
3749b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
375a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
376a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
3779b94acceSBarry Smith 
3789b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
379a86d99e1SLois Curfman McInnes 
380a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
3819b94acceSBarry Smith @*/
3829b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
3839b94acceSBarry Smith {
3843a40ed3dSBarry Smith   PetscFunctionBegin;
38577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
38674679c65SBarry Smith   PetscValidScalarPointer(gnorm);
38774679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
388d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
38974679c65SBarry Smith   }
3909b94acceSBarry Smith   *gnorm = snes->norm;
3913a40ed3dSBarry Smith   PetscFunctionReturn(0);
3929b94acceSBarry Smith }
39374679c65SBarry Smith 
3945615d1e5SSatish Balay #undef __FUNC__
395d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps"
3969b94acceSBarry Smith /*@
3979b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
3989b94acceSBarry Smith    attempted by the nonlinear solver.
3999b94acceSBarry Smith 
4009b94acceSBarry Smith    Input Parameter:
4019b94acceSBarry Smith .  snes - SNES context
4029b94acceSBarry Smith 
4039b94acceSBarry Smith    Output Parameter:
4049b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4059b94acceSBarry Smith 
406c96a6f78SLois Curfman McInnes    Notes:
407c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
408c96a6f78SLois Curfman McInnes 
4099b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4109b94acceSBarry Smith @*/
4119b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4129b94acceSBarry Smith {
4133a40ed3dSBarry Smith   PetscFunctionBegin;
41477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
41574679c65SBarry Smith   PetscValidIntPointer(nfails);
4169b94acceSBarry Smith   *nfails = snes->nfailures;
4173a40ed3dSBarry Smith   PetscFunctionReturn(0);
4189b94acceSBarry Smith }
419a847f771SSatish Balay 
4205615d1e5SSatish Balay #undef __FUNC__
421d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations"
422c96a6f78SLois Curfman McInnes /*@
423c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
424c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
425c96a6f78SLois Curfman McInnes 
426c96a6f78SLois Curfman McInnes    Input Parameter:
427c96a6f78SLois Curfman McInnes .  snes - SNES context
428c96a6f78SLois Curfman McInnes 
429c96a6f78SLois Curfman McInnes    Output Parameter:
430c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
431c96a6f78SLois Curfman McInnes 
432c96a6f78SLois Curfman McInnes    Notes:
433c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
434c96a6f78SLois Curfman McInnes 
435c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
436c96a6f78SLois Curfman McInnes @*/
43786bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
438c96a6f78SLois Curfman McInnes {
4393a40ed3dSBarry Smith   PetscFunctionBegin;
440c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
441c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
442c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
4433a40ed3dSBarry Smith   PetscFunctionReturn(0);
444c96a6f78SLois Curfman McInnes }
445c96a6f78SLois Curfman McInnes 
446c96a6f78SLois Curfman McInnes #undef __FUNC__
447d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES"
4489b94acceSBarry Smith /*@C
4499b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4509b94acceSBarry Smith 
4519b94acceSBarry Smith    Input Parameter:
4529b94acceSBarry Smith .  snes - the SNES context
4539b94acceSBarry Smith 
4549b94acceSBarry Smith    Output Parameter:
4559b94acceSBarry Smith .  sles - the SLES context
4569b94acceSBarry Smith 
4579b94acceSBarry Smith    Notes:
4589b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4599b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4609b94acceSBarry Smith    KSP and PC contexts as well.
4619b94acceSBarry Smith 
4629b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4639b94acceSBarry Smith 
4649b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4659b94acceSBarry Smith @*/
4669b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4679b94acceSBarry Smith {
4683a40ed3dSBarry Smith   PetscFunctionBegin;
46977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4709b94acceSBarry Smith   *sles = snes->sles;
4713a40ed3dSBarry Smith   PetscFunctionReturn(0);
4729b94acceSBarry Smith }
47382bf6240SBarry Smith 
4749b94acceSBarry Smith /* -----------------------------------------------------------*/
4755615d1e5SSatish Balay #undef __FUNC__
4765615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
4779b94acceSBarry Smith /*@C
4789b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
4799b94acceSBarry Smith 
4809b94acceSBarry Smith    Input Parameter:
4819b94acceSBarry Smith .  comm - MPI communicator
4829b94acceSBarry Smith .  type - type of method, one of
4839b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
4849b94acceSBarry Smith $      (for systems of nonlinear equations)
4859b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
4869b94acceSBarry Smith $      (for unconstrained minimization)
4879b94acceSBarry Smith 
4889b94acceSBarry Smith    Output Parameter:
4899b94acceSBarry Smith .  outsnes - the new SNES context
4909b94acceSBarry Smith 
491c1f60f51SBarry Smith    Options Database Key:
49219bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
49319bd6747SLois Curfman McInnes $              and no preconditioning matrix
49419bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
49519bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
49619bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
49719bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
498c1f60f51SBarry Smith 
4999b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5009b94acceSBarry Smith 
50163a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5029b94acceSBarry Smith @*/
5034b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5049b94acceSBarry Smith {
5059b94acceSBarry Smith   int                 ierr;
5069b94acceSBarry Smith   SNES                snes;
5079b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
50837fcc0dbSBarry Smith 
5093a40ed3dSBarry Smith   PetscFunctionBegin;
5109b94acceSBarry Smith   *outsnes = 0;
511d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
512d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
513d64ed03dSBarry Smith   }
514f830108cSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,comm,SNESDestroy,SNESView);
5159b94acceSBarry Smith   PLogObjectCreate(snes);
5169b94acceSBarry Smith   snes->max_its           = 50;
5179750a799SBarry Smith   snes->max_funcs	  = 10000;
5189b94acceSBarry Smith   snes->norm		  = 0.0;
5195a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
520b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
521b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5229b94acceSBarry Smith     snes->atol		  = 1.e-10;
523d64ed03dSBarry Smith   } else {
524b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
525b4874afaSBarry Smith     snes->ttol            = 0.0;
526b4874afaSBarry Smith     snes->atol		  = 1.e-50;
527b4874afaSBarry Smith   }
5289b94acceSBarry Smith   snes->xtol		  = 1.e-8;
529b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5309b94acceSBarry Smith   snes->nfuncs            = 0;
5319b94acceSBarry Smith   snes->nfailures         = 0;
5327a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
533639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5349b94acceSBarry Smith   snes->data              = 0;
5359b94acceSBarry Smith   snes->view              = 0;
5369b94acceSBarry Smith   snes->computeumfunction = 0;
5379b94acceSBarry Smith   snes->umfunP            = 0;
5389b94acceSBarry Smith   snes->fc                = 0;
5399b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5409b94acceSBarry Smith   snes->fmin              = -1.e30;
5419b94acceSBarry Smith   snes->method_class      = type;
5429b94acceSBarry Smith   snes->set_method_called = 0;
54382bf6240SBarry Smith   snes->setupcalled      = 0;
5449b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5456f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5466f24a144SLois Curfman McInnes   snes->vwork             = 0;
5476f24a144SLois Curfman McInnes   snes->nwork             = 0;
548c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
549c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
550c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
5519b94acceSBarry Smith 
5529b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5530452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
554eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
5559b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5569b94acceSBarry Smith   kctx->version     = 2;
5579b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5589b94acceSBarry Smith                              this was too large for some test cases */
5599b94acceSBarry Smith   kctx->rtol_last   = 0;
5609b94acceSBarry Smith   kctx->rtol_max    = .9;
5619b94acceSBarry Smith   kctx->gamma       = 1.0;
5629b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5639b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5649b94acceSBarry Smith   kctx->threshold   = .1;
5659b94acceSBarry Smith   kctx->lresid_last = 0;
5669b94acceSBarry Smith   kctx->norm_last   = 0;
5679b94acceSBarry Smith 
5689b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5699b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5705334005bSBarry Smith 
5719b94acceSBarry Smith   *outsnes = snes;
5723a40ed3dSBarry Smith   PetscFunctionReturn(0);
5739b94acceSBarry Smith }
5749b94acceSBarry Smith 
5759b94acceSBarry Smith /* --------------------------------------------------------------- */
5765615d1e5SSatish Balay #undef __FUNC__
577d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction"
5789b94acceSBarry Smith /*@C
5799b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
5809b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
5819b94acceSBarry Smith    equations.
5829b94acceSBarry Smith 
5839b94acceSBarry Smith    Input Parameters:
5849b94acceSBarry Smith .  snes - the SNES context
5859b94acceSBarry Smith .  func - function evaluation routine
5869b94acceSBarry Smith .  r - vector to store function value
5872cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
5882cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
5899b94acceSBarry Smith 
5909b94acceSBarry Smith    Calling sequence of func:
591313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
5929b94acceSBarry Smith 
5939b94acceSBarry Smith .  x - input vector
594313e4042SLois Curfman McInnes .  f - function vector
5952cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
5969b94acceSBarry Smith 
5979b94acceSBarry Smith    Notes:
5989b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
5999b94acceSBarry Smith $      f'(x) x = -f(x),
6009b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
6019b94acceSBarry Smith 
6029b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6039b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6049b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6059b94acceSBarry Smith 
6069b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6079b94acceSBarry Smith 
608a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6099b94acceSBarry Smith @*/
6105334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6119b94acceSBarry Smith {
6123a40ed3dSBarry Smith   PetscFunctionBegin;
61377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
614a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
615a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
616a8c6a408SBarry Smith   }
6179b94acceSBarry Smith   snes->computefunction     = func;
6189b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6199b94acceSBarry Smith   snes->funP                = ctx;
6203a40ed3dSBarry Smith   PetscFunctionReturn(0);
6219b94acceSBarry Smith }
6229b94acceSBarry Smith 
6235615d1e5SSatish Balay #undef __FUNC__
6245615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
6259b94acceSBarry Smith /*@
6269b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6279b94acceSBarry Smith    SNESSetFunction().
6289b94acceSBarry Smith 
6299b94acceSBarry Smith    Input Parameters:
6309b94acceSBarry Smith .  snes - the SNES context
6319b94acceSBarry Smith .  x - input vector
6329b94acceSBarry Smith 
6339b94acceSBarry Smith    Output Parameter:
6343638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6359b94acceSBarry Smith 
6361bffabb2SLois Curfman McInnes    Notes:
6379b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6389b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6399b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6409b94acceSBarry Smith 
6419b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6429b94acceSBarry Smith 
643a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6449b94acceSBarry Smith @*/
6459b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6469b94acceSBarry Smith {
6479b94acceSBarry Smith   int    ierr;
6489b94acceSBarry Smith 
6493a40ed3dSBarry Smith   PetscFunctionBegin;
650d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
651a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
652d4bb536fSBarry Smith   }
6539b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
654d64ed03dSBarry Smith   PetscStackPush("SNES user function");
6559b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
656d64ed03dSBarry Smith   PetscStackPop;
657ae3c334cSLois Curfman McInnes   snes->nfuncs++;
6589b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6593a40ed3dSBarry Smith   PetscFunctionReturn(0);
6609b94acceSBarry Smith }
6619b94acceSBarry Smith 
6625615d1e5SSatish Balay #undef __FUNC__
663d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction"
6649b94acceSBarry Smith /*@C
6659b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6669b94acceSBarry Smith    unconstrained minimization.
6679b94acceSBarry Smith 
6689b94acceSBarry Smith    Input Parameters:
6699b94acceSBarry Smith .  snes - the SNES context
6709b94acceSBarry Smith .  func - function evaluation routine
671044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
672044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6739b94acceSBarry Smith 
6749b94acceSBarry Smith    Calling sequence of func:
6759b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
6769b94acceSBarry Smith 
6779b94acceSBarry Smith .  x - input vector
6789b94acceSBarry Smith .  f - function
679044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
6809b94acceSBarry Smith 
6819b94acceSBarry Smith    Notes:
6829b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
6839b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
6849b94acceSBarry Smith    SNESSetFunction().
6859b94acceSBarry Smith 
6869b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
6879b94acceSBarry Smith 
688a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
689a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
6909b94acceSBarry Smith @*/
6919b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
6929b94acceSBarry Smith                       void *ctx)
6939b94acceSBarry Smith {
6943a40ed3dSBarry Smith   PetscFunctionBegin;
69577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
696a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
697a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
698a8c6a408SBarry Smith   }
6999b94acceSBarry Smith   snes->computeumfunction   = func;
7009b94acceSBarry Smith   snes->umfunP              = ctx;
7013a40ed3dSBarry Smith   PetscFunctionReturn(0);
7029b94acceSBarry Smith }
7039b94acceSBarry Smith 
7045615d1e5SSatish Balay #undef __FUNC__
7055615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7069b94acceSBarry Smith /*@
7079b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7089b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7099b94acceSBarry Smith 
7109b94acceSBarry Smith    Input Parameters:
7119b94acceSBarry Smith .  snes - the SNES context
7129b94acceSBarry Smith .  x - input vector
7139b94acceSBarry Smith 
7149b94acceSBarry Smith    Output Parameter:
7159b94acceSBarry Smith .  y - function value
7169b94acceSBarry Smith 
7179b94acceSBarry Smith    Notes:
7189b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7199b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7209b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
721a86d99e1SLois Curfman McInnes 
722a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
723a86d99e1SLois Curfman McInnes 
724a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
725a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7269b94acceSBarry Smith @*/
7279b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7289b94acceSBarry Smith {
7299b94acceSBarry Smith   int    ierr;
7303a40ed3dSBarry Smith 
7313a40ed3dSBarry Smith   PetscFunctionBegin;
732a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
733a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
734a8c6a408SBarry Smith   }
7359b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
736d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
7379b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
738d64ed03dSBarry Smith   PetscStackPop;
739ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7409b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7413a40ed3dSBarry Smith   PetscFunctionReturn(0);
7429b94acceSBarry Smith }
7439b94acceSBarry Smith 
7445615d1e5SSatish Balay #undef __FUNC__
745d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient"
7469b94acceSBarry Smith /*@C
7479b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7489b94acceSBarry Smith    vector for use by the SNES routines.
7499b94acceSBarry Smith 
7509b94acceSBarry Smith    Input Parameters:
7519b94acceSBarry Smith .  snes - the SNES context
7529b94acceSBarry Smith .  func - function evaluation routine
753044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
754044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7559b94acceSBarry Smith .  r - vector to store gradient value
7569b94acceSBarry Smith 
7579b94acceSBarry Smith    Calling sequence of func:
7589b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7599b94acceSBarry Smith 
7609b94acceSBarry Smith .  x - input vector
7619b94acceSBarry Smith .  g - gradient vector
762044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7639b94acceSBarry Smith 
7649b94acceSBarry Smith    Notes:
7659b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7669b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7679b94acceSBarry Smith    SNESSetFunction().
7689b94acceSBarry Smith 
7699b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7709b94acceSBarry Smith 
771a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
772a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
7739b94acceSBarry Smith @*/
77474679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
7759b94acceSBarry Smith {
7763a40ed3dSBarry Smith   PetscFunctionBegin;
77777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
778a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
779a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
780a8c6a408SBarry Smith   }
7819b94acceSBarry Smith   snes->computefunction     = func;
7829b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7839b94acceSBarry Smith   snes->funP                = ctx;
7843a40ed3dSBarry Smith   PetscFunctionReturn(0);
7859b94acceSBarry Smith }
7869b94acceSBarry Smith 
7875615d1e5SSatish Balay #undef __FUNC__
7885615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
7899b94acceSBarry Smith /*@
790a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
791a86d99e1SLois Curfman McInnes    SNESSetGradient().
7929b94acceSBarry Smith 
7939b94acceSBarry Smith    Input Parameters:
7949b94acceSBarry Smith .  snes - the SNES context
7959b94acceSBarry Smith .  x - input vector
7969b94acceSBarry Smith 
7979b94acceSBarry Smith    Output Parameter:
7989b94acceSBarry Smith .  y - gradient vector
7999b94acceSBarry Smith 
8009b94acceSBarry Smith    Notes:
8019b94acceSBarry Smith    SNESComputeGradient() is valid only for
8029b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8039b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
804a86d99e1SLois Curfman McInnes 
805a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
806a86d99e1SLois Curfman McInnes 
807a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
808a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8099b94acceSBarry Smith @*/
8109b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8119b94acceSBarry Smith {
8129b94acceSBarry Smith   int    ierr;
8133a40ed3dSBarry Smith 
8143a40ed3dSBarry Smith   PetscFunctionBegin;
8153a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
816a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8173a40ed3dSBarry Smith   }
8189b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
819d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
8209b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
821d64ed03dSBarry Smith   PetscStackPop;
8229b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8233a40ed3dSBarry Smith   PetscFunctionReturn(0);
8249b94acceSBarry Smith }
8259b94acceSBarry Smith 
8265615d1e5SSatish Balay #undef __FUNC__
8275615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
82862fef451SLois Curfman McInnes /*@
82962fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
83062fef451SLois Curfman McInnes    set with SNESSetJacobian().
83162fef451SLois Curfman McInnes 
83262fef451SLois Curfman McInnes    Input Parameters:
83362fef451SLois Curfman McInnes .  snes - the SNES context
83462fef451SLois Curfman McInnes .  x - input vector
83562fef451SLois Curfman McInnes 
83662fef451SLois Curfman McInnes    Output Parameters:
83762fef451SLois Curfman McInnes .  A - Jacobian matrix
83862fef451SLois Curfman McInnes .  B - optional preconditioning matrix
83962fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
84062fef451SLois Curfman McInnes 
84162fef451SLois Curfman McInnes    Notes:
84262fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
84362fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
84462fef451SLois Curfman McInnes 
845dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
846dc5a77f8SLois Curfman McInnes    flag parameter.
84762fef451SLois Curfman McInnes 
84862fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
84962fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
850005c665bSBarry Smith    methods is SNESComputeHessian().
85162fef451SLois Curfman McInnes 
85262fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
85362fef451SLois Curfman McInnes 
85462fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
85562fef451SLois Curfman McInnes @*/
8569b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8579b94acceSBarry Smith {
8589b94acceSBarry Smith   int    ierr;
8593a40ed3dSBarry Smith 
8603a40ed3dSBarry Smith   PetscFunctionBegin;
8613a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
862a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
8633a40ed3dSBarry Smith   }
8643a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
8659b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
866c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
867d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
8689b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
869d64ed03dSBarry Smith   PetscStackPop;
8709b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8716d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
87277c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
87377c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8743a40ed3dSBarry Smith   PetscFunctionReturn(0);
8759b94acceSBarry Smith }
8769b94acceSBarry Smith 
8775615d1e5SSatish Balay #undef __FUNC__
8785615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
87962fef451SLois Curfman McInnes /*@
88062fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
88162fef451SLois Curfman McInnes    set with SNESSetHessian().
88262fef451SLois Curfman McInnes 
88362fef451SLois Curfman McInnes    Input Parameters:
88462fef451SLois Curfman McInnes .  snes - the SNES context
88562fef451SLois Curfman McInnes .  x - input vector
88662fef451SLois Curfman McInnes 
88762fef451SLois Curfman McInnes    Output Parameters:
88862fef451SLois Curfman McInnes .  A - Hessian matrix
88962fef451SLois Curfman McInnes .  B - optional preconditioning matrix
89062fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
89162fef451SLois Curfman McInnes 
89262fef451SLois Curfman McInnes    Notes:
89362fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
89462fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
89562fef451SLois Curfman McInnes 
896dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
897dc5a77f8SLois Curfman McInnes    flag parameter.
89862fef451SLois Curfman McInnes 
89962fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
90062fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
90162fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
90262fef451SLois Curfman McInnes 
90362fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
90462fef451SLois Curfman McInnes 
905a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
906a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
90762fef451SLois Curfman McInnes @*/
90862fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
9099b94acceSBarry Smith {
9109b94acceSBarry Smith   int    ierr;
9113a40ed3dSBarry Smith 
9123a40ed3dSBarry Smith   PetscFunctionBegin;
9133a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
914a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9153a40ed3dSBarry Smith   }
9163a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
91762fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9180f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
919d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
92062fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
921d64ed03dSBarry Smith   PetscStackPop;
92262fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9230f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
92477c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
92577c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9263a40ed3dSBarry Smith   PetscFunctionReturn(0);
9279b94acceSBarry Smith }
9289b94acceSBarry Smith 
9295615d1e5SSatish Balay #undef __FUNC__
930d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian"
9319b94acceSBarry Smith /*@C
9329b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
933044dda88SLois Curfman McInnes    location to store the matrix.
9349b94acceSBarry Smith 
9359b94acceSBarry Smith    Input Parameters:
9369b94acceSBarry Smith .  snes - the SNES context
9379b94acceSBarry Smith .  A - Jacobian matrix
9389b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9399b94acceSBarry Smith .  func - Jacobian evaluation routine
9402cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9412cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9429b94acceSBarry Smith 
9439b94acceSBarry Smith    Calling sequence of func:
944313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9459b94acceSBarry Smith 
9469b94acceSBarry Smith .  x - input vector
9479b94acceSBarry Smith .  A - Jacobian matrix
9489b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
949ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
950ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9512cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
9529b94acceSBarry Smith 
9539b94acceSBarry Smith    Notes:
954dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9552cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
956ac21db08SLois Curfman McInnes 
957ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9589b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9599b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9609b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9619b94acceSBarry Smith    throughout the global iterations.
9629b94acceSBarry Smith 
9639b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9649b94acceSBarry Smith 
965ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
9669b94acceSBarry Smith @*/
9679b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9689b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9699b94acceSBarry Smith {
9703a40ed3dSBarry Smith   PetscFunctionBegin;
97177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
972a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
973a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
974a8c6a408SBarry Smith   }
9759b94acceSBarry Smith   snes->computejacobian = func;
9769b94acceSBarry Smith   snes->jacP            = ctx;
9779b94acceSBarry Smith   snes->jacobian        = A;
9789b94acceSBarry Smith   snes->jacobian_pre    = B;
9793a40ed3dSBarry Smith   PetscFunctionReturn(0);
9809b94acceSBarry Smith }
98162fef451SLois Curfman McInnes 
9825615d1e5SSatish Balay #undef __FUNC__
983d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian"
984b4fd4287SBarry Smith /*@
985b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
986b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
987b4fd4287SBarry Smith 
988b4fd4287SBarry Smith    Input Parameter:
989b4fd4287SBarry Smith .  snes - the nonlinear solver context
990b4fd4287SBarry Smith 
991b4fd4287SBarry Smith    Output Parameters:
992b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
993b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
994b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
995b4fd4287SBarry Smith 
996b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
997b4fd4287SBarry Smith @*/
998b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
999b4fd4287SBarry Smith {
10003a40ed3dSBarry Smith   PetscFunctionBegin;
10013a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1002a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
10033a40ed3dSBarry Smith   }
1004b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
1005b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
1006b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
10073a40ed3dSBarry Smith   PetscFunctionReturn(0);
1008b4fd4287SBarry Smith }
1009b4fd4287SBarry Smith 
10105615d1e5SSatish Balay #undef __FUNC__
1011d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian"
10129b94acceSBarry Smith /*@C
10139b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1014044dda88SLois Curfman McInnes    location to store the matrix.
10159b94acceSBarry Smith 
10169b94acceSBarry Smith    Input Parameters:
10179b94acceSBarry Smith .  snes - the SNES context
10189b94acceSBarry Smith .  A - Hessian matrix
10199b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
10209b94acceSBarry Smith .  func - Jacobian evaluation routine
1021313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1022313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10239b94acceSBarry Smith 
10249b94acceSBarry Smith    Calling sequence of func:
1025313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10269b94acceSBarry Smith 
10279b94acceSBarry Smith .  x - input vector
10289b94acceSBarry Smith .  A - Hessian matrix
10299b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1030ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1031ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10322cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10339b94acceSBarry Smith 
10349b94acceSBarry Smith    Notes:
1035dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10362cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1037ac21db08SLois Curfman McInnes 
10389b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10399b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10409b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10419b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10429b94acceSBarry Smith    throughout the global iterations.
10439b94acceSBarry Smith 
10449b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10459b94acceSBarry Smith 
1046ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10479b94acceSBarry Smith @*/
10489b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10499b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10509b94acceSBarry Smith {
10513a40ed3dSBarry Smith   PetscFunctionBegin;
105277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1053d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1054a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1055d4bb536fSBarry Smith   }
10569b94acceSBarry Smith   snes->computejacobian = func;
10579b94acceSBarry Smith   snes->jacP            = ctx;
10589b94acceSBarry Smith   snes->jacobian        = A;
10599b94acceSBarry Smith   snes->jacobian_pre    = B;
10603a40ed3dSBarry Smith   PetscFunctionReturn(0);
10619b94acceSBarry Smith }
10629b94acceSBarry Smith 
10635615d1e5SSatish Balay #undef __FUNC__
1064d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian"
106562fef451SLois Curfman McInnes /*@
106662fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
106762fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
106862fef451SLois Curfman McInnes 
106962fef451SLois Curfman McInnes    Input Parameter:
107062fef451SLois Curfman McInnes .  snes - the nonlinear solver context
107162fef451SLois Curfman McInnes 
107262fef451SLois Curfman McInnes    Output Parameters:
107362fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
107462fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
107562fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
107662fef451SLois Curfman McInnes 
107762fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
107862fef451SLois Curfman McInnes @*/
107962fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
108062fef451SLois Curfman McInnes {
10813a40ed3dSBarry Smith   PetscFunctionBegin;
10823a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1083a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10843a40ed3dSBarry Smith   }
108562fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
108662fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
108762fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
10883a40ed3dSBarry Smith   PetscFunctionReturn(0);
108962fef451SLois Curfman McInnes }
109062fef451SLois Curfman McInnes 
10919b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10929b94acceSBarry Smith 
10935615d1e5SSatish Balay #undef __FUNC__
10945615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
10959b94acceSBarry Smith /*@
10969b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1097272ac6f2SLois Curfman McInnes    of a nonlinear solver.
10989b94acceSBarry Smith 
10999b94acceSBarry Smith    Input Parameter:
11009b94acceSBarry Smith .  snes - the SNES context
11018ddd3da0SLois Curfman McInnes .  x - the solution vector
11029b94acceSBarry Smith 
1103272ac6f2SLois Curfman McInnes    Notes:
1104272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1105272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1106272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1107272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1108272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1109272ac6f2SLois Curfman McInnes 
11109b94acceSBarry Smith .keywords: SNES, nonlinear, setup
11119b94acceSBarry Smith 
11129b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
11139b94acceSBarry Smith @*/
11148ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
11159b94acceSBarry Smith {
1116272ac6f2SLois Curfman McInnes   int ierr, flg;
11173a40ed3dSBarry Smith 
11183a40ed3dSBarry Smith   PetscFunctionBegin;
111977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
112077c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
11218ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
11229b94acceSBarry Smith 
1123c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1124c1f60f51SBarry Smith   /*
1125c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1126dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1127c1f60f51SBarry Smith   */
1128c1f60f51SBarry Smith   if (flg) {
1129c1f60f51SBarry Smith     Mat J;
1130c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1131c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1132c1f60f51SBarry Smith     snes->mfshell = J;
1133c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1134c1f60f51SBarry Smith       snes->jacobian = J;
1135c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1136d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1137c1f60f51SBarry Smith       snes->jacobian = J;
1138c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1139d4bb536fSBarry Smith     } else {
1140a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option");
1141d4bb536fSBarry Smith     }
1142c1f60f51SBarry Smith   }
1143272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1144c1f60f51SBarry Smith   /*
1145dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1146c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1147c1f60f51SBarry Smith    */
114831615d04SLois Curfman McInnes   if (flg) {
1149272ac6f2SLois Curfman McInnes     Mat J;
1150272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1151272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1152272ac6f2SLois Curfman McInnes     snes->mfshell = J;
115383e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
115483e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
115594a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1156d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
115783e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
115894a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1159d4bb536fSBarry Smith     } else {
1160a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option");
1161d4bb536fSBarry Smith     }
1162272ac6f2SLois Curfman McInnes   }
11639b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1164a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1165a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1166a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first");
1167a8c6a408SBarry Smith     if (snes->vec_func == snes->vec_sol) {
1168a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector");
1169a8c6a408SBarry Smith     }
1170a703fe33SLois Curfman McInnes 
1171a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
117282bf6240SBarry Smith     if (snes->ksp_ewconv && PetscStrcmp(snes->type_name,SNES_EQ_TR)) {
1173a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1174a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1175a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1176a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1177a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1178a703fe33SLois Curfman McInnes     }
11799b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1180a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1181a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1182a8c6a408SBarry Smith     if (!snes->computeumfunction) {
1183a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first");
1184a8c6a408SBarry Smith     }
1185a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first");
1186d4bb536fSBarry Smith   } else {
1187a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
1188d4bb536fSBarry Smith   }
1189a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
119082bf6240SBarry Smith   snes->setupcalled = 1;
11913a40ed3dSBarry Smith   PetscFunctionReturn(0);
11929b94acceSBarry Smith }
11939b94acceSBarry Smith 
11945615d1e5SSatish Balay #undef __FUNC__
1195d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy"
11969b94acceSBarry Smith /*@C
11979b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11989b94acceSBarry Smith    with SNESCreate().
11999b94acceSBarry Smith 
12009b94acceSBarry Smith    Input Parameter:
12019b94acceSBarry Smith .  snes - the SNES context
12029b94acceSBarry Smith 
12039b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
12049b94acceSBarry Smith 
120563a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
12069b94acceSBarry Smith @*/
12079b94acceSBarry Smith int SNESDestroy(SNES snes)
12089b94acceSBarry Smith {
12099b94acceSBarry Smith   int ierr;
12103a40ed3dSBarry Smith 
12113a40ed3dSBarry Smith   PetscFunctionBegin;
121277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12133a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1214d4bb536fSBarry Smith 
12159750a799SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);}
12160452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
1217522c5e43SBarry Smith   if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);}
12189b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
1219522c5e43SBarry Smith   if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);}
1220522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
12219b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
12220452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
12233a40ed3dSBarry Smith   PetscFunctionReturn(0);
12249b94acceSBarry Smith }
12259b94acceSBarry Smith 
12269b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12279b94acceSBarry Smith 
12285615d1e5SSatish Balay #undef __FUNC__
12295615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
12309b94acceSBarry Smith /*@
1231d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12329b94acceSBarry Smith 
12339b94acceSBarry Smith    Input Parameters:
12349b94acceSBarry Smith .  snes - the SNES context
123533174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
123633174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
123733174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
123833174efeSLois Curfman McInnes            of the change in the solution between steps
123933174efeSLois Curfman McInnes .  maxit - maximum number of iterations
124033174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
12419b94acceSBarry Smith 
124233174efeSLois Curfman McInnes    Options Database Keys:
124333174efeSLois Curfman McInnes $    -snes_atol <atol>
124433174efeSLois Curfman McInnes $    -snes_rtol <rtol>
124533174efeSLois Curfman McInnes $    -snes_stol <stol>
124633174efeSLois Curfman McInnes $    -snes_max_it <maxit>
124733174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
12489b94acceSBarry Smith 
1249d7a720efSLois Curfman McInnes    Notes:
12509b94acceSBarry Smith    The default maximum number of iterations is 50.
12519b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
12529b94acceSBarry Smith 
125333174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
12549b94acceSBarry Smith 
1255d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
12569b94acceSBarry Smith @*/
1257d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
12589b94acceSBarry Smith {
12593a40ed3dSBarry Smith   PetscFunctionBegin;
126077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1261d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1262d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1263d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1264d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1265d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
12663a40ed3dSBarry Smith   PetscFunctionReturn(0);
12679b94acceSBarry Smith }
12689b94acceSBarry Smith 
12695615d1e5SSatish Balay #undef __FUNC__
12705615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
12719b94acceSBarry Smith /*@
127233174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
127333174efeSLois Curfman McInnes 
127433174efeSLois Curfman McInnes    Input Parameters:
127533174efeSLois Curfman McInnes .  snes - the SNES context
127633174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
127733174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
127833174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
127933174efeSLois Curfman McInnes            of the change in the solution between steps
128033174efeSLois Curfman McInnes .  maxit - maximum number of iterations
128133174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
128233174efeSLois Curfman McInnes 
128333174efeSLois Curfman McInnes    Notes:
128433174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
128533174efeSLois Curfman McInnes 
128633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
128733174efeSLois Curfman McInnes 
128833174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
128933174efeSLois Curfman McInnes @*/
129033174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
129133174efeSLois Curfman McInnes {
12923a40ed3dSBarry Smith   PetscFunctionBegin;
129333174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
129433174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
129533174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
129633174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
129733174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
129833174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
12993a40ed3dSBarry Smith   PetscFunctionReturn(0);
130033174efeSLois Curfman McInnes }
130133174efeSLois Curfman McInnes 
13025615d1e5SSatish Balay #undef __FUNC__
13035615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
130433174efeSLois Curfman McInnes /*@
13059b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
13069b94acceSBarry Smith 
13079b94acceSBarry Smith    Input Parameters:
13089b94acceSBarry Smith .  snes - the SNES context
13099b94acceSBarry Smith .  tol - tolerance
13109b94acceSBarry Smith 
13119b94acceSBarry Smith    Options Database Key:
1312d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
13139b94acceSBarry Smith 
13149b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
13159b94acceSBarry Smith 
1316d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
13179b94acceSBarry Smith @*/
13189b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
13199b94acceSBarry Smith {
13203a40ed3dSBarry Smith   PetscFunctionBegin;
132177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13229b94acceSBarry Smith   snes->deltatol = tol;
13233a40ed3dSBarry Smith   PetscFunctionReturn(0);
13249b94acceSBarry Smith }
13259b94acceSBarry Smith 
13265615d1e5SSatish Balay #undef __FUNC__
13275615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
13289b94acceSBarry Smith /*@
132977c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
13309b94acceSBarry Smith    for unconstrained minimization solvers.
13319b94acceSBarry Smith 
13329b94acceSBarry Smith    Input Parameters:
13339b94acceSBarry Smith .  snes - the SNES context
13349b94acceSBarry Smith .  ftol - minimum function tolerance
13359b94acceSBarry Smith 
13369b94acceSBarry Smith    Options Database Key:
1337d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
13389b94acceSBarry Smith 
13399b94acceSBarry Smith    Note:
134077c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
13419b94acceSBarry Smith    methods only.
13429b94acceSBarry Smith 
13439b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
13449b94acceSBarry Smith 
1345d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
13469b94acceSBarry Smith @*/
134777c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
13489b94acceSBarry Smith {
13493a40ed3dSBarry Smith   PetscFunctionBegin;
135077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13519b94acceSBarry Smith   snes->fmin = ftol;
13523a40ed3dSBarry Smith   PetscFunctionReturn(0);
13539b94acceSBarry Smith }
13549b94acceSBarry Smith 
13559b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
13569b94acceSBarry Smith 
13575615d1e5SSatish Balay #undef __FUNC__
1358d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor"
13599b94acceSBarry Smith /*@C
13609b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
13619b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
13629b94acceSBarry Smith    progress.
13639b94acceSBarry Smith 
13649b94acceSBarry Smith    Input Parameters:
13659b94acceSBarry Smith .  snes - the SNES context
13669b94acceSBarry Smith .  func - monitoring routine
1367044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1368044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
13699b94acceSBarry Smith 
13709b94acceSBarry Smith    Calling sequence of func:
1371682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
13729b94acceSBarry Smith 
13739b94acceSBarry Smith $    snes - the SNES context
13749b94acceSBarry Smith $    its - iteration number
1375044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
13769b94acceSBarry Smith $
13779b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13789b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
13799b94acceSBarry Smith $
13809b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13819b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
13829b94acceSBarry Smith 
13839665c990SLois Curfman McInnes    Options Database Keys:
13849665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
13859665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
13869665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
13879665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
13889665c990SLois Curfman McInnes $                           been hardwired into a code by
13899665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
13909665c990SLois Curfman McInnes $                           does not cancel those set via
13919665c990SLois Curfman McInnes $                           the options database.
13929665c990SLois Curfman McInnes 
13939665c990SLois Curfman McInnes 
1394639f9d9dSBarry Smith    Notes:
13956bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13966bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13976bc08f3fSLois Curfman McInnes    order in which they were set.
1398639f9d9dSBarry Smith 
13999b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
14009b94acceSBarry Smith 
14019b94acceSBarry Smith .seealso: SNESDefaultMonitor()
14029b94acceSBarry Smith @*/
140374679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
14049b94acceSBarry Smith {
14053a40ed3dSBarry Smith   PetscFunctionBegin;
1406639f9d9dSBarry Smith   if (!func) {
1407639f9d9dSBarry Smith     snes->numbermonitors = 0;
14083a40ed3dSBarry Smith     PetscFunctionReturn(0);
1409639f9d9dSBarry Smith   }
1410639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1411a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
1412639f9d9dSBarry Smith   }
1413639f9d9dSBarry Smith 
1414639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1415639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
14163a40ed3dSBarry Smith   PetscFunctionReturn(0);
14179b94acceSBarry Smith }
14189b94acceSBarry Smith 
14195615d1e5SSatish Balay #undef __FUNC__
1420d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest"
14219b94acceSBarry Smith /*@C
14229b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
14239b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
14249b94acceSBarry Smith 
14259b94acceSBarry Smith    Input Parameters:
14269b94acceSBarry Smith .  snes - the SNES context
14279b94acceSBarry Smith .  func - routine to test for convergence
1428044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1429044dda88SLois Curfman McInnes           (may be PETSC_NULL)
14309b94acceSBarry Smith 
14319b94acceSBarry Smith    Calling sequence of func:
14329b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
14339b94acceSBarry Smith              double f,void *cctx)
14349b94acceSBarry Smith 
14359b94acceSBarry Smith $    snes - the SNES context
1436044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
14379b94acceSBarry Smith $    xnorm - 2-norm of current iterate
14389b94acceSBarry Smith $
14399b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14409b94acceSBarry Smith $    gnorm - 2-norm of current step
14419b94acceSBarry Smith $    f - 2-norm of function
14429b94acceSBarry Smith $
14439b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14449b94acceSBarry Smith $    gnorm - 2-norm of current gradient
14459b94acceSBarry Smith $    f - function value
14469b94acceSBarry Smith 
14479b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14489b94acceSBarry Smith 
144940191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
145040191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
14519b94acceSBarry Smith @*/
145274679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
14539b94acceSBarry Smith {
14543a40ed3dSBarry Smith   PetscFunctionBegin;
14559b94acceSBarry Smith   (snes)->converged = func;
14569b94acceSBarry Smith   (snes)->cnvP      = cctx;
14573a40ed3dSBarry Smith   PetscFunctionReturn(0);
14589b94acceSBarry Smith }
14599b94acceSBarry Smith 
14605615d1e5SSatish Balay #undef __FUNC__
1461d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory"
1462c9005455SLois Curfman McInnes /*@
1463c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1464c9005455SLois Curfman McInnes 
1465c9005455SLois Curfman McInnes    Input Parameters:
1466c9005455SLois Curfman McInnes .  snes - iterative context obtained from SNESCreate()
1467c9005455SLois Curfman McInnes .  a   - array to hold history
1468c9005455SLois Curfman McInnes .  na  - size of a
1469c9005455SLois Curfman McInnes 
1470c9005455SLois Curfman McInnes    Notes:
1471c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1472c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1473c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1474c9005455SLois Curfman McInnes    at each step.
1475c9005455SLois Curfman McInnes 
1476c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1477c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1478c9005455SLois Curfman McInnes    during the section of code that is being timed.
1479c9005455SLois Curfman McInnes 
1480c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1481c9005455SLois Curfman McInnes @*/
1482c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1483c9005455SLois Curfman McInnes {
14843a40ed3dSBarry Smith   PetscFunctionBegin;
1485c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1486c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1487c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1488c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
14893a40ed3dSBarry Smith   PetscFunctionReturn(0);
1490c9005455SLois Curfman McInnes }
1491c9005455SLois Curfman McInnes 
1492c9005455SLois Curfman McInnes #undef __FUNC__
14935615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
14949b94acceSBarry Smith /*
14959b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
14969b94acceSBarry Smith    positive parameter delta.
14979b94acceSBarry Smith 
14989b94acceSBarry Smith     Input Parameters:
14999b94acceSBarry Smith .   snes - the SNES context
15009b94acceSBarry Smith .   y - approximate solution of linear system
15019b94acceSBarry Smith .   fnorm - 2-norm of current function
15029b94acceSBarry Smith .   delta - trust region size
15039b94acceSBarry Smith 
15049b94acceSBarry Smith     Output Parameters:
15059b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
15069b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
15079b94acceSBarry Smith     region, and exceeds zero otherwise.
15089b94acceSBarry Smith .   ynorm - 2-norm of the step
15099b94acceSBarry Smith 
15109b94acceSBarry Smith     Note:
151140191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
15129b94acceSBarry Smith     is set to be the maximum allowable step size.
15139b94acceSBarry Smith 
15149b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
15159b94acceSBarry Smith */
15169b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
15179b94acceSBarry Smith                           double *gpnorm,double *ynorm)
15189b94acceSBarry Smith {
15199b94acceSBarry Smith   double norm;
15209b94acceSBarry Smith   Scalar cnorm;
15213a40ed3dSBarry Smith   int    ierr;
15223a40ed3dSBarry Smith 
15233a40ed3dSBarry Smith   PetscFunctionBegin;
15243a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
15259b94acceSBarry Smith   if (norm > *delta) {
15269b94acceSBarry Smith      norm = *delta/norm;
15279b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
15289b94acceSBarry Smith      cnorm = norm;
15299b94acceSBarry Smith      VecScale( &cnorm, y );
15309b94acceSBarry Smith      *ynorm = *delta;
15319b94acceSBarry Smith   } else {
15329b94acceSBarry Smith      *gpnorm = 0.0;
15339b94acceSBarry Smith      *ynorm = norm;
15349b94acceSBarry Smith   }
15353a40ed3dSBarry Smith   PetscFunctionReturn(0);
15369b94acceSBarry Smith }
15379b94acceSBarry Smith 
15385615d1e5SSatish Balay #undef __FUNC__
15395615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
15409b94acceSBarry Smith /*@
15419b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
154263a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
15439b94acceSBarry Smith 
15449b94acceSBarry Smith    Input Parameter:
15459b94acceSBarry Smith .  snes - the SNES context
15468ddd3da0SLois Curfman McInnes .  x - the solution vector
15479b94acceSBarry Smith 
15489b94acceSBarry Smith    Output Parameter:
15499b94acceSBarry Smith    its - number of iterations until termination
15509b94acceSBarry Smith 
15518ddd3da0SLois Curfman McInnes    Note:
15528ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
15538ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
15548ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
15558ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
15568ddd3da0SLois Curfman McInnes 
15579b94acceSBarry Smith .keywords: SNES, nonlinear, solve
15589b94acceSBarry Smith 
155963a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
15609b94acceSBarry Smith @*/
15618ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
15629b94acceSBarry Smith {
15633c7409f5SSatish Balay   int ierr, flg;
1564052efed2SBarry Smith 
15653a40ed3dSBarry Smith   PetscFunctionBegin;
156677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
156774679c65SBarry Smith   PetscValidIntPointer(its);
156882bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1569c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
15709b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1571c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
15729b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
15739b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
15743c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
15756d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
15763a40ed3dSBarry Smith   PetscFunctionReturn(0);
15779b94acceSBarry Smith }
15789b94acceSBarry Smith 
15799b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
15809b94acceSBarry Smith 
15815615d1e5SSatish Balay #undef __FUNC__
15825615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
158382bf6240SBarry Smith /*@C
15844b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
15859b94acceSBarry Smith 
15869b94acceSBarry Smith    Input Parameters:
15879b94acceSBarry Smith .  snes - the SNES context
15889b94acceSBarry Smith .  method - a known method
15899b94acceSBarry Smith 
1590ae12b187SLois Curfman McInnes   Options Database Command:
1591ae12b187SLois Curfman McInnes $ -snes_type  <method>
1592ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1593ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1594ae12b187SLois Curfman McInnes 
15959b94acceSBarry Smith    Notes:
15969b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
15979b94acceSBarry Smith $  Systems of nonlinear equations:
159840191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
159940191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
16009b94acceSBarry Smith $  Unconstrained minimization:
160140191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
160240191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
16039b94acceSBarry Smith 
1604ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1605ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1606ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1607ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1608ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1609ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1610ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1611ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1612ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1613ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1614a703fe33SLois Curfman McInnes 
1615f536c699SSatish Balay .keywords: SNES, set, method
16169b94acceSBarry Smith @*/
16174b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
16189b94acceSBarry Smith {
161984cb2905SBarry Smith   int ierr;
16209b94acceSBarry Smith   int (*r)(SNES);
16213a40ed3dSBarry Smith 
16223a40ed3dSBarry Smith   PetscFunctionBegin;
162377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
162482bf6240SBarry Smith 
162582bf6240SBarry Smith   if (!PetscStrcmp(snes->type_name,method)) PetscFunctionReturn(0);
162682bf6240SBarry Smith 
162782bf6240SBarry Smith   if (snes->setupcalled) {
162882bf6240SBarry Smith     ierr       = (*(snes)->destroy)((PetscObject)snes);CHKERRQ(ierr);
162982bf6240SBarry Smith     snes->data = 0;
163082bf6240SBarry Smith   }
16317d1a2b2bSBarry Smith 
16329b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
163382bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);}
163482bf6240SBarry Smith 
1635*ecf371e4SBarry Smith   ierr =  DLRegisterFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr);
163682bf6240SBarry Smith 
16370452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
163882bf6240SBarry Smith   snes->data = 0;
16393a40ed3dSBarry Smith   ierr = (*r)(snes); CHKERRQ(ierr);
164082bf6240SBarry Smith 
164182bf6240SBarry Smith   if (snes->type_name) PetscFree(snes->type_name);
164282bf6240SBarry Smith   snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name);
164382bf6240SBarry Smith   PetscStrcpy(snes->type_name,method);
164482bf6240SBarry Smith   snes->set_method_called = 1;
164582bf6240SBarry Smith 
16463a40ed3dSBarry Smith   PetscFunctionReturn(0);
16479b94acceSBarry Smith }
16489b94acceSBarry Smith 
1649a847f771SSatish Balay 
16509b94acceSBarry Smith /* --------------------------------------------------------------------- */
16515615d1e5SSatish Balay #undef __FUNC__
1652d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy"
16539b94acceSBarry Smith /*@C
16549b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
16559b94acceSBarry Smith    registered by SNESRegister().
16569b94acceSBarry Smith 
16579b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
16589b94acceSBarry Smith 
16599b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
16609b94acceSBarry Smith @*/
16619b94acceSBarry Smith int SNESRegisterDestroy()
16629b94acceSBarry Smith {
166382bf6240SBarry Smith   int ierr;
166482bf6240SBarry Smith 
16653a40ed3dSBarry Smith   PetscFunctionBegin;
166682bf6240SBarry Smith   if (SNESList) {
166782bf6240SBarry Smith     ierr = DLRegisterDestroy( SNESList );CHKERRQ(ierr);
166882bf6240SBarry Smith     SNESList = 0;
16699b94acceSBarry Smith   }
167084cb2905SBarry Smith   SNESRegisterAllCalled = 0;
16713a40ed3dSBarry Smith   PetscFunctionReturn(0);
16729b94acceSBarry Smith }
16739b94acceSBarry Smith 
16745615d1e5SSatish Balay #undef __FUNC__
1675d4bb536fSBarry Smith #define __FUNC__ "SNESGetType"
16769b94acceSBarry Smith /*@C
16779a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
16789b94acceSBarry Smith 
16799b94acceSBarry Smith    Input Parameter:
16804b0e389bSBarry Smith .  snes - nonlinear solver context
16819b94acceSBarry Smith 
16829b94acceSBarry Smith    Output Parameter:
168382bf6240SBarry Smith .  method - SNES method (a charactor string)
16849b94acceSBarry Smith 
16859b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
16869b94acceSBarry Smith @*/
168782bf6240SBarry Smith int SNESGetType(SNES snes, SNESType *method)
16889b94acceSBarry Smith {
16893a40ed3dSBarry Smith   PetscFunctionBegin;
169082bf6240SBarry Smith   *method = snes->type_name;
16913a40ed3dSBarry Smith   PetscFunctionReturn(0);
16929b94acceSBarry Smith }
16939b94acceSBarry Smith 
16945615d1e5SSatish Balay #undef __FUNC__
1695d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution"
16969b94acceSBarry Smith /*@C
16979b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
16989b94acceSBarry Smith    stored.
16999b94acceSBarry Smith 
17009b94acceSBarry Smith    Input Parameter:
17019b94acceSBarry Smith .  snes - the SNES context
17029b94acceSBarry Smith 
17039b94acceSBarry Smith    Output Parameter:
17049b94acceSBarry Smith .  x - the solution
17059b94acceSBarry Smith 
17069b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
17079b94acceSBarry Smith 
1708a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
17099b94acceSBarry Smith @*/
17109b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
17119b94acceSBarry Smith {
17123a40ed3dSBarry Smith   PetscFunctionBegin;
171377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17149b94acceSBarry Smith   *x = snes->vec_sol_always;
17153a40ed3dSBarry Smith   PetscFunctionReturn(0);
17169b94acceSBarry Smith }
17179b94acceSBarry Smith 
17185615d1e5SSatish Balay #undef __FUNC__
1719d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate"
17209b94acceSBarry Smith /*@C
17219b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
17229b94acceSBarry Smith    stored.
17239b94acceSBarry Smith 
17249b94acceSBarry Smith    Input Parameter:
17259b94acceSBarry Smith .  snes - the SNES context
17269b94acceSBarry Smith 
17279b94acceSBarry Smith    Output Parameter:
17289b94acceSBarry Smith .  x - the solution update
17299b94acceSBarry Smith 
17309b94acceSBarry Smith    Notes:
17319b94acceSBarry Smith    This vector is implementation dependent.
17329b94acceSBarry Smith 
17339b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
17349b94acceSBarry Smith 
17359b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
17369b94acceSBarry Smith @*/
17379b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
17389b94acceSBarry Smith {
17393a40ed3dSBarry Smith   PetscFunctionBegin;
174077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17419b94acceSBarry Smith   *x = snes->vec_sol_update_always;
17423a40ed3dSBarry Smith   PetscFunctionReturn(0);
17439b94acceSBarry Smith }
17449b94acceSBarry Smith 
17455615d1e5SSatish Balay #undef __FUNC__
1746d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction"
17479b94acceSBarry Smith /*@C
17483638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
17499b94acceSBarry Smith 
17509b94acceSBarry Smith    Input Parameter:
17519b94acceSBarry Smith .  snes - the SNES context
17529b94acceSBarry Smith 
17539b94acceSBarry Smith    Output Parameter:
17543638b69dSLois Curfman McInnes .  r - the function
17559b94acceSBarry Smith 
17569b94acceSBarry Smith    Notes:
17579b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
17589b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
17599b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
17609b94acceSBarry Smith 
1761a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
17629b94acceSBarry Smith 
176331615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
176431615d04SLois Curfman McInnes           SNESGetGradient()
17659b94acceSBarry Smith @*/
17669b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
17679b94acceSBarry Smith {
17683a40ed3dSBarry Smith   PetscFunctionBegin;
176977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1770a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1771a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1772a8c6a408SBarry Smith   }
17739b94acceSBarry Smith   *r = snes->vec_func_always;
17743a40ed3dSBarry Smith   PetscFunctionReturn(0);
17759b94acceSBarry Smith }
17769b94acceSBarry Smith 
17775615d1e5SSatish Balay #undef __FUNC__
1778d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient"
17799b94acceSBarry Smith /*@C
17803638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
17819b94acceSBarry Smith 
17829b94acceSBarry Smith    Input Parameter:
17839b94acceSBarry Smith .  snes - the SNES context
17849b94acceSBarry Smith 
17859b94acceSBarry Smith    Output Parameter:
17869b94acceSBarry Smith .  r - the gradient
17879b94acceSBarry Smith 
17889b94acceSBarry Smith    Notes:
17899b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
17909b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
17919b94acceSBarry Smith    SNESGetFunction().
17929b94acceSBarry Smith 
17939b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
17949b94acceSBarry Smith 
179531615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
17969b94acceSBarry Smith @*/
17979b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
17989b94acceSBarry Smith {
17993a40ed3dSBarry Smith   PetscFunctionBegin;
180077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18013a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1802a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
18033a40ed3dSBarry Smith   }
18049b94acceSBarry Smith   *r = snes->vec_func_always;
18053a40ed3dSBarry Smith   PetscFunctionReturn(0);
18069b94acceSBarry Smith }
18079b94acceSBarry Smith 
18085615d1e5SSatish Balay #undef __FUNC__
1809d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction"
18109b94acceSBarry Smith /*@
18119b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
18129b94acceSBarry Smith    unconstrained minimization problems.
18139b94acceSBarry Smith 
18149b94acceSBarry Smith    Input Parameter:
18159b94acceSBarry Smith .  snes - the SNES context
18169b94acceSBarry Smith 
18179b94acceSBarry Smith    Output Parameter:
18189b94acceSBarry Smith .  r - the function
18199b94acceSBarry Smith 
18209b94acceSBarry Smith    Notes:
18219b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
18229b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18239b94acceSBarry Smith    SNESGetFunction().
18249b94acceSBarry Smith 
18259b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
18269b94acceSBarry Smith 
182731615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
18289b94acceSBarry Smith @*/
18299b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
18309b94acceSBarry Smith {
18313a40ed3dSBarry Smith   PetscFunctionBegin;
183277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
183374679c65SBarry Smith   PetscValidScalarPointer(r);
18343a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1835a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
18363a40ed3dSBarry Smith   }
18379b94acceSBarry Smith   *r = snes->fc;
18383a40ed3dSBarry Smith   PetscFunctionReturn(0);
18399b94acceSBarry Smith }
18409b94acceSBarry Smith 
18415615d1e5SSatish Balay #undef __FUNC__
1842d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix"
18433c7409f5SSatish Balay /*@C
18443c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1845d850072dSLois Curfman McInnes    SNES options in the database.
18463c7409f5SSatish Balay 
18473c7409f5SSatish Balay    Input Parameter:
18483c7409f5SSatish Balay .  snes - the SNES context
18493c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18503c7409f5SSatish Balay 
1851d850072dSLois Curfman McInnes    Notes:
1852a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1853a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1854a83b1b31SSatish Balay    hyphen.
1855d850072dSLois Curfman McInnes 
18563c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1857a86d99e1SLois Curfman McInnes 
1858a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
18593c7409f5SSatish Balay @*/
18603c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
18613c7409f5SSatish Balay {
18623c7409f5SSatish Balay   int ierr;
18633c7409f5SSatish Balay 
18643a40ed3dSBarry Smith   PetscFunctionBegin;
186577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1866639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18673c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
18683a40ed3dSBarry Smith   PetscFunctionReturn(0);
18693c7409f5SSatish Balay }
18703c7409f5SSatish Balay 
18715615d1e5SSatish Balay #undef __FUNC__
1872d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix"
18733c7409f5SSatish Balay /*@C
1874f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1875d850072dSLois Curfman McInnes    SNES options in the database.
18763c7409f5SSatish Balay 
18773c7409f5SSatish Balay    Input Parameter:
18783c7409f5SSatish Balay .  snes - the SNES context
18793c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18803c7409f5SSatish Balay 
1881d850072dSLois Curfman McInnes    Notes:
1882a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1883a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1884a83b1b31SSatish Balay    hyphen.
1885d850072dSLois Curfman McInnes 
18863c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1887a86d99e1SLois Curfman McInnes 
1888a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
18893c7409f5SSatish Balay @*/
18903c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
18913c7409f5SSatish Balay {
18923c7409f5SSatish Balay   int ierr;
18933c7409f5SSatish Balay 
18943a40ed3dSBarry Smith   PetscFunctionBegin;
189577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1896639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18973c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
18983a40ed3dSBarry Smith   PetscFunctionReturn(0);
18993c7409f5SSatish Balay }
19003c7409f5SSatish Balay 
19015615d1e5SSatish Balay #undef __FUNC__
1902d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix"
1903c83590e2SSatish Balay /*@
19043c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
19053c7409f5SSatish Balay    SNES options in the database.
19063c7409f5SSatish Balay 
19073c7409f5SSatish Balay    Input Parameter:
19083c7409f5SSatish Balay .  snes - the SNES context
19093c7409f5SSatish Balay 
19103c7409f5SSatish Balay    Output Parameter:
19113c7409f5SSatish Balay .  prefix - pointer to the prefix string used
19123c7409f5SSatish Balay 
19133c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1914a86d99e1SLois Curfman McInnes 
1915a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
19163c7409f5SSatish Balay @*/
19173c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
19183c7409f5SSatish Balay {
19193c7409f5SSatish Balay   int ierr;
19203c7409f5SSatish Balay 
19213a40ed3dSBarry Smith   PetscFunctionBegin;
192277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1923639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19243a40ed3dSBarry Smith   PetscFunctionReturn(0);
19253c7409f5SSatish Balay }
19263c7409f5SSatish Balay 
1927ca161407SBarry Smith #undef __FUNC__
1928ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp"
1929ca161407SBarry Smith /*@
1930ca161407SBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
1931ca161407SBarry Smith 
1932ca161407SBarry Smith    Input Parameter:
1933ca161407SBarry Smith .  snes - the SNES context
1934ca161407SBarry Smith 
1935ca161407SBarry Smith    Options Database Keys:
1936ca161407SBarry Smith $  -help, -h
1937ca161407SBarry Smith 
1938ca161407SBarry Smith .keywords: SNES, nonlinear, help
1939ca161407SBarry Smith 
1940ca161407SBarry Smith .seealso: SNESSetFromOptions()
1941ca161407SBarry Smith @*/
1942ca161407SBarry Smith int SNESPrintHelp(SNES snes)
1943ca161407SBarry Smith {
1944ca161407SBarry Smith   char                p[64];
1945ca161407SBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
1946ca161407SBarry Smith   int                 ierr;
1947ca161407SBarry Smith 
1948ca161407SBarry Smith   PetscFunctionBegin;
1949ca161407SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1950ca161407SBarry Smith 
1951ca161407SBarry Smith   PetscStrcpy(p,"-");
1952ca161407SBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
1953ca161407SBarry Smith 
1954ca161407SBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
1955ca161407SBarry Smith 
195676be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");
195782bf6240SBarry Smith   ierr = DLRegisterPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr);
195876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
195976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
196076be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
196176be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
196276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
196376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
196476be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
196576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");
196676be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
196776be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
1968ca161407SBarry Smith     residual norm at each iteration.\n",p);
196976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
1970ca161407SBarry Smith     residual norm for small residual norms. This is useful to conceal\n\
1971ca161407SBarry Smith     meaningless digits that may be different on different machines.\n",p);
197276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
1973ca161407SBarry Smith   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
197476be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1975ca161407SBarry Smith      " Options for solving systems of nonlinear equations only:\n");
197676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
197776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
197876be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
197976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
198076be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
198176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1982ca161407SBarry Smith      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
198376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1984ca161407SBarry Smith      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
198576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1986ca161407SBarry Smith      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
198776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1988ca161407SBarry Smith      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
198976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1990ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
199176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1992ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
199376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1994ca161407SBarry Smith      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
1995ca161407SBarry Smith   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
199676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1997ca161407SBarry Smith      " Options for solving unconstrained minimization problems only:\n");
199876be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
199976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
200076be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
2001ca161407SBarry Smith   }
200276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
200376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"a particular method\n");
2004ca161407SBarry Smith   if (snes->printhelp) {
2005ca161407SBarry Smith     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2006ca161407SBarry Smith   }
2007ca161407SBarry Smith   PetscFunctionReturn(0);
2008ca161407SBarry Smith }
20093c7409f5SSatish Balay 
20109b94acceSBarry Smith 
20119b94acceSBarry Smith 
20129b94acceSBarry Smith 
20133a40ed3dSBarry Smith 
2014