xref: /petsc/src/snes/interface/snes.c (revision ef1dfb1141b5dd2eae9a1667e2446c2015e70e16)
19b94acceSBarry Smith #ifndef lint
2*ef1dfb11SLois Curfman McInnes static char vcid[] = "$Id: snes.c,v 1.84 1996/08/26 02:35:13 bsmith Exp curfman $";
39b94acceSBarry Smith #endif
49b94acceSBarry Smith 
59b94acceSBarry Smith #include "draw.h"          /*I "draw.h"  I*/
670f55243SBarry Smith #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
7f5eb4b81SSatish Balay #include "src/sys/nreg.h"
89b94acceSBarry Smith #include "pinclude/pviewer.h"
99b94acceSBarry Smith #include <math.h>
109b94acceSBarry Smith 
11052efed2SBarry Smith extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*);
1233455573SSatish Balay extern int SNESPrintTypes_Private(MPI_Comm,char*,char*);
139b94acceSBarry Smith 
149b94acceSBarry Smith /*@
159b94acceSBarry Smith    SNESView - Prints the SNES data structure.
169b94acceSBarry Smith 
179b94acceSBarry Smith    Input Parameters:
189b94acceSBarry Smith .  SNES - the SNES context
199b94acceSBarry Smith .  viewer - visualization context
209b94acceSBarry Smith 
219b94acceSBarry Smith    Options Database Key:
229b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
239b94acceSBarry Smith 
249b94acceSBarry Smith    Notes:
259b94acceSBarry Smith    The available visualization contexts include
266d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
276d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
289b94acceSBarry Smith $       output where only the first processor opens
299b94acceSBarry Smith $       the file.  All other processors send their
309b94acceSBarry Smith $       data to the first processor to print.
319b94acceSBarry Smith 
329b94acceSBarry Smith    The user can open alternative vistualization contexts with
33dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
349b94acceSBarry Smith 
359b94acceSBarry Smith .keywords: SNES, view
369b94acceSBarry Smith 
37dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
389b94acceSBarry Smith @*/
399b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
409b94acceSBarry Smith {
419b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
429b94acceSBarry Smith   FILE                *fd;
439b94acceSBarry Smith   int                 ierr;
449b94acceSBarry Smith   SLES                sles;
459b94acceSBarry Smith   char                *method;
4619bcc07fSBarry Smith   ViewerType          vtype;
479b94acceSBarry Smith 
4874679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4974679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5074679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5174679c65SBarry Smith 
5219bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5319bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5490ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5577c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
564b0e389bSBarry Smith     SNESGetType(snes,PETSC_NULL,&method);
5777c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
589b94acceSBarry Smith     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
5977c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
609b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
619b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6277c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
639b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
649b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
659b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
6677c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
679b94acceSBarry Smith     if (snes->ksp_ewconv) {
689b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
699b94acceSBarry Smith       if (kctx) {
7077c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
719b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
729b94acceSBarry Smith         kctx->version);
7377c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
749b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
759b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
7677c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
779b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
789b94acceSBarry Smith       }
799b94acceSBarry Smith     }
80c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
81c30f285eSLois Curfman McInnes     SNESGetType(snes,PETSC_NULL,&method);
82c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
8319bcc07fSBarry Smith   }
849b94acceSBarry Smith   SNESGetSLES(snes,&sles);
859b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
869b94acceSBarry Smith   return 0;
879b94acceSBarry Smith }
889b94acceSBarry Smith 
899b94acceSBarry Smith /*@
90682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
919b94acceSBarry Smith 
929b94acceSBarry Smith    Input Parameter:
939b94acceSBarry Smith .  snes - the SNES context
949b94acceSBarry Smith 
959b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
969b94acceSBarry Smith 
97a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
989b94acceSBarry Smith @*/
999b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1009b94acceSBarry Smith {
1014b0e389bSBarry Smith   SNESType method;
1029b94acceSBarry Smith   double   tmp;
1039b94acceSBarry Smith   SLES     sles;
104d31a3109SLois Curfman McInnes   int      ierr, flg, mset = 0;
1053c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1069b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1079b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1089b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1099b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1109b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1119b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1129b94acceSBarry Smith 
11377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
11474679c65SBarry Smith   if (snes->setup_called) SETERRQ(1,"SNESSetFromOptions:Must call prior to SNESSetUp");
115052efed2SBarry Smith   ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr);
116052efed2SBarry Smith   if (flg) {
117052efed2SBarry Smith     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
1189b94acceSBarry Smith   }
1195334005bSBarry Smith   else if (!snes->set_method_called) {
1205334005bSBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
12140191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
1225334005bSBarry Smith     }
1235334005bSBarry Smith     else {
12440191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
1255334005bSBarry Smith     }
1265334005bSBarry Smith   }
1275334005bSBarry Smith 
1283c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1293c7409f5SSatish Balay   if (flg) { SNESPrintHelp(snes); }
1303c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
131d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);}
1323c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
133d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);}
1343c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
135d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); }
1363c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1373c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
138d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
139d7a720efSLois Curfman McInnes   if (flg) { SNESSetTrustRegionTolerance(snes,tmp); }
140d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
141d7a720efSLois Curfman McInnes   if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); }
1423c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1433c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
144b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
145b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
146b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
147b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
148b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
149b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
150b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
1519b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
1529b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
1533c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
154d31a3109SLois Curfman McInnes   if (flg) {SNESSetMonitor(snes,SNESDefaultMonitor,0); mset = 1;}
1553c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
156d31a3109SLois Curfman McInnes   if (flg) {
157d31a3109SLois Curfman McInnes    if (mset) SETERRQ(1,"SNESSetFromOptions:Monitor for SNES is already set");
158d31a3109SLois Curfman McInnes    SNESSetMonitor(snes,SNESDefaultSMonitor,0); mset = 1;
159d31a3109SLois Curfman McInnes   }
1603c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_xmonitor",&flg); CHKERRQ(ierr);
1613c7409f5SSatish Balay   if (flg) {
16217699dbbSLois Curfman McInnes     int    rank = 0;
163d7e8b826SBarry Smith     DrawLG lg;
164d31a3109SLois Curfman McInnes     if (mset) SETERRQ(1,"SNESSetFromOptions:Monitor for SNES is already set");
16517699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
16617699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
16717699dbbSLois Curfman McInnes     if (!rank) {
1689b94acceSBarry Smith       ierr = SNESLGMonitorCreate(0,0,0,0,300,300,&lg); CHKERRQ(ierr);
1699b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
170f8c826e1SBarry Smith       snes->xmonitor = lg;
1719b94acceSBarry Smith       PLogObjectParent(snes,lg);
1729b94acceSBarry Smith     }
173d31a3109SLois Curfman McInnes     mset = 1;
1749b94acceSBarry Smith   }
1753c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
1763c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1779b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
1789b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
17994a424c1SBarry Smith     PLogInfo(snes,
18031615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
18131615d04SLois Curfman McInnes   }
18231615d04SLois Curfman McInnes   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
18331615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
18431615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
18594a424c1SBarry Smith     PLogInfo(snes,
18631615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
1879b94acceSBarry Smith   }
1889b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1899b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
1909b94acceSBarry Smith   if (!snes->setfromoptions) return 0;
1919b94acceSBarry Smith   return (*snes->setfromoptions)(snes);
1929b94acceSBarry Smith }
1939b94acceSBarry Smith 
1949b94acceSBarry Smith /*@
1959b94acceSBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
1969b94acceSBarry Smith 
1979b94acceSBarry Smith    Input Parameter:
1989b94acceSBarry Smith .  snes - the SNES context
1999b94acceSBarry Smith 
200a703fe33SLois Curfman McInnes    Options Database Keys:
201a703fe33SLois Curfman McInnes $  -help, -h
202a703fe33SLois Curfman McInnes 
2039b94acceSBarry Smith .keywords: SNES, nonlinear, help
2049b94acceSBarry Smith 
205682d7d0cSBarry Smith .seealso: SNESSetFromOptions()
2069b94acceSBarry Smith @*/
2079b94acceSBarry Smith int SNESPrintHelp(SNES snes)
2089b94acceSBarry Smith {
2096daaf66cSBarry Smith   char                p[64];
2106daaf66cSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2116daaf66cSBarry Smith 
21277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2136daaf66cSBarry Smith 
2146daaf66cSBarry Smith   PetscStrcpy(p,"-");
2156daaf66cSBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2166daaf66cSBarry Smith 
2176daaf66cSBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2186daaf66cSBarry Smith 
219d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
22033455573SSatish Balay   SNESPrintTypes_Private(snes->comm,p,"snes_type");
22177c4ece6SBarry Smith   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
222b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
223b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
224b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
225b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
226b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
227b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
228d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm," SNES Monitoring Options: Choose 1 of the following\n");
229d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor\n",p);
230d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
231b18e04deSLois Curfman McInnes   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
23277c4ece6SBarry Smith     PetscPrintf(snes->comm,
233d31a3109SLois Curfman McInnes      " Options for solving systems of nonlinear equations only:\n");
23477c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
23577c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
236*ef1dfb11SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
2371650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p);
2381650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p);
23977c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
24077c4ece6SBarry Smith     PetscPrintf(snes->comm,
241b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
24277c4ece6SBarry Smith     PetscPrintf(snes->comm,
243b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
24477c4ece6SBarry Smith     PetscPrintf(snes->comm,
245b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
24677c4ece6SBarry Smith     PetscPrintf(snes->comm,
247b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
24877c4ece6SBarry Smith     PetscPrintf(snes->comm,
249b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
25077c4ece6SBarry Smith     PetscPrintf(snes->comm,
251b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
25277c4ece6SBarry Smith     PetscPrintf(snes->comm,
253b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
254b18e04deSLois Curfman McInnes   }
255b18e04deSLois Curfman McInnes   else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
25677c4ece6SBarry Smith     PetscPrintf(snes->comm,
257d31a3109SLois Curfman McInnes      " Options for solving unconstrained minimization problems only:\n");
258b18e04deSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
25977c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
26077c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
261d31a3109SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p);
262d31a3109SLois Curfman McInnes      PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p);
263b18e04deSLois Curfman McInnes   }
2644537a946SLois Curfman McInnes   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
26577c4ece6SBarry Smith   PetscPrintf(snes->comm,"a particular method\n");
2666daaf66cSBarry Smith   if (snes->printhelp) (*snes->printhelp)(snes,p);
2679b94acceSBarry Smith   return 0;
2689b94acceSBarry Smith }
2699b94acceSBarry Smith /*@
2709b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2719b94acceSBarry Smith    the nonlinear solvers.
2729b94acceSBarry Smith 
2739b94acceSBarry Smith    Input Parameters:
2749b94acceSBarry Smith .  snes - the SNES context
2759b94acceSBarry Smith .  usrP - optional user context
2769b94acceSBarry Smith 
2779b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2789b94acceSBarry Smith 
2799b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2809b94acceSBarry Smith @*/
2819b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
2829b94acceSBarry Smith {
28377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2849b94acceSBarry Smith   snes->user		= usrP;
2859b94acceSBarry Smith   return 0;
2869b94acceSBarry Smith }
28774679c65SBarry Smith 
2889b94acceSBarry Smith /*@C
2899b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
2909b94acceSBarry Smith    nonlinear solvers.
2919b94acceSBarry Smith 
2929b94acceSBarry Smith    Input Parameter:
2939b94acceSBarry Smith .  snes - SNES context
2949b94acceSBarry Smith 
2959b94acceSBarry Smith    Output Parameter:
2969b94acceSBarry Smith .  usrP - user context
2979b94acceSBarry Smith 
2989b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
2999b94acceSBarry Smith 
3009b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3019b94acceSBarry Smith @*/
3029b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3039b94acceSBarry Smith {
30477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3059b94acceSBarry Smith   *usrP = snes->user;
3069b94acceSBarry Smith   return 0;
3079b94acceSBarry Smith }
30874679c65SBarry Smith 
3099b94acceSBarry Smith /*@
3109b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3119b94acceSBarry Smith    nonlinear solver.
3129b94acceSBarry Smith 
3139b94acceSBarry Smith    Input Parameter:
3149b94acceSBarry Smith .  snes - SNES context
3159b94acceSBarry Smith 
3169b94acceSBarry Smith    Output Parameter:
3179b94acceSBarry Smith .  iter - iteration number
3189b94acceSBarry Smith 
3199b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3209b94acceSBarry Smith @*/
3219b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3229b94acceSBarry Smith {
32377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
32474679c65SBarry Smith   PetscValidIntPointer(iter);
3259b94acceSBarry Smith   *iter = snes->iter;
3269b94acceSBarry Smith   return 0;
3279b94acceSBarry Smith }
32874679c65SBarry Smith 
3299b94acceSBarry Smith /*@
3309b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3319b94acceSBarry Smith    with SNESSSetFunction().
3329b94acceSBarry Smith 
3339b94acceSBarry Smith    Input Parameter:
3349b94acceSBarry Smith .  snes - SNES context
3359b94acceSBarry Smith 
3369b94acceSBarry Smith    Output Parameter:
3379b94acceSBarry Smith .  fnorm - 2-norm of function
3389b94acceSBarry Smith 
3399b94acceSBarry Smith    Note:
3409b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
341a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
342a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3439b94acceSBarry Smith 
3449b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
345a86d99e1SLois Curfman McInnes 
346a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3479b94acceSBarry Smith @*/
3489b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3499b94acceSBarry Smith {
35077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
35174679c65SBarry Smith   PetscValidScalarPointer(fnorm);
35274679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
35374679c65SBarry Smith     SETERRQ(1,"SNESGetFunctionNorm:For SNES_NONLINEAR_EQUATIONS only");
35474679c65SBarry Smith   }
3559b94acceSBarry Smith   *fnorm = snes->norm;
3569b94acceSBarry Smith   return 0;
3579b94acceSBarry Smith }
35874679c65SBarry Smith 
3599b94acceSBarry Smith /*@
3609b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
3619b94acceSBarry Smith    with SNESSSetGradient().
3629b94acceSBarry Smith 
3639b94acceSBarry Smith    Input Parameter:
3649b94acceSBarry Smith .  snes - SNES context
3659b94acceSBarry Smith 
3669b94acceSBarry Smith    Output Parameter:
3679b94acceSBarry Smith .  fnorm - 2-norm of gradient
3689b94acceSBarry Smith 
3699b94acceSBarry Smith    Note:
3709b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
371a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
372a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
3739b94acceSBarry Smith 
3749b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
375a86d99e1SLois Curfman McInnes 
376a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
3779b94acceSBarry Smith @*/
3789b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
3799b94acceSBarry Smith {
38077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
38174679c65SBarry Smith   PetscValidScalarPointer(gnorm);
38274679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
38374679c65SBarry Smith     SETERRQ(1,"SNESGetGradientNorm:For SNES_UNCONSTRAINED_MINIMIZATION only");
38474679c65SBarry Smith   }
3859b94acceSBarry Smith   *gnorm = snes->norm;
3869b94acceSBarry Smith   return 0;
3879b94acceSBarry Smith }
38874679c65SBarry Smith 
3899b94acceSBarry Smith /*@
3909b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
3919b94acceSBarry Smith    attempted by the nonlinear solver.
3929b94acceSBarry Smith 
3939b94acceSBarry Smith    Input Parameter:
3949b94acceSBarry Smith .  snes - SNES context
3959b94acceSBarry Smith 
3969b94acceSBarry Smith    Output Parameter:
3979b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
3989b94acceSBarry Smith 
3999b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4009b94acceSBarry Smith @*/
4019b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4029b94acceSBarry Smith {
40377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
40474679c65SBarry Smith   PetscValidIntPointer(nfails);
4059b94acceSBarry Smith   *nfails = snes->nfailures;
4069b94acceSBarry Smith   return 0;
4079b94acceSBarry Smith }
4089b94acceSBarry Smith /*@C
4099b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4109b94acceSBarry Smith 
4119b94acceSBarry Smith    Input Parameter:
4129b94acceSBarry Smith .  snes - the SNES context
4139b94acceSBarry Smith 
4149b94acceSBarry Smith    Output Parameter:
4159b94acceSBarry Smith .  sles - the SLES context
4169b94acceSBarry Smith 
4179b94acceSBarry Smith    Notes:
4189b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4199b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4209b94acceSBarry Smith    KSP and PC contexts as well.
4219b94acceSBarry Smith 
4229b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4239b94acceSBarry Smith 
4249b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4259b94acceSBarry Smith @*/
4269b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4279b94acceSBarry Smith {
42877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4299b94acceSBarry Smith   *sles = snes->sles;
4309b94acceSBarry Smith   return 0;
4319b94acceSBarry Smith }
4329b94acceSBarry Smith /* -----------------------------------------------------------*/
4339b94acceSBarry Smith /*@C
4349b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
4359b94acceSBarry Smith 
4369b94acceSBarry Smith    Input Parameter:
4379b94acceSBarry Smith .  comm - MPI communicator
4389b94acceSBarry Smith .  type - type of method, one of
4399b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
4409b94acceSBarry Smith $      (for systems of nonlinear equations)
4419b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
4429b94acceSBarry Smith $      (for unconstrained minimization)
4439b94acceSBarry Smith 
4449b94acceSBarry Smith    Output Parameter:
4459b94acceSBarry Smith .  outsnes - the new SNES context
4469b94acceSBarry Smith 
447c1f60f51SBarry Smith    Options Database Key:
448c1f60f51SBarry Smith $   -snes_mf - use default matrix free for Jacobian and preconditioner
449c1f60f51SBarry Smith $   -snes_mf_operator - use default matrix free for Jacobian but not for preconditioner
450c1f60f51SBarry Smith $   -snes_fd - use slow finite differences to compute Jacobian
451c1f60f51SBarry Smith 
4529b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
4539b94acceSBarry Smith 
45463a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
4559b94acceSBarry Smith @*/
4564b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
4579b94acceSBarry Smith {
4589b94acceSBarry Smith   int                 ierr;
4599b94acceSBarry Smith   SNES                snes;
4609b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
46137fcc0dbSBarry Smith 
4629b94acceSBarry Smith   *outsnes = 0;
4632a0acf01SLois Curfman McInnes   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS)
4642a0acf01SLois Curfman McInnes     SETERRQ(1,"SNESCreate:incorrect method type");
4650452661fSBarry Smith   PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm);
4669b94acceSBarry Smith   PLogObjectCreate(snes);
4679b94acceSBarry Smith   snes->max_its           = 50;
4689b94acceSBarry Smith   snes->max_funcs	  = 1000;
4699b94acceSBarry Smith   snes->norm		  = 0.0;
4705a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
471b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
472b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
4739b94acceSBarry Smith     snes->atol		  = 1.e-10;
4745a2d0531SBarry Smith   }
475b4874afaSBarry Smith   else {
476b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
477b4874afaSBarry Smith     snes->ttol            = 0.0;
478b4874afaSBarry Smith     snes->atol		  = 1.e-50;
479b4874afaSBarry Smith   }
4809b94acceSBarry Smith   snes->xtol		  = 1.e-8;
481b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
4829b94acceSBarry Smith   snes->nfuncs            = 0;
4839b94acceSBarry Smith   snes->nfailures         = 0;
4849b94acceSBarry Smith   snes->monitor           = 0;
4859b94acceSBarry Smith   snes->data              = 0;
4869b94acceSBarry Smith   snes->view              = 0;
4879b94acceSBarry Smith   snes->computeumfunction = 0;
4889b94acceSBarry Smith   snes->umfunP            = 0;
4899b94acceSBarry Smith   snes->fc                = 0;
4909b94acceSBarry Smith   snes->deltatol          = 1.e-12;
4919b94acceSBarry Smith   snes->fmin              = -1.e30;
4929b94acceSBarry Smith   snes->method_class      = type;
4939b94acceSBarry Smith   snes->set_method_called = 0;
494a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
4959b94acceSBarry Smith   snes->ksp_ewconv        = 0;
4967d1a2b2bSBarry Smith   snes->type              = -1;
4976f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
4986f24a144SLois Curfman McInnes   snes->vwork             = 0;
4996f24a144SLois Curfman McInnes   snes->nwork             = 0;
5009b94acceSBarry Smith 
5019b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5020452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
5039b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5049b94acceSBarry Smith   kctx->version     = 2;
5059b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5069b94acceSBarry Smith                              this was too large for some test cases */
5079b94acceSBarry Smith   kctx->rtol_last   = 0;
5089b94acceSBarry Smith   kctx->rtol_max    = .9;
5099b94acceSBarry Smith   kctx->gamma       = 1.0;
5109b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5119b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5129b94acceSBarry Smith   kctx->threshold   = .1;
5139b94acceSBarry Smith   kctx->lresid_last = 0;
5149b94acceSBarry Smith   kctx->norm_last   = 0;
5159b94acceSBarry Smith 
5169b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5179b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5185334005bSBarry Smith 
5199b94acceSBarry Smith   *outsnes = snes;
5209b94acceSBarry Smith   return 0;
5219b94acceSBarry Smith }
5229b94acceSBarry Smith 
5239b94acceSBarry Smith /* --------------------------------------------------------------- */
5249b94acceSBarry Smith /*@C
5259b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
5269b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
5279b94acceSBarry Smith    equations.
5289b94acceSBarry Smith 
5299b94acceSBarry Smith    Input Parameters:
5309b94acceSBarry Smith .  snes - the SNES context
5319b94acceSBarry Smith .  func - function evaluation routine
5329b94acceSBarry Smith .  ctx - optional user-defined function context
5339b94acceSBarry Smith .  r - vector to store function value
5349b94acceSBarry Smith 
5359b94acceSBarry Smith    Calling sequence of func:
5369b94acceSBarry Smith .  func (SNES, Vec x, Vec f, void *ctx);
5379b94acceSBarry Smith 
5389b94acceSBarry Smith .  x - input vector
5395334005bSBarry Smith .  f - vector function
5409b94acceSBarry Smith .  ctx - optional user-defined context for private data for the
5419b94acceSBarry Smith          function evaluation routine (may be null)
5429b94acceSBarry Smith 
5439b94acceSBarry Smith    Notes:
5449b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
5459b94acceSBarry Smith $      f'(x) x = -f(x),
5469b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
5479b94acceSBarry Smith 
5489b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
5499b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
5509b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
5519b94acceSBarry Smith 
5529b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
5539b94acceSBarry Smith 
554a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
5559b94acceSBarry Smith @*/
5565334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
5579b94acceSBarry Smith {
55877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5599b94acceSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
56048d91487SBarry Smith     "SNESSetFunction:For SNES_NONLINEAR_EQUATIONS only");
5619b94acceSBarry Smith   snes->computefunction     = func;
5629b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
5639b94acceSBarry Smith   snes->funP                = ctx;
5649b94acceSBarry Smith   return 0;
5659b94acceSBarry Smith }
5669b94acceSBarry Smith 
5679b94acceSBarry Smith /*@
5689b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
5699b94acceSBarry Smith    SNESSetFunction().
5709b94acceSBarry Smith 
5719b94acceSBarry Smith    Input Parameters:
5729b94acceSBarry Smith .  snes - the SNES context
5739b94acceSBarry Smith .  x - input vector
5749b94acceSBarry Smith 
5759b94acceSBarry Smith    Output Parameter:
5763638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
5779b94acceSBarry Smith 
5781bffabb2SLois Curfman McInnes    Notes:
5799b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
5809b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
5819b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
5829b94acceSBarry Smith 
5839b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
5849b94acceSBarry Smith 
585a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
5869b94acceSBarry Smith @*/
5879b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
5889b94acceSBarry Smith {
5899b94acceSBarry Smith   int    ierr;
5909b94acceSBarry Smith 
59174679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
59274679c65SBarry Smith     SETERRQ(1,"SNESComputeFunction: For SNES_NONLINEAR_EQUATIONS only");
5939b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
5949b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
5959b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
5969b94acceSBarry Smith   return 0;
5979b94acceSBarry Smith }
5989b94acceSBarry Smith 
5999b94acceSBarry Smith /*@C
6009b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6019b94acceSBarry Smith    unconstrained minimization.
6029b94acceSBarry Smith 
6039b94acceSBarry Smith    Input Parameters:
6049b94acceSBarry Smith .  snes - the SNES context
6059b94acceSBarry Smith .  func - function evaluation routine
6069b94acceSBarry Smith .  ctx - optional user-defined function context
6079b94acceSBarry Smith 
6089b94acceSBarry Smith    Calling sequence of func:
6099b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
6109b94acceSBarry Smith 
6119b94acceSBarry Smith .  x - input vector
6129b94acceSBarry Smith .  f - function
6139b94acceSBarry Smith .  ctx - optional user-defined context for private data for the
6149b94acceSBarry Smith          function evaluation routine (may be null)
6159b94acceSBarry Smith 
6169b94acceSBarry Smith    Notes:
6179b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
6189b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
6199b94acceSBarry Smith    SNESSetFunction().
6209b94acceSBarry Smith 
6219b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
6229b94acceSBarry Smith 
623a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
624a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
6259b94acceSBarry Smith @*/
6269b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
6279b94acceSBarry Smith                       void *ctx)
6289b94acceSBarry Smith {
62977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
6309b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
63148d91487SBarry Smith     "SNESSetMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION");
6329b94acceSBarry Smith   snes->computeumfunction   = func;
6339b94acceSBarry Smith   snes->umfunP              = ctx;
6349b94acceSBarry Smith   return 0;
6359b94acceSBarry Smith }
6369b94acceSBarry Smith 
6379b94acceSBarry Smith /*@
6389b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
6399b94acceSBarry Smith    set with SNESSetMinimizationFunction().
6409b94acceSBarry Smith 
6419b94acceSBarry Smith    Input Parameters:
6429b94acceSBarry Smith .  snes - the SNES context
6439b94acceSBarry Smith .  x - input vector
6449b94acceSBarry Smith 
6459b94acceSBarry Smith    Output Parameter:
6469b94acceSBarry Smith .  y - function value
6479b94acceSBarry Smith 
6489b94acceSBarry Smith    Notes:
6499b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
6509b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
6519b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
652a86d99e1SLois Curfman McInnes 
653a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
654a86d99e1SLois Curfman McInnes 
655a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
656a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
6579b94acceSBarry Smith @*/
6589b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
6599b94acceSBarry Smith {
6609b94acceSBarry Smith   int    ierr;
6619b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
66248d91487SBarry Smith     "SNESComputeMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION");
6639b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
6649b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
6659b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
6669b94acceSBarry Smith   return 0;
6679b94acceSBarry Smith }
6689b94acceSBarry Smith 
6699b94acceSBarry Smith /*@C
6709b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
6719b94acceSBarry Smith    vector for use by the SNES routines.
6729b94acceSBarry Smith 
6739b94acceSBarry Smith    Input Parameters:
6749b94acceSBarry Smith .  snes - the SNES context
6759b94acceSBarry Smith .  func - function evaluation routine
6769b94acceSBarry Smith .  ctx - optional user-defined function context
6779b94acceSBarry Smith .  r - vector to store gradient value
6789b94acceSBarry Smith 
6799b94acceSBarry Smith    Calling sequence of func:
6809b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
6819b94acceSBarry Smith 
6829b94acceSBarry Smith .  x - input vector
6839b94acceSBarry Smith .  g - gradient vector
6849b94acceSBarry Smith .  ctx - optional user-defined context for private data for the
6859b94acceSBarry Smith          function evaluation routine (may be null)
6869b94acceSBarry Smith 
6879b94acceSBarry Smith    Notes:
6889b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
6899b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
6909b94acceSBarry Smith    SNESSetFunction().
6919b94acceSBarry Smith 
6929b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6939b94acceSBarry Smith 
694a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
695a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
6969b94acceSBarry Smith @*/
69774679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
6989b94acceSBarry Smith {
69977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
7009b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
70148d91487SBarry Smith     "SNESSetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
7029b94acceSBarry Smith   snes->computefunction     = func;
7039b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7049b94acceSBarry Smith   snes->funP                = ctx;
7059b94acceSBarry Smith   return 0;
7069b94acceSBarry Smith }
7079b94acceSBarry Smith 
7089b94acceSBarry Smith /*@
709a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
710a86d99e1SLois Curfman McInnes    SNESSetGradient().
7119b94acceSBarry Smith 
7129b94acceSBarry Smith    Input Parameters:
7139b94acceSBarry Smith .  snes - the SNES context
7149b94acceSBarry Smith .  x - input vector
7159b94acceSBarry Smith 
7169b94acceSBarry Smith    Output Parameter:
7179b94acceSBarry Smith .  y - gradient vector
7189b94acceSBarry Smith 
7199b94acceSBarry Smith    Notes:
7209b94acceSBarry Smith    SNESComputeGradient() is valid only for
7219b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7229b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
723a86d99e1SLois Curfman McInnes 
724a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
725a86d99e1SLois Curfman McInnes 
726a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
727a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
7289b94acceSBarry Smith @*/
7299b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
7309b94acceSBarry Smith {
7319b94acceSBarry Smith   int    ierr;
73274679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
73374679c65SBarry Smith     SETERRQ(1,"SNESComputeGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
7349b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
7359b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
7369b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
7379b94acceSBarry Smith   return 0;
7389b94acceSBarry Smith }
7399b94acceSBarry Smith 
74062fef451SLois Curfman McInnes /*@
74162fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
74262fef451SLois Curfman McInnes    set with SNESSetJacobian().
74362fef451SLois Curfman McInnes 
74462fef451SLois Curfman McInnes    Input Parameters:
74562fef451SLois Curfman McInnes .  snes - the SNES context
74662fef451SLois Curfman McInnes .  x - input vector
74762fef451SLois Curfman McInnes 
74862fef451SLois Curfman McInnes    Output Parameters:
74962fef451SLois Curfman McInnes .  A - Jacobian matrix
75062fef451SLois Curfman McInnes .  B - optional preconditioning matrix
75162fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
75262fef451SLois Curfman McInnes 
75362fef451SLois Curfman McInnes    Notes:
75462fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
75562fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
75662fef451SLois Curfman McInnes 
75762fef451SLois Curfman McInnes    See SLESSetOperators() for information about setting the flag parameter.
75862fef451SLois Curfman McInnes 
75962fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
76062fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
76162fef451SLois Curfman McInnes    methods is SNESComputeJacobian().
76262fef451SLois Curfman McInnes 
76362fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
76462fef451SLois Curfman McInnes 
76562fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
76662fef451SLois Curfman McInnes @*/
7679b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
7689b94acceSBarry Smith {
7699b94acceSBarry Smith   int    ierr;
77074679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
77174679c65SBarry Smith     SETERRQ(1,"SNESComputeJacobian: For SNES_NONLINEAR_EQUATIONS only");
7729b94acceSBarry Smith   if (!snes->computejacobian) return 0;
7739b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
774c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
7759b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
7769b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
7776d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
77877c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
77977c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
7809b94acceSBarry Smith   return 0;
7819b94acceSBarry Smith }
7829b94acceSBarry Smith 
78362fef451SLois Curfman McInnes /*@
78462fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
78562fef451SLois Curfman McInnes    set with SNESSetHessian().
78662fef451SLois Curfman McInnes 
78762fef451SLois Curfman McInnes    Input Parameters:
78862fef451SLois Curfman McInnes .  snes - the SNES context
78962fef451SLois Curfman McInnes .  x - input vector
79062fef451SLois Curfman McInnes 
79162fef451SLois Curfman McInnes    Output Parameters:
79262fef451SLois Curfman McInnes .  A - Hessian matrix
79362fef451SLois Curfman McInnes .  B - optional preconditioning matrix
79462fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
79562fef451SLois Curfman McInnes 
79662fef451SLois Curfman McInnes    Notes:
79762fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
79862fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
79962fef451SLois Curfman McInnes 
80062fef451SLois Curfman McInnes    See SLESSetOperators() for information about setting the flag parameter.
80162fef451SLois Curfman McInnes 
80262fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
80362fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
80462fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
80562fef451SLois Curfman McInnes 
80662fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
80762fef451SLois Curfman McInnes 
808a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
809a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
81062fef451SLois Curfman McInnes @*/
81162fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
8129b94acceSBarry Smith {
8139b94acceSBarry Smith   int    ierr;
81474679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
81574679c65SBarry Smith     SETERRQ(1,"SNESComputeHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
8169b94acceSBarry Smith   if (!snes->computejacobian) return 0;
81762fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
8180f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
81962fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
82062fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
8210f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
82277c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
82377c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8249b94acceSBarry Smith   return 0;
8259b94acceSBarry Smith }
8269b94acceSBarry Smith 
8279b94acceSBarry Smith /*@C
8289b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
8299b94acceSBarry Smith    location to store it.
8309b94acceSBarry Smith 
8319b94acceSBarry Smith    Input Parameters:
8329b94acceSBarry Smith .  snes - the SNES context
8339b94acceSBarry Smith .  A - Jacobian matrix
8349b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
8359b94acceSBarry Smith .  func - Jacobian evaluation routine
8369b94acceSBarry Smith .  ctx - optional user-defined context for private data for the
8379b94acceSBarry Smith          Jacobian evaluation routine (may be null)
8389b94acceSBarry Smith 
8399b94acceSBarry Smith    Calling sequence of func:
8409b94acceSBarry Smith .  func (SNES,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
8419b94acceSBarry Smith 
8429b94acceSBarry Smith .  x - input vector
8439b94acceSBarry Smith .  A - Jacobian matrix
8449b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
845ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
846ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
8479b94acceSBarry Smith .  ctx - optional user-defined Jacobian context
8489b94acceSBarry Smith 
8499b94acceSBarry Smith    Notes:
850ac21db08SLois Curfman McInnes    See SLESSetOperators() for information about setting the flag input
851ac21db08SLois Curfman McInnes    parameter in the routine func().  Be sure to read this information!
852ac21db08SLois Curfman McInnes 
853ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
8549b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
8559b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
8569b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
8579b94acceSBarry Smith    throughout the global iterations.
8589b94acceSBarry Smith 
8599b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
8609b94acceSBarry Smith 
861ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
8629b94acceSBarry Smith @*/
8639b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
8649b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
8659b94acceSBarry Smith {
86677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
86774679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
86874679c65SBarry Smith     SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only");
8699b94acceSBarry Smith   snes->computejacobian = func;
8709b94acceSBarry Smith   snes->jacP            = ctx;
8719b94acceSBarry Smith   snes->jacobian        = A;
8729b94acceSBarry Smith   snes->jacobian_pre    = B;
8739b94acceSBarry Smith   return 0;
8749b94acceSBarry Smith }
87562fef451SLois Curfman McInnes 
876b4fd4287SBarry Smith /*@
877b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
878b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
879b4fd4287SBarry Smith 
880b4fd4287SBarry Smith    Input Parameter:
881b4fd4287SBarry Smith .  snes - the nonlinear solver context
882b4fd4287SBarry Smith 
883b4fd4287SBarry Smith    Output Parameters:
884b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
885b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
886b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
887b4fd4287SBarry Smith 
888b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
889b4fd4287SBarry Smith @*/
890b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
891b4fd4287SBarry Smith {
89274679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
89374679c65SBarry Smith     SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only");
894b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
895b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
896b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
897b4fd4287SBarry Smith   return 0;
898b4fd4287SBarry Smith }
899b4fd4287SBarry Smith 
9009b94acceSBarry Smith /*@C
9019b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
9029b94acceSBarry Smith    location to store it.
9039b94acceSBarry Smith 
9049b94acceSBarry Smith    Input Parameters:
9059b94acceSBarry Smith .  snes - the SNES context
9069b94acceSBarry Smith .  A - Hessian matrix
9079b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
9089b94acceSBarry Smith .  func - Jacobian evaluation routine
9099b94acceSBarry Smith .  ctx - optional user-defined context for private data for the
9109b94acceSBarry Smith          Hessian evaluation routine (may be null)
9119b94acceSBarry Smith 
9129b94acceSBarry Smith    Calling sequence of func:
9139b94acceSBarry Smith .  func (SNES,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9149b94acceSBarry Smith 
9159b94acceSBarry Smith .  x - input vector
9169b94acceSBarry Smith .  A - Hessian matrix
9179b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
918ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
919ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9209b94acceSBarry Smith .  ctx - optional user-defined Hessian context
9219b94acceSBarry Smith 
9229b94acceSBarry Smith    Notes:
923ac21db08SLois Curfman McInnes    See SLESSetOperators() for information about setting the flag input
924ac21db08SLois Curfman McInnes    parameter in the routine func().  Be sure to read this information!
925ac21db08SLois Curfman McInnes 
9269b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
9279b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
9289b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9299b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9309b94acceSBarry Smith    throughout the global iterations.
9319b94acceSBarry Smith 
9329b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
9339b94acceSBarry Smith 
934ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
9359b94acceSBarry Smith @*/
9369b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9379b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9389b94acceSBarry Smith {
93977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
94074679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
94174679c65SBarry Smith     SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
9429b94acceSBarry Smith   snes->computejacobian = func;
9439b94acceSBarry Smith   snes->jacP            = ctx;
9449b94acceSBarry Smith   snes->jacobian        = A;
9459b94acceSBarry Smith   snes->jacobian_pre    = B;
9469b94acceSBarry Smith   return 0;
9479b94acceSBarry Smith }
9489b94acceSBarry Smith 
94962fef451SLois Curfman McInnes /*@
95062fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
95162fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
95262fef451SLois Curfman McInnes 
95362fef451SLois Curfman McInnes    Input Parameter:
95462fef451SLois Curfman McInnes .  snes - the nonlinear solver context
95562fef451SLois Curfman McInnes 
95662fef451SLois Curfman McInnes    Output Parameters:
95762fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
95862fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
95962fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
96062fef451SLois Curfman McInnes 
96162fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
96262fef451SLois Curfman McInnes @*/
96362fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
96462fef451SLois Curfman McInnes {
96574679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
96674679c65SBarry Smith     SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
96762fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
96862fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
96962fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
97062fef451SLois Curfman McInnes   return 0;
97162fef451SLois Curfman McInnes }
97262fef451SLois Curfman McInnes 
9739b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
9749b94acceSBarry Smith 
9759b94acceSBarry Smith /*@
9769b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
977272ac6f2SLois Curfman McInnes    of a nonlinear solver.
9789b94acceSBarry Smith 
9799b94acceSBarry Smith    Input Parameter:
9809b94acceSBarry Smith .  snes - the SNES context
9818ddd3da0SLois Curfman McInnes .  x - the solution vector
9829b94acceSBarry Smith 
983272ac6f2SLois Curfman McInnes    Notes:
984272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
985272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
986272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
987272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
988272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
989272ac6f2SLois Curfman McInnes 
9909b94acceSBarry Smith .keywords: SNES, nonlinear, setup
9919b94acceSBarry Smith 
9929b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
9939b94acceSBarry Smith @*/
9948ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
9959b94acceSBarry Smith {
996272ac6f2SLois Curfman McInnes   int ierr, flg;
99777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
99877c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
9998ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
10009b94acceSBarry Smith 
1001c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1002c1f60f51SBarry Smith   /*
1003c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1004c1f60f51SBarry Smith       matrix free version but still uses the user provided preconditioner matrix
1005c1f60f51SBarry Smith   */
1006c1f60f51SBarry Smith   if (flg) {
1007c1f60f51SBarry Smith     Mat J;
1008c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1009c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1010c1f60f51SBarry Smith     snes->mfshell = J;
1011c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1012c1f60f51SBarry Smith       snes->jacobian = J;
1013c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1014c1f60f51SBarry Smith     }
1015c1f60f51SBarry Smith     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1016c1f60f51SBarry Smith       snes->jacobian = J;
1017c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1018c1f60f51SBarry Smith     } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free operator option");
1019c1f60f51SBarry Smith   }
1020272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1021c1f60f51SBarry Smith   /*
1022c1f60f51SBarry Smith       This version replaces both the user provided Jacobian and the user
1023c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1024c1f60f51SBarry Smith    */
102531615d04SLois Curfman McInnes   if (flg) {
1026272ac6f2SLois Curfman McInnes     Mat J;
1027272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1028272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1029272ac6f2SLois Curfman McInnes     snes->mfshell = J;
103083e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
103183e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
103294a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
103383e56afdSLois Curfman McInnes     }
103483e56afdSLois Curfman McInnes     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
103583e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
103694a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
103783e56afdSLois Curfman McInnes     } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free option");
1038272ac6f2SLois Curfman McInnes   }
10399b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
104037fcc0dbSBarry Smith     if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first");
104137fcc0dbSBarry Smith     if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first");
104237fcc0dbSBarry Smith     if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetJacobian() first");
1043d804570eSBarry Smith     if (snes->vec_func == snes->vec_sol) SETERRQ(1,"SNESSetUp:Solution vector cannot be function vector");
1044a703fe33SLois Curfman McInnes 
1045a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
104640191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1047a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1048a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1049a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1050a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1051a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1052a703fe33SLois Curfman McInnes     }
10539b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
105437fcc0dbSBarry Smith     if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first");
105537fcc0dbSBarry Smith     if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first");
105637fcc0dbSBarry Smith     if (!snes->computeumfunction)
105737fcc0dbSBarry Smith       SETERRQ(1,"SNESSetUp:Must call SNESSetMinimizationFunction() first");
105837fcc0dbSBarry Smith     if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetHessian() first");
10599b94acceSBarry Smith   } else SETERRQ(1,"SNESSetUp:Unknown method class");
1060a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1061a703fe33SLois Curfman McInnes   snes->setup_called = 1;
1062a703fe33SLois Curfman McInnes   return 0;
10639b94acceSBarry Smith }
10649b94acceSBarry Smith 
10659b94acceSBarry Smith /*@C
10669b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
10679b94acceSBarry Smith    with SNESCreate().
10689b94acceSBarry Smith 
10699b94acceSBarry Smith    Input Parameter:
10709b94acceSBarry Smith .  snes - the SNES context
10719b94acceSBarry Smith 
10729b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
10739b94acceSBarry Smith 
107463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
10759b94acceSBarry Smith @*/
10769b94acceSBarry Smith int SNESDestroy(SNES snes)
10779b94acceSBarry Smith {
10789b94acceSBarry Smith   int ierr;
107977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
10809b94acceSBarry Smith   ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);
10810452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
10829b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
10839b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
10843e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
10856f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
10869b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
10870452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
10889b94acceSBarry Smith   return 0;
10899b94acceSBarry Smith }
10909b94acceSBarry Smith 
10919b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
10929b94acceSBarry Smith 
10939b94acceSBarry Smith /*@
1094d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
10959b94acceSBarry Smith 
10969b94acceSBarry Smith    Input Parameters:
10979b94acceSBarry Smith .  snes - the SNES context
109833174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
109933174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
110033174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
110133174efeSLois Curfman McInnes            of the change in the solution between steps
110233174efeSLois Curfman McInnes .  maxit - maximum number of iterations
110333174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
11049b94acceSBarry Smith 
110533174efeSLois Curfman McInnes    Options Database Keys:
110633174efeSLois Curfman McInnes $    -snes_atol <atol>
110733174efeSLois Curfman McInnes $    -snes_rtol <rtol>
110833174efeSLois Curfman McInnes $    -snes_stol <stol>
110933174efeSLois Curfman McInnes $    -snes_max_it <maxit>
111033174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
11119b94acceSBarry Smith 
1112d7a720efSLois Curfman McInnes    Notes:
11139b94acceSBarry Smith    The default maximum number of iterations is 50.
11149b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
11159b94acceSBarry Smith 
111633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
11179b94acceSBarry Smith 
1118d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
11199b94acceSBarry Smith @*/
1120d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
11219b94acceSBarry Smith {
112277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1123d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1124d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1125d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1126d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1127d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
11289b94acceSBarry Smith   return 0;
11299b94acceSBarry Smith }
11309b94acceSBarry Smith 
11319b94acceSBarry Smith /*@
113233174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
113333174efeSLois Curfman McInnes 
113433174efeSLois Curfman McInnes    Input Parameters:
113533174efeSLois Curfman McInnes .  snes - the SNES context
113633174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
113733174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
113833174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
113933174efeSLois Curfman McInnes            of the change in the solution between steps
114033174efeSLois Curfman McInnes .  maxit - maximum number of iterations
114133174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
114233174efeSLois Curfman McInnes 
114333174efeSLois Curfman McInnes    Notes:
114433174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
114533174efeSLois Curfman McInnes 
114633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
114733174efeSLois Curfman McInnes 
114833174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
114933174efeSLois Curfman McInnes @*/
115033174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
115133174efeSLois Curfman McInnes {
115233174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
115333174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
115433174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
115533174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
115633174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
115733174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
115833174efeSLois Curfman McInnes   return 0;
115933174efeSLois Curfman McInnes }
116033174efeSLois Curfman McInnes 
116133174efeSLois Curfman McInnes /*@
11629b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
11639b94acceSBarry Smith 
11649b94acceSBarry Smith    Input Parameters:
11659b94acceSBarry Smith .  snes - the SNES context
11669b94acceSBarry Smith .  tol - tolerance
11679b94acceSBarry Smith 
11689b94acceSBarry Smith    Options Database Key:
1169d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
11709b94acceSBarry Smith 
11719b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
11729b94acceSBarry Smith 
1173d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
11749b94acceSBarry Smith @*/
11759b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
11769b94acceSBarry Smith {
117777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
11789b94acceSBarry Smith   snes->deltatol = tol;
11799b94acceSBarry Smith   return 0;
11809b94acceSBarry Smith }
11819b94acceSBarry Smith 
11829b94acceSBarry Smith /*@
118377c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
11849b94acceSBarry Smith    for unconstrained minimization solvers.
11859b94acceSBarry Smith 
11869b94acceSBarry Smith    Input Parameters:
11879b94acceSBarry Smith .  snes - the SNES context
11889b94acceSBarry Smith .  ftol - minimum function tolerance
11899b94acceSBarry Smith 
11909b94acceSBarry Smith    Options Database Key:
1191d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
11929b94acceSBarry Smith 
11939b94acceSBarry Smith    Note:
119477c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
11959b94acceSBarry Smith    methods only.
11969b94acceSBarry Smith 
11979b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
11989b94acceSBarry Smith 
1199d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
12009b94acceSBarry Smith @*/
120177c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
12029b94acceSBarry Smith {
120377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12049b94acceSBarry Smith   snes->fmin = ftol;
12059b94acceSBarry Smith   return 0;
12069b94acceSBarry Smith }
12079b94acceSBarry Smith 
12089b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
12099b94acceSBarry Smith 
12109b94acceSBarry Smith /*@C
12119b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
12129b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
12139b94acceSBarry Smith    progress.
12149b94acceSBarry Smith 
12159b94acceSBarry Smith    Input Parameters:
12169b94acceSBarry Smith .  snes - the SNES context
12179b94acceSBarry Smith .  func - monitoring routine
12189b94acceSBarry Smith .  mctx - optional user-defined context for private data for the
12199b94acceSBarry Smith           monitor routine (may be null)
12209b94acceSBarry Smith 
12219b94acceSBarry Smith    Calling sequence of func:
1222682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
12239b94acceSBarry Smith 
12249b94acceSBarry Smith $    snes - the SNES context
12259b94acceSBarry Smith $    its - iteration number
12269b94acceSBarry Smith $    mctx - optional monitoring context
12279b94acceSBarry Smith $
12289b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
12299b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
12309b94acceSBarry Smith $
12319b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
12329b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
12339b94acceSBarry Smith 
12349b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
12359b94acceSBarry Smith 
12369b94acceSBarry Smith .seealso: SNESDefaultMonitor()
12379b94acceSBarry Smith @*/
123874679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
12399b94acceSBarry Smith {
12409b94acceSBarry Smith   snes->monitor = func;
124174679c65SBarry Smith   snes->monP    = mctx;
12429b94acceSBarry Smith   return 0;
12439b94acceSBarry Smith }
12449b94acceSBarry Smith 
12459b94acceSBarry Smith /*@C
12469b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
12479b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
12489b94acceSBarry Smith 
12499b94acceSBarry Smith    Input Parameters:
12509b94acceSBarry Smith .  snes - the SNES context
12519b94acceSBarry Smith .  func - routine to test for convergence
12529b94acceSBarry Smith .  cctx - optional context for private data for the convergence routine
12539b94acceSBarry Smith           (may be null)
12549b94acceSBarry Smith 
12559b94acceSBarry Smith    Calling sequence of func:
12569b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
12579b94acceSBarry Smith              double f,void *cctx)
12589b94acceSBarry Smith 
12599b94acceSBarry Smith $    snes - the SNES context
12609b94acceSBarry Smith $    cctx - optional convergence context
12619b94acceSBarry Smith $    xnorm - 2-norm of current iterate
12629b94acceSBarry Smith $
12639b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
12649b94acceSBarry Smith $    gnorm - 2-norm of current step
12659b94acceSBarry Smith $    f - 2-norm of function
12669b94acceSBarry Smith $
12679b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
12689b94acceSBarry Smith $    gnorm - 2-norm of current gradient
12699b94acceSBarry Smith $    f - function value
12709b94acceSBarry Smith 
12719b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
12729b94acceSBarry Smith 
127340191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
127440191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
12759b94acceSBarry Smith @*/
127674679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
12779b94acceSBarry Smith {
12789b94acceSBarry Smith   (snes)->converged = func;
12799b94acceSBarry Smith   (snes)->cnvP      = cctx;
12809b94acceSBarry Smith   return 0;
12819b94acceSBarry Smith }
12829b94acceSBarry Smith 
12839b94acceSBarry Smith /*
12849b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
12859b94acceSBarry Smith    positive parameter delta.
12869b94acceSBarry Smith 
12879b94acceSBarry Smith     Input Parameters:
12889b94acceSBarry Smith .   snes - the SNES context
12899b94acceSBarry Smith .   y - approximate solution of linear system
12909b94acceSBarry Smith .   fnorm - 2-norm of current function
12919b94acceSBarry Smith .   delta - trust region size
12929b94acceSBarry Smith 
12939b94acceSBarry Smith     Output Parameters:
12949b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
12959b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
12969b94acceSBarry Smith     region, and exceeds zero otherwise.
12979b94acceSBarry Smith .   ynorm - 2-norm of the step
12989b94acceSBarry Smith 
12999b94acceSBarry Smith     Note:
130040191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
13019b94acceSBarry Smith     is set to be the maximum allowable step size.
13029b94acceSBarry Smith 
13039b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
13049b94acceSBarry Smith */
13059b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
13069b94acceSBarry Smith                           double *gpnorm,double *ynorm)
13079b94acceSBarry Smith {
13089b94acceSBarry Smith   double norm;
13099b94acceSBarry Smith   Scalar cnorm;
1310cddf8d76SBarry Smith   VecNorm(y,NORM_2, &norm );
13119b94acceSBarry Smith   if (norm > *delta) {
13129b94acceSBarry Smith      norm = *delta/norm;
13139b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
13149b94acceSBarry Smith      cnorm = norm;
13159b94acceSBarry Smith      VecScale( &cnorm, y );
13169b94acceSBarry Smith      *ynorm = *delta;
13179b94acceSBarry Smith   } else {
13189b94acceSBarry Smith      *gpnorm = 0.0;
13199b94acceSBarry Smith      *ynorm = norm;
13209b94acceSBarry Smith   }
13219b94acceSBarry Smith   return 0;
13229b94acceSBarry Smith }
13239b94acceSBarry Smith 
13249b94acceSBarry Smith /*@
13259b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
132663a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
13279b94acceSBarry Smith 
13289b94acceSBarry Smith    Input Parameter:
13299b94acceSBarry Smith .  snes - the SNES context
13308ddd3da0SLois Curfman McInnes .  x - the solution vector
13319b94acceSBarry Smith 
13329b94acceSBarry Smith    Output Parameter:
13339b94acceSBarry Smith    its - number of iterations until termination
13349b94acceSBarry Smith 
13358ddd3da0SLois Curfman McInnes    Note:
13368ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
13378ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
13388ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
13398ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
13408ddd3da0SLois Curfman McInnes 
13419b94acceSBarry Smith .keywords: SNES, nonlinear, solve
13429b94acceSBarry Smith 
134363a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
13449b94acceSBarry Smith @*/
13458ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
13469b94acceSBarry Smith {
13473c7409f5SSatish Balay   int ierr, flg;
1348052efed2SBarry Smith 
134977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
135074679c65SBarry Smith   PetscValidIntPointer(its);
1351c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1352c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
13539b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
13549b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
13559b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
13563c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
13576d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
13589b94acceSBarry Smith   return 0;
13599b94acceSBarry Smith }
13609b94acceSBarry Smith 
13619b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
136237fcc0dbSBarry Smith static NRList *__SNESList = 0;
13639b94acceSBarry Smith 
13649b94acceSBarry Smith /*@
13654b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
13669b94acceSBarry Smith 
13679b94acceSBarry Smith    Input Parameters:
13689b94acceSBarry Smith .  snes - the SNES context
13699b94acceSBarry Smith .  method - a known method
13709b94acceSBarry Smith 
13719b94acceSBarry Smith    Notes:
13729b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
13739b94acceSBarry Smith $  Systems of nonlinear equations:
137440191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
137540191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
13769b94acceSBarry Smith $  Unconstrained minimization:
137740191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
137840191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
13799b94acceSBarry Smith 
13809b94acceSBarry Smith   Options Database Command:
13814b0e389bSBarry Smith $ -snes_type  <method>
13829b94acceSBarry Smith $    Use -help for a list of available methods
13839b94acceSBarry Smith $    (for instance, ls or tr)
1384a703fe33SLois Curfman McInnes 
1385a703fe33SLois Curfman McInnes .keysords: SNES, set, method
13869b94acceSBarry Smith @*/
13874b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
13889b94acceSBarry Smith {
13899b94acceSBarry Smith   int (*r)(SNES);
139077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13917d1a2b2bSBarry Smith   if (snes->type == (int) method) return 0;
13927d1a2b2bSBarry Smith 
13939b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
139437fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
139537fcc0dbSBarry Smith   if (!__SNESList) {SETERRQ(1,"SNESSetType:Could not get methods");}
139637fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
13974b0e389bSBarry Smith   if (!r) {SETERRQ(1,"SNESSetType:Unknown method");}
13980452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
13999b94acceSBarry Smith   snes->set_method_called = 1;
14009b94acceSBarry Smith   return (*r)(snes);
14019b94acceSBarry Smith }
14029b94acceSBarry Smith 
14039b94acceSBarry Smith /* --------------------------------------------------------------------- */
14049b94acceSBarry Smith /*@C
14059b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
14064b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
14079b94acceSBarry Smith 
14089b94acceSBarry Smith    Input Parameters:
140940191667SLois Curfman McInnes .  name - for instance SNES_EQ_LS, SNES_EQ_TR, ...
141040191667SLois Curfman McInnes .  sname - corresponding string for name
14119b94acceSBarry Smith .  create - routine to create method context
14129b94acceSBarry Smith 
14139b94acceSBarry Smith .keywords: SNES, nonlinear, register
14149b94acceSBarry Smith 
14159b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
14169b94acceSBarry Smith @*/
14179b94acceSBarry Smith int SNESRegister(int name, char *sname, int (*create)(SNES))
14189b94acceSBarry Smith {
14199b94acceSBarry Smith   int ierr;
142037fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
142137fcc0dbSBarry Smith   NRRegister( __SNESList, name, sname, (int (*)(void*))create );
14229b94acceSBarry Smith   return 0;
14239b94acceSBarry Smith }
14249b94acceSBarry Smith /* --------------------------------------------------------------------- */
14259b94acceSBarry Smith /*@C
14269b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
14279b94acceSBarry Smith    registered by SNESRegister().
14289b94acceSBarry Smith 
14299b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
14309b94acceSBarry Smith 
14319b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
14329b94acceSBarry Smith @*/
14339b94acceSBarry Smith int SNESRegisterDestroy()
14349b94acceSBarry Smith {
143537fcc0dbSBarry Smith   if (__SNESList) {
143637fcc0dbSBarry Smith     NRDestroy( __SNESList );
143737fcc0dbSBarry Smith     __SNESList = 0;
14389b94acceSBarry Smith   }
14399b94acceSBarry Smith   return 0;
14409b94acceSBarry Smith }
14419b94acceSBarry Smith 
14429b94acceSBarry Smith /*
14434b0e389bSBarry Smith    SNESGetTypeFromOptions_Private - Sets the selected method from the
14449b94acceSBarry Smith    options database.
14459b94acceSBarry Smith 
14469b94acceSBarry Smith    Input Parameter:
14479b94acceSBarry Smith .  ctx - the SNES context
14489b94acceSBarry Smith 
14499b94acceSBarry Smith    Output Parameter:
14509b94acceSBarry Smith .  method -  solver method
14519b94acceSBarry Smith 
14529b94acceSBarry Smith    Returns:
14539b94acceSBarry Smith    Returns 1 if the method is found; 0 otherwise.
14549b94acceSBarry Smith 
14559b94acceSBarry Smith    Options Database Key:
14564b0e389bSBarry Smith $  -snes_type  method
14579b94acceSBarry Smith */
1458052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
14599b94acceSBarry Smith {
1460052efed2SBarry Smith   int  ierr;
14619b94acceSBarry Smith   char sbuf[50];
14625baf8537SBarry Smith 
1463052efed2SBarry Smith   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1464052efed2SBarry Smith   if (*flg) {
1465052efed2SBarry Smith     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
146637fcc0dbSBarry Smith     *method = (SNESType)NRFindID( __SNESList, sbuf );
14679b94acceSBarry Smith   }
14689b94acceSBarry Smith   return 0;
14699b94acceSBarry Smith }
14709b94acceSBarry Smith 
14719b94acceSBarry Smith /*@C
14729a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
14739b94acceSBarry Smith 
14749b94acceSBarry Smith    Input Parameter:
14754b0e389bSBarry Smith .  snes - nonlinear solver context
14769b94acceSBarry Smith 
14779b94acceSBarry Smith    Output Parameter:
14789a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
14799a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
14809b94acceSBarry Smith 
14819b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
14829b94acceSBarry Smith @*/
14834b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
14849b94acceSBarry Smith {
148506a9b0adSLois Curfman McInnes   int ierr;
148637fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
14874b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
148837fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
14899b94acceSBarry Smith   return 0;
14909b94acceSBarry Smith }
14919b94acceSBarry Smith 
14929b94acceSBarry Smith #include <stdio.h>
14939b94acceSBarry Smith /*
14944b0e389bSBarry Smith    SNESPrintTypes_Private - Prints the SNES methods available from the
14959b94acceSBarry Smith    options database.
14969b94acceSBarry Smith 
14979b94acceSBarry Smith    Input Parameters:
149833455573SSatish Balay .  comm   - communicator (usually MPI_COMM_WORLD)
14999b94acceSBarry Smith .  prefix - prefix (usually "-")
15004b0e389bSBarry Smith .  name   - the options database name (by default "snes_type")
15019b94acceSBarry Smith */
150233455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
15039b94acceSBarry Smith {
15049b94acceSBarry Smith   FuncList *entry;
150537fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
150637fcc0dbSBarry Smith   entry = __SNESList->head;
150777c4ece6SBarry Smith   PetscPrintf(comm," %s%s (one of)",prefix,name);
15089b94acceSBarry Smith   while (entry) {
150977c4ece6SBarry Smith     PetscPrintf(comm," %s",entry->name);
15109b94acceSBarry Smith     entry = entry->next;
15119b94acceSBarry Smith   }
151277c4ece6SBarry Smith   PetscPrintf(comm,"\n");
15139b94acceSBarry Smith   return 0;
15149b94acceSBarry Smith }
15159b94acceSBarry Smith 
15169b94acceSBarry Smith /*@C
15179b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
15189b94acceSBarry Smith    stored.
15199b94acceSBarry Smith 
15209b94acceSBarry Smith    Input Parameter:
15219b94acceSBarry Smith .  snes - the SNES context
15229b94acceSBarry Smith 
15239b94acceSBarry Smith    Output Parameter:
15249b94acceSBarry Smith .  x - the solution
15259b94acceSBarry Smith 
15269b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
15279b94acceSBarry Smith 
1528a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
15299b94acceSBarry Smith @*/
15309b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
15319b94acceSBarry Smith {
153277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15339b94acceSBarry Smith   *x = snes->vec_sol_always;
15349b94acceSBarry Smith   return 0;
15359b94acceSBarry Smith }
15369b94acceSBarry Smith 
15379b94acceSBarry Smith /*@C
15389b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
15399b94acceSBarry Smith    stored.
15409b94acceSBarry Smith 
15419b94acceSBarry Smith    Input Parameter:
15429b94acceSBarry Smith .  snes - the SNES context
15439b94acceSBarry Smith 
15449b94acceSBarry Smith    Output Parameter:
15459b94acceSBarry Smith .  x - the solution update
15469b94acceSBarry Smith 
15479b94acceSBarry Smith    Notes:
15489b94acceSBarry Smith    This vector is implementation dependent.
15499b94acceSBarry Smith 
15509b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
15519b94acceSBarry Smith 
15529b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
15539b94acceSBarry Smith @*/
15549b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
15559b94acceSBarry Smith {
155677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15579b94acceSBarry Smith   *x = snes->vec_sol_update_always;
15589b94acceSBarry Smith   return 0;
15599b94acceSBarry Smith }
15609b94acceSBarry Smith 
15619b94acceSBarry Smith /*@C
15623638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
15639b94acceSBarry Smith 
15649b94acceSBarry Smith    Input Parameter:
15659b94acceSBarry Smith .  snes - the SNES context
15669b94acceSBarry Smith 
15679b94acceSBarry Smith    Output Parameter:
15683638b69dSLois Curfman McInnes .  r - the function
15699b94acceSBarry Smith 
15709b94acceSBarry Smith    Notes:
15719b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
15729b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
15739b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
15749b94acceSBarry Smith 
1575a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
15769b94acceSBarry Smith 
157731615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
157831615d04SLois Curfman McInnes           SNESGetGradient()
15799b94acceSBarry Smith @*/
15809b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
15819b94acceSBarry Smith {
158277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15839b94acceSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
158448d91487SBarry Smith     "SNESGetFunction:For SNES_NONLINEAR_EQUATIONS only");
15859b94acceSBarry Smith   *r = snes->vec_func_always;
15869b94acceSBarry Smith   return 0;
15879b94acceSBarry Smith }
15889b94acceSBarry Smith 
15899b94acceSBarry Smith /*@C
15903638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
15919b94acceSBarry Smith 
15929b94acceSBarry Smith    Input Parameter:
15939b94acceSBarry Smith .  snes - the SNES context
15949b94acceSBarry Smith 
15959b94acceSBarry Smith    Output Parameter:
15969b94acceSBarry Smith .  r - the gradient
15979b94acceSBarry Smith 
15989b94acceSBarry Smith    Notes:
15999b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
16009b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
16019b94acceSBarry Smith    SNESGetFunction().
16029b94acceSBarry Smith 
16039b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
16049b94acceSBarry Smith 
160531615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
16069b94acceSBarry Smith @*/
16079b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
16089b94acceSBarry Smith {
160977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16109b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
161148d91487SBarry Smith     "SNESGetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
16129b94acceSBarry Smith   *r = snes->vec_func_always;
16139b94acceSBarry Smith   return 0;
16149b94acceSBarry Smith }
16159b94acceSBarry Smith 
16169b94acceSBarry Smith /*@
16179b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
16189b94acceSBarry Smith    unconstrained minimization problems.
16199b94acceSBarry Smith 
16209b94acceSBarry Smith    Input Parameter:
16219b94acceSBarry Smith .  snes - the SNES context
16229b94acceSBarry Smith 
16239b94acceSBarry Smith    Output Parameter:
16249b94acceSBarry Smith .  r - the function
16259b94acceSBarry Smith 
16269b94acceSBarry Smith    Notes:
16279b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
16289b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
16299b94acceSBarry Smith    SNESGetFunction().
16309b94acceSBarry Smith 
16319b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
16329b94acceSBarry Smith 
163331615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
16349b94acceSBarry Smith @*/
16359b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
16369b94acceSBarry Smith {
163777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
163874679c65SBarry Smith   PetscValidScalarPointer(r);
16399b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
164048d91487SBarry Smith     "SNESGetMinimizationFunction:For SNES_UNCONSTRAINED_MINIMIZATION only");
16419b94acceSBarry Smith   *r = snes->fc;
16429b94acceSBarry Smith   return 0;
16439b94acceSBarry Smith }
16449b94acceSBarry Smith 
16453c7409f5SSatish Balay /*@C
16463c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
16473c7409f5SSatish Balay    SNES options in the database.
16483c7409f5SSatish Balay 
16493c7409f5SSatish Balay    Input Parameter:
16503c7409f5SSatish Balay .  snes - the SNES context
16513c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
16523c7409f5SSatish Balay 
16533c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1654a86d99e1SLois Curfman McInnes 
1655a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
16563c7409f5SSatish Balay @*/
16573c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
16583c7409f5SSatish Balay {
16593c7409f5SSatish Balay   int ierr;
16603c7409f5SSatish Balay 
166177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16623c7409f5SSatish Balay   ierr = PetscObjectSetPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
16633c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
16643c7409f5SSatish Balay   return 0;
16653c7409f5SSatish Balay }
16663c7409f5SSatish Balay 
16673c7409f5SSatish Balay /*@C
16683c7409f5SSatish Balay    SNESAppendOptionsPrefix - Append to the prefix used for searching for all
16693c7409f5SSatish Balay    SNES options in the database.
16703c7409f5SSatish Balay 
16713c7409f5SSatish Balay    Input Parameter:
16723c7409f5SSatish Balay .  snes - the SNES context
16733c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
16743c7409f5SSatish Balay 
16753c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1676a86d99e1SLois Curfman McInnes 
1677a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
16783c7409f5SSatish Balay @*/
16793c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
16803c7409f5SSatish Balay {
16813c7409f5SSatish Balay   int ierr;
16823c7409f5SSatish Balay 
168377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16843c7409f5SSatish Balay   ierr = PetscObjectAppendPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
16853c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
16863c7409f5SSatish Balay   return 0;
16873c7409f5SSatish Balay }
16883c7409f5SSatish Balay 
1689c83590e2SSatish Balay /*@
16903c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
16913c7409f5SSatish Balay    SNES options in the database.
16923c7409f5SSatish Balay 
16933c7409f5SSatish Balay    Input Parameter:
16943c7409f5SSatish Balay .  snes - the SNES context
16953c7409f5SSatish Balay 
16963c7409f5SSatish Balay    Output Parameter:
16973c7409f5SSatish Balay .  prefix - pointer to the prefix string used
16983c7409f5SSatish Balay 
16993c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1700a86d99e1SLois Curfman McInnes 
1701a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
17023c7409f5SSatish Balay @*/
17033c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
17043c7409f5SSatish Balay {
17053c7409f5SSatish Balay   int ierr;
17063c7409f5SSatish Balay 
170777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17083c7409f5SSatish Balay   ierr = PetscObjectGetPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
17093c7409f5SSatish Balay   return 0;
17103c7409f5SSatish Balay }
17113c7409f5SSatish Balay 
17123c7409f5SSatish Balay 
17139b94acceSBarry Smith 
17149b94acceSBarry Smith 
17159b94acceSBarry Smith 
1716