xref: /petsc/src/snes/interface/snes.c (revision c9005455e5a2a46f3e06d732e71c5ab89b366105)
19b94acceSBarry Smith #ifndef lint
2*c9005455SLois Curfman McInnes static char vcid[] = "$Id: snes.c,v 1.111 1997/01/21 19:13:51 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);
69ae3c334cSLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
7068632166SLois Curfman McInnes      "  total number of function evaluations=%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 
458c96a6f78SLois Curfman McInnes    Notes:
459c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
460c96a6f78SLois Curfman McInnes 
4619b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4629b94acceSBarry Smith @*/
4639b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4649b94acceSBarry Smith {
46577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
46674679c65SBarry Smith   PetscValidIntPointer(nfails);
4679b94acceSBarry Smith   *nfails = snes->nfailures;
4689b94acceSBarry Smith   return 0;
4699b94acceSBarry Smith }
470a847f771SSatish Balay 
4715615d1e5SSatish Balay #undef __FUNC__
472c96a6f78SLois Curfman McInnes #define __FUNC__ "SNESGetNumberLinearIterations"
473c96a6f78SLois Curfman McInnes /*@
474c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
475c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
476c96a6f78SLois Curfman McInnes 
477c96a6f78SLois Curfman McInnes    Input Parameter:
478c96a6f78SLois Curfman McInnes .  snes - SNES context
479c96a6f78SLois Curfman McInnes 
480c96a6f78SLois Curfman McInnes    Output Parameter:
481c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
482c96a6f78SLois Curfman McInnes 
483c96a6f78SLois Curfman McInnes    Notes:
484c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
485c96a6f78SLois Curfman McInnes 
486c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
487c96a6f78SLois Curfman McInnes @*/
48886bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
489c96a6f78SLois Curfman McInnes {
490c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
491c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
492c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
493c96a6f78SLois Curfman McInnes   return 0;
494c96a6f78SLois Curfman McInnes }
495c96a6f78SLois Curfman McInnes 
496c96a6f78SLois Curfman McInnes #undef __FUNC__
4975615d1e5SSatish Balay #define __FUNC__ "SNESGetSLES"
4989b94acceSBarry Smith /*@C
4999b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
5009b94acceSBarry Smith 
5019b94acceSBarry Smith    Input Parameter:
5029b94acceSBarry Smith .  snes - the SNES context
5039b94acceSBarry Smith 
5049b94acceSBarry Smith    Output Parameter:
5059b94acceSBarry Smith .  sles - the SLES context
5069b94acceSBarry Smith 
5079b94acceSBarry Smith    Notes:
5089b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
5099b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5109b94acceSBarry Smith    KSP and PC contexts as well.
5119b94acceSBarry Smith 
5129b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
5139b94acceSBarry Smith 
5149b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
5159b94acceSBarry Smith @*/
5169b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
5179b94acceSBarry Smith {
51877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5199b94acceSBarry Smith   *sles = snes->sles;
5209b94acceSBarry Smith   return 0;
5219b94acceSBarry Smith }
5229b94acceSBarry Smith /* -----------------------------------------------------------*/
5235615d1e5SSatish Balay #undef __FUNC__
5245615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
5259b94acceSBarry Smith /*@C
5269b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5279b94acceSBarry Smith 
5289b94acceSBarry Smith    Input Parameter:
5299b94acceSBarry Smith .  comm - MPI communicator
5309b94acceSBarry Smith .  type - type of method, one of
5319b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
5329b94acceSBarry Smith $      (for systems of nonlinear equations)
5339b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
5349b94acceSBarry Smith $      (for unconstrained minimization)
5359b94acceSBarry Smith 
5369b94acceSBarry Smith    Output Parameter:
5379b94acceSBarry Smith .  outsnes - the new SNES context
5389b94acceSBarry Smith 
539c1f60f51SBarry Smith    Options Database Key:
54019bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
54119bd6747SLois Curfman McInnes $              and no preconditioning matrix
54219bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
54319bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
54419bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
54519bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
546c1f60f51SBarry Smith 
5479b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5489b94acceSBarry Smith 
54963a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5509b94acceSBarry Smith @*/
5514b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5529b94acceSBarry Smith {
5539b94acceSBarry Smith   int                 ierr;
5549b94acceSBarry Smith   SNES                snes;
5559b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
55637fcc0dbSBarry Smith 
5579b94acceSBarry Smith   *outsnes = 0;
5582a0acf01SLois Curfman McInnes   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS)
559e3372554SBarry Smith     SETERRQ(1,0,"incorrect method type");
5600452661fSBarry Smith   PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm);
5619b94acceSBarry Smith   PLogObjectCreate(snes);
5629b94acceSBarry Smith   snes->max_its           = 50;
5639b94acceSBarry Smith   snes->max_funcs	  = 1000;
5649b94acceSBarry Smith   snes->norm		  = 0.0;
5655a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
566b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
567b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5689b94acceSBarry Smith     snes->atol		  = 1.e-10;
5695a2d0531SBarry Smith   }
570b4874afaSBarry Smith   else {
571b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
572b4874afaSBarry Smith     snes->ttol            = 0.0;
573b4874afaSBarry Smith     snes->atol		  = 1.e-50;
574b4874afaSBarry Smith   }
5759b94acceSBarry Smith   snes->xtol		  = 1.e-8;
576b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5779b94acceSBarry Smith   snes->nfuncs            = 0;
5789b94acceSBarry Smith   snes->nfailures         = 0;
5797a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
580639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5819b94acceSBarry Smith   snes->data              = 0;
5829b94acceSBarry Smith   snes->view              = 0;
5839b94acceSBarry Smith   snes->computeumfunction = 0;
5849b94acceSBarry Smith   snes->umfunP            = 0;
5859b94acceSBarry Smith   snes->fc                = 0;
5869b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5879b94acceSBarry Smith   snes->fmin              = -1.e30;
5889b94acceSBarry Smith   snes->method_class      = type;
5899b94acceSBarry Smith   snes->set_method_called = 0;
590a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
5919b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5927d1a2b2bSBarry Smith   snes->type              = -1;
5936f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5946f24a144SLois Curfman McInnes   snes->vwork             = 0;
5956f24a144SLois Curfman McInnes   snes->nwork             = 0;
596*c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
597*c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
598*c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
5999b94acceSBarry Smith 
6009b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
6010452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
6029b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6039b94acceSBarry Smith   kctx->version     = 2;
6049b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6059b94acceSBarry Smith                              this was too large for some test cases */
6069b94acceSBarry Smith   kctx->rtol_last   = 0;
6079b94acceSBarry Smith   kctx->rtol_max    = .9;
6089b94acceSBarry Smith   kctx->gamma       = 1.0;
6099b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6109b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6119b94acceSBarry Smith   kctx->threshold   = .1;
6129b94acceSBarry Smith   kctx->lresid_last = 0;
6139b94acceSBarry Smith   kctx->norm_last   = 0;
6149b94acceSBarry Smith 
6159b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
6169b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
6175334005bSBarry Smith 
6189b94acceSBarry Smith   *outsnes = snes;
6199b94acceSBarry Smith   return 0;
6209b94acceSBarry Smith }
6219b94acceSBarry Smith 
6229b94acceSBarry Smith /* --------------------------------------------------------------- */
6235615d1e5SSatish Balay #undef __FUNC__
6245615d1e5SSatish Balay #define __FUNC__ "SNESSetFunction"
6259b94acceSBarry Smith /*@C
6269b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6279b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6289b94acceSBarry Smith    equations.
6299b94acceSBarry Smith 
6309b94acceSBarry Smith    Input Parameters:
6319b94acceSBarry Smith .  snes - the SNES context
6329b94acceSBarry Smith .  func - function evaluation routine
6339b94acceSBarry Smith .  r - vector to store function value
6342cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
6352cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6369b94acceSBarry Smith 
6379b94acceSBarry Smith    Calling sequence of func:
638313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
6399b94acceSBarry Smith 
6409b94acceSBarry Smith .  x - input vector
641313e4042SLois Curfman McInnes .  f - function vector
6422cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
6439b94acceSBarry Smith 
6449b94acceSBarry Smith    Notes:
6459b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6469b94acceSBarry Smith $      f'(x) x = -f(x),
6479b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
6489b94acceSBarry Smith 
6499b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6509b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6519b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6529b94acceSBarry Smith 
6539b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6549b94acceSBarry Smith 
655a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6569b94acceSBarry Smith @*/
6575334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6589b94acceSBarry Smith {
65977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
660e3372554SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
6619b94acceSBarry Smith   snes->computefunction     = func;
6629b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6639b94acceSBarry Smith   snes->funP                = ctx;
6649b94acceSBarry Smith   return 0;
6659b94acceSBarry Smith }
6669b94acceSBarry Smith 
6675615d1e5SSatish Balay #undef __FUNC__
6685615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
6699b94acceSBarry Smith /*@
6709b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6719b94acceSBarry Smith    SNESSetFunction().
6729b94acceSBarry Smith 
6739b94acceSBarry Smith    Input Parameters:
6749b94acceSBarry Smith .  snes - the SNES context
6759b94acceSBarry Smith .  x - input vector
6769b94acceSBarry Smith 
6779b94acceSBarry Smith    Output Parameter:
6783638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6799b94acceSBarry Smith 
6801bffabb2SLois Curfman McInnes    Notes:
6819b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6829b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6839b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6849b94acceSBarry Smith 
6859b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6869b94acceSBarry Smith 
687a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6889b94acceSBarry Smith @*/
6899b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6909b94acceSBarry Smith {
6919b94acceSBarry Smith   int    ierr;
6929b94acceSBarry Smith 
69374679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
694e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
6959b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
6969b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
697ae3c334cSLois Curfman McInnes   snes->nfuncs++;
6989b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6999b94acceSBarry Smith   return 0;
7009b94acceSBarry Smith }
7019b94acceSBarry Smith 
7025615d1e5SSatish Balay #undef __FUNC__
7035615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunction"
7049b94acceSBarry Smith /*@C
7059b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
7069b94acceSBarry Smith    unconstrained minimization.
7079b94acceSBarry Smith 
7089b94acceSBarry Smith    Input Parameters:
7099b94acceSBarry Smith .  snes - the SNES context
7109b94acceSBarry Smith .  func - function evaluation routine
711044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
712044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7139b94acceSBarry Smith 
7149b94acceSBarry Smith    Calling sequence of func:
7159b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
7169b94acceSBarry Smith 
7179b94acceSBarry Smith .  x - input vector
7189b94acceSBarry Smith .  f - function
719044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
7209b94acceSBarry Smith 
7219b94acceSBarry Smith    Notes:
7229b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7239b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7249b94acceSBarry Smith    SNESSetFunction().
7259b94acceSBarry Smith 
7269b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
7279b94acceSBarry Smith 
728a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
729a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
7309b94acceSBarry Smith @*/
7319b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
7329b94acceSBarry Smith                       void *ctx)
7339b94acceSBarry Smith {
73477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
735e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
7369b94acceSBarry Smith   snes->computeumfunction   = func;
7379b94acceSBarry Smith   snes->umfunP              = ctx;
7389b94acceSBarry Smith   return 0;
7399b94acceSBarry Smith }
7409b94acceSBarry Smith 
7415615d1e5SSatish Balay #undef __FUNC__
7425615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7439b94acceSBarry Smith /*@
7449b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7459b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7469b94acceSBarry Smith 
7479b94acceSBarry Smith    Input Parameters:
7489b94acceSBarry Smith .  snes - the SNES context
7499b94acceSBarry Smith .  x - input vector
7509b94acceSBarry Smith 
7519b94acceSBarry Smith    Output Parameter:
7529b94acceSBarry Smith .  y - function value
7539b94acceSBarry Smith 
7549b94acceSBarry Smith    Notes:
7559b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7569b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7579b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
758a86d99e1SLois Curfman McInnes 
759a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
760a86d99e1SLois Curfman McInnes 
761a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
762a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7639b94acceSBarry Smith @*/
7649b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7659b94acceSBarry Smith {
7669b94acceSBarry Smith   int    ierr;
767e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
7689b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
7699b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
770ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7719b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7729b94acceSBarry Smith   return 0;
7739b94acceSBarry Smith }
7749b94acceSBarry Smith 
7755615d1e5SSatish Balay #undef __FUNC__
7765615d1e5SSatish Balay #define __FUNC__ "SNESSetGradient"
7779b94acceSBarry Smith /*@C
7789b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7799b94acceSBarry Smith    vector for use by the SNES routines.
7809b94acceSBarry Smith 
7819b94acceSBarry Smith    Input Parameters:
7829b94acceSBarry Smith .  snes - the SNES context
7839b94acceSBarry Smith .  func - function evaluation routine
784044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
785044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7869b94acceSBarry Smith .  r - vector to store gradient value
7879b94acceSBarry Smith 
7889b94acceSBarry Smith    Calling sequence of func:
7899b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7909b94acceSBarry Smith 
7919b94acceSBarry Smith .  x - input vector
7929b94acceSBarry Smith .  g - gradient vector
793044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7949b94acceSBarry Smith 
7959b94acceSBarry Smith    Notes:
7969b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7979b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7989b94acceSBarry Smith    SNESSetFunction().
7999b94acceSBarry Smith 
8009b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
8019b94acceSBarry Smith 
802a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
803a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
8049b94acceSBarry Smith @*/
80574679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
8069b94acceSBarry Smith {
80777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
808e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8099b94acceSBarry Smith   snes->computefunction     = func;
8109b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
8119b94acceSBarry Smith   snes->funP                = ctx;
8129b94acceSBarry Smith   return 0;
8139b94acceSBarry Smith }
8149b94acceSBarry Smith 
8155615d1e5SSatish Balay #undef __FUNC__
8165615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
8179b94acceSBarry Smith /*@
818a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
819a86d99e1SLois Curfman McInnes    SNESSetGradient().
8209b94acceSBarry Smith 
8219b94acceSBarry Smith    Input Parameters:
8229b94acceSBarry Smith .  snes - the SNES context
8239b94acceSBarry Smith .  x - input vector
8249b94acceSBarry Smith 
8259b94acceSBarry Smith    Output Parameter:
8269b94acceSBarry Smith .  y - gradient vector
8279b94acceSBarry Smith 
8289b94acceSBarry Smith    Notes:
8299b94acceSBarry Smith    SNESComputeGradient() is valid only for
8309b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8319b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
832a86d99e1SLois Curfman McInnes 
833a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
834a86d99e1SLois Curfman McInnes 
835a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
836a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8379b94acceSBarry Smith @*/
8389b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8399b94acceSBarry Smith {
8409b94acceSBarry Smith   int    ierr;
84174679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
842e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8439b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
8449b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
8459b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8469b94acceSBarry Smith   return 0;
8479b94acceSBarry Smith }
8489b94acceSBarry Smith 
8495615d1e5SSatish Balay #undef __FUNC__
8505615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
85162fef451SLois Curfman McInnes /*@
85262fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
85362fef451SLois Curfman McInnes    set with SNESSetJacobian().
85462fef451SLois Curfman McInnes 
85562fef451SLois Curfman McInnes    Input Parameters:
85662fef451SLois Curfman McInnes .  snes - the SNES context
85762fef451SLois Curfman McInnes .  x - input vector
85862fef451SLois Curfman McInnes 
85962fef451SLois Curfman McInnes    Output Parameters:
86062fef451SLois Curfman McInnes .  A - Jacobian matrix
86162fef451SLois Curfman McInnes .  B - optional preconditioning matrix
86262fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
86362fef451SLois Curfman McInnes 
86462fef451SLois Curfman McInnes    Notes:
86562fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
86662fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
86762fef451SLois Curfman McInnes 
868dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
869dc5a77f8SLois Curfman McInnes    flag parameter.
87062fef451SLois Curfman McInnes 
87162fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
87262fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
87362fef451SLois Curfman McInnes    methods is SNESComputeJacobian().
87462fef451SLois Curfman McInnes 
87562fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
87662fef451SLois Curfman McInnes 
87762fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
87862fef451SLois Curfman McInnes @*/
8799b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8809b94acceSBarry Smith {
8819b94acceSBarry Smith   int    ierr;
88274679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
883e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
8849b94acceSBarry Smith   if (!snes->computejacobian) return 0;
8859b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
886c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
8879b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
8889b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8896d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
89077c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
89177c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8929b94acceSBarry Smith   return 0;
8939b94acceSBarry Smith }
8949b94acceSBarry Smith 
8955615d1e5SSatish Balay #undef __FUNC__
8965615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
89762fef451SLois Curfman McInnes /*@
89862fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
89962fef451SLois Curfman McInnes    set with SNESSetHessian().
90062fef451SLois Curfman McInnes 
90162fef451SLois Curfman McInnes    Input Parameters:
90262fef451SLois Curfman McInnes .  snes - the SNES context
90362fef451SLois Curfman McInnes .  x - input vector
90462fef451SLois Curfman McInnes 
90562fef451SLois Curfman McInnes    Output Parameters:
90662fef451SLois Curfman McInnes .  A - Hessian matrix
90762fef451SLois Curfman McInnes .  B - optional preconditioning matrix
90862fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
90962fef451SLois Curfman McInnes 
91062fef451SLois Curfman McInnes    Notes:
91162fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
91262fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
91362fef451SLois Curfman McInnes 
914dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
915dc5a77f8SLois Curfman McInnes    flag parameter.
91662fef451SLois Curfman McInnes 
91762fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
91862fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
91962fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
92062fef451SLois Curfman McInnes 
92162fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
92262fef451SLois Curfman McInnes 
923a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
924a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
92562fef451SLois Curfman McInnes @*/
92662fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
9279b94acceSBarry Smith {
9289b94acceSBarry Smith   int    ierr;
92974679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
930e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9319b94acceSBarry Smith   if (!snes->computejacobian) return 0;
93262fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9330f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
93462fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
93562fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9360f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
93777c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
93877c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9399b94acceSBarry Smith   return 0;
9409b94acceSBarry Smith }
9419b94acceSBarry Smith 
9425615d1e5SSatish Balay #undef __FUNC__
9435615d1e5SSatish Balay #define __FUNC__ "SNESSetJacobian"
9449b94acceSBarry Smith /*@C
9459b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
946044dda88SLois Curfman McInnes    location to store the matrix.
9479b94acceSBarry Smith 
9489b94acceSBarry Smith    Input Parameters:
9499b94acceSBarry Smith .  snes - the SNES context
9509b94acceSBarry Smith .  A - Jacobian matrix
9519b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9529b94acceSBarry Smith .  func - Jacobian evaluation routine
9532cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9542cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9559b94acceSBarry Smith 
9569b94acceSBarry Smith    Calling sequence of func:
957313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9589b94acceSBarry Smith 
9599b94acceSBarry Smith .  x - input vector
9609b94acceSBarry Smith .  A - Jacobian matrix
9619b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
962ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
963ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9642cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
9659b94acceSBarry Smith 
9669b94acceSBarry Smith    Notes:
967dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9682cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
969ac21db08SLois Curfman McInnes 
970ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9719b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9729b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9739b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9749b94acceSBarry Smith    throughout the global iterations.
9759b94acceSBarry Smith 
9769b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9779b94acceSBarry Smith 
978ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
9799b94acceSBarry Smith @*/
9809b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9819b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9829b94acceSBarry Smith {
98377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
98474679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
985e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
9869b94acceSBarry Smith   snes->computejacobian = func;
9879b94acceSBarry Smith   snes->jacP            = ctx;
9889b94acceSBarry Smith   snes->jacobian        = A;
9899b94acceSBarry Smith   snes->jacobian_pre    = B;
9909b94acceSBarry Smith   return 0;
9919b94acceSBarry Smith }
99262fef451SLois Curfman McInnes 
9935615d1e5SSatish Balay #undef __FUNC__
9945615d1e5SSatish Balay #define __FUNC__ "SNESGetJacobian"
995b4fd4287SBarry Smith /*@
996b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
997b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
998b4fd4287SBarry Smith 
999b4fd4287SBarry Smith    Input Parameter:
1000b4fd4287SBarry Smith .  snes - the nonlinear solver context
1001b4fd4287SBarry Smith 
1002b4fd4287SBarry Smith    Output Parameters:
1003b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
1004b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
1005b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1006b4fd4287SBarry Smith 
1007b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1008b4fd4287SBarry Smith @*/
1009b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
1010b4fd4287SBarry Smith {
101174679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
1012e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
1013b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
1014b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
1015b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
1016b4fd4287SBarry Smith   return 0;
1017b4fd4287SBarry Smith }
1018b4fd4287SBarry Smith 
10195615d1e5SSatish Balay #undef __FUNC__
10205615d1e5SSatish Balay #define __FUNC__ "SNESSetHessian"
10219b94acceSBarry Smith /*@C
10229b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1023044dda88SLois Curfman McInnes    location to store the matrix.
10249b94acceSBarry Smith 
10259b94acceSBarry Smith    Input Parameters:
10269b94acceSBarry Smith .  snes - the SNES context
10279b94acceSBarry Smith .  A - Hessian matrix
10289b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
10299b94acceSBarry Smith .  func - Jacobian evaluation routine
1030313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1031313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10329b94acceSBarry Smith 
10339b94acceSBarry Smith    Calling sequence of func:
1034313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10359b94acceSBarry Smith 
10369b94acceSBarry Smith .  x - input vector
10379b94acceSBarry Smith .  A - Hessian matrix
10389b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1039ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1040ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10412cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10429b94acceSBarry Smith 
10439b94acceSBarry Smith    Notes:
1044dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10452cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1046ac21db08SLois Curfman McInnes 
10479b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10489b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10499b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10509b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10519b94acceSBarry Smith    throughout the global iterations.
10529b94acceSBarry Smith 
10539b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10549b94acceSBarry Smith 
1055ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10569b94acceSBarry Smith @*/
10579b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10589b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10599b94acceSBarry Smith {
106077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
106174679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
1062e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10639b94acceSBarry Smith   snes->computejacobian = func;
10649b94acceSBarry Smith   snes->jacP            = ctx;
10659b94acceSBarry Smith   snes->jacobian        = A;
10669b94acceSBarry Smith   snes->jacobian_pre    = B;
10679b94acceSBarry Smith   return 0;
10689b94acceSBarry Smith }
10699b94acceSBarry Smith 
10705615d1e5SSatish Balay #undef __FUNC__
10715615d1e5SSatish Balay #define __FUNC__ "SNESGetHessian"
107262fef451SLois Curfman McInnes /*@
107362fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
107462fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
107562fef451SLois Curfman McInnes 
107662fef451SLois Curfman McInnes    Input Parameter:
107762fef451SLois Curfman McInnes .  snes - the nonlinear solver context
107862fef451SLois Curfman McInnes 
107962fef451SLois Curfman McInnes    Output Parameters:
108062fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
108162fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
108262fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
108362fef451SLois Curfman McInnes 
108462fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
108562fef451SLois Curfman McInnes @*/
108662fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
108762fef451SLois Curfman McInnes {
108874679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
1089e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
109062fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
109162fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
109262fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
109362fef451SLois Curfman McInnes   return 0;
109462fef451SLois Curfman McInnes }
109562fef451SLois Curfman McInnes 
10969b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10979b94acceSBarry Smith 
10985615d1e5SSatish Balay #undef __FUNC__
10995615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
11009b94acceSBarry Smith /*@
11019b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1102272ac6f2SLois Curfman McInnes    of a nonlinear solver.
11039b94acceSBarry Smith 
11049b94acceSBarry Smith    Input Parameter:
11059b94acceSBarry Smith .  snes - the SNES context
11068ddd3da0SLois Curfman McInnes .  x - the solution vector
11079b94acceSBarry Smith 
1108272ac6f2SLois Curfman McInnes    Notes:
1109272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1110272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1111272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1112272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1113272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1114272ac6f2SLois Curfman McInnes 
11159b94acceSBarry Smith .keywords: SNES, nonlinear, setup
11169b94acceSBarry Smith 
11179b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
11189b94acceSBarry Smith @*/
11198ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
11209b94acceSBarry Smith {
1121272ac6f2SLois Curfman McInnes   int ierr, flg;
112277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
112377c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
11248ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
11259b94acceSBarry Smith 
1126c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1127c1f60f51SBarry Smith   /*
1128c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1129dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1130c1f60f51SBarry Smith   */
1131c1f60f51SBarry Smith   if (flg) {
1132c1f60f51SBarry Smith     Mat J;
1133c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1134c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1135c1f60f51SBarry Smith     snes->mfshell = J;
1136c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1137c1f60f51SBarry Smith       snes->jacobian = J;
1138c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1139c1f60f51SBarry Smith     }
1140c1f60f51SBarry Smith     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1141c1f60f51SBarry Smith       snes->jacobian = J;
1142c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1143e3372554SBarry Smith     } else SETERRQ(1,0,"Method class doesn't support matrix-free operator option");
1144c1f60f51SBarry Smith   }
1145272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1146c1f60f51SBarry Smith   /*
1147dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1148c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1149c1f60f51SBarry Smith    */
115031615d04SLois Curfman McInnes   if (flg) {
1151272ac6f2SLois Curfman McInnes     Mat J;
1152272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1153272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1154272ac6f2SLois Curfman McInnes     snes->mfshell = J;
115583e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
115683e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
115794a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
115883e56afdSLois Curfman McInnes     }
115983e56afdSLois Curfman McInnes     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
116083e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
116194a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1162e3372554SBarry Smith     } else SETERRQ(1,0,"Method class doesn't support matrix-free option");
1163272ac6f2SLois Curfman McInnes   }
11649b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1165e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first");
1166e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first");
1167e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first");
1168e3372554SBarry Smith     if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector");
1169a703fe33SLois Curfman McInnes 
1170a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
117140191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1172a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1173a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1174a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1175a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1176a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1177a703fe33SLois Curfman McInnes     }
11789b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1179e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first");
1180e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first");
118137fcc0dbSBarry Smith     if (!snes->computeumfunction)
1182e3372554SBarry Smith       SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first");
1183e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first");
1184e3372554SBarry Smith   } else SETERRQ(1,0,"Unknown method class");
1185a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1186a703fe33SLois Curfman McInnes   snes->setup_called = 1;
1187a703fe33SLois Curfman McInnes   return 0;
11889b94acceSBarry Smith }
11899b94acceSBarry Smith 
11905615d1e5SSatish Balay #undef __FUNC__
11915615d1e5SSatish Balay #define __FUNC__ "SNESDestroy"
11929b94acceSBarry Smith /*@C
11939b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11949b94acceSBarry Smith    with SNESCreate().
11959b94acceSBarry Smith 
11969b94acceSBarry Smith    Input Parameter:
11979b94acceSBarry Smith .  snes - the SNES context
11989b94acceSBarry Smith 
11999b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
12009b94acceSBarry Smith 
120163a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
12029b94acceSBarry Smith @*/
12039b94acceSBarry Smith int SNESDestroy(SNES snes)
12049b94acceSBarry Smith {
12059b94acceSBarry Smith   int ierr;
120677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12079b94acceSBarry Smith   ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);
12080452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
12099b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
12109b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
12113e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
12126f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
12139b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
12140452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
12159b94acceSBarry Smith   return 0;
12169b94acceSBarry Smith }
12179b94acceSBarry Smith 
12189b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12199b94acceSBarry Smith 
12205615d1e5SSatish Balay #undef __FUNC__
12215615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
12229b94acceSBarry Smith /*@
1223d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12249b94acceSBarry Smith 
12259b94acceSBarry Smith    Input Parameters:
12269b94acceSBarry Smith .  snes - the SNES context
122733174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
122833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
122933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
123033174efeSLois Curfman McInnes            of the change in the solution between steps
123133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
123233174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
12339b94acceSBarry Smith 
123433174efeSLois Curfman McInnes    Options Database Keys:
123533174efeSLois Curfman McInnes $    -snes_atol <atol>
123633174efeSLois Curfman McInnes $    -snes_rtol <rtol>
123733174efeSLois Curfman McInnes $    -snes_stol <stol>
123833174efeSLois Curfman McInnes $    -snes_max_it <maxit>
123933174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
12409b94acceSBarry Smith 
1241d7a720efSLois Curfman McInnes    Notes:
12429b94acceSBarry Smith    The default maximum number of iterations is 50.
12439b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
12449b94acceSBarry Smith 
124533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
12469b94acceSBarry Smith 
1247d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
12489b94acceSBarry Smith @*/
1249d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
12509b94acceSBarry Smith {
125177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1252d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1253d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1254d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1255d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1256d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
12579b94acceSBarry Smith   return 0;
12589b94acceSBarry Smith }
12599b94acceSBarry Smith 
12605615d1e5SSatish Balay #undef __FUNC__
12615615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
12629b94acceSBarry Smith /*@
126333174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
126433174efeSLois Curfman McInnes 
126533174efeSLois Curfman McInnes    Input Parameters:
126633174efeSLois Curfman McInnes .  snes - the SNES context
126733174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
126833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
126933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
127033174efeSLois Curfman McInnes            of the change in the solution between steps
127133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
127233174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
127333174efeSLois Curfman McInnes 
127433174efeSLois Curfman McInnes    Notes:
127533174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
127633174efeSLois Curfman McInnes 
127733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
127833174efeSLois Curfman McInnes 
127933174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
128033174efeSLois Curfman McInnes @*/
128133174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
128233174efeSLois Curfman McInnes {
128333174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
128433174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
128533174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
128633174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
128733174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
128833174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
128933174efeSLois Curfman McInnes   return 0;
129033174efeSLois Curfman McInnes }
129133174efeSLois Curfman McInnes 
12925615d1e5SSatish Balay #undef __FUNC__
12935615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
129433174efeSLois Curfman McInnes /*@
12959b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12969b94acceSBarry Smith 
12979b94acceSBarry Smith    Input Parameters:
12989b94acceSBarry Smith .  snes - the SNES context
12999b94acceSBarry Smith .  tol - tolerance
13009b94acceSBarry Smith 
13019b94acceSBarry Smith    Options Database Key:
1302d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
13039b94acceSBarry Smith 
13049b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
13059b94acceSBarry Smith 
1306d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
13079b94acceSBarry Smith @*/
13089b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
13099b94acceSBarry Smith {
131077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13119b94acceSBarry Smith   snes->deltatol = tol;
13129b94acceSBarry Smith   return 0;
13139b94acceSBarry Smith }
13149b94acceSBarry Smith 
13155615d1e5SSatish Balay #undef __FUNC__
13165615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
13179b94acceSBarry Smith /*@
131877c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
13199b94acceSBarry Smith    for unconstrained minimization solvers.
13209b94acceSBarry Smith 
13219b94acceSBarry Smith    Input Parameters:
13229b94acceSBarry Smith .  snes - the SNES context
13239b94acceSBarry Smith .  ftol - minimum function tolerance
13249b94acceSBarry Smith 
13259b94acceSBarry Smith    Options Database Key:
1326d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
13279b94acceSBarry Smith 
13289b94acceSBarry Smith    Note:
132977c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
13309b94acceSBarry Smith    methods only.
13319b94acceSBarry Smith 
13329b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
13339b94acceSBarry Smith 
1334d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
13359b94acceSBarry Smith @*/
133677c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
13379b94acceSBarry Smith {
133877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13399b94acceSBarry Smith   snes->fmin = ftol;
13409b94acceSBarry Smith   return 0;
13419b94acceSBarry Smith }
13429b94acceSBarry Smith 
13439b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
13449b94acceSBarry Smith 
13455615d1e5SSatish Balay #undef __FUNC__
13465615d1e5SSatish Balay #define __FUNC__ "SNESSetMonitor"
13479b94acceSBarry Smith /*@C
13489b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
13499b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
13509b94acceSBarry Smith    progress.
13519b94acceSBarry Smith 
13529b94acceSBarry Smith    Input Parameters:
13539b94acceSBarry Smith .  snes - the SNES context
13549b94acceSBarry Smith .  func - monitoring routine
1355044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1356044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
13579b94acceSBarry Smith 
13589b94acceSBarry Smith    Calling sequence of func:
1359682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
13609b94acceSBarry Smith 
13619b94acceSBarry Smith $    snes - the SNES context
13629b94acceSBarry Smith $    its - iteration number
1363044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
13649b94acceSBarry Smith $
13659b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13669b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
13679b94acceSBarry Smith $
13689b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13699b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
13709b94acceSBarry Smith 
13719665c990SLois Curfman McInnes    Options Database Keys:
13729665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
13739665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
13749665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
13759665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
13769665c990SLois Curfman McInnes $                           been hardwired into a code by
13779665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
13789665c990SLois Curfman McInnes $                           does not cancel those set via
13799665c990SLois Curfman McInnes $                           the options database.
13809665c990SLois Curfman McInnes 
13819665c990SLois Curfman McInnes 
1382639f9d9dSBarry Smith    Notes:
13836bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13846bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13856bc08f3fSLois Curfman McInnes    order in which they were set.
1386639f9d9dSBarry Smith 
13879b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13889b94acceSBarry Smith 
13899b94acceSBarry Smith .seealso: SNESDefaultMonitor()
13909b94acceSBarry Smith @*/
139174679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
13929b94acceSBarry Smith {
1393639f9d9dSBarry Smith   if (!func) {
1394639f9d9dSBarry Smith     snes->numbermonitors = 0;
1395639f9d9dSBarry Smith     return 0;
1396639f9d9dSBarry Smith   }
1397639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1398e3372554SBarry Smith     SETERRQ(1,0,"Too many monitors set");
1399639f9d9dSBarry Smith   }
1400639f9d9dSBarry Smith 
1401639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1402639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
14039b94acceSBarry Smith   return 0;
14049b94acceSBarry Smith }
14059b94acceSBarry Smith 
14065615d1e5SSatish Balay #undef __FUNC__
14075615d1e5SSatish Balay #define __FUNC__ "SNESSetConvergenceTest"
14089b94acceSBarry Smith /*@C
14099b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
14109b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
14119b94acceSBarry Smith 
14129b94acceSBarry Smith    Input Parameters:
14139b94acceSBarry Smith .  snes - the SNES context
14149b94acceSBarry Smith .  func - routine to test for convergence
1415044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1416044dda88SLois Curfman McInnes           (may be PETSC_NULL)
14179b94acceSBarry Smith 
14189b94acceSBarry Smith    Calling sequence of func:
14199b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
14209b94acceSBarry Smith              double f,void *cctx)
14219b94acceSBarry Smith 
14229b94acceSBarry Smith $    snes - the SNES context
1423044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
14249b94acceSBarry Smith $    xnorm - 2-norm of current iterate
14259b94acceSBarry Smith $
14269b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14279b94acceSBarry Smith $    gnorm - 2-norm of current step
14289b94acceSBarry Smith $    f - 2-norm of function
14299b94acceSBarry Smith $
14309b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14319b94acceSBarry Smith $    gnorm - 2-norm of current gradient
14329b94acceSBarry Smith $    f - function value
14339b94acceSBarry Smith 
14349b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14359b94acceSBarry Smith 
143640191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
143740191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
14389b94acceSBarry Smith @*/
143974679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
14409b94acceSBarry Smith {
14419b94acceSBarry Smith   (snes)->converged = func;
14429b94acceSBarry Smith   (snes)->cnvP      = cctx;
14439b94acceSBarry Smith   return 0;
14449b94acceSBarry Smith }
14459b94acceSBarry Smith 
14465615d1e5SSatish Balay #undef __FUNC__
1447*c9005455SLois Curfman McInnes #define __FUNC__ "SNESSetConvergenceHistory"
1448*c9005455SLois Curfman McInnes /*@
1449*c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1450*c9005455SLois Curfman McInnes 
1451*c9005455SLois Curfman McInnes    Input Parameters:
1452*c9005455SLois Curfman McInnes .  snes - iterative context obtained from SNESCreate()
1453*c9005455SLois Curfman McInnes .  a   - array to hold history
1454*c9005455SLois Curfman McInnes .  na  - size of a
1455*c9005455SLois Curfman McInnes 
1456*c9005455SLois Curfman McInnes    Notes:
1457*c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1458*c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1459*c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1460*c9005455SLois Curfman McInnes    at each step.
1461*c9005455SLois Curfman McInnes 
1462*c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1463*c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1464*c9005455SLois Curfman McInnes    during the section of code that is being timed.
1465*c9005455SLois Curfman McInnes 
1466*c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1467*c9005455SLois Curfman McInnes @*/
1468*c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1469*c9005455SLois Curfman McInnes {
1470*c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1471*c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1472*c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1473*c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
1474*c9005455SLois Curfman McInnes   return 0;
1475*c9005455SLois Curfman McInnes }
1476*c9005455SLois Curfman McInnes 
1477*c9005455SLois Curfman McInnes #undef __FUNC__
14785615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
14799b94acceSBarry Smith /*
14809b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
14819b94acceSBarry Smith    positive parameter delta.
14829b94acceSBarry Smith 
14839b94acceSBarry Smith     Input Parameters:
14849b94acceSBarry Smith .   snes - the SNES context
14859b94acceSBarry Smith .   y - approximate solution of linear system
14869b94acceSBarry Smith .   fnorm - 2-norm of current function
14879b94acceSBarry Smith .   delta - trust region size
14889b94acceSBarry Smith 
14899b94acceSBarry Smith     Output Parameters:
14909b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
14919b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
14929b94acceSBarry Smith     region, and exceeds zero otherwise.
14939b94acceSBarry Smith .   ynorm - 2-norm of the step
14949b94acceSBarry Smith 
14959b94acceSBarry Smith     Note:
149640191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
14979b94acceSBarry Smith     is set to be the maximum allowable step size.
14989b94acceSBarry Smith 
14999b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
15009b94acceSBarry Smith */
15019b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
15029b94acceSBarry Smith                           double *gpnorm,double *ynorm)
15039b94acceSBarry Smith {
15049b94acceSBarry Smith   double norm;
15059b94acceSBarry Smith   Scalar cnorm;
1506cddf8d76SBarry Smith   VecNorm(y,NORM_2, &norm );
15079b94acceSBarry Smith   if (norm > *delta) {
15089b94acceSBarry Smith      norm = *delta/norm;
15099b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
15109b94acceSBarry Smith      cnorm = norm;
15119b94acceSBarry Smith      VecScale( &cnorm, y );
15129b94acceSBarry Smith      *ynorm = *delta;
15139b94acceSBarry Smith   } else {
15149b94acceSBarry Smith      *gpnorm = 0.0;
15159b94acceSBarry Smith      *ynorm = norm;
15169b94acceSBarry Smith   }
15179b94acceSBarry Smith   return 0;
15189b94acceSBarry Smith }
15199b94acceSBarry Smith 
15205615d1e5SSatish Balay #undef __FUNC__
15215615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
15229b94acceSBarry Smith /*@
15239b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
152463a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
15259b94acceSBarry Smith 
15269b94acceSBarry Smith    Input Parameter:
15279b94acceSBarry Smith .  snes - the SNES context
15288ddd3da0SLois Curfman McInnes .  x - the solution vector
15299b94acceSBarry Smith 
15309b94acceSBarry Smith    Output Parameter:
15319b94acceSBarry Smith    its - number of iterations until termination
15329b94acceSBarry Smith 
15338ddd3da0SLois Curfman McInnes    Note:
15348ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
15358ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
15368ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
15378ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
15388ddd3da0SLois Curfman McInnes 
15399b94acceSBarry Smith .keywords: SNES, nonlinear, solve
15409b94acceSBarry Smith 
154163a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
15429b94acceSBarry Smith @*/
15438ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
15449b94acceSBarry Smith {
15453c7409f5SSatish Balay   int ierr, flg;
1546052efed2SBarry Smith 
154777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
154874679c65SBarry Smith   PetscValidIntPointer(its);
1549c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1550c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
15519b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1552c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
15539b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
15549b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
15553c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
15566d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
15579b94acceSBarry Smith   return 0;
15589b94acceSBarry Smith }
15599b94acceSBarry Smith 
15609b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
156137fcc0dbSBarry Smith static NRList *__SNESList = 0;
15629b94acceSBarry Smith 
15635615d1e5SSatish Balay #undef __FUNC__
15645615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
15659b94acceSBarry Smith /*@
15664b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
15679b94acceSBarry Smith 
15689b94acceSBarry Smith    Input Parameters:
15699b94acceSBarry Smith .  snes - the SNES context
15709b94acceSBarry Smith .  method - a known method
15719b94acceSBarry Smith 
1572ae12b187SLois Curfman McInnes   Options Database Command:
1573ae12b187SLois Curfman McInnes $ -snes_type  <method>
1574ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1575ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1576ae12b187SLois Curfman McInnes 
15779b94acceSBarry Smith    Notes:
15789b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
15799b94acceSBarry Smith $  Systems of nonlinear equations:
158040191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
158140191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
15829b94acceSBarry Smith $  Unconstrained minimization:
158340191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
158440191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
15859b94acceSBarry Smith 
1586ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1587ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1588ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1589ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1590ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1591ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1592ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1593ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1594ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1595ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1596a703fe33SLois Curfman McInnes 
1597f536c699SSatish Balay .keywords: SNES, set, method
15989b94acceSBarry Smith @*/
15994b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
16009b94acceSBarry Smith {
16019b94acceSBarry Smith   int (*r)(SNES);
160277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16037d1a2b2bSBarry Smith   if (snes->type == (int) method) return 0;
16047d1a2b2bSBarry Smith 
16059b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
160637fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
1607e3372554SBarry Smith   if (!__SNESList) {SETERRQ(1,0,"Could not get methods");}
160837fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1609e3372554SBarry Smith   if (!r) {SETERRQ(1,0,"Unknown method");}
16100452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
16119b94acceSBarry Smith   snes->set_method_called = 1;
16129b94acceSBarry Smith   return (*r)(snes);
16139b94acceSBarry Smith }
16149b94acceSBarry Smith 
16159b94acceSBarry Smith /* --------------------------------------------------------------------- */
16165615d1e5SSatish Balay #undef __FUNC__
16175615d1e5SSatish Balay #define __FUNC__ "SNESRegister"
16189b94acceSBarry Smith /*@C
16199b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
16204b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
16219b94acceSBarry Smith 
16229b94acceSBarry Smith    Input Parameters:
162340191667SLois Curfman McInnes .  name - for instance SNES_EQ_LS, SNES_EQ_TR, ...
162440191667SLois Curfman McInnes .  sname - corresponding string for name
16259b94acceSBarry Smith .  create - routine to create method context
16269b94acceSBarry Smith 
16279b94acceSBarry Smith .keywords: SNES, nonlinear, register
16289b94acceSBarry Smith 
16299b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
16309b94acceSBarry Smith @*/
16319b94acceSBarry Smith int SNESRegister(int name, char *sname, int (*create)(SNES))
16329b94acceSBarry Smith {
16339b94acceSBarry Smith   int ierr;
163437fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
163537fcc0dbSBarry Smith   NRRegister( __SNESList, name, sname, (int (*)(void*))create );
16369b94acceSBarry Smith   return 0;
16379b94acceSBarry Smith }
1638a847f771SSatish Balay 
16399b94acceSBarry Smith /* --------------------------------------------------------------------- */
16405615d1e5SSatish Balay #undef __FUNC__
16415615d1e5SSatish Balay #define __FUNC__ "SNESRegisterDestroy"
16429b94acceSBarry Smith /*@C
16439b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
16449b94acceSBarry Smith    registered by SNESRegister().
16459b94acceSBarry Smith 
16469b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
16479b94acceSBarry Smith 
16489b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
16499b94acceSBarry Smith @*/
16509b94acceSBarry Smith int SNESRegisterDestroy()
16519b94acceSBarry Smith {
165237fcc0dbSBarry Smith   if (__SNESList) {
165337fcc0dbSBarry Smith     NRDestroy( __SNESList );
165437fcc0dbSBarry Smith     __SNESList = 0;
16559b94acceSBarry Smith   }
16569b94acceSBarry Smith   return 0;
16579b94acceSBarry Smith }
16589b94acceSBarry Smith 
16595615d1e5SSatish Balay #undef __FUNC__
16605615d1e5SSatish Balay #define __FUNC__ "SNESGetTypeFromOptions_Private"
16619b94acceSBarry Smith /*
16624b0e389bSBarry Smith    SNESGetTypeFromOptions_Private - Sets the selected method from the
16639b94acceSBarry Smith    options database.
16649b94acceSBarry Smith 
16659b94acceSBarry Smith    Input Parameter:
16669b94acceSBarry Smith .  ctx - the SNES context
16679b94acceSBarry Smith 
16689b94acceSBarry Smith    Output Parameter:
16699b94acceSBarry Smith .  method -  solver method
16709b94acceSBarry Smith 
16719b94acceSBarry Smith    Returns:
16729b94acceSBarry Smith    Returns 1 if the method is found; 0 otherwise.
16739b94acceSBarry Smith 
16749b94acceSBarry Smith    Options Database Key:
16754b0e389bSBarry Smith $  -snes_type  method
16769b94acceSBarry Smith */
1677052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
16789b94acceSBarry Smith {
1679052efed2SBarry Smith   int  ierr;
16809b94acceSBarry Smith   char sbuf[50];
16815baf8537SBarry Smith 
1682052efed2SBarry Smith   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1683052efed2SBarry Smith   if (*flg) {
1684052efed2SBarry Smith     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
168537fcc0dbSBarry Smith     *method = (SNESType)NRFindID( __SNESList, sbuf );
16869b94acceSBarry Smith   }
16879b94acceSBarry Smith   return 0;
16889b94acceSBarry Smith }
16899b94acceSBarry Smith 
16905615d1e5SSatish Balay #undef __FUNC__
16915615d1e5SSatish Balay #define __FUNC__ "SNESGetType"
16929b94acceSBarry Smith /*@C
16939a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
16949b94acceSBarry Smith 
16959b94acceSBarry Smith    Input Parameter:
16964b0e389bSBarry Smith .  snes - nonlinear solver context
16979b94acceSBarry Smith 
16989b94acceSBarry Smith    Output Parameter:
16999a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
17009a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
17019b94acceSBarry Smith 
17029b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
17039b94acceSBarry Smith @*/
17044b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
17059b94acceSBarry Smith {
170606a9b0adSLois Curfman McInnes   int ierr;
170737fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
17084b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
170937fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
17109b94acceSBarry Smith   return 0;
17119b94acceSBarry Smith }
17129b94acceSBarry Smith 
17139b94acceSBarry Smith #include <stdio.h>
17145615d1e5SSatish Balay #undef __FUNC__
17155615d1e5SSatish Balay #define __FUNC__ "SNESPrintTypes_Private"
17169b94acceSBarry Smith /*
17174b0e389bSBarry Smith    SNESPrintTypes_Private - Prints the SNES methods available from the
17189b94acceSBarry Smith    options database.
17199b94acceSBarry Smith 
17209b94acceSBarry Smith    Input Parameters:
172133455573SSatish Balay .  comm   - communicator (usually MPI_COMM_WORLD)
17229b94acceSBarry Smith .  prefix - prefix (usually "-")
17234b0e389bSBarry Smith .  name   - the options database name (by default "snes_type")
17249b94acceSBarry Smith */
172533455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
17269b94acceSBarry Smith {
17279b94acceSBarry Smith   FuncList *entry;
172837fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
172937fcc0dbSBarry Smith   entry = __SNESList->head;
173077c4ece6SBarry Smith   PetscPrintf(comm," %s%s (one of)",prefix,name);
17319b94acceSBarry Smith   while (entry) {
173277c4ece6SBarry Smith     PetscPrintf(comm," %s",entry->name);
17339b94acceSBarry Smith     entry = entry->next;
17349b94acceSBarry Smith   }
173577c4ece6SBarry Smith   PetscPrintf(comm,"\n");
17369b94acceSBarry Smith   return 0;
17379b94acceSBarry Smith }
17389b94acceSBarry Smith 
17395615d1e5SSatish Balay #undef __FUNC__
17405615d1e5SSatish Balay #define __FUNC__ "SNESGetSolution"
17419b94acceSBarry Smith /*@C
17429b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
17439b94acceSBarry Smith    stored.
17449b94acceSBarry Smith 
17459b94acceSBarry Smith    Input Parameter:
17469b94acceSBarry Smith .  snes - the SNES context
17479b94acceSBarry Smith 
17489b94acceSBarry Smith    Output Parameter:
17499b94acceSBarry Smith .  x - the solution
17509b94acceSBarry Smith 
17519b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
17529b94acceSBarry Smith 
1753a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
17549b94acceSBarry Smith @*/
17559b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
17569b94acceSBarry Smith {
175777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17589b94acceSBarry Smith   *x = snes->vec_sol_always;
17599b94acceSBarry Smith   return 0;
17609b94acceSBarry Smith }
17619b94acceSBarry Smith 
17625615d1e5SSatish Balay #undef __FUNC__
17635615d1e5SSatish Balay #define __FUNC__ "SNESGetSolutionUpdate"
17649b94acceSBarry Smith /*@C
17659b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
17669b94acceSBarry Smith    stored.
17679b94acceSBarry Smith 
17689b94acceSBarry Smith    Input Parameter:
17699b94acceSBarry Smith .  snes - the SNES context
17709b94acceSBarry Smith 
17719b94acceSBarry Smith    Output Parameter:
17729b94acceSBarry Smith .  x - the solution update
17739b94acceSBarry Smith 
17749b94acceSBarry Smith    Notes:
17759b94acceSBarry Smith    This vector is implementation dependent.
17769b94acceSBarry Smith 
17779b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
17789b94acceSBarry Smith 
17799b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
17809b94acceSBarry Smith @*/
17819b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
17829b94acceSBarry Smith {
178377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17849b94acceSBarry Smith   *x = snes->vec_sol_update_always;
17859b94acceSBarry Smith   return 0;
17869b94acceSBarry Smith }
17879b94acceSBarry Smith 
17885615d1e5SSatish Balay #undef __FUNC__
17895615d1e5SSatish Balay #define __FUNC__ "SNESGetFunction"
17909b94acceSBarry Smith /*@C
17913638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
17929b94acceSBarry Smith 
17939b94acceSBarry Smith    Input Parameter:
17949b94acceSBarry Smith .  snes - the SNES context
17959b94acceSBarry Smith 
17969b94acceSBarry Smith    Output Parameter:
17973638b69dSLois Curfman McInnes .  r - the function
17989b94acceSBarry Smith 
17999b94acceSBarry Smith    Notes:
18009b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
18019b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
18029b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
18039b94acceSBarry Smith 
1804a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
18059b94acceSBarry Smith 
180631615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
180731615d04SLois Curfman McInnes           SNESGetGradient()
18089b94acceSBarry Smith @*/
18099b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
18109b94acceSBarry Smith {
181177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1812e3372554SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,
18137c792091SSatish Balay     "For SNES_NONLINEAR_EQUATIONS only");
18149b94acceSBarry Smith   *r = snes->vec_func_always;
18159b94acceSBarry Smith   return 0;
18169b94acceSBarry Smith }
18179b94acceSBarry Smith 
18185615d1e5SSatish Balay #undef __FUNC__
18195615d1e5SSatish Balay #define __FUNC__ "SNESGetGradient"
18209b94acceSBarry Smith /*@C
18213638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
18229b94acceSBarry Smith 
18239b94acceSBarry Smith    Input Parameter:
18249b94acceSBarry Smith .  snes - the SNES context
18259b94acceSBarry Smith 
18269b94acceSBarry Smith    Output Parameter:
18279b94acceSBarry Smith .  r - the gradient
18289b94acceSBarry Smith 
18299b94acceSBarry Smith    Notes:
18309b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
18319b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18329b94acceSBarry Smith    SNESGetFunction().
18339b94acceSBarry Smith 
18349b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
18359b94acceSBarry Smith 
183631615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
18379b94acceSBarry Smith @*/
18389b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
18399b94acceSBarry Smith {
184077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1841e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,
184263c41f6aSSatish Balay     "For SNES_UNCONSTRAINED_MINIMIZATION only");
18439b94acceSBarry Smith   *r = snes->vec_func_always;
18449b94acceSBarry Smith   return 0;
18459b94acceSBarry Smith }
18469b94acceSBarry Smith 
18475615d1e5SSatish Balay #undef __FUNC__
18485615d1e5SSatish Balay #define __FUNC__ "SNESGetMinimizationFunction"
18499b94acceSBarry Smith /*@
18509b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
18519b94acceSBarry Smith    unconstrained minimization problems.
18529b94acceSBarry Smith 
18539b94acceSBarry Smith    Input Parameter:
18549b94acceSBarry Smith .  snes - the SNES context
18559b94acceSBarry Smith 
18569b94acceSBarry Smith    Output Parameter:
18579b94acceSBarry Smith .  r - the function
18589b94acceSBarry Smith 
18599b94acceSBarry Smith    Notes:
18609b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
18619b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18629b94acceSBarry Smith    SNESGetFunction().
18639b94acceSBarry Smith 
18649b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
18659b94acceSBarry Smith 
186631615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
18679b94acceSBarry Smith @*/
18689b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
18699b94acceSBarry Smith {
187077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
187174679c65SBarry Smith   PetscValidScalarPointer(r);
1872e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,
187363c41f6aSSatish Balay     "For SNES_UNCONSTRAINED_MINIMIZATION only");
18749b94acceSBarry Smith   *r = snes->fc;
18759b94acceSBarry Smith   return 0;
18769b94acceSBarry Smith }
18779b94acceSBarry Smith 
18785615d1e5SSatish Balay #undef __FUNC__
18795615d1e5SSatish Balay #define __FUNC__ "SNESSetOptionsPrefix"
18803c7409f5SSatish Balay /*@C
18813c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1882639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1883639f9d9dSBarry Smith    the prefix name.
18843c7409f5SSatish Balay 
18853c7409f5SSatish Balay    Input Parameter:
18863c7409f5SSatish Balay .  snes - the SNES context
18873c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18883c7409f5SSatish Balay 
18893c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1890a86d99e1SLois Curfman McInnes 
1891a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
18923c7409f5SSatish Balay @*/
18933c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
18943c7409f5SSatish Balay {
18953c7409f5SSatish Balay   int ierr;
18963c7409f5SSatish Balay 
189777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1898639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18993c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19003c7409f5SSatish Balay   return 0;
19013c7409f5SSatish Balay }
19023c7409f5SSatish Balay 
19035615d1e5SSatish Balay #undef __FUNC__
19045615d1e5SSatish Balay #define __FUNC__ "SNESAppendOptionsPrefix"
19053c7409f5SSatish Balay /*@C
1906f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1907639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1908639f9d9dSBarry Smith    the prefix name.
19093c7409f5SSatish Balay 
19103c7409f5SSatish Balay    Input Parameter:
19113c7409f5SSatish Balay .  snes - the SNES context
19123c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
19133c7409f5SSatish Balay 
19143c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1915a86d99e1SLois Curfman McInnes 
1916a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
19173c7409f5SSatish Balay @*/
19183c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
19193c7409f5SSatish Balay {
19203c7409f5SSatish Balay   int ierr;
19213c7409f5SSatish Balay 
192277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1923639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19243c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19253c7409f5SSatish Balay   return 0;
19263c7409f5SSatish Balay }
19273c7409f5SSatish Balay 
19285615d1e5SSatish Balay #undef __FUNC__
19295615d1e5SSatish Balay #define __FUNC__ "SNESGetOptionsPrefix"
1930c83590e2SSatish Balay /*@
19313c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
19323c7409f5SSatish Balay    SNES options in the database.
19333c7409f5SSatish Balay 
19343c7409f5SSatish Balay    Input Parameter:
19353c7409f5SSatish Balay .  snes - the SNES context
19363c7409f5SSatish Balay 
19373c7409f5SSatish Balay    Output Parameter:
19383c7409f5SSatish Balay .  prefix - pointer to the prefix string used
19393c7409f5SSatish Balay 
19403c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1941a86d99e1SLois Curfman McInnes 
1942a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
19433c7409f5SSatish Balay @*/
19443c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
19453c7409f5SSatish Balay {
19463c7409f5SSatish Balay   int ierr;
19473c7409f5SSatish Balay 
194877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1949639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19503c7409f5SSatish Balay   return 0;
19513c7409f5SSatish Balay }
19523c7409f5SSatish Balay 
19533c7409f5SSatish Balay 
19549b94acceSBarry Smith 
19559b94acceSBarry Smith 
19569b94acceSBarry Smith 
1957