xref: /petsc/src/snes/interface/snes.c (revision ae3c334c5d8e947687ef6a015c77d2112d2bf88a)
19b94acceSBarry Smith #ifndef lint
2*ae3c334cSLois Curfman McInnes static char vcid[] = "$Id: snes.c,v 1.106 1997/01/14 22:57:02 curfman 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 
145615d1e5SSatish Balay #undef __FUNC__
155615d1e5SSatish Balay #define __FUNC__ "SNESView"
169b94acceSBarry Smith /*@
179b94acceSBarry Smith    SNESView - Prints the SNES data structure.
189b94acceSBarry Smith 
199b94acceSBarry Smith    Input Parameters:
209b94acceSBarry Smith .  SNES - the SNES context
219b94acceSBarry Smith .  viewer - visualization context
229b94acceSBarry Smith 
239b94acceSBarry Smith    Options Database Key:
249b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
259b94acceSBarry Smith 
269b94acceSBarry Smith    Notes:
279b94acceSBarry Smith    The available visualization contexts include
286d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
296d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
309b94acceSBarry Smith $       output where only the first processor opens
319b94acceSBarry Smith $       the file.  All other processors send their
329b94acceSBarry Smith $       data to the first processor to print.
339b94acceSBarry Smith 
349b94acceSBarry Smith    The user can open alternative vistualization contexts with
35dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
369b94acceSBarry Smith 
379b94acceSBarry Smith .keywords: SNES, view
389b94acceSBarry Smith 
39dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
409b94acceSBarry Smith @*/
419b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
429b94acceSBarry Smith {
439b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
449b94acceSBarry Smith   FILE                *fd;
459b94acceSBarry Smith   int                 ierr;
469b94acceSBarry Smith   SLES                sles;
479b94acceSBarry Smith   char                *method;
4819bcc07fSBarry Smith   ViewerType          vtype;
499b94acceSBarry Smith 
5074679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5174679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5274679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5374679c65SBarry Smith 
5419bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5519bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5690ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5777c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
584b0e389bSBarry Smith     SNESGetType(snes,PETSC_NULL,&method);
5977c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
609b94acceSBarry Smith     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
6177c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
629b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
639b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6477c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
659b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
669b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
677a00f4a9SLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
687a00f4a9SLois Curfman McInnes     "  total number of linear solver iterations=%d\n",snes->linear_its);
69*ae3c334cSLois Curfman McInnes      PetscFPrintf(snes->comm,fd,
70*ae3c334cSLois Curfman McInnes      "  total number of function evaluation=%d\n",snes->nfuncs);
719b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
7277c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
739b94acceSBarry Smith     if (snes->ksp_ewconv) {
749b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
759b94acceSBarry Smith       if (kctx) {
7677c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
779b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
789b94acceSBarry Smith         kctx->version);
7977c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
809b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
819b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
8277c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
839b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
849b94acceSBarry Smith       }
859b94acceSBarry Smith     }
86c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
87c30f285eSLois Curfman McInnes     SNESGetType(snes,PETSC_NULL,&method);
88c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
8919bcc07fSBarry Smith   }
909b94acceSBarry Smith   SNESGetSLES(snes,&sles);
919b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
929b94acceSBarry Smith   return 0;
939b94acceSBarry Smith }
949b94acceSBarry Smith 
95639f9d9dSBarry Smith /*
96639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
97639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
98639f9d9dSBarry Smith */
99639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
100639f9d9dSBarry Smith static int numberofsetfromoptions;
101639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
102639f9d9dSBarry Smith 
1035615d1e5SSatish Balay #undef __FUNC__
1045615d1e5SSatish Balay #define __FUNC__ "SNESAddOptionsChecker"
105639f9d9dSBarry Smith /*@
106639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
107639f9d9dSBarry Smith 
108639f9d9dSBarry Smith   Input Parameter:
109639f9d9dSBarry Smith .   snescheck - function that checks for options
110639f9d9dSBarry Smith 
111639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
112639f9d9dSBarry Smith @*/
113639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
114639f9d9dSBarry Smith {
115639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
116e3372554SBarry Smith     SETERRQ(1,0,"Too many options checkers, only 5 allowed");
117639f9d9dSBarry Smith   }
118639f9d9dSBarry Smith 
119639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
120639f9d9dSBarry Smith   return 0;
121639f9d9dSBarry Smith }
122639f9d9dSBarry Smith 
1235615d1e5SSatish Balay #undef __FUNC__
1245615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1259b94acceSBarry Smith /*@
126682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1279b94acceSBarry Smith 
1289b94acceSBarry Smith    Input Parameter:
1299b94acceSBarry Smith .  snes - the SNES context
1309b94acceSBarry Smith 
1319b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1329b94acceSBarry Smith 
133a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1349b94acceSBarry Smith @*/
1359b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1369b94acceSBarry Smith {
1374b0e389bSBarry Smith   SNESType method;
1389b94acceSBarry Smith   double   tmp;
1399b94acceSBarry Smith   SLES     sles;
14081f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1413c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1429b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1439b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1449b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1459b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1469b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1479b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1489b94acceSBarry Smith 
14981f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
15081f57ec7SBarry Smith 
15181f57ec7SBarry Smith 
15277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
153e3372554SBarry Smith   if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp");
154052efed2SBarry Smith   ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr);
155052efed2SBarry Smith   if (flg) {
156052efed2SBarry Smith     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
1579b94acceSBarry Smith   }
1585334005bSBarry Smith   else if (!snes->set_method_called) {
1595334005bSBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
16040191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
1615334005bSBarry Smith     }
1625334005bSBarry Smith     else {
16340191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
1645334005bSBarry Smith     }
1655334005bSBarry Smith   }
1665334005bSBarry Smith 
1673c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1683c7409f5SSatish Balay   if (flg) { SNESPrintHelp(snes); }
1693c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
170d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);}
1713c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
172d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);}
1733c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
174d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); }
1753c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1763c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
177d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
178d7a720efSLois Curfman McInnes   if (flg) { SNESSetTrustRegionTolerance(snes,tmp); }
179d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
180d7a720efSLois Curfman McInnes   if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); }
1813c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1823c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
183b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
184b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
185b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
186b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
187b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
188b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
189b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
1909b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
1919b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
1925bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
1935bbad29bSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);}
1943c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
195639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
1963c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
197d31a3109SLois Curfman McInnes   if (flg) {
198639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
199d31a3109SLois Curfman McInnes   }
20081f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
2013c7409f5SSatish Balay   if (flg) {
20217699dbbSLois Curfman McInnes     int    rank = 0;
203d7e8b826SBarry Smith     DrawLG lg;
20417699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
20517699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
20617699dbbSLois Curfman McInnes     if (!rank) {
20781f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2089b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
209f8c826e1SBarry Smith       snes->xmonitor = lg;
2109b94acceSBarry Smith       PLogObjectParent(snes,lg);
2119b94acceSBarry Smith     }
2129b94acceSBarry Smith   }
2133c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2143c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2159b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2169b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
21794a424c1SBarry Smith     PLogInfo(snes,
21831615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
21931615d04SLois Curfman McInnes   }
22031615d04SLois Curfman McInnes   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
22131615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
22231615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
22394a424c1SBarry Smith     PLogInfo(snes,
22431615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2259b94acceSBarry Smith   }
226639f9d9dSBarry Smith 
227639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
228639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
229639f9d9dSBarry Smith   }
230639f9d9dSBarry Smith 
2319b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2329b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
2339b94acceSBarry Smith   if (!snes->setfromoptions) return 0;
2349b94acceSBarry Smith   return (*snes->setfromoptions)(snes);
2359b94acceSBarry Smith }
2369b94acceSBarry Smith 
2375615d1e5SSatish Balay #undef __FUNC__
2385615d1e5SSatish Balay #define __FUNC__ "SNESPrintHelp"
2399b94acceSBarry Smith /*@
2409b94acceSBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2419b94acceSBarry Smith 
2429b94acceSBarry Smith    Input Parameter:
2439b94acceSBarry Smith .  snes - the SNES context
2449b94acceSBarry Smith 
245a703fe33SLois Curfman McInnes    Options Database Keys:
246a703fe33SLois Curfman McInnes $  -help, -h
247a703fe33SLois Curfman McInnes 
2489b94acceSBarry Smith .keywords: SNES, nonlinear, help
2499b94acceSBarry Smith 
250682d7d0cSBarry Smith .seealso: SNESSetFromOptions()
2519b94acceSBarry Smith @*/
2529b94acceSBarry Smith int SNESPrintHelp(SNES snes)
2539b94acceSBarry Smith {
2546daaf66cSBarry Smith   char                p[64];
2556daaf66cSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2566daaf66cSBarry Smith 
25777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2586daaf66cSBarry Smith 
2596daaf66cSBarry Smith   PetscStrcpy(p,"-");
2606daaf66cSBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2616daaf66cSBarry Smith 
2626daaf66cSBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2636daaf66cSBarry Smith 
264d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
26533455573SSatish Balay   SNESPrintTypes_Private(snes->comm,p,"snes_type");
26677c4ece6SBarry Smith   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
267b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
268b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
269b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
270b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
271b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
272b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
2735bbad29bSBarry Smith   PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n");
2745bbad29bSBarry Smith   PetscPrintf(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
275d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor\n",p);
276d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
277b18e04deSLois Curfman McInnes   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
27877c4ece6SBarry Smith     PetscPrintf(snes->comm,
279d31a3109SLois Curfman McInnes      " Options for solving systems of nonlinear equations only:\n");
28077c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
28177c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
282ef1dfb11SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
2831650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p);
2841650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p);
28577c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
28677c4ece6SBarry Smith     PetscPrintf(snes->comm,
287b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
28877c4ece6SBarry Smith     PetscPrintf(snes->comm,
289b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
29077c4ece6SBarry Smith     PetscPrintf(snes->comm,
291b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
29277c4ece6SBarry Smith     PetscPrintf(snes->comm,
293b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
29477c4ece6SBarry Smith     PetscPrintf(snes->comm,
295b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
29677c4ece6SBarry Smith     PetscPrintf(snes->comm,
297b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
29877c4ece6SBarry Smith     PetscPrintf(snes->comm,
299b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
300b18e04deSLois Curfman McInnes   }
301b18e04deSLois Curfman McInnes   else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
30277c4ece6SBarry Smith     PetscPrintf(snes->comm,
303d31a3109SLois Curfman McInnes      " Options for solving unconstrained minimization problems only:\n");
304b18e04deSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
30577c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
30677c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
307d31a3109SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p);
308d31a3109SLois Curfman McInnes      PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p);
309b18e04deSLois Curfman McInnes   }
3104537a946SLois Curfman McInnes   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
31177c4ece6SBarry Smith   PetscPrintf(snes->comm,"a particular method\n");
3126daaf66cSBarry Smith   if (snes->printhelp) (*snes->printhelp)(snes,p);
3139b94acceSBarry Smith   return 0;
3149b94acceSBarry Smith }
315a847f771SSatish Balay 
3165615d1e5SSatish Balay #undef __FUNC__
3175615d1e5SSatish Balay #define __FUNC__ "SNESSetApplicationContext"
3189b94acceSBarry Smith /*@
3199b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3209b94acceSBarry Smith    the nonlinear solvers.
3219b94acceSBarry Smith 
3229b94acceSBarry Smith    Input Parameters:
3239b94acceSBarry Smith .  snes - the SNES context
3249b94acceSBarry Smith .  usrP - optional user context
3259b94acceSBarry Smith 
3269b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3279b94acceSBarry Smith 
3289b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3299b94acceSBarry Smith @*/
3309b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3319b94acceSBarry Smith {
33277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3339b94acceSBarry Smith   snes->user		= usrP;
3349b94acceSBarry Smith   return 0;
3359b94acceSBarry Smith }
33674679c65SBarry Smith 
3375615d1e5SSatish Balay #undef __FUNC__
3385615d1e5SSatish Balay #define __FUNC__ "SNESGetApplicationContext"
3399b94acceSBarry Smith /*@C
3409b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3419b94acceSBarry Smith    nonlinear solvers.
3429b94acceSBarry Smith 
3439b94acceSBarry Smith    Input Parameter:
3449b94acceSBarry Smith .  snes - SNES context
3459b94acceSBarry Smith 
3469b94acceSBarry Smith    Output Parameter:
3479b94acceSBarry Smith .  usrP - user context
3489b94acceSBarry Smith 
3499b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3509b94acceSBarry Smith 
3519b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3529b94acceSBarry Smith @*/
3539b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3549b94acceSBarry Smith {
35577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3569b94acceSBarry Smith   *usrP = snes->user;
3579b94acceSBarry Smith   return 0;
3589b94acceSBarry Smith }
35974679c65SBarry Smith 
3605615d1e5SSatish Balay #undef __FUNC__
3615615d1e5SSatish Balay #define __FUNC__ "SNESGetIterationNumber"
3629b94acceSBarry Smith /*@
3639b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3649b94acceSBarry Smith    nonlinear solver.
3659b94acceSBarry Smith 
3669b94acceSBarry Smith    Input Parameter:
3679b94acceSBarry Smith .  snes - SNES context
3689b94acceSBarry Smith 
3699b94acceSBarry Smith    Output Parameter:
3709b94acceSBarry Smith .  iter - iteration number
3719b94acceSBarry Smith 
3729b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3739b94acceSBarry Smith @*/
3749b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3759b94acceSBarry Smith {
37677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
37774679c65SBarry Smith   PetscValidIntPointer(iter);
3789b94acceSBarry Smith   *iter = snes->iter;
3799b94acceSBarry Smith   return 0;
3809b94acceSBarry Smith }
38174679c65SBarry Smith 
3825615d1e5SSatish Balay #undef __FUNC__
3835615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
3849b94acceSBarry Smith /*@
3859b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3869b94acceSBarry Smith    with SNESSSetFunction().
3879b94acceSBarry Smith 
3889b94acceSBarry Smith    Input Parameter:
3899b94acceSBarry Smith .  snes - SNES context
3909b94acceSBarry Smith 
3919b94acceSBarry Smith    Output Parameter:
3929b94acceSBarry Smith .  fnorm - 2-norm of function
3939b94acceSBarry Smith 
3949b94acceSBarry Smith    Note:
3959b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
396a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
397a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3989b94acceSBarry Smith 
3999b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
400a86d99e1SLois Curfman McInnes 
401a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
4029b94acceSBarry Smith @*/
4039b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
4049b94acceSBarry Smith {
40577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
40674679c65SBarry Smith   PetscValidScalarPointer(fnorm);
40774679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
408e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
40974679c65SBarry Smith   }
4109b94acceSBarry Smith   *fnorm = snes->norm;
4119b94acceSBarry Smith   return 0;
4129b94acceSBarry Smith }
41374679c65SBarry Smith 
4145615d1e5SSatish Balay #undef __FUNC__
4155615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
4169b94acceSBarry Smith /*@
4179b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
4189b94acceSBarry Smith    with SNESSSetGradient().
4199b94acceSBarry Smith 
4209b94acceSBarry Smith    Input Parameter:
4219b94acceSBarry Smith .  snes - SNES context
4229b94acceSBarry Smith 
4239b94acceSBarry Smith    Output Parameter:
4249b94acceSBarry Smith .  fnorm - 2-norm of gradient
4259b94acceSBarry Smith 
4269b94acceSBarry Smith    Note:
4279b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
428a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
429a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4309b94acceSBarry Smith 
4319b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
432a86d99e1SLois Curfman McInnes 
433a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4349b94acceSBarry Smith @*/
4359b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
4369b94acceSBarry Smith {
43777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43874679c65SBarry Smith   PetscValidScalarPointer(gnorm);
43974679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
440e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
44174679c65SBarry Smith   }
4429b94acceSBarry Smith   *gnorm = snes->norm;
4439b94acceSBarry Smith   return 0;
4449b94acceSBarry Smith }
44574679c65SBarry Smith 
4465615d1e5SSatish Balay #undef __FUNC__
4475615d1e5SSatish Balay #define __FUNC__ "SNESGetNumberUnsuccessfulSteps"
4489b94acceSBarry Smith /*@
4499b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4509b94acceSBarry Smith    attempted by the nonlinear solver.
4519b94acceSBarry Smith 
4529b94acceSBarry Smith    Input Parameter:
4539b94acceSBarry Smith .  snes - SNES context
4549b94acceSBarry Smith 
4559b94acceSBarry Smith    Output Parameter:
4569b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4579b94acceSBarry Smith 
4589b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4599b94acceSBarry Smith @*/
4609b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4619b94acceSBarry Smith {
46277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
46374679c65SBarry Smith   PetscValidIntPointer(nfails);
4649b94acceSBarry Smith   *nfails = snes->nfailures;
4659b94acceSBarry Smith   return 0;
4669b94acceSBarry Smith }
467a847f771SSatish Balay 
4685615d1e5SSatish Balay #undef __FUNC__
4695615d1e5SSatish Balay #define __FUNC__ "SNESGetSLES"
4709b94acceSBarry Smith /*@C
4719b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4729b94acceSBarry Smith 
4739b94acceSBarry Smith    Input Parameter:
4749b94acceSBarry Smith .  snes - the SNES context
4759b94acceSBarry Smith 
4769b94acceSBarry Smith    Output Parameter:
4779b94acceSBarry Smith .  sles - the SLES context
4789b94acceSBarry Smith 
4799b94acceSBarry Smith    Notes:
4809b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4819b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4829b94acceSBarry Smith    KSP and PC contexts as well.
4839b94acceSBarry Smith 
4849b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4859b94acceSBarry Smith 
4869b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4879b94acceSBarry Smith @*/
4889b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4899b94acceSBarry Smith {
49077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4919b94acceSBarry Smith   *sles = snes->sles;
4929b94acceSBarry Smith   return 0;
4939b94acceSBarry Smith }
4949b94acceSBarry Smith /* -----------------------------------------------------------*/
4955615d1e5SSatish Balay #undef __FUNC__
4965615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
4979b94acceSBarry Smith /*@C
4989b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
4999b94acceSBarry Smith 
5009b94acceSBarry Smith    Input Parameter:
5019b94acceSBarry Smith .  comm - MPI communicator
5029b94acceSBarry Smith .  type - type of method, one of
5039b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
5049b94acceSBarry Smith $      (for systems of nonlinear equations)
5059b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
5069b94acceSBarry Smith $      (for unconstrained minimization)
5079b94acceSBarry Smith 
5089b94acceSBarry Smith    Output Parameter:
5099b94acceSBarry Smith .  outsnes - the new SNES context
5109b94acceSBarry Smith 
511c1f60f51SBarry Smith    Options Database Key:
51219bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
51319bd6747SLois Curfman McInnes $              and no preconditioning matrix
51419bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
51519bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
51619bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
51719bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
518c1f60f51SBarry Smith 
5199b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5209b94acceSBarry Smith 
52163a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5229b94acceSBarry Smith @*/
5234b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5249b94acceSBarry Smith {
5259b94acceSBarry Smith   int                 ierr;
5269b94acceSBarry Smith   SNES                snes;
5279b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
52837fcc0dbSBarry Smith 
5299b94acceSBarry Smith   *outsnes = 0;
5302a0acf01SLois Curfman McInnes   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS)
531e3372554SBarry Smith     SETERRQ(1,0,"incorrect method type");
5320452661fSBarry Smith   PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm);
5339b94acceSBarry Smith   PLogObjectCreate(snes);
5349b94acceSBarry Smith   snes->max_its           = 50;
5359b94acceSBarry Smith   snes->max_funcs	  = 1000;
5369b94acceSBarry Smith   snes->norm		  = 0.0;
5375a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
538b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
539b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5409b94acceSBarry Smith     snes->atol		  = 1.e-10;
5415a2d0531SBarry Smith   }
542b4874afaSBarry Smith   else {
543b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
544b4874afaSBarry Smith     snes->ttol            = 0.0;
545b4874afaSBarry Smith     snes->atol		  = 1.e-50;
546b4874afaSBarry Smith   }
5479b94acceSBarry Smith   snes->xtol		  = 1.e-8;
548b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5499b94acceSBarry Smith   snes->nfuncs            = 0;
5509b94acceSBarry Smith   snes->nfailures         = 0;
5517a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
552639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5539b94acceSBarry Smith   snes->data              = 0;
5549b94acceSBarry Smith   snes->view              = 0;
5559b94acceSBarry Smith   snes->computeumfunction = 0;
5569b94acceSBarry Smith   snes->umfunP            = 0;
5579b94acceSBarry Smith   snes->fc                = 0;
5589b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5599b94acceSBarry Smith   snes->fmin              = -1.e30;
5609b94acceSBarry Smith   snes->method_class      = type;
5619b94acceSBarry Smith   snes->set_method_called = 0;
562a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
5639b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5647d1a2b2bSBarry Smith   snes->type              = -1;
5656f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5666f24a144SLois Curfman McInnes   snes->vwork             = 0;
5676f24a144SLois Curfman McInnes   snes->nwork             = 0;
5689b94acceSBarry Smith 
5699b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5700452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
5719b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5729b94acceSBarry Smith   kctx->version     = 2;
5739b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5749b94acceSBarry Smith                              this was too large for some test cases */
5759b94acceSBarry Smith   kctx->rtol_last   = 0;
5769b94acceSBarry Smith   kctx->rtol_max    = .9;
5779b94acceSBarry Smith   kctx->gamma       = 1.0;
5789b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5799b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5809b94acceSBarry Smith   kctx->threshold   = .1;
5819b94acceSBarry Smith   kctx->lresid_last = 0;
5829b94acceSBarry Smith   kctx->norm_last   = 0;
5839b94acceSBarry Smith 
5849b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5859b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5865334005bSBarry Smith 
5879b94acceSBarry Smith   *outsnes = snes;
5889b94acceSBarry Smith   return 0;
5899b94acceSBarry Smith }
5909b94acceSBarry Smith 
5919b94acceSBarry Smith /* --------------------------------------------------------------- */
5925615d1e5SSatish Balay #undef __FUNC__
5935615d1e5SSatish Balay #define __FUNC__ "SNESSetFunction"
5949b94acceSBarry Smith /*@C
5959b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
5969b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
5979b94acceSBarry Smith    equations.
5989b94acceSBarry Smith 
5999b94acceSBarry Smith    Input Parameters:
6009b94acceSBarry Smith .  snes - the SNES context
6019b94acceSBarry Smith .  func - function evaluation routine
6029b94acceSBarry Smith .  r - vector to store function value
6032cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
6042cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6059b94acceSBarry Smith 
6069b94acceSBarry Smith    Calling sequence of func:
607313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
6089b94acceSBarry Smith 
6099b94acceSBarry Smith .  x - input vector
610313e4042SLois Curfman McInnes .  f - function vector
6112cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
6129b94acceSBarry Smith 
6139b94acceSBarry Smith    Notes:
6149b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6159b94acceSBarry Smith $      f'(x) x = -f(x),
6169b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
6179b94acceSBarry Smith 
6189b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6199b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6209b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6219b94acceSBarry Smith 
6229b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6239b94acceSBarry Smith 
624a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6259b94acceSBarry Smith @*/
6265334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6279b94acceSBarry Smith {
62877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
629e3372554SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
6309b94acceSBarry Smith   snes->computefunction     = func;
6319b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6329b94acceSBarry Smith   snes->funP                = ctx;
6339b94acceSBarry Smith   return 0;
6349b94acceSBarry Smith }
6359b94acceSBarry Smith 
6365615d1e5SSatish Balay #undef __FUNC__
6375615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
6389b94acceSBarry Smith /*@
6399b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6409b94acceSBarry Smith    SNESSetFunction().
6419b94acceSBarry Smith 
6429b94acceSBarry Smith    Input Parameters:
6439b94acceSBarry Smith .  snes - the SNES context
6449b94acceSBarry Smith .  x - input vector
6459b94acceSBarry Smith 
6469b94acceSBarry Smith    Output Parameter:
6473638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6489b94acceSBarry Smith 
6491bffabb2SLois Curfman McInnes    Notes:
6509b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6519b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6529b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6539b94acceSBarry Smith 
6549b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6559b94acceSBarry Smith 
656a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6579b94acceSBarry Smith @*/
6589b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6599b94acceSBarry Smith {
6609b94acceSBarry Smith   int    ierr;
6619b94acceSBarry Smith 
66274679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
663e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
6649b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
6659b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
666*ae3c334cSLois Curfman McInnes   snes->nfuncs++;
6679b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6689b94acceSBarry Smith   return 0;
6699b94acceSBarry Smith }
6709b94acceSBarry Smith 
6715615d1e5SSatish Balay #undef __FUNC__
6725615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunction"
6739b94acceSBarry Smith /*@C
6749b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6759b94acceSBarry Smith    unconstrained minimization.
6769b94acceSBarry Smith 
6779b94acceSBarry Smith    Input Parameters:
6789b94acceSBarry Smith .  snes - the SNES context
6799b94acceSBarry Smith .  func - function evaluation routine
680044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
681044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6829b94acceSBarry Smith 
6839b94acceSBarry Smith    Calling sequence of func:
6849b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
6859b94acceSBarry Smith 
6869b94acceSBarry Smith .  x - input vector
6879b94acceSBarry Smith .  f - function
688044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
6899b94acceSBarry Smith 
6909b94acceSBarry Smith    Notes:
6919b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
6929b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
6939b94acceSBarry Smith    SNESSetFunction().
6949b94acceSBarry Smith 
6959b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
6969b94acceSBarry Smith 
697a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
698a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
6999b94acceSBarry Smith @*/
7009b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
7019b94acceSBarry Smith                       void *ctx)
7029b94acceSBarry Smith {
70377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
704e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
7059b94acceSBarry Smith   snes->computeumfunction   = func;
7069b94acceSBarry Smith   snes->umfunP              = ctx;
7079b94acceSBarry Smith   return 0;
7089b94acceSBarry Smith }
7099b94acceSBarry Smith 
7105615d1e5SSatish Balay #undef __FUNC__
7115615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7129b94acceSBarry Smith /*@
7139b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7149b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7159b94acceSBarry Smith 
7169b94acceSBarry Smith    Input Parameters:
7179b94acceSBarry Smith .  snes - the SNES context
7189b94acceSBarry Smith .  x - input vector
7199b94acceSBarry Smith 
7209b94acceSBarry Smith    Output Parameter:
7219b94acceSBarry Smith .  y - function value
7229b94acceSBarry Smith 
7239b94acceSBarry Smith    Notes:
7249b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7259b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7269b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
727a86d99e1SLois Curfman McInnes 
728a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
729a86d99e1SLois Curfman McInnes 
730a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
731a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7329b94acceSBarry Smith @*/
7339b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7349b94acceSBarry Smith {
7359b94acceSBarry Smith   int    ierr;
736e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
7379b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
7389b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
739*ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7409b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7419b94acceSBarry Smith   return 0;
7429b94acceSBarry Smith }
7439b94acceSBarry Smith 
7445615d1e5SSatish Balay #undef __FUNC__
7455615d1e5SSatish Balay #define __FUNC__ "SNESSetGradient"
7469b94acceSBarry Smith /*@C
7479b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7489b94acceSBarry Smith    vector for use by the SNES routines.
7499b94acceSBarry Smith 
7509b94acceSBarry Smith    Input Parameters:
7519b94acceSBarry Smith .  snes - the SNES context
7529b94acceSBarry Smith .  func - function evaluation routine
753044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
754044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7559b94acceSBarry Smith .  r - vector to store gradient value
7569b94acceSBarry Smith 
7579b94acceSBarry Smith    Calling sequence of func:
7589b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7599b94acceSBarry Smith 
7609b94acceSBarry Smith .  x - input vector
7619b94acceSBarry Smith .  g - gradient vector
762044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7639b94acceSBarry Smith 
7649b94acceSBarry Smith    Notes:
7659b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7669b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7679b94acceSBarry Smith    SNESSetFunction().
7689b94acceSBarry Smith 
7699b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7709b94acceSBarry Smith 
771a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
772a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
7739b94acceSBarry Smith @*/
77474679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
7759b94acceSBarry Smith {
77677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
777e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
7789b94acceSBarry Smith   snes->computefunction     = func;
7799b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7809b94acceSBarry Smith   snes->funP                = ctx;
7819b94acceSBarry Smith   return 0;
7829b94acceSBarry Smith }
7839b94acceSBarry Smith 
7845615d1e5SSatish Balay #undef __FUNC__
7855615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
7869b94acceSBarry Smith /*@
787a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
788a86d99e1SLois Curfman McInnes    SNESSetGradient().
7899b94acceSBarry Smith 
7909b94acceSBarry Smith    Input Parameters:
7919b94acceSBarry Smith .  snes - the SNES context
7929b94acceSBarry Smith .  x - input vector
7939b94acceSBarry Smith 
7949b94acceSBarry Smith    Output Parameter:
7959b94acceSBarry Smith .  y - gradient vector
7969b94acceSBarry Smith 
7979b94acceSBarry Smith    Notes:
7989b94acceSBarry Smith    SNESComputeGradient() is valid only for
7999b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8009b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
801a86d99e1SLois Curfman McInnes 
802a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
803a86d99e1SLois Curfman McInnes 
804a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
805a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8069b94acceSBarry Smith @*/
8079b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8089b94acceSBarry Smith {
8099b94acceSBarry Smith   int    ierr;
81074679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
811e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8129b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
8139b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
8149b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8159b94acceSBarry Smith   return 0;
8169b94acceSBarry Smith }
8179b94acceSBarry Smith 
8185615d1e5SSatish Balay #undef __FUNC__
8195615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
82062fef451SLois Curfman McInnes /*@
82162fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
82262fef451SLois Curfman McInnes    set with SNESSetJacobian().
82362fef451SLois Curfman McInnes 
82462fef451SLois Curfman McInnes    Input Parameters:
82562fef451SLois Curfman McInnes .  snes - the SNES context
82662fef451SLois Curfman McInnes .  x - input vector
82762fef451SLois Curfman McInnes 
82862fef451SLois Curfman McInnes    Output Parameters:
82962fef451SLois Curfman McInnes .  A - Jacobian matrix
83062fef451SLois Curfman McInnes .  B - optional preconditioning matrix
83162fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
83262fef451SLois Curfman McInnes 
83362fef451SLois Curfman McInnes    Notes:
83462fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
83562fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
83662fef451SLois Curfman McInnes 
837dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
838dc5a77f8SLois Curfman McInnes    flag parameter.
83962fef451SLois Curfman McInnes 
84062fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
84162fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
84262fef451SLois Curfman McInnes    methods is SNESComputeJacobian().
84362fef451SLois Curfman McInnes 
84462fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
84562fef451SLois Curfman McInnes 
84662fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
84762fef451SLois Curfman McInnes @*/
8489b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8499b94acceSBarry Smith {
8509b94acceSBarry Smith   int    ierr;
85174679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
852e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
8539b94acceSBarry Smith   if (!snes->computejacobian) return 0;
8549b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
855c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
8569b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
8579b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8586d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
85977c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
86077c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8619b94acceSBarry Smith   return 0;
8629b94acceSBarry Smith }
8639b94acceSBarry Smith 
8645615d1e5SSatish Balay #undef __FUNC__
8655615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
86662fef451SLois Curfman McInnes /*@
86762fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
86862fef451SLois Curfman McInnes    set with SNESSetHessian().
86962fef451SLois Curfman McInnes 
87062fef451SLois Curfman McInnes    Input Parameters:
87162fef451SLois Curfman McInnes .  snes - the SNES context
87262fef451SLois Curfman McInnes .  x - input vector
87362fef451SLois Curfman McInnes 
87462fef451SLois Curfman McInnes    Output Parameters:
87562fef451SLois Curfman McInnes .  A - Hessian matrix
87662fef451SLois Curfman McInnes .  B - optional preconditioning matrix
87762fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
87862fef451SLois Curfman McInnes 
87962fef451SLois Curfman McInnes    Notes:
88062fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
88162fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
88262fef451SLois Curfman McInnes 
883dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
884dc5a77f8SLois Curfman McInnes    flag parameter.
88562fef451SLois Curfman McInnes 
88662fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
88762fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
88862fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
88962fef451SLois Curfman McInnes 
89062fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
89162fef451SLois Curfman McInnes 
892a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
893a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
89462fef451SLois Curfman McInnes @*/
89562fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
8969b94acceSBarry Smith {
8979b94acceSBarry Smith   int    ierr;
89874679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
899e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9009b94acceSBarry Smith   if (!snes->computejacobian) return 0;
90162fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9020f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
90362fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
90462fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9050f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
90677c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
90777c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9089b94acceSBarry Smith   return 0;
9099b94acceSBarry Smith }
9109b94acceSBarry Smith 
9115615d1e5SSatish Balay #undef __FUNC__
9125615d1e5SSatish Balay #define __FUNC__ "SNESSetJacobian"
9139b94acceSBarry Smith /*@C
9149b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
915044dda88SLois Curfman McInnes    location to store the matrix.
9169b94acceSBarry Smith 
9179b94acceSBarry Smith    Input Parameters:
9189b94acceSBarry Smith .  snes - the SNES context
9199b94acceSBarry Smith .  A - Jacobian matrix
9209b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9219b94acceSBarry Smith .  func - Jacobian evaluation routine
9222cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9232cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9249b94acceSBarry Smith 
9259b94acceSBarry Smith    Calling sequence of func:
926313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9279b94acceSBarry Smith 
9289b94acceSBarry Smith .  x - input vector
9299b94acceSBarry Smith .  A - Jacobian matrix
9309b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
931ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
932ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9332cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
9349b94acceSBarry Smith 
9359b94acceSBarry Smith    Notes:
936dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9372cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
938ac21db08SLois Curfman McInnes 
939ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9409b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9419b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9429b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9439b94acceSBarry Smith    throughout the global iterations.
9449b94acceSBarry Smith 
9459b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9469b94acceSBarry Smith 
947ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
9489b94acceSBarry Smith @*/
9499b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9509b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9519b94acceSBarry Smith {
95277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
95374679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
954e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
9559b94acceSBarry Smith   snes->computejacobian = func;
9569b94acceSBarry Smith   snes->jacP            = ctx;
9579b94acceSBarry Smith   snes->jacobian        = A;
9589b94acceSBarry Smith   snes->jacobian_pre    = B;
9599b94acceSBarry Smith   return 0;
9609b94acceSBarry Smith }
96162fef451SLois Curfman McInnes 
9625615d1e5SSatish Balay #undef __FUNC__
9635615d1e5SSatish Balay #define __FUNC__ "SNESGetJacobian"
964b4fd4287SBarry Smith /*@
965b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
966b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
967b4fd4287SBarry Smith 
968b4fd4287SBarry Smith    Input Parameter:
969b4fd4287SBarry Smith .  snes - the nonlinear solver context
970b4fd4287SBarry Smith 
971b4fd4287SBarry Smith    Output Parameters:
972b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
973b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
974b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
975b4fd4287SBarry Smith 
976b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
977b4fd4287SBarry Smith @*/
978b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
979b4fd4287SBarry Smith {
98074679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
981e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
982b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
983b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
984b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
985b4fd4287SBarry Smith   return 0;
986b4fd4287SBarry Smith }
987b4fd4287SBarry Smith 
9885615d1e5SSatish Balay #undef __FUNC__
9895615d1e5SSatish Balay #define __FUNC__ "SNESSetHessian"
9909b94acceSBarry Smith /*@C
9919b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
992044dda88SLois Curfman McInnes    location to store the matrix.
9939b94acceSBarry Smith 
9949b94acceSBarry Smith    Input Parameters:
9959b94acceSBarry Smith .  snes - the SNES context
9969b94acceSBarry Smith .  A - Hessian matrix
9979b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
9989b94acceSBarry Smith .  func - Jacobian evaluation routine
999313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1000313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10019b94acceSBarry Smith 
10029b94acceSBarry Smith    Calling sequence of func:
1003313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10049b94acceSBarry Smith 
10059b94acceSBarry Smith .  x - input vector
10069b94acceSBarry Smith .  A - Hessian matrix
10079b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1008ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1009ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10102cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10119b94acceSBarry Smith 
10129b94acceSBarry Smith    Notes:
1013dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10142cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1015ac21db08SLois Curfman McInnes 
10169b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10179b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10189b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10199b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10209b94acceSBarry Smith    throughout the global iterations.
10219b94acceSBarry Smith 
10229b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10239b94acceSBarry Smith 
1024ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10259b94acceSBarry Smith @*/
10269b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10279b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10289b94acceSBarry Smith {
102977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
103074679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
1031e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10329b94acceSBarry Smith   snes->computejacobian = func;
10339b94acceSBarry Smith   snes->jacP            = ctx;
10349b94acceSBarry Smith   snes->jacobian        = A;
10359b94acceSBarry Smith   snes->jacobian_pre    = B;
10369b94acceSBarry Smith   return 0;
10379b94acceSBarry Smith }
10389b94acceSBarry Smith 
10395615d1e5SSatish Balay #undef __FUNC__
10405615d1e5SSatish Balay #define __FUNC__ "SNESGetHessian"
104162fef451SLois Curfman McInnes /*@
104262fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
104362fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
104462fef451SLois Curfman McInnes 
104562fef451SLois Curfman McInnes    Input Parameter:
104662fef451SLois Curfman McInnes .  snes - the nonlinear solver context
104762fef451SLois Curfman McInnes 
104862fef451SLois Curfman McInnes    Output Parameters:
104962fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
105062fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
105162fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
105262fef451SLois Curfman McInnes 
105362fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
105462fef451SLois Curfman McInnes @*/
105562fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
105662fef451SLois Curfman McInnes {
105774679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
1058e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
105962fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
106062fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
106162fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
106262fef451SLois Curfman McInnes   return 0;
106362fef451SLois Curfman McInnes }
106462fef451SLois Curfman McInnes 
10659b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10669b94acceSBarry Smith 
10675615d1e5SSatish Balay #undef __FUNC__
10685615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
10699b94acceSBarry Smith /*@
10709b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1071272ac6f2SLois Curfman McInnes    of a nonlinear solver.
10729b94acceSBarry Smith 
10739b94acceSBarry Smith    Input Parameter:
10749b94acceSBarry Smith .  snes - the SNES context
10758ddd3da0SLois Curfman McInnes .  x - the solution vector
10769b94acceSBarry Smith 
1077272ac6f2SLois Curfman McInnes    Notes:
1078272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1079272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1080272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1081272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1082272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1083272ac6f2SLois Curfman McInnes 
10849b94acceSBarry Smith .keywords: SNES, nonlinear, setup
10859b94acceSBarry Smith 
10869b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
10879b94acceSBarry Smith @*/
10888ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
10899b94acceSBarry Smith {
1090272ac6f2SLois Curfman McInnes   int ierr, flg;
109177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
109277c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
10938ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
10949b94acceSBarry Smith 
1095c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1096c1f60f51SBarry Smith   /*
1097c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1098dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1099c1f60f51SBarry Smith   */
1100c1f60f51SBarry Smith   if (flg) {
1101c1f60f51SBarry Smith     Mat J;
1102c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1103c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1104c1f60f51SBarry Smith     snes->mfshell = J;
1105c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1106c1f60f51SBarry Smith       snes->jacobian = J;
1107c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1108c1f60f51SBarry Smith     }
1109c1f60f51SBarry Smith     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1110c1f60f51SBarry Smith       snes->jacobian = J;
1111c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1112e3372554SBarry Smith     } else SETERRQ(1,0,"Method class doesn't support matrix-free operator option");
1113c1f60f51SBarry Smith   }
1114272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1115c1f60f51SBarry Smith   /*
1116dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1117c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1118c1f60f51SBarry Smith    */
111931615d04SLois Curfman McInnes   if (flg) {
1120272ac6f2SLois Curfman McInnes     Mat J;
1121272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1122272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1123272ac6f2SLois Curfman McInnes     snes->mfshell = J;
112483e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
112583e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
112694a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
112783e56afdSLois Curfman McInnes     }
112883e56afdSLois Curfman McInnes     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
112983e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
113094a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1131e3372554SBarry Smith     } else SETERRQ(1,0,"Method class doesn't support matrix-free option");
1132272ac6f2SLois Curfman McInnes   }
11339b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1134e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first");
1135e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first");
1136e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first");
1137e3372554SBarry Smith     if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector");
1138a703fe33SLois Curfman McInnes 
1139a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
114040191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1141a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1142a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1143a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1144a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1145a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1146a703fe33SLois Curfman McInnes     }
11479b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1148e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first");
1149e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first");
115037fcc0dbSBarry Smith     if (!snes->computeumfunction)
1151e3372554SBarry Smith       SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first");
1152e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first");
1153e3372554SBarry Smith   } else SETERRQ(1,0,"Unknown method class");
1154a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1155a703fe33SLois Curfman McInnes   snes->setup_called = 1;
1156a703fe33SLois Curfman McInnes   return 0;
11579b94acceSBarry Smith }
11589b94acceSBarry Smith 
11595615d1e5SSatish Balay #undef __FUNC__
11605615d1e5SSatish Balay #define __FUNC__ "SNESDestroy"
11619b94acceSBarry Smith /*@C
11629b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11639b94acceSBarry Smith    with SNESCreate().
11649b94acceSBarry Smith 
11659b94acceSBarry Smith    Input Parameter:
11669b94acceSBarry Smith .  snes - the SNES context
11679b94acceSBarry Smith 
11689b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
11699b94acceSBarry Smith 
117063a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
11719b94acceSBarry Smith @*/
11729b94acceSBarry Smith int SNESDestroy(SNES snes)
11739b94acceSBarry Smith {
11749b94acceSBarry Smith   int ierr;
117577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
11769b94acceSBarry Smith   ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);
11770452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
11789b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
11799b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
11803e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
11816f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
11829b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
11830452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
11849b94acceSBarry Smith   return 0;
11859b94acceSBarry Smith }
11869b94acceSBarry Smith 
11879b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11889b94acceSBarry Smith 
11895615d1e5SSatish Balay #undef __FUNC__
11905615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
11919b94acceSBarry Smith /*@
1192d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11939b94acceSBarry Smith 
11949b94acceSBarry Smith    Input Parameters:
11959b94acceSBarry Smith .  snes - the SNES context
119633174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
119733174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
119833174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
119933174efeSLois Curfman McInnes            of the change in the solution between steps
120033174efeSLois Curfman McInnes .  maxit - maximum number of iterations
120133174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
12029b94acceSBarry Smith 
120333174efeSLois Curfman McInnes    Options Database Keys:
120433174efeSLois Curfman McInnes $    -snes_atol <atol>
120533174efeSLois Curfman McInnes $    -snes_rtol <rtol>
120633174efeSLois Curfman McInnes $    -snes_stol <stol>
120733174efeSLois Curfman McInnes $    -snes_max_it <maxit>
120833174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
12099b94acceSBarry Smith 
1210d7a720efSLois Curfman McInnes    Notes:
12119b94acceSBarry Smith    The default maximum number of iterations is 50.
12129b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
12139b94acceSBarry Smith 
121433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
12159b94acceSBarry Smith 
1216d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
12179b94acceSBarry Smith @*/
1218d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
12199b94acceSBarry Smith {
122077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1221d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1222d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1223d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1224d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1225d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
12269b94acceSBarry Smith   return 0;
12279b94acceSBarry Smith }
12289b94acceSBarry Smith 
12295615d1e5SSatish Balay #undef __FUNC__
12305615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
12319b94acceSBarry Smith /*@
123233174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
123333174efeSLois Curfman McInnes 
123433174efeSLois Curfman McInnes    Input Parameters:
123533174efeSLois Curfman McInnes .  snes - the SNES context
123633174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
123733174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
123833174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
123933174efeSLois Curfman McInnes            of the change in the solution between steps
124033174efeSLois Curfman McInnes .  maxit - maximum number of iterations
124133174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
124233174efeSLois Curfman McInnes 
124333174efeSLois Curfman McInnes    Notes:
124433174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
124533174efeSLois Curfman McInnes 
124633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
124733174efeSLois Curfman McInnes 
124833174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
124933174efeSLois Curfman McInnes @*/
125033174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
125133174efeSLois Curfman McInnes {
125233174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
125333174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
125433174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
125533174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
125633174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
125733174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
125833174efeSLois Curfman McInnes   return 0;
125933174efeSLois Curfman McInnes }
126033174efeSLois Curfman McInnes 
12615615d1e5SSatish Balay #undef __FUNC__
12625615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
126333174efeSLois Curfman McInnes /*@
12649b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12659b94acceSBarry Smith 
12669b94acceSBarry Smith    Input Parameters:
12679b94acceSBarry Smith .  snes - the SNES context
12689b94acceSBarry Smith .  tol - tolerance
12699b94acceSBarry Smith 
12709b94acceSBarry Smith    Options Database Key:
1271d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
12729b94acceSBarry Smith 
12739b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12749b94acceSBarry Smith 
1275d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
12769b94acceSBarry Smith @*/
12779b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
12789b94acceSBarry Smith {
127977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12809b94acceSBarry Smith   snes->deltatol = tol;
12819b94acceSBarry Smith   return 0;
12829b94acceSBarry Smith }
12839b94acceSBarry Smith 
12845615d1e5SSatish Balay #undef __FUNC__
12855615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
12869b94acceSBarry Smith /*@
128777c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
12889b94acceSBarry Smith    for unconstrained minimization solvers.
12899b94acceSBarry Smith 
12909b94acceSBarry Smith    Input Parameters:
12919b94acceSBarry Smith .  snes - the SNES context
12929b94acceSBarry Smith .  ftol - minimum function tolerance
12939b94acceSBarry Smith 
12949b94acceSBarry Smith    Options Database Key:
1295d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
12969b94acceSBarry Smith 
12979b94acceSBarry Smith    Note:
129877c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
12999b94acceSBarry Smith    methods only.
13009b94acceSBarry Smith 
13019b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
13029b94acceSBarry Smith 
1303d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
13049b94acceSBarry Smith @*/
130577c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
13069b94acceSBarry Smith {
130777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13089b94acceSBarry Smith   snes->fmin = ftol;
13099b94acceSBarry Smith   return 0;
13109b94acceSBarry Smith }
13119b94acceSBarry Smith 
13129b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
13139b94acceSBarry Smith 
13145615d1e5SSatish Balay #undef __FUNC__
13155615d1e5SSatish Balay #define __FUNC__ "SNESSetMonitor"
13169b94acceSBarry Smith /*@C
13179b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
13189b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
13199b94acceSBarry Smith    progress.
13209b94acceSBarry Smith 
13219b94acceSBarry Smith    Input Parameters:
13229b94acceSBarry Smith .  snes - the SNES context
13239b94acceSBarry Smith .  func - monitoring routine
1324044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1325044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
13269b94acceSBarry Smith 
13279b94acceSBarry Smith    Calling sequence of func:
1328682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
13299b94acceSBarry Smith 
13309b94acceSBarry Smith $    snes - the SNES context
13319b94acceSBarry Smith $    its - iteration number
1332044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
13339b94acceSBarry Smith $
13349b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13359b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
13369b94acceSBarry Smith $
13379b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13389b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
13399b94acceSBarry Smith 
13409665c990SLois Curfman McInnes    Options Database Keys:
13419665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
13429665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
13439665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
13449665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
13459665c990SLois Curfman McInnes $                           been hardwired into a code by
13469665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
13479665c990SLois Curfman McInnes $                           does not cancel those set via
13489665c990SLois Curfman McInnes $                           the options database.
13499665c990SLois Curfman McInnes 
13509665c990SLois Curfman McInnes 
1351639f9d9dSBarry Smith    Notes:
13526bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13536bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13546bc08f3fSLois Curfman McInnes    order in which they were set.
1355639f9d9dSBarry Smith 
13569b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13579b94acceSBarry Smith 
13589b94acceSBarry Smith .seealso: SNESDefaultMonitor()
13599b94acceSBarry Smith @*/
136074679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
13619b94acceSBarry Smith {
1362639f9d9dSBarry Smith   if (!func) {
1363639f9d9dSBarry Smith     snes->numbermonitors = 0;
1364639f9d9dSBarry Smith     return 0;
1365639f9d9dSBarry Smith   }
1366639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1367e3372554SBarry Smith     SETERRQ(1,0,"Too many monitors set");
1368639f9d9dSBarry Smith   }
1369639f9d9dSBarry Smith 
1370639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1371639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
13729b94acceSBarry Smith   return 0;
13739b94acceSBarry Smith }
13749b94acceSBarry Smith 
13755615d1e5SSatish Balay #undef __FUNC__
13765615d1e5SSatish Balay #define __FUNC__ "SNESSetConvergenceTest"
13779b94acceSBarry Smith /*@C
13789b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
13799b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
13809b94acceSBarry Smith 
13819b94acceSBarry Smith    Input Parameters:
13829b94acceSBarry Smith .  snes - the SNES context
13839b94acceSBarry Smith .  func - routine to test for convergence
1384044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1385044dda88SLois Curfman McInnes           (may be PETSC_NULL)
13869b94acceSBarry Smith 
13879b94acceSBarry Smith    Calling sequence of func:
13889b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
13899b94acceSBarry Smith              double f,void *cctx)
13909b94acceSBarry Smith 
13919b94acceSBarry Smith $    snes - the SNES context
1392044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
13939b94acceSBarry Smith $    xnorm - 2-norm of current iterate
13949b94acceSBarry Smith $
13959b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13969b94acceSBarry Smith $    gnorm - 2-norm of current step
13979b94acceSBarry Smith $    f - 2-norm of function
13989b94acceSBarry Smith $
13999b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14009b94acceSBarry Smith $    gnorm - 2-norm of current gradient
14019b94acceSBarry Smith $    f - function value
14029b94acceSBarry Smith 
14039b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14049b94acceSBarry Smith 
140540191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
140640191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
14079b94acceSBarry Smith @*/
140874679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
14099b94acceSBarry Smith {
14109b94acceSBarry Smith   (snes)->converged = func;
14119b94acceSBarry Smith   (snes)->cnvP      = cctx;
14129b94acceSBarry Smith   return 0;
14139b94acceSBarry Smith }
14149b94acceSBarry Smith 
14155615d1e5SSatish Balay #undef __FUNC__
14165615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
14179b94acceSBarry Smith /*
14189b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
14199b94acceSBarry Smith    positive parameter delta.
14209b94acceSBarry Smith 
14219b94acceSBarry Smith     Input Parameters:
14229b94acceSBarry Smith .   snes - the SNES context
14239b94acceSBarry Smith .   y - approximate solution of linear system
14249b94acceSBarry Smith .   fnorm - 2-norm of current function
14259b94acceSBarry Smith .   delta - trust region size
14269b94acceSBarry Smith 
14279b94acceSBarry Smith     Output Parameters:
14289b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
14299b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
14309b94acceSBarry Smith     region, and exceeds zero otherwise.
14319b94acceSBarry Smith .   ynorm - 2-norm of the step
14329b94acceSBarry Smith 
14339b94acceSBarry Smith     Note:
143440191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
14359b94acceSBarry Smith     is set to be the maximum allowable step size.
14369b94acceSBarry Smith 
14379b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
14389b94acceSBarry Smith */
14399b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
14409b94acceSBarry Smith                           double *gpnorm,double *ynorm)
14419b94acceSBarry Smith {
14429b94acceSBarry Smith   double norm;
14439b94acceSBarry Smith   Scalar cnorm;
1444cddf8d76SBarry Smith   VecNorm(y,NORM_2, &norm );
14459b94acceSBarry Smith   if (norm > *delta) {
14469b94acceSBarry Smith      norm = *delta/norm;
14479b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
14489b94acceSBarry Smith      cnorm = norm;
14499b94acceSBarry Smith      VecScale( &cnorm, y );
14509b94acceSBarry Smith      *ynorm = *delta;
14519b94acceSBarry Smith   } else {
14529b94acceSBarry Smith      *gpnorm = 0.0;
14539b94acceSBarry Smith      *ynorm = norm;
14549b94acceSBarry Smith   }
14559b94acceSBarry Smith   return 0;
14569b94acceSBarry Smith }
14579b94acceSBarry Smith 
14585615d1e5SSatish Balay #undef __FUNC__
14595615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
14609b94acceSBarry Smith /*@
14619b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
146263a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
14639b94acceSBarry Smith 
14649b94acceSBarry Smith    Input Parameter:
14659b94acceSBarry Smith .  snes - the SNES context
14668ddd3da0SLois Curfman McInnes .  x - the solution vector
14679b94acceSBarry Smith 
14689b94acceSBarry Smith    Output Parameter:
14699b94acceSBarry Smith    its - number of iterations until termination
14709b94acceSBarry Smith 
14718ddd3da0SLois Curfman McInnes    Note:
14728ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
14738ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
14748ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
14758ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
14768ddd3da0SLois Curfman McInnes 
14779b94acceSBarry Smith .keywords: SNES, nonlinear, solve
14789b94acceSBarry Smith 
147963a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
14809b94acceSBarry Smith @*/
14818ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
14829b94acceSBarry Smith {
14833c7409f5SSatish Balay   int ierr, flg;
1484052efed2SBarry Smith 
148577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
148674679c65SBarry Smith   PetscValidIntPointer(its);
1487c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1488c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
14899b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
14909b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
14919b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
14923c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
14936d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
14949b94acceSBarry Smith   return 0;
14959b94acceSBarry Smith }
14969b94acceSBarry Smith 
14979b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
149837fcc0dbSBarry Smith static NRList *__SNESList = 0;
14999b94acceSBarry Smith 
15005615d1e5SSatish Balay #undef __FUNC__
15015615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
15029b94acceSBarry Smith /*@
15034b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
15049b94acceSBarry Smith 
15059b94acceSBarry Smith    Input Parameters:
15069b94acceSBarry Smith .  snes - the SNES context
15079b94acceSBarry Smith .  method - a known method
15089b94acceSBarry Smith 
1509ae12b187SLois Curfman McInnes   Options Database Command:
1510ae12b187SLois Curfman McInnes $ -snes_type  <method>
1511ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1512ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1513ae12b187SLois Curfman McInnes 
15149b94acceSBarry Smith    Notes:
15159b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
15169b94acceSBarry Smith $  Systems of nonlinear equations:
151740191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
151840191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
15199b94acceSBarry Smith $  Unconstrained minimization:
152040191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
152140191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
15229b94acceSBarry Smith 
1523ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1524ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1525ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1526ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1527ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1528ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1529ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1530ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1531ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1532ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1533a703fe33SLois Curfman McInnes 
1534f536c699SSatish Balay .keywords: SNES, set, method
15359b94acceSBarry Smith @*/
15364b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
15379b94acceSBarry Smith {
15389b94acceSBarry Smith   int (*r)(SNES);
153977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15407d1a2b2bSBarry Smith   if (snes->type == (int) method) return 0;
15417d1a2b2bSBarry Smith 
15429b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
154337fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
1544e3372554SBarry Smith   if (!__SNESList) {SETERRQ(1,0,"Could not get methods");}
154537fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1546e3372554SBarry Smith   if (!r) {SETERRQ(1,0,"Unknown method");}
15470452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
15489b94acceSBarry Smith   snes->set_method_called = 1;
15499b94acceSBarry Smith   return (*r)(snes);
15509b94acceSBarry Smith }
15519b94acceSBarry Smith 
15529b94acceSBarry Smith /* --------------------------------------------------------------------- */
15535615d1e5SSatish Balay #undef __FUNC__
15545615d1e5SSatish Balay #define __FUNC__ "SNESRegister"
15559b94acceSBarry Smith /*@C
15569b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
15574b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
15589b94acceSBarry Smith 
15599b94acceSBarry Smith    Input Parameters:
156040191667SLois Curfman McInnes .  name - for instance SNES_EQ_LS, SNES_EQ_TR, ...
156140191667SLois Curfman McInnes .  sname - corresponding string for name
15629b94acceSBarry Smith .  create - routine to create method context
15639b94acceSBarry Smith 
15649b94acceSBarry Smith .keywords: SNES, nonlinear, register
15659b94acceSBarry Smith 
15669b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
15679b94acceSBarry Smith @*/
15689b94acceSBarry Smith int SNESRegister(int name, char *sname, int (*create)(SNES))
15699b94acceSBarry Smith {
15709b94acceSBarry Smith   int ierr;
157137fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
157237fcc0dbSBarry Smith   NRRegister( __SNESList, name, sname, (int (*)(void*))create );
15739b94acceSBarry Smith   return 0;
15749b94acceSBarry Smith }
1575a847f771SSatish Balay 
15769b94acceSBarry Smith /* --------------------------------------------------------------------- */
15775615d1e5SSatish Balay #undef __FUNC__
15785615d1e5SSatish Balay #define __FUNC__ "SNESRegisterDestroy"
15799b94acceSBarry Smith /*@C
15809b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
15819b94acceSBarry Smith    registered by SNESRegister().
15829b94acceSBarry Smith 
15839b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
15849b94acceSBarry Smith 
15859b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
15869b94acceSBarry Smith @*/
15879b94acceSBarry Smith int SNESRegisterDestroy()
15889b94acceSBarry Smith {
158937fcc0dbSBarry Smith   if (__SNESList) {
159037fcc0dbSBarry Smith     NRDestroy( __SNESList );
159137fcc0dbSBarry Smith     __SNESList = 0;
15929b94acceSBarry Smith   }
15939b94acceSBarry Smith   return 0;
15949b94acceSBarry Smith }
15959b94acceSBarry Smith 
15965615d1e5SSatish Balay #undef __FUNC__
15975615d1e5SSatish Balay #define __FUNC__ "SNESGetTypeFromOptions_Private"
15989b94acceSBarry Smith /*
15994b0e389bSBarry Smith    SNESGetTypeFromOptions_Private - Sets the selected method from the
16009b94acceSBarry Smith    options database.
16019b94acceSBarry Smith 
16029b94acceSBarry Smith    Input Parameter:
16039b94acceSBarry Smith .  ctx - the SNES context
16049b94acceSBarry Smith 
16059b94acceSBarry Smith    Output Parameter:
16069b94acceSBarry Smith .  method -  solver method
16079b94acceSBarry Smith 
16089b94acceSBarry Smith    Returns:
16099b94acceSBarry Smith    Returns 1 if the method is found; 0 otherwise.
16109b94acceSBarry Smith 
16119b94acceSBarry Smith    Options Database Key:
16124b0e389bSBarry Smith $  -snes_type  method
16139b94acceSBarry Smith */
1614052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
16159b94acceSBarry Smith {
1616052efed2SBarry Smith   int  ierr;
16179b94acceSBarry Smith   char sbuf[50];
16185baf8537SBarry Smith 
1619052efed2SBarry Smith   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1620052efed2SBarry Smith   if (*flg) {
1621052efed2SBarry Smith     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
162237fcc0dbSBarry Smith     *method = (SNESType)NRFindID( __SNESList, sbuf );
16239b94acceSBarry Smith   }
16249b94acceSBarry Smith   return 0;
16259b94acceSBarry Smith }
16269b94acceSBarry Smith 
16275615d1e5SSatish Balay #undef __FUNC__
16285615d1e5SSatish Balay #define __FUNC__ "SNESGetType"
16299b94acceSBarry Smith /*@C
16309a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
16319b94acceSBarry Smith 
16329b94acceSBarry Smith    Input Parameter:
16334b0e389bSBarry Smith .  snes - nonlinear solver context
16349b94acceSBarry Smith 
16359b94acceSBarry Smith    Output Parameter:
16369a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
16379a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
16389b94acceSBarry Smith 
16399b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
16409b94acceSBarry Smith @*/
16414b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
16429b94acceSBarry Smith {
164306a9b0adSLois Curfman McInnes   int ierr;
164437fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
16454b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
164637fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
16479b94acceSBarry Smith   return 0;
16489b94acceSBarry Smith }
16499b94acceSBarry Smith 
16509b94acceSBarry Smith #include <stdio.h>
16515615d1e5SSatish Balay #undef __FUNC__
16525615d1e5SSatish Balay #define __FUNC__ "SNESPrintTypes_Private"
16539b94acceSBarry Smith /*
16544b0e389bSBarry Smith    SNESPrintTypes_Private - Prints the SNES methods available from the
16559b94acceSBarry Smith    options database.
16569b94acceSBarry Smith 
16579b94acceSBarry Smith    Input Parameters:
165833455573SSatish Balay .  comm   - communicator (usually MPI_COMM_WORLD)
16599b94acceSBarry Smith .  prefix - prefix (usually "-")
16604b0e389bSBarry Smith .  name   - the options database name (by default "snes_type")
16619b94acceSBarry Smith */
166233455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
16639b94acceSBarry Smith {
16649b94acceSBarry Smith   FuncList *entry;
166537fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
166637fcc0dbSBarry Smith   entry = __SNESList->head;
166777c4ece6SBarry Smith   PetscPrintf(comm," %s%s (one of)",prefix,name);
16689b94acceSBarry Smith   while (entry) {
166977c4ece6SBarry Smith     PetscPrintf(comm," %s",entry->name);
16709b94acceSBarry Smith     entry = entry->next;
16719b94acceSBarry Smith   }
167277c4ece6SBarry Smith   PetscPrintf(comm,"\n");
16739b94acceSBarry Smith   return 0;
16749b94acceSBarry Smith }
16759b94acceSBarry Smith 
16765615d1e5SSatish Balay #undef __FUNC__
16775615d1e5SSatish Balay #define __FUNC__ "SNESGetSolution"
16789b94acceSBarry Smith /*@C
16799b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
16809b94acceSBarry Smith    stored.
16819b94acceSBarry Smith 
16829b94acceSBarry Smith    Input Parameter:
16839b94acceSBarry Smith .  snes - the SNES context
16849b94acceSBarry Smith 
16859b94acceSBarry Smith    Output Parameter:
16869b94acceSBarry Smith .  x - the solution
16879b94acceSBarry Smith 
16889b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
16899b94acceSBarry Smith 
1690a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
16919b94acceSBarry Smith @*/
16929b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
16939b94acceSBarry Smith {
169477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16959b94acceSBarry Smith   *x = snes->vec_sol_always;
16969b94acceSBarry Smith   return 0;
16979b94acceSBarry Smith }
16989b94acceSBarry Smith 
16995615d1e5SSatish Balay #undef __FUNC__
17005615d1e5SSatish Balay #define __FUNC__ "SNESGetSolutionUpdate"
17019b94acceSBarry Smith /*@C
17029b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
17039b94acceSBarry Smith    stored.
17049b94acceSBarry Smith 
17059b94acceSBarry Smith    Input Parameter:
17069b94acceSBarry Smith .  snes - the SNES context
17079b94acceSBarry Smith 
17089b94acceSBarry Smith    Output Parameter:
17099b94acceSBarry Smith .  x - the solution update
17109b94acceSBarry Smith 
17119b94acceSBarry Smith    Notes:
17129b94acceSBarry Smith    This vector is implementation dependent.
17139b94acceSBarry Smith 
17149b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
17159b94acceSBarry Smith 
17169b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
17179b94acceSBarry Smith @*/
17189b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
17199b94acceSBarry Smith {
172077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17219b94acceSBarry Smith   *x = snes->vec_sol_update_always;
17229b94acceSBarry Smith   return 0;
17239b94acceSBarry Smith }
17249b94acceSBarry Smith 
17255615d1e5SSatish Balay #undef __FUNC__
17265615d1e5SSatish Balay #define __FUNC__ "SNESGetFunction"
17279b94acceSBarry Smith /*@C
17283638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
17299b94acceSBarry Smith 
17309b94acceSBarry Smith    Input Parameter:
17319b94acceSBarry Smith .  snes - the SNES context
17329b94acceSBarry Smith 
17339b94acceSBarry Smith    Output Parameter:
17343638b69dSLois Curfman McInnes .  r - the function
17359b94acceSBarry Smith 
17369b94acceSBarry Smith    Notes:
17379b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
17389b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
17399b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
17409b94acceSBarry Smith 
1741a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
17429b94acceSBarry Smith 
174331615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
174431615d04SLois Curfman McInnes           SNESGetGradient()
17459b94acceSBarry Smith @*/
17469b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
17479b94acceSBarry Smith {
174877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1749e3372554SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,
17507c792091SSatish Balay     "For SNES_NONLINEAR_EQUATIONS only");
17519b94acceSBarry Smith   *r = snes->vec_func_always;
17529b94acceSBarry Smith   return 0;
17539b94acceSBarry Smith }
17549b94acceSBarry Smith 
17555615d1e5SSatish Balay #undef __FUNC__
17565615d1e5SSatish Balay #define __FUNC__ "SNESGetGradient"
17579b94acceSBarry Smith /*@C
17583638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
17599b94acceSBarry Smith 
17609b94acceSBarry Smith    Input Parameter:
17619b94acceSBarry Smith .  snes - the SNES context
17629b94acceSBarry Smith 
17639b94acceSBarry Smith    Output Parameter:
17649b94acceSBarry Smith .  r - the gradient
17659b94acceSBarry Smith 
17669b94acceSBarry Smith    Notes:
17679b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
17689b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
17699b94acceSBarry Smith    SNESGetFunction().
17709b94acceSBarry Smith 
17719b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
17729b94acceSBarry Smith 
177331615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
17749b94acceSBarry Smith @*/
17759b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
17769b94acceSBarry Smith {
177777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1778e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,
177963c41f6aSSatish Balay     "For SNES_UNCONSTRAINED_MINIMIZATION only");
17809b94acceSBarry Smith   *r = snes->vec_func_always;
17819b94acceSBarry Smith   return 0;
17829b94acceSBarry Smith }
17839b94acceSBarry Smith 
17845615d1e5SSatish Balay #undef __FUNC__
17855615d1e5SSatish Balay #define __FUNC__ "SNESGetMinimizationFunction"
17869b94acceSBarry Smith /*@
17879b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
17889b94acceSBarry Smith    unconstrained minimization problems.
17899b94acceSBarry Smith 
17909b94acceSBarry Smith    Input Parameter:
17919b94acceSBarry Smith .  snes - the SNES context
17929b94acceSBarry Smith 
17939b94acceSBarry Smith    Output Parameter:
17949b94acceSBarry Smith .  r - the function
17959b94acceSBarry Smith 
17969b94acceSBarry Smith    Notes:
17979b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
17989b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
17999b94acceSBarry Smith    SNESGetFunction().
18009b94acceSBarry Smith 
18019b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
18029b94acceSBarry Smith 
180331615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
18049b94acceSBarry Smith @*/
18059b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
18069b94acceSBarry Smith {
180777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
180874679c65SBarry Smith   PetscValidScalarPointer(r);
1809e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,
181063c41f6aSSatish Balay     "For SNES_UNCONSTRAINED_MINIMIZATION only");
18119b94acceSBarry Smith   *r = snes->fc;
18129b94acceSBarry Smith   return 0;
18139b94acceSBarry Smith }
18149b94acceSBarry Smith 
18155615d1e5SSatish Balay #undef __FUNC__
18165615d1e5SSatish Balay #define __FUNC__ "SNESSetOptionsPrefix"
18173c7409f5SSatish Balay /*@C
18183c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1819639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1820639f9d9dSBarry Smith    the prefix name.
18213c7409f5SSatish Balay 
18223c7409f5SSatish Balay    Input Parameter:
18233c7409f5SSatish Balay .  snes - the SNES context
18243c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18253c7409f5SSatish Balay 
18263c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1827a86d99e1SLois Curfman McInnes 
1828a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
18293c7409f5SSatish Balay @*/
18303c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
18313c7409f5SSatish Balay {
18323c7409f5SSatish Balay   int ierr;
18333c7409f5SSatish Balay 
183477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1835639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18363c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
18373c7409f5SSatish Balay   return 0;
18383c7409f5SSatish Balay }
18393c7409f5SSatish Balay 
18405615d1e5SSatish Balay #undef __FUNC__
18415615d1e5SSatish Balay #define __FUNC__ "SNESAppendOptionsPrefix"
18423c7409f5SSatish Balay /*@C
1843f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1844639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1845639f9d9dSBarry Smith    the prefix name.
18463c7409f5SSatish Balay 
18473c7409f5SSatish Balay    Input Parameter:
18483c7409f5SSatish Balay .  snes - the SNES context
18493c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18503c7409f5SSatish Balay 
18513c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1852a86d99e1SLois Curfman McInnes 
1853a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
18543c7409f5SSatish Balay @*/
18553c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
18563c7409f5SSatish Balay {
18573c7409f5SSatish Balay   int ierr;
18583c7409f5SSatish Balay 
185977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1860639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18613c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
18623c7409f5SSatish Balay   return 0;
18633c7409f5SSatish Balay }
18643c7409f5SSatish Balay 
18655615d1e5SSatish Balay #undef __FUNC__
18665615d1e5SSatish Balay #define __FUNC__ "SNESGetOptionsPrefix"
1867c83590e2SSatish Balay /*@
18683c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
18693c7409f5SSatish Balay    SNES options in the database.
18703c7409f5SSatish Balay 
18713c7409f5SSatish Balay    Input Parameter:
18723c7409f5SSatish Balay .  snes - the SNES context
18733c7409f5SSatish Balay 
18743c7409f5SSatish Balay    Output Parameter:
18753c7409f5SSatish Balay .  prefix - pointer to the prefix string used
18763c7409f5SSatish Balay 
18773c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1878a86d99e1SLois Curfman McInnes 
1879a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
18803c7409f5SSatish Balay @*/
18813c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
18823c7409f5SSatish Balay {
18833c7409f5SSatish Balay   int ierr;
18843c7409f5SSatish Balay 
188577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1886639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18873c7409f5SSatish Balay   return 0;
18883c7409f5SSatish Balay }
18893c7409f5SSatish Balay 
18903c7409f5SSatish Balay 
18919b94acceSBarry Smith 
18929b94acceSBarry Smith 
18939b94acceSBarry Smith 
1894