xref: /petsc/src/snes/interface/snes.c (revision 3a40ed3dce77c081171d005ae1a6ff4bb9d13b6f)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*3a40ed3dSBarry Smith static char vcid[] = "$Id: snes.c,v 1.130 1997/09/26 02:20:57 bsmith Exp bsmith $";
39b94acceSBarry Smith #endif
49b94acceSBarry Smith 
570f55243SBarry Smith #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
6f5eb4b81SSatish Balay #include "src/sys/nreg.h"
79b94acceSBarry Smith #include "pinclude/pviewer.h"
89b94acceSBarry Smith #include <math.h>
99b94acceSBarry Smith 
10052efed2SBarry Smith extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*);
1133455573SSatish Balay extern int SNESPrintTypes_Private(MPI_Comm,char*,char*);
129b94acceSBarry Smith 
1384cb2905SBarry Smith int SNESRegisterAllCalled = 0;
1484cb2905SBarry Smith 
155615d1e5SSatish Balay #undef __FUNC__
16d4bb536fSBarry Smith #define __FUNC__ "SNESView"
179b94acceSBarry Smith /*@
189b94acceSBarry Smith    SNESView - Prints the SNES data structure.
199b94acceSBarry Smith 
209b94acceSBarry Smith    Input Parameters:
219b94acceSBarry Smith .  SNES - the SNES context
229b94acceSBarry Smith .  viewer - visualization context
239b94acceSBarry Smith 
249b94acceSBarry Smith    Options Database Key:
259b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
269b94acceSBarry Smith 
279b94acceSBarry Smith    Notes:
289b94acceSBarry Smith    The available visualization contexts include
296d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
306d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
319b94acceSBarry Smith $       output where only the first processor opens
329b94acceSBarry Smith $       the file.  All other processors send their
339b94acceSBarry Smith $       data to the first processor to print.
349b94acceSBarry Smith 
359b94acceSBarry Smith    The user can open alternative vistualization contexts with
36dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
379b94acceSBarry Smith 
389b94acceSBarry Smith .keywords: SNES, view
399b94acceSBarry Smith 
40dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
419b94acceSBarry Smith @*/
429b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
439b94acceSBarry Smith {
449b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
459b94acceSBarry Smith   FILE                *fd;
469b94acceSBarry Smith   int                 ierr;
479b94acceSBarry Smith   SLES                sles;
489b94acceSBarry Smith   char                *method;
4919bcc07fSBarry Smith   ViewerType          vtype;
509b94acceSBarry Smith 
51*3a40ed3dSBarry Smith   PetscFunctionBegin;
5274679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5374679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5474679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5574679c65SBarry Smith 
5619bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5719bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5890ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5977c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
604b0e389bSBarry Smith     SNESGetType(snes,PETSC_NULL,&method);
6177c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
629b94acceSBarry Smith     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
6377c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
649b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
659b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6677c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
679b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
689b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
697a00f4a9SLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
707a00f4a9SLois Curfman McInnes     "  total number of linear solver iterations=%d\n",snes->linear_its);
71ae3c334cSLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
7268632166SLois Curfman McInnes      "  total number of function evaluations=%d\n",snes->nfuncs);
739b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
7477c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
759b94acceSBarry Smith     if (snes->ksp_ewconv) {
769b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
779b94acceSBarry Smith       if (kctx) {
7877c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
799b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
809b94acceSBarry Smith         kctx->version);
8177c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
829b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
839b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
8477c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
859b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
869b94acceSBarry Smith       }
879b94acceSBarry Smith     }
88c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
89c30f285eSLois Curfman McInnes     SNESGetType(snes,PETSC_NULL,&method);
90c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
9119bcc07fSBarry Smith   }
929b94acceSBarry Smith   SNESGetSLES(snes,&sles);
939b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
94*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
959b94acceSBarry Smith }
969b94acceSBarry Smith 
97639f9d9dSBarry Smith /*
98639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
99639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
100639f9d9dSBarry Smith */
101639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
102639f9d9dSBarry Smith static int numberofsetfromoptions;
103639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
104639f9d9dSBarry Smith 
1055615d1e5SSatish Balay #undef __FUNC__
106d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker"
107639f9d9dSBarry Smith /*@
108639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
109639f9d9dSBarry Smith 
110639f9d9dSBarry Smith     Input Parameter:
111639f9d9dSBarry Smith .   snescheck - function that checks for options
112639f9d9dSBarry Smith 
113639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
114639f9d9dSBarry Smith @*/
115639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
116639f9d9dSBarry Smith {
117*3a40ed3dSBarry Smith   PetscFunctionBegin;
118639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
119e3372554SBarry Smith     SETERRQ(1,0,"Too many options checkers, only 5 allowed");
120639f9d9dSBarry Smith   }
121639f9d9dSBarry Smith 
122639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
123*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
124639f9d9dSBarry Smith }
125639f9d9dSBarry Smith 
1265615d1e5SSatish Balay #undef __FUNC__
1275615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1289b94acceSBarry Smith /*@
129682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1309b94acceSBarry Smith 
1319b94acceSBarry Smith    Input Parameter:
1329b94acceSBarry Smith .  snes - the SNES context
1339b94acceSBarry Smith 
13411ca99fdSLois Curfman McInnes    Notes:
13511ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
13611ca99fdSLois Curfman McInnes    the users manual.
13783e2fdc7SBarry Smith 
1389b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1399b94acceSBarry Smith 
140a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1419b94acceSBarry Smith @*/
1429b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1439b94acceSBarry Smith {
1444b0e389bSBarry Smith   SNESType method;
1459b94acceSBarry Smith   double   tmp;
1469b94acceSBarry Smith   SLES     sles;
14781f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1483c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1499b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1509b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1519b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1529b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1539b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1549b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1559b94acceSBarry Smith 
156*3a40ed3dSBarry Smith   PetscFunctionBegin;
15781f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
15881f57ec7SBarry Smith 
15981f57ec7SBarry Smith 
16077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
161e3372554SBarry Smith   if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp");
162052efed2SBarry Smith   ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr);
163052efed2SBarry Smith   if (flg) {
164052efed2SBarry Smith     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
1659b94acceSBarry Smith   }
1665334005bSBarry Smith   else if (!snes->set_method_called) {
1675334005bSBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
16840191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
1695334005bSBarry Smith     }
1705334005bSBarry Smith     else {
17140191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
1725334005bSBarry Smith     }
1735334005bSBarry Smith   }
1745334005bSBarry Smith 
1753c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1763c7409f5SSatish Balay   if (flg) { SNESPrintHelp(snes); }
1773c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
178d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);}
1793c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
180d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);}
1813c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
182d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); }
1833c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1843c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
185d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
186d7a720efSLois Curfman McInnes   if (flg) { SNESSetTrustRegionTolerance(snes,tmp); }
187d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
188d7a720efSLois Curfman McInnes   if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); }
1893c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1903c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
191b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
192b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
193b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
194b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
195b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
196b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
197b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
19893c39befSBarry Smith 
19993c39befSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr);
20093c39befSBarry Smith   if (flg) {snes->converged = 0;}
20193c39befSBarry Smith 
2029b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
2039b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
2045bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
2055bbad29bSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);}
2063c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
207639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
2083c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
209d31a3109SLois Curfman McInnes   if (flg) {
210639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
211d31a3109SLois Curfman McInnes   }
21281f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
2133c7409f5SSatish Balay   if (flg) {
21417699dbbSLois Curfman McInnes     int    rank = 0;
215d7e8b826SBarry Smith     DrawLG lg;
21617699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
21717699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
21817699dbbSLois Curfman McInnes     if (!rank) {
21981f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2209b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
221f8c826e1SBarry Smith       snes->xmonitor = lg;
2229b94acceSBarry Smith       PLogObjectParent(snes,lg);
2239b94acceSBarry Smith     }
2249b94acceSBarry Smith   }
2253c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2263c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2279b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2289b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
22994a424c1SBarry Smith     PLogInfo(snes,
23031615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
23131615d04SLois Curfman McInnes   }
23231615d04SLois Curfman McInnes   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
23331615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
23431615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
23594a424c1SBarry Smith     PLogInfo(snes,
23631615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2379b94acceSBarry Smith   }
238639f9d9dSBarry Smith 
239639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
240639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
241639f9d9dSBarry Smith   }
242639f9d9dSBarry Smith 
2439b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2449b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
245*3a40ed3dSBarry Smith   if (snes->setfromoptions) {
246*3a40ed3dSBarry Smith     ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
247*3a40ed3dSBarry Smith   }
248*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
2499b94acceSBarry Smith }
2509b94acceSBarry Smith 
2515615d1e5SSatish Balay #undef __FUNC__
252d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp"
2539b94acceSBarry Smith /*@
2549b94acceSBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2559b94acceSBarry Smith 
2569b94acceSBarry Smith    Input Parameter:
2579b94acceSBarry Smith .  snes - the SNES context
2589b94acceSBarry Smith 
259a703fe33SLois Curfman McInnes    Options Database Keys:
260a703fe33SLois Curfman McInnes $  -help, -h
261a703fe33SLois Curfman McInnes 
2629b94acceSBarry Smith .keywords: SNES, nonlinear, help
2639b94acceSBarry Smith 
264682d7d0cSBarry Smith .seealso: SNESSetFromOptions()
2659b94acceSBarry Smith @*/
2669b94acceSBarry Smith int SNESPrintHelp(SNES snes)
2679b94acceSBarry Smith {
2686daaf66cSBarry Smith   char                p[64];
2696daaf66cSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2706daaf66cSBarry Smith 
271*3a40ed3dSBarry Smith   PetscFunctionBegin;
27277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2736daaf66cSBarry Smith 
2746daaf66cSBarry Smith   PetscStrcpy(p,"-");
2756daaf66cSBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2766daaf66cSBarry Smith 
2776daaf66cSBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2786daaf66cSBarry Smith 
279d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
28033455573SSatish Balay   SNESPrintTypes_Private(snes->comm,p,"snes_type");
28177c4ece6SBarry Smith   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
282b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
283b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
284b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
285b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
286b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
287b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
2885bbad29bSBarry Smith   PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n");
2895bbad29bSBarry Smith   PetscPrintf(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
29083e2fdc7SBarry Smith   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
29183e2fdc7SBarry Smith     residual norm at each iteration.\n",p);
29283e2fdc7SBarry Smith   PetscPrintf(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
29383e2fdc7SBarry Smith     residual norm for small residual norms. This is useful to conceal\n\
29483e2fdc7SBarry Smith     meaningless digits that may be different on different machines.\n",p);
295d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
296b18e04deSLois Curfman McInnes   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
29777c4ece6SBarry Smith     PetscPrintf(snes->comm,
298d31a3109SLois Curfman McInnes      " Options for solving systems of nonlinear equations only:\n");
29977c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
30077c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
301ef1dfb11SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
30293c39befSBarry Smith     PetscPrintf(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
30377c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
30477c4ece6SBarry Smith     PetscPrintf(snes->comm,
305b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
30677c4ece6SBarry Smith     PetscPrintf(snes->comm,
307b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
30877c4ece6SBarry Smith     PetscPrintf(snes->comm,
309b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
31077c4ece6SBarry Smith     PetscPrintf(snes->comm,
311b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
31277c4ece6SBarry Smith     PetscPrintf(snes->comm,
313b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
31477c4ece6SBarry Smith     PetscPrintf(snes->comm,
315b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
31677c4ece6SBarry Smith     PetscPrintf(snes->comm,
317b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
318b18e04deSLois Curfman McInnes   }
319b18e04deSLois Curfman McInnes   else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
32077c4ece6SBarry Smith     PetscPrintf(snes->comm,
321d31a3109SLois Curfman McInnes      " Options for solving unconstrained minimization problems only:\n");
322b18e04deSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
32377c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
32477c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
325b18e04deSLois Curfman McInnes   }
3264537a946SLois Curfman McInnes   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
32777c4ece6SBarry Smith   PetscPrintf(snes->comm,"a particular method\n");
3286daaf66cSBarry Smith   if (snes->printhelp) (*snes->printhelp)(snes,p);
329*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
3309b94acceSBarry Smith }
331a847f771SSatish Balay 
3325615d1e5SSatish Balay #undef __FUNC__
333d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext"
3349b94acceSBarry Smith /*@
3359b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3369b94acceSBarry Smith    the nonlinear solvers.
3379b94acceSBarry Smith 
3389b94acceSBarry Smith    Input Parameters:
3399b94acceSBarry Smith .  snes - the SNES context
3409b94acceSBarry Smith .  usrP - optional user context
3419b94acceSBarry Smith 
3429b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3439b94acceSBarry Smith 
3449b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3459b94acceSBarry Smith @*/
3469b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3479b94acceSBarry Smith {
348*3a40ed3dSBarry Smith   PetscFunctionBegin;
34977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3509b94acceSBarry Smith   snes->user		= usrP;
351*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
3529b94acceSBarry Smith }
35374679c65SBarry Smith 
3545615d1e5SSatish Balay #undef __FUNC__
355d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext"
3569b94acceSBarry Smith /*@C
3579b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3589b94acceSBarry Smith    nonlinear solvers.
3599b94acceSBarry Smith 
3609b94acceSBarry Smith    Input Parameter:
3619b94acceSBarry Smith .  snes - SNES context
3629b94acceSBarry Smith 
3639b94acceSBarry Smith    Output Parameter:
3649b94acceSBarry Smith .  usrP - user context
3659b94acceSBarry Smith 
3669b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3679b94acceSBarry Smith 
3689b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3699b94acceSBarry Smith @*/
3709b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3719b94acceSBarry Smith {
372*3a40ed3dSBarry Smith   PetscFunctionBegin;
37377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3749b94acceSBarry Smith   *usrP = snes->user;
375*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
3769b94acceSBarry Smith }
37774679c65SBarry Smith 
3785615d1e5SSatish Balay #undef __FUNC__
379d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber"
3809b94acceSBarry Smith /*@
3819b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3829b94acceSBarry Smith    nonlinear solver.
3839b94acceSBarry Smith 
3849b94acceSBarry Smith    Input Parameter:
3859b94acceSBarry Smith .  snes - SNES context
3869b94acceSBarry Smith 
3879b94acceSBarry Smith    Output Parameter:
3889b94acceSBarry Smith .  iter - iteration number
3899b94acceSBarry Smith 
3909b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3919b94acceSBarry Smith @*/
3929b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3939b94acceSBarry Smith {
394*3a40ed3dSBarry Smith   PetscFunctionBegin;
39577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
39674679c65SBarry Smith   PetscValidIntPointer(iter);
3979b94acceSBarry Smith   *iter = snes->iter;
398*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
3999b94acceSBarry Smith }
40074679c65SBarry Smith 
4015615d1e5SSatish Balay #undef __FUNC__
4025615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
4039b94acceSBarry Smith /*@
4049b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
4059b94acceSBarry Smith    with SNESSSetFunction().
4069b94acceSBarry Smith 
4079b94acceSBarry Smith    Input Parameter:
4089b94acceSBarry Smith .  snes - SNES context
4099b94acceSBarry Smith 
4109b94acceSBarry Smith    Output Parameter:
4119b94acceSBarry Smith .  fnorm - 2-norm of function
4129b94acceSBarry Smith 
4139b94acceSBarry Smith    Note:
4149b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
415a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
416a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
4179b94acceSBarry Smith 
4189b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
419a86d99e1SLois Curfman McInnes 
420a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
4219b94acceSBarry Smith @*/
4229b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
4239b94acceSBarry Smith {
424*3a40ed3dSBarry Smith   PetscFunctionBegin;
42577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
42674679c65SBarry Smith   PetscValidScalarPointer(fnorm);
42774679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
428d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only");
42974679c65SBarry Smith   }
4309b94acceSBarry Smith   *fnorm = snes->norm;
431*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
4329b94acceSBarry Smith }
43374679c65SBarry Smith 
4345615d1e5SSatish Balay #undef __FUNC__
4355615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
4369b94acceSBarry Smith /*@
4379b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
4389b94acceSBarry Smith    with SNESSSetGradient().
4399b94acceSBarry Smith 
4409b94acceSBarry Smith    Input Parameter:
4419b94acceSBarry Smith .  snes - SNES context
4429b94acceSBarry Smith 
4439b94acceSBarry Smith    Output Parameter:
4449b94acceSBarry Smith .  fnorm - 2-norm of gradient
4459b94acceSBarry Smith 
4469b94acceSBarry Smith    Note:
4479b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
448a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
449a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4509b94acceSBarry Smith 
4519b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
452a86d99e1SLois Curfman McInnes 
453a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4549b94acceSBarry Smith @*/
4559b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
4569b94acceSBarry Smith {
457*3a40ed3dSBarry Smith   PetscFunctionBegin;
45877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
45974679c65SBarry Smith   PetscValidScalarPointer(gnorm);
46074679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
461d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
46274679c65SBarry Smith   }
4639b94acceSBarry Smith   *gnorm = snes->norm;
464*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
4659b94acceSBarry Smith }
46674679c65SBarry Smith 
4675615d1e5SSatish Balay #undef __FUNC__
468d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps"
4699b94acceSBarry Smith /*@
4709b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4719b94acceSBarry Smith    attempted by the nonlinear solver.
4729b94acceSBarry Smith 
4739b94acceSBarry Smith    Input Parameter:
4749b94acceSBarry Smith .  snes - SNES context
4759b94acceSBarry Smith 
4769b94acceSBarry Smith    Output Parameter:
4779b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4789b94acceSBarry Smith 
479c96a6f78SLois Curfman McInnes    Notes:
480c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
481c96a6f78SLois Curfman McInnes 
4829b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4839b94acceSBarry Smith @*/
4849b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4859b94acceSBarry Smith {
486*3a40ed3dSBarry Smith   PetscFunctionBegin;
48777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
48874679c65SBarry Smith   PetscValidIntPointer(nfails);
4899b94acceSBarry Smith   *nfails = snes->nfailures;
490*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
4919b94acceSBarry Smith }
492a847f771SSatish Balay 
4935615d1e5SSatish Balay #undef __FUNC__
494d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations"
495c96a6f78SLois Curfman McInnes /*@
496c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
497c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
498c96a6f78SLois Curfman McInnes 
499c96a6f78SLois Curfman McInnes    Input Parameter:
500c96a6f78SLois Curfman McInnes .  snes - SNES context
501c96a6f78SLois Curfman McInnes 
502c96a6f78SLois Curfman McInnes    Output Parameter:
503c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
504c96a6f78SLois Curfman McInnes 
505c96a6f78SLois Curfman McInnes    Notes:
506c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
507c96a6f78SLois Curfman McInnes 
508c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
509c96a6f78SLois Curfman McInnes @*/
51086bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
511c96a6f78SLois Curfman McInnes {
512*3a40ed3dSBarry Smith   PetscFunctionBegin;
513c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
514c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
515c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
516*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
517c96a6f78SLois Curfman McInnes }
518c96a6f78SLois Curfman McInnes 
519c96a6f78SLois Curfman McInnes #undef __FUNC__
520d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES"
5219b94acceSBarry Smith /*@C
5229b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
5239b94acceSBarry Smith 
5249b94acceSBarry Smith    Input Parameter:
5259b94acceSBarry Smith .  snes - the SNES context
5269b94acceSBarry Smith 
5279b94acceSBarry Smith    Output Parameter:
5289b94acceSBarry Smith .  sles - the SLES context
5299b94acceSBarry Smith 
5309b94acceSBarry Smith    Notes:
5319b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
5329b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5339b94acceSBarry Smith    KSP and PC contexts as well.
5349b94acceSBarry Smith 
5359b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
5369b94acceSBarry Smith 
5379b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
5389b94acceSBarry Smith @*/
5399b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
5409b94acceSBarry Smith {
541*3a40ed3dSBarry Smith   PetscFunctionBegin;
54277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5439b94acceSBarry Smith   *sles = snes->sles;
544*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
5459b94acceSBarry Smith }
5469b94acceSBarry Smith /* -----------------------------------------------------------*/
5475615d1e5SSatish Balay #undef __FUNC__
5485615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
5499b94acceSBarry Smith /*@C
5509b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5519b94acceSBarry Smith 
5529b94acceSBarry Smith    Input Parameter:
5539b94acceSBarry Smith .  comm - MPI communicator
5549b94acceSBarry Smith .  type - type of method, one of
5559b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
5569b94acceSBarry Smith $      (for systems of nonlinear equations)
5579b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
5589b94acceSBarry Smith $      (for unconstrained minimization)
5599b94acceSBarry Smith 
5609b94acceSBarry Smith    Output Parameter:
5619b94acceSBarry Smith .  outsnes - the new SNES context
5629b94acceSBarry Smith 
563c1f60f51SBarry Smith    Options Database Key:
56419bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
56519bd6747SLois Curfman McInnes $              and no preconditioning matrix
56619bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
56719bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
56819bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
56919bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
570c1f60f51SBarry Smith 
5719b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5729b94acceSBarry Smith 
57363a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5749b94acceSBarry Smith @*/
5754b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5769b94acceSBarry Smith {
5779b94acceSBarry Smith   int                 ierr;
5789b94acceSBarry Smith   SNES                snes;
5799b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
58037fcc0dbSBarry Smith 
581*3a40ed3dSBarry Smith   PetscFunctionBegin;
5829b94acceSBarry Smith   *outsnes = 0;
5832a0acf01SLois Curfman McInnes   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS)
584d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
585d4bb536fSBarry Smith   PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm,SNESDestroy,SNESView);
5869b94acceSBarry Smith   PLogObjectCreate(snes);
5879b94acceSBarry Smith   snes->max_its           = 50;
5889750a799SBarry Smith   snes->max_funcs	  = 10000;
5899b94acceSBarry Smith   snes->norm		  = 0.0;
5905a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
591b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
592b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5939b94acceSBarry Smith     snes->atol		  = 1.e-10;
5945a2d0531SBarry Smith   }
595b4874afaSBarry Smith   else {
596b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
597b4874afaSBarry Smith     snes->ttol            = 0.0;
598b4874afaSBarry Smith     snes->atol		  = 1.e-50;
599b4874afaSBarry Smith   }
6009b94acceSBarry Smith   snes->xtol		  = 1.e-8;
601b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
6029b94acceSBarry Smith   snes->nfuncs            = 0;
6039b94acceSBarry Smith   snes->nfailures         = 0;
6047a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
605639f9d9dSBarry Smith   snes->numbermonitors    = 0;
6069b94acceSBarry Smith   snes->data              = 0;
6079b94acceSBarry Smith   snes->view              = 0;
6089b94acceSBarry Smith   snes->computeumfunction = 0;
6099b94acceSBarry Smith   snes->umfunP            = 0;
6109b94acceSBarry Smith   snes->fc                = 0;
6119b94acceSBarry Smith   snes->deltatol          = 1.e-12;
6129b94acceSBarry Smith   snes->fmin              = -1.e30;
6139b94acceSBarry Smith   snes->method_class      = type;
6149b94acceSBarry Smith   snes->set_method_called = 0;
615a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
6169b94acceSBarry Smith   snes->ksp_ewconv        = 0;
6177d1a2b2bSBarry Smith   snes->type              = -1;
6186f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
6196f24a144SLois Curfman McInnes   snes->vwork             = 0;
6206f24a144SLois Curfman McInnes   snes->nwork             = 0;
621c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
622c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
623c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
6249b94acceSBarry Smith 
6259b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
6260452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
627eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
6289b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6299b94acceSBarry Smith   kctx->version     = 2;
6309b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6319b94acceSBarry Smith                              this was too large for some test cases */
6329b94acceSBarry Smith   kctx->rtol_last   = 0;
6339b94acceSBarry Smith   kctx->rtol_max    = .9;
6349b94acceSBarry Smith   kctx->gamma       = 1.0;
6359b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6369b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6379b94acceSBarry Smith   kctx->threshold   = .1;
6389b94acceSBarry Smith   kctx->lresid_last = 0;
6399b94acceSBarry Smith   kctx->norm_last   = 0;
6409b94acceSBarry Smith 
6419b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
6429b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
6435334005bSBarry Smith 
6449b94acceSBarry Smith   *outsnes = snes;
645*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
6469b94acceSBarry Smith }
6479b94acceSBarry Smith 
6489b94acceSBarry Smith /* --------------------------------------------------------------- */
6495615d1e5SSatish Balay #undef __FUNC__
650d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction"
6519b94acceSBarry Smith /*@C
6529b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6539b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6549b94acceSBarry Smith    equations.
6559b94acceSBarry Smith 
6569b94acceSBarry Smith    Input Parameters:
6579b94acceSBarry Smith .  snes - the SNES context
6589b94acceSBarry Smith .  func - function evaluation routine
6599b94acceSBarry Smith .  r - vector to store function value
6602cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
6612cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6629b94acceSBarry Smith 
6639b94acceSBarry Smith    Calling sequence of func:
664313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
6659b94acceSBarry Smith 
6669b94acceSBarry Smith .  x - input vector
667313e4042SLois Curfman McInnes .  f - function vector
6682cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
6699b94acceSBarry Smith 
6709b94acceSBarry Smith    Notes:
6719b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6729b94acceSBarry Smith $      f'(x) x = -f(x),
6739b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
6749b94acceSBarry Smith 
6759b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6769b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6779b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6789b94acceSBarry Smith 
6799b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6809b94acceSBarry Smith 
681a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6829b94acceSBarry Smith @*/
6835334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6849b94acceSBarry Smith {
685*3a40ed3dSBarry Smith   PetscFunctionBegin;
68677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
687e3372554SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
6889b94acceSBarry Smith   snes->computefunction     = func;
6899b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6909b94acceSBarry Smith   snes->funP                = ctx;
691*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
6929b94acceSBarry Smith }
6939b94acceSBarry Smith 
6945615d1e5SSatish Balay #undef __FUNC__
6955615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
6969b94acceSBarry Smith /*@
6979b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6989b94acceSBarry Smith    SNESSetFunction().
6999b94acceSBarry Smith 
7009b94acceSBarry Smith    Input Parameters:
7019b94acceSBarry Smith .  snes - the SNES context
7029b94acceSBarry Smith .  x - input vector
7039b94acceSBarry Smith 
7049b94acceSBarry Smith    Output Parameter:
7053638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
7069b94acceSBarry Smith 
7071bffabb2SLois Curfman McInnes    Notes:
7089b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
7099b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
7109b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
7119b94acceSBarry Smith 
7129b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7139b94acceSBarry Smith 
714a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7159b94acceSBarry Smith @*/
7169b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
7179b94acceSBarry Smith {
7189b94acceSBarry Smith   int    ierr;
7199b94acceSBarry Smith 
720*3a40ed3dSBarry Smith   PetscFunctionBegin;
721d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
722e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
723d4bb536fSBarry Smith   }
7249b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
7259b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
726ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7279b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
728*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
7299b94acceSBarry Smith }
7309b94acceSBarry Smith 
7315615d1e5SSatish Balay #undef __FUNC__
732d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction"
7339b94acceSBarry Smith /*@C
7349b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
7359b94acceSBarry Smith    unconstrained minimization.
7369b94acceSBarry Smith 
7379b94acceSBarry Smith    Input Parameters:
7389b94acceSBarry Smith .  snes - the SNES context
7399b94acceSBarry Smith .  func - function evaluation routine
740044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
741044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7429b94acceSBarry Smith 
7439b94acceSBarry Smith    Calling sequence of func:
7449b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
7459b94acceSBarry Smith 
7469b94acceSBarry Smith .  x - input vector
7479b94acceSBarry Smith .  f - function
748044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
7499b94acceSBarry Smith 
7509b94acceSBarry Smith    Notes:
7519b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7529b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7539b94acceSBarry Smith    SNESSetFunction().
7549b94acceSBarry Smith 
7559b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
7569b94acceSBarry Smith 
757a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
758a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
7599b94acceSBarry Smith @*/
7609b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
7619b94acceSBarry Smith                       void *ctx)
7629b94acceSBarry Smith {
763*3a40ed3dSBarry Smith   PetscFunctionBegin;
76477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
765e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
7669b94acceSBarry Smith   snes->computeumfunction   = func;
7679b94acceSBarry Smith   snes->umfunP              = ctx;
768*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
7699b94acceSBarry Smith }
7709b94acceSBarry Smith 
7715615d1e5SSatish Balay #undef __FUNC__
7725615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7739b94acceSBarry Smith /*@
7749b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7759b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7769b94acceSBarry Smith 
7779b94acceSBarry Smith    Input Parameters:
7789b94acceSBarry Smith .  snes - the SNES context
7799b94acceSBarry Smith .  x - input vector
7809b94acceSBarry Smith 
7819b94acceSBarry Smith    Output Parameter:
7829b94acceSBarry Smith .  y - function value
7839b94acceSBarry Smith 
7849b94acceSBarry Smith    Notes:
7859b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7869b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7879b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
788a86d99e1SLois Curfman McInnes 
789a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
790a86d99e1SLois Curfman McInnes 
791a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
792a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7939b94acceSBarry Smith @*/
7949b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7959b94acceSBarry Smith {
7969b94acceSBarry Smith   int    ierr;
797*3a40ed3dSBarry Smith 
798*3a40ed3dSBarry Smith   PetscFunctionBegin;
799e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
8009b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
8019b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
802ae3c334cSLois Curfman McInnes   snes->nfuncs++;
8039b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
804*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
8059b94acceSBarry Smith }
8069b94acceSBarry Smith 
8075615d1e5SSatish Balay #undef __FUNC__
808d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient"
8099b94acceSBarry Smith /*@C
8109b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
8119b94acceSBarry Smith    vector for use by the SNES routines.
8129b94acceSBarry Smith 
8139b94acceSBarry Smith    Input Parameters:
8149b94acceSBarry Smith .  snes - the SNES context
8159b94acceSBarry Smith .  func - function evaluation routine
816044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
817044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
8189b94acceSBarry Smith .  r - vector to store gradient value
8199b94acceSBarry Smith 
8209b94acceSBarry Smith    Calling sequence of func:
8219b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
8229b94acceSBarry Smith 
8239b94acceSBarry Smith .  x - input vector
8249b94acceSBarry Smith .  g - gradient vector
825044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
8269b94acceSBarry Smith 
8279b94acceSBarry Smith    Notes:
8289b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
8299b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
8309b94acceSBarry Smith    SNESSetFunction().
8319b94acceSBarry Smith 
8329b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
8339b94acceSBarry Smith 
834a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
835a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
8369b94acceSBarry Smith @*/
83774679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
8389b94acceSBarry Smith {
839*3a40ed3dSBarry Smith   PetscFunctionBegin;
84077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
841e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8429b94acceSBarry Smith   snes->computefunction     = func;
8439b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
8449b94acceSBarry Smith   snes->funP                = ctx;
845*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
8469b94acceSBarry Smith }
8479b94acceSBarry Smith 
8485615d1e5SSatish Balay #undef __FUNC__
8495615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
8509b94acceSBarry Smith /*@
851a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
852a86d99e1SLois Curfman McInnes    SNESSetGradient().
8539b94acceSBarry Smith 
8549b94acceSBarry Smith    Input Parameters:
8559b94acceSBarry Smith .  snes - the SNES context
8569b94acceSBarry Smith .  x - input vector
8579b94acceSBarry Smith 
8589b94acceSBarry Smith    Output Parameter:
8599b94acceSBarry Smith .  y - gradient vector
8609b94acceSBarry Smith 
8619b94acceSBarry Smith    Notes:
8629b94acceSBarry Smith    SNESComputeGradient() is valid only for
8639b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8649b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
865a86d99e1SLois Curfman McInnes 
866a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
867a86d99e1SLois Curfman McInnes 
868a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
869a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8709b94acceSBarry Smith @*/
8719b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8729b94acceSBarry Smith {
8739b94acceSBarry Smith   int    ierr;
874*3a40ed3dSBarry Smith 
875*3a40ed3dSBarry Smith   PetscFunctionBegin;
876*3a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
877e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
878*3a40ed3dSBarry Smith   }
8799b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
8809b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
8819b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
882*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
8839b94acceSBarry Smith }
8849b94acceSBarry Smith 
8855615d1e5SSatish Balay #undef __FUNC__
8865615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
88762fef451SLois Curfman McInnes /*@
88862fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
88962fef451SLois Curfman McInnes    set with SNESSetJacobian().
89062fef451SLois Curfman McInnes 
89162fef451SLois Curfman McInnes    Input Parameters:
89262fef451SLois Curfman McInnes .  snes - the SNES context
89362fef451SLois Curfman McInnes .  x - input vector
89462fef451SLois Curfman McInnes 
89562fef451SLois Curfman McInnes    Output Parameters:
89662fef451SLois Curfman McInnes .  A - Jacobian matrix
89762fef451SLois Curfman McInnes .  B - optional preconditioning matrix
89862fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
89962fef451SLois Curfman McInnes 
90062fef451SLois Curfman McInnes    Notes:
90162fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
90262fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
90362fef451SLois Curfman McInnes 
904dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
905dc5a77f8SLois Curfman McInnes    flag parameter.
90662fef451SLois Curfman McInnes 
90762fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
90862fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
909005c665bSBarry Smith    methods is SNESComputeHessian().
91062fef451SLois Curfman McInnes 
91162fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
91262fef451SLois Curfman McInnes 
91362fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
91462fef451SLois Curfman McInnes @*/
9159b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
9169b94acceSBarry Smith {
9179b94acceSBarry Smith   int    ierr;
918*3a40ed3dSBarry Smith 
919*3a40ed3dSBarry Smith   PetscFunctionBegin;
920*3a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
921e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
922*3a40ed3dSBarry Smith   }
923*3a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
9249b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
925c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
9269b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
9279b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
9286d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
92977c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
93077c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
931*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
9329b94acceSBarry Smith }
9339b94acceSBarry Smith 
9345615d1e5SSatish Balay #undef __FUNC__
9355615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
93662fef451SLois Curfman McInnes /*@
93762fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
93862fef451SLois Curfman McInnes    set with SNESSetHessian().
93962fef451SLois Curfman McInnes 
94062fef451SLois Curfman McInnes    Input Parameters:
94162fef451SLois Curfman McInnes .  snes - the SNES context
94262fef451SLois Curfman McInnes .  x - input vector
94362fef451SLois Curfman McInnes 
94462fef451SLois Curfman McInnes    Output Parameters:
94562fef451SLois Curfman McInnes .  A - Hessian matrix
94662fef451SLois Curfman McInnes .  B - optional preconditioning matrix
94762fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
94862fef451SLois Curfman McInnes 
94962fef451SLois Curfman McInnes    Notes:
95062fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
95162fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
95262fef451SLois Curfman McInnes 
953dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
954dc5a77f8SLois Curfman McInnes    flag parameter.
95562fef451SLois Curfman McInnes 
95662fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
95762fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
95862fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
95962fef451SLois Curfman McInnes 
96062fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
96162fef451SLois Curfman McInnes 
962a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
963a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
96462fef451SLois Curfman McInnes @*/
96562fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
9669b94acceSBarry Smith {
9679b94acceSBarry Smith   int    ierr;
968*3a40ed3dSBarry Smith 
969*3a40ed3dSBarry Smith   PetscFunctionBegin;
970*3a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
971e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
972*3a40ed3dSBarry Smith   }
973*3a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
97462fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9750f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
97662fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
97762fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9780f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
97977c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
98077c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
981*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
9829b94acceSBarry Smith }
9839b94acceSBarry Smith 
9845615d1e5SSatish Balay #undef __FUNC__
985d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian"
9869b94acceSBarry Smith /*@C
9879b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
988044dda88SLois Curfman McInnes    location to store the matrix.
9899b94acceSBarry Smith 
9909b94acceSBarry Smith    Input Parameters:
9919b94acceSBarry Smith .  snes - the SNES context
9929b94acceSBarry Smith .  A - Jacobian matrix
9939b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9949b94acceSBarry Smith .  func - Jacobian evaluation routine
9952cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9962cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9979b94acceSBarry Smith 
9989b94acceSBarry Smith    Calling sequence of func:
999313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10009b94acceSBarry Smith 
10019b94acceSBarry Smith .  x - input vector
10029b94acceSBarry Smith .  A - Jacobian matrix
10039b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1004ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1005ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10062cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
10079b94acceSBarry Smith 
10089b94acceSBarry Smith    Notes:
1009dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10102cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1011ac21db08SLois Curfman McInnes 
1012ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
10139b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
10149b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10159b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10169b94acceSBarry Smith    throughout the global iterations.
10179b94acceSBarry Smith 
10189b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
10199b94acceSBarry Smith 
1020ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
10219b94acceSBarry Smith @*/
10229b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10239b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10249b94acceSBarry Smith {
1025*3a40ed3dSBarry Smith   PetscFunctionBegin;
102677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1027d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
10289b94acceSBarry Smith   snes->computejacobian = func;
10299b94acceSBarry Smith   snes->jacP            = ctx;
10309b94acceSBarry Smith   snes->jacobian        = A;
10319b94acceSBarry Smith   snes->jacobian_pre    = B;
1032*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
10339b94acceSBarry Smith }
103462fef451SLois Curfman McInnes 
10355615d1e5SSatish Balay #undef __FUNC__
1036d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian"
1037b4fd4287SBarry Smith /*@
1038b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1039b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1040b4fd4287SBarry Smith 
1041b4fd4287SBarry Smith    Input Parameter:
1042b4fd4287SBarry Smith .  snes - the nonlinear solver context
1043b4fd4287SBarry Smith 
1044b4fd4287SBarry Smith    Output Parameters:
1045b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
1046b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
1047b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1048b4fd4287SBarry Smith 
1049b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1050b4fd4287SBarry Smith @*/
1051b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
1052b4fd4287SBarry Smith {
1053*3a40ed3dSBarry Smith   PetscFunctionBegin;
1054*3a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1055e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
1056*3a40ed3dSBarry Smith   }
1057b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
1058b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
1059b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
1060*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
1061b4fd4287SBarry Smith }
1062b4fd4287SBarry Smith 
10635615d1e5SSatish Balay #undef __FUNC__
1064d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian"
10659b94acceSBarry Smith /*@C
10669b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1067044dda88SLois Curfman McInnes    location to store the matrix.
10689b94acceSBarry Smith 
10699b94acceSBarry Smith    Input Parameters:
10709b94acceSBarry Smith .  snes - the SNES context
10719b94acceSBarry Smith .  A - Hessian matrix
10729b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
10739b94acceSBarry Smith .  func - Jacobian evaluation routine
1074313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1075313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10769b94acceSBarry Smith 
10779b94acceSBarry Smith    Calling sequence of func:
1078313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10799b94acceSBarry Smith 
10809b94acceSBarry Smith .  x - input vector
10819b94acceSBarry Smith .  A - Hessian matrix
10829b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1083ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1084ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10852cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10869b94acceSBarry Smith 
10879b94acceSBarry Smith    Notes:
1088dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10892cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1090ac21db08SLois Curfman McInnes 
10919b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10929b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10939b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10949b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10959b94acceSBarry Smith    throughout the global iterations.
10969b94acceSBarry Smith 
10979b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10989b94acceSBarry Smith 
1099ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
11009b94acceSBarry Smith @*/
11019b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
11029b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
11039b94acceSBarry Smith {
1104*3a40ed3dSBarry Smith   PetscFunctionBegin;
110577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1106d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1107e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1108d4bb536fSBarry Smith   }
11099b94acceSBarry Smith   snes->computejacobian = func;
11109b94acceSBarry Smith   snes->jacP            = ctx;
11119b94acceSBarry Smith   snes->jacobian        = A;
11129b94acceSBarry Smith   snes->jacobian_pre    = B;
1113*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
11149b94acceSBarry Smith }
11159b94acceSBarry Smith 
11165615d1e5SSatish Balay #undef __FUNC__
1117d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian"
111862fef451SLois Curfman McInnes /*@
111962fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
112062fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
112162fef451SLois Curfman McInnes 
112262fef451SLois Curfman McInnes    Input Parameter:
112362fef451SLois Curfman McInnes .  snes - the nonlinear solver context
112462fef451SLois Curfman McInnes 
112562fef451SLois Curfman McInnes    Output Parameters:
112662fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
112762fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
112862fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
112962fef451SLois Curfman McInnes 
113062fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
113162fef451SLois Curfman McInnes @*/
113262fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
113362fef451SLois Curfman McInnes {
1134*3a40ed3dSBarry Smith   PetscFunctionBegin;
1135*3a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1136e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1137*3a40ed3dSBarry Smith   }
113862fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
113962fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
114062fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
1141*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
114262fef451SLois Curfman McInnes }
114362fef451SLois Curfman McInnes 
11449b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
11459b94acceSBarry Smith 
11465615d1e5SSatish Balay #undef __FUNC__
11475615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
11489b94acceSBarry Smith /*@
11499b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1150272ac6f2SLois Curfman McInnes    of a nonlinear solver.
11519b94acceSBarry Smith 
11529b94acceSBarry Smith    Input Parameter:
11539b94acceSBarry Smith .  snes - the SNES context
11548ddd3da0SLois Curfman McInnes .  x - the solution vector
11559b94acceSBarry Smith 
1156272ac6f2SLois Curfman McInnes    Notes:
1157272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1158272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1159272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1160272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1161272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1162272ac6f2SLois Curfman McInnes 
11639b94acceSBarry Smith .keywords: SNES, nonlinear, setup
11649b94acceSBarry Smith 
11659b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
11669b94acceSBarry Smith @*/
11678ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
11689b94acceSBarry Smith {
1169272ac6f2SLois Curfman McInnes   int ierr, flg;
1170*3a40ed3dSBarry Smith 
1171*3a40ed3dSBarry Smith   PetscFunctionBegin;
117277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
117377c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
11748ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
11759b94acceSBarry Smith 
1176c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1177c1f60f51SBarry Smith   /*
1178c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1179dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1180c1f60f51SBarry Smith   */
1181c1f60f51SBarry Smith   if (flg) {
1182c1f60f51SBarry Smith     Mat J;
1183c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1184c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1185c1f60f51SBarry Smith     snes->mfshell = J;
1186c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1187c1f60f51SBarry Smith       snes->jacobian = J;
1188c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1189c1f60f51SBarry Smith     }
1190c1f60f51SBarry Smith     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1191c1f60f51SBarry Smith       snes->jacobian = J;
1192c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1193d4bb536fSBarry Smith     } else {
1194d4bb536fSBarry Smith       SETERRQ(1,0,"Method class doesn't support matrix-free operator option");
1195d4bb536fSBarry Smith     }
1196c1f60f51SBarry Smith   }
1197272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1198c1f60f51SBarry Smith   /*
1199dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1200c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1201c1f60f51SBarry Smith    */
120231615d04SLois Curfman McInnes   if (flg) {
1203272ac6f2SLois Curfman McInnes     Mat J;
1204272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1205272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1206272ac6f2SLois Curfman McInnes     snes->mfshell = J;
120783e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
120883e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
120994a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
121083e56afdSLois Curfman McInnes     }
121183e56afdSLois Curfman McInnes     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
121283e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
121394a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1214d4bb536fSBarry Smith     } else {
1215d4bb536fSBarry Smith       SETERRQ(1,0,"Method class doesn't support matrix-free option");
1216d4bb536fSBarry Smith     }
1217272ac6f2SLois Curfman McInnes   }
12189b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1219e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first");
1220e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first");
1221e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first");
1222e3372554SBarry Smith     if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector");
1223a703fe33SLois Curfman McInnes 
1224a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
122540191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1226a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1227a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1228a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1229a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1230a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1231a703fe33SLois Curfman McInnes     }
12329b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1233e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first");
1234e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first");
1235d4bb536fSBarry Smith     if (!snes->computeumfunction) SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first");
1236e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first");
1237d4bb536fSBarry Smith   } else {
1238d4bb536fSBarry Smith     SETERRQ(1,0,"Unknown method class");
1239d4bb536fSBarry Smith   }
1240a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1241a703fe33SLois Curfman McInnes   snes->setup_called = 1;
1242*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
12439b94acceSBarry Smith }
12449b94acceSBarry Smith 
12455615d1e5SSatish Balay #undef __FUNC__
1246d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy"
12479b94acceSBarry Smith /*@C
12489b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
12499b94acceSBarry Smith    with SNESCreate().
12509b94acceSBarry Smith 
12519b94acceSBarry Smith    Input Parameter:
12529b94acceSBarry Smith .  snes - the SNES context
12539b94acceSBarry Smith 
12549b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
12559b94acceSBarry Smith 
125663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
12579b94acceSBarry Smith @*/
12589b94acceSBarry Smith int SNESDestroy(SNES snes)
12599b94acceSBarry Smith {
12609b94acceSBarry Smith   int ierr;
1261*3a40ed3dSBarry Smith 
1262*3a40ed3dSBarry Smith   PetscFunctionBegin;
126377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1264*3a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1265d4bb536fSBarry Smith 
12669750a799SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);}
12670452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
12689b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
12699b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
12703e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
12716f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
12729b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
12730452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
1274*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
12759b94acceSBarry Smith }
12769b94acceSBarry Smith 
12779b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12789b94acceSBarry Smith 
12795615d1e5SSatish Balay #undef __FUNC__
12805615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
12819b94acceSBarry Smith /*@
1282d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12839b94acceSBarry Smith 
12849b94acceSBarry Smith    Input Parameters:
12859b94acceSBarry Smith .  snes - the SNES context
128633174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
128733174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
128833174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
128933174efeSLois Curfman McInnes            of the change in the solution between steps
129033174efeSLois Curfman McInnes .  maxit - maximum number of iterations
129133174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
12929b94acceSBarry Smith 
129333174efeSLois Curfman McInnes    Options Database Keys:
129433174efeSLois Curfman McInnes $    -snes_atol <atol>
129533174efeSLois Curfman McInnes $    -snes_rtol <rtol>
129633174efeSLois Curfman McInnes $    -snes_stol <stol>
129733174efeSLois Curfman McInnes $    -snes_max_it <maxit>
129833174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
12999b94acceSBarry Smith 
1300d7a720efSLois Curfman McInnes    Notes:
13019b94acceSBarry Smith    The default maximum number of iterations is 50.
13029b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
13039b94acceSBarry Smith 
130433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
13059b94acceSBarry Smith 
1306d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
13079b94acceSBarry Smith @*/
1308d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
13099b94acceSBarry Smith {
1310*3a40ed3dSBarry Smith   PetscFunctionBegin;
131177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1312d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1313d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1314d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1315d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1316d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1317*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
13189b94acceSBarry Smith }
13199b94acceSBarry Smith 
13205615d1e5SSatish Balay #undef __FUNC__
13215615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
13229b94acceSBarry Smith /*@
132333174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
132433174efeSLois Curfman McInnes 
132533174efeSLois Curfman McInnes    Input Parameters:
132633174efeSLois Curfman McInnes .  snes - the SNES context
132733174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
132833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
132933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
133033174efeSLois Curfman McInnes            of the change in the solution between steps
133133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
133233174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
133333174efeSLois Curfman McInnes 
133433174efeSLois Curfman McInnes    Notes:
133533174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
133633174efeSLois Curfman McInnes 
133733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
133833174efeSLois Curfman McInnes 
133933174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
134033174efeSLois Curfman McInnes @*/
134133174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
134233174efeSLois Curfman McInnes {
1343*3a40ed3dSBarry Smith   PetscFunctionBegin;
134433174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
134533174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
134633174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
134733174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
134833174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
134933174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
1350*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
135133174efeSLois Curfman McInnes }
135233174efeSLois Curfman McInnes 
13535615d1e5SSatish Balay #undef __FUNC__
13545615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
135533174efeSLois Curfman McInnes /*@
13569b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
13579b94acceSBarry Smith 
13589b94acceSBarry Smith    Input Parameters:
13599b94acceSBarry Smith .  snes - the SNES context
13609b94acceSBarry Smith .  tol - tolerance
13619b94acceSBarry Smith 
13629b94acceSBarry Smith    Options Database Key:
1363d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
13649b94acceSBarry Smith 
13659b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
13669b94acceSBarry Smith 
1367d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
13689b94acceSBarry Smith @*/
13699b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
13709b94acceSBarry Smith {
1371*3a40ed3dSBarry Smith   PetscFunctionBegin;
137277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13739b94acceSBarry Smith   snes->deltatol = tol;
1374*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
13759b94acceSBarry Smith }
13769b94acceSBarry Smith 
13775615d1e5SSatish Balay #undef __FUNC__
13785615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
13799b94acceSBarry Smith /*@
138077c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
13819b94acceSBarry Smith    for unconstrained minimization solvers.
13829b94acceSBarry Smith 
13839b94acceSBarry Smith    Input Parameters:
13849b94acceSBarry Smith .  snes - the SNES context
13859b94acceSBarry Smith .  ftol - minimum function tolerance
13869b94acceSBarry Smith 
13879b94acceSBarry Smith    Options Database Key:
1388d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
13899b94acceSBarry Smith 
13909b94acceSBarry Smith    Note:
139177c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
13929b94acceSBarry Smith    methods only.
13939b94acceSBarry Smith 
13949b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
13959b94acceSBarry Smith 
1396d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
13979b94acceSBarry Smith @*/
139877c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
13999b94acceSBarry Smith {
1400*3a40ed3dSBarry Smith   PetscFunctionBegin;
140177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
14029b94acceSBarry Smith   snes->fmin = ftol;
1403*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
14049b94acceSBarry Smith }
14059b94acceSBarry Smith 
14069b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
14079b94acceSBarry Smith 
14085615d1e5SSatish Balay #undef __FUNC__
1409d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor"
14109b94acceSBarry Smith /*@C
14119b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
14129b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
14139b94acceSBarry Smith    progress.
14149b94acceSBarry Smith 
14159b94acceSBarry Smith    Input Parameters:
14169b94acceSBarry Smith .  snes - the SNES context
14179b94acceSBarry Smith .  func - monitoring routine
1418044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1419044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
14209b94acceSBarry Smith 
14219b94acceSBarry Smith    Calling sequence of func:
1422682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
14239b94acceSBarry Smith 
14249b94acceSBarry Smith $    snes - the SNES context
14259b94acceSBarry Smith $    its - iteration number
1426044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
14279b94acceSBarry Smith $
14289b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14299b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
14309b94acceSBarry Smith $
14319b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14329b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
14339b94acceSBarry Smith 
14349665c990SLois Curfman McInnes    Options Database Keys:
14359665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
14369665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
14379665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
14389665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
14399665c990SLois Curfman McInnes $                           been hardwired into a code by
14409665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
14419665c990SLois Curfman McInnes $                           does not cancel those set via
14429665c990SLois Curfman McInnes $                           the options database.
14439665c990SLois Curfman McInnes 
14449665c990SLois Curfman McInnes 
1445639f9d9dSBarry Smith    Notes:
14466bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
14476bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
14486bc08f3fSLois Curfman McInnes    order in which they were set.
1449639f9d9dSBarry Smith 
14509b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
14519b94acceSBarry Smith 
14529b94acceSBarry Smith .seealso: SNESDefaultMonitor()
14539b94acceSBarry Smith @*/
145474679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
14559b94acceSBarry Smith {
1456*3a40ed3dSBarry Smith   PetscFunctionBegin;
1457639f9d9dSBarry Smith   if (!func) {
1458639f9d9dSBarry Smith     snes->numbermonitors = 0;
1459*3a40ed3dSBarry Smith     PetscFunctionReturn(0);
1460639f9d9dSBarry Smith   }
1461639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1462e3372554SBarry Smith     SETERRQ(1,0,"Too many monitors set");
1463639f9d9dSBarry Smith   }
1464639f9d9dSBarry Smith 
1465639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1466639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1467*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
14689b94acceSBarry Smith }
14699b94acceSBarry Smith 
14705615d1e5SSatish Balay #undef __FUNC__
1471d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest"
14729b94acceSBarry Smith /*@C
14739b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
14749b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
14759b94acceSBarry Smith 
14769b94acceSBarry Smith    Input Parameters:
14779b94acceSBarry Smith .  snes - the SNES context
14789b94acceSBarry Smith .  func - routine to test for convergence
1479044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1480044dda88SLois Curfman McInnes           (may be PETSC_NULL)
14819b94acceSBarry Smith 
14829b94acceSBarry Smith    Calling sequence of func:
14839b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
14849b94acceSBarry Smith              double f,void *cctx)
14859b94acceSBarry Smith 
14869b94acceSBarry Smith $    snes - the SNES context
1487044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
14889b94acceSBarry Smith $    xnorm - 2-norm of current iterate
14899b94acceSBarry Smith $
14909b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14919b94acceSBarry Smith $    gnorm - 2-norm of current step
14929b94acceSBarry Smith $    f - 2-norm of function
14939b94acceSBarry Smith $
14949b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14959b94acceSBarry Smith $    gnorm - 2-norm of current gradient
14969b94acceSBarry Smith $    f - function value
14979b94acceSBarry Smith 
14989b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14999b94acceSBarry Smith 
150040191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
150140191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
15029b94acceSBarry Smith @*/
150374679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
15049b94acceSBarry Smith {
1505*3a40ed3dSBarry Smith   PetscFunctionBegin;
15069b94acceSBarry Smith   (snes)->converged = func;
15079b94acceSBarry Smith   (snes)->cnvP      = cctx;
1508*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
15099b94acceSBarry Smith }
15109b94acceSBarry Smith 
15115615d1e5SSatish Balay #undef __FUNC__
1512d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory"
1513c9005455SLois Curfman McInnes /*@
1514c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1515c9005455SLois Curfman McInnes 
1516c9005455SLois Curfman McInnes    Input Parameters:
1517c9005455SLois Curfman McInnes .  snes - iterative context obtained from SNESCreate()
1518c9005455SLois Curfman McInnes .  a   - array to hold history
1519c9005455SLois Curfman McInnes .  na  - size of a
1520c9005455SLois Curfman McInnes 
1521c9005455SLois Curfman McInnes    Notes:
1522c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1523c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1524c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1525c9005455SLois Curfman McInnes    at each step.
1526c9005455SLois Curfman McInnes 
1527c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1528c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1529c9005455SLois Curfman McInnes    during the section of code that is being timed.
1530c9005455SLois Curfman McInnes 
1531c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1532c9005455SLois Curfman McInnes @*/
1533c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1534c9005455SLois Curfman McInnes {
1535*3a40ed3dSBarry Smith   PetscFunctionBegin;
1536c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1537c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1538c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1539c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
1540*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
1541c9005455SLois Curfman McInnes }
1542c9005455SLois Curfman McInnes 
1543c9005455SLois Curfman McInnes #undef __FUNC__
15445615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
15459b94acceSBarry Smith /*
15469b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
15479b94acceSBarry Smith    positive parameter delta.
15489b94acceSBarry Smith 
15499b94acceSBarry Smith     Input Parameters:
15509b94acceSBarry Smith .   snes - the SNES context
15519b94acceSBarry Smith .   y - approximate solution of linear system
15529b94acceSBarry Smith .   fnorm - 2-norm of current function
15539b94acceSBarry Smith .   delta - trust region size
15549b94acceSBarry Smith 
15559b94acceSBarry Smith     Output Parameters:
15569b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
15579b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
15589b94acceSBarry Smith     region, and exceeds zero otherwise.
15599b94acceSBarry Smith .   ynorm - 2-norm of the step
15609b94acceSBarry Smith 
15619b94acceSBarry Smith     Note:
156240191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
15639b94acceSBarry Smith     is set to be the maximum allowable step size.
15649b94acceSBarry Smith 
15659b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
15669b94acceSBarry Smith */
15679b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
15689b94acceSBarry Smith                           double *gpnorm,double *ynorm)
15699b94acceSBarry Smith {
15709b94acceSBarry Smith   double norm;
15719b94acceSBarry Smith   Scalar cnorm;
1572*3a40ed3dSBarry Smith   int    ierr;
1573*3a40ed3dSBarry Smith 
1574*3a40ed3dSBarry Smith   PetscFunctionBegin;
1575*3a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
15769b94acceSBarry Smith   if (norm > *delta) {
15779b94acceSBarry Smith      norm = *delta/norm;
15789b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
15799b94acceSBarry Smith      cnorm = norm;
15809b94acceSBarry Smith      VecScale( &cnorm, y );
15819b94acceSBarry Smith      *ynorm = *delta;
15829b94acceSBarry Smith   } else {
15839b94acceSBarry Smith      *gpnorm = 0.0;
15849b94acceSBarry Smith      *ynorm = norm;
15859b94acceSBarry Smith   }
1586*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
15879b94acceSBarry Smith }
15889b94acceSBarry Smith 
15895615d1e5SSatish Balay #undef __FUNC__
15905615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
15919b94acceSBarry Smith /*@
15929b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
159363a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
15949b94acceSBarry Smith 
15959b94acceSBarry Smith    Input Parameter:
15969b94acceSBarry Smith .  snes - the SNES context
15978ddd3da0SLois Curfman McInnes .  x - the solution vector
15989b94acceSBarry Smith 
15999b94acceSBarry Smith    Output Parameter:
16009b94acceSBarry Smith    its - number of iterations until termination
16019b94acceSBarry Smith 
16028ddd3da0SLois Curfman McInnes    Note:
16038ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
16048ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
16058ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
16068ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
16078ddd3da0SLois Curfman McInnes 
16089b94acceSBarry Smith .keywords: SNES, nonlinear, solve
16099b94acceSBarry Smith 
161063a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
16119b94acceSBarry Smith @*/
16128ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
16139b94acceSBarry Smith {
16143c7409f5SSatish Balay   int ierr, flg;
1615052efed2SBarry Smith 
1616*3a40ed3dSBarry Smith   PetscFunctionBegin;
161777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
161874679c65SBarry Smith   PetscValidIntPointer(its);
1619c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1620c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
16219b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1622c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
16239b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
16249b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
16253c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
16266d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
1627*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
16289b94acceSBarry Smith }
16299b94acceSBarry Smith 
16309b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
163137fcc0dbSBarry Smith static NRList *__SNESList = 0;
16329b94acceSBarry Smith 
16335615d1e5SSatish Balay #undef __FUNC__
16345615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
16359b94acceSBarry Smith /*@
16364b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
16379b94acceSBarry Smith 
16389b94acceSBarry Smith    Input Parameters:
16399b94acceSBarry Smith .  snes - the SNES context
16409b94acceSBarry Smith .  method - a known method
16419b94acceSBarry Smith 
1642ae12b187SLois Curfman McInnes   Options Database Command:
1643ae12b187SLois Curfman McInnes $ -snes_type  <method>
1644ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1645ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1646ae12b187SLois Curfman McInnes 
16479b94acceSBarry Smith    Notes:
16489b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
16499b94acceSBarry Smith $  Systems of nonlinear equations:
165040191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
165140191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
16529b94acceSBarry Smith $  Unconstrained minimization:
165340191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
165440191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
16559b94acceSBarry Smith 
1656ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1657ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1658ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1659ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1660ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1661ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1662ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1663ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1664ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1665ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1666a703fe33SLois Curfman McInnes 
1667f536c699SSatish Balay .keywords: SNES, set, method
16689b94acceSBarry Smith @*/
16694b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
16709b94acceSBarry Smith {
167184cb2905SBarry Smith   int ierr;
16729b94acceSBarry Smith   int (*r)(SNES);
1673*3a40ed3dSBarry Smith 
1674*3a40ed3dSBarry Smith   PetscFunctionBegin;
167577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1676*3a40ed3dSBarry Smith   if (snes->type == (int) method) PetscFunctionReturn(0);
16777d1a2b2bSBarry Smith 
16789b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
167984cb2905SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
1680e3372554SBarry Smith   if (!__SNESList) {SETERRQ(1,0,"Could not get methods");}
168137fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1682e3372554SBarry Smith   if (!r) {SETERRQ(1,0,"Unknown method");}
16830452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
16849b94acceSBarry Smith   snes->set_method_called = 1;
1685*3a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
1686*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
16879b94acceSBarry Smith }
16889b94acceSBarry Smith 
16899b94acceSBarry Smith /* --------------------------------------------------------------------- */
16905615d1e5SSatish Balay #undef __FUNC__
1691d4bb536fSBarry Smith #define __FUNC__ "SNESRegister"
16929b94acceSBarry Smith /*@C
16939b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
16944b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
16959b94acceSBarry Smith 
16969b94acceSBarry Smith    Input Parameters:
16972d872ea7SLois Curfman McInnes .  name - either a predefined name such as SNES_EQ_LS, or SNES_NEW
16982d872ea7SLois Curfman McInnes           to indicate a new user-defined solver
169940191667SLois Curfman McInnes .  sname - corresponding string for name
17009b94acceSBarry Smith .  create - routine to create method context
17019b94acceSBarry Smith 
170284cb2905SBarry Smith    Output Parameter:
170384cb2905SBarry Smith .  oname - type associated with this new solver
170484cb2905SBarry Smith 
17052d872ea7SLois Curfman McInnes    Notes:
17062d872ea7SLois Curfman McInnes    Multiple user-defined nonlinear solvers can be added by calling
17072d872ea7SLois Curfman McInnes    SNESRegister() with the input parameter "name" set to be SNES_NEW;
17082d872ea7SLois Curfman McInnes    each call will return a unique solver type in the output
17092d872ea7SLois Curfman McInnes    parameter "oname".
17102d872ea7SLois Curfman McInnes 
17119b94acceSBarry Smith .keywords: SNES, nonlinear, register
17129b94acceSBarry Smith 
17139b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
17149b94acceSBarry Smith @*/
171584cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES))
17169b94acceSBarry Smith {
17179b94acceSBarry Smith   int ierr;
171884cb2905SBarry Smith   static int numberregistered = 0;
171984cb2905SBarry Smith 
1720*3a40ed3dSBarry Smith   PetscFunctionBegin;
1721d252947aSBarry Smith   if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++);
172284cb2905SBarry Smith 
172384cb2905SBarry Smith   if (oname) *oname = name;
172437fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
172584cb2905SBarry Smith   NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create );
1726*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
17279b94acceSBarry Smith }
1728a847f771SSatish Balay 
17299b94acceSBarry Smith /* --------------------------------------------------------------------- */
17305615d1e5SSatish Balay #undef __FUNC__
1731d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy"
17329b94acceSBarry Smith /*@C
17339b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
17349b94acceSBarry Smith    registered by SNESRegister().
17359b94acceSBarry Smith 
17369b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
17379b94acceSBarry Smith 
17389b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
17399b94acceSBarry Smith @*/
17409b94acceSBarry Smith int SNESRegisterDestroy()
17419b94acceSBarry Smith {
1742*3a40ed3dSBarry Smith   PetscFunctionBegin;
174337fcc0dbSBarry Smith   if (__SNESList) {
174437fcc0dbSBarry Smith     NRDestroy( __SNESList );
174537fcc0dbSBarry Smith     __SNESList = 0;
17469b94acceSBarry Smith   }
174784cb2905SBarry Smith   SNESRegisterAllCalled = 0;
1748*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
17499b94acceSBarry Smith }
17509b94acceSBarry Smith 
17515615d1e5SSatish Balay #undef __FUNC__
1752d4bb536fSBarry Smith #define __FUNC__ "SNESGetTypeFromOptions_Private"
17539b94acceSBarry Smith /*
17544b0e389bSBarry Smith    SNESGetTypeFromOptions_Private - Sets the selected method from the
17559b94acceSBarry Smith    options database.
17569b94acceSBarry Smith 
17579b94acceSBarry Smith    Input Parameter:
17589b94acceSBarry Smith .  ctx - the SNES context
17599b94acceSBarry Smith 
17609b94acceSBarry Smith    Output Parameter:
17619b94acceSBarry Smith .  method -  solver method
17629b94acceSBarry Smith 
17639b94acceSBarry Smith    Returns:
17649b94acceSBarry Smith    Returns 1 if the method is found; 0 otherwise.
17659b94acceSBarry Smith 
17669b94acceSBarry Smith    Options Database Key:
17674b0e389bSBarry Smith $  -snes_type  method
17689b94acceSBarry Smith */
1769052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
17709b94acceSBarry Smith {
1771052efed2SBarry Smith   int  ierr;
17729b94acceSBarry Smith   char sbuf[50];
17735baf8537SBarry Smith 
1774*3a40ed3dSBarry Smith   PetscFunctionBegin;
1775052efed2SBarry Smith   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1776052efed2SBarry Smith   if (*flg) {
1777052efed2SBarry Smith     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
177837fcc0dbSBarry Smith     *method = (SNESType)NRFindID( __SNESList, sbuf );
17799b94acceSBarry Smith   }
1780*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
17819b94acceSBarry Smith }
17829b94acceSBarry Smith 
17835615d1e5SSatish Balay #undef __FUNC__
1784d4bb536fSBarry Smith #define __FUNC__ "SNESGetType"
17859b94acceSBarry Smith /*@C
17869a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
17879b94acceSBarry Smith 
17889b94acceSBarry Smith    Input Parameter:
17894b0e389bSBarry Smith .  snes - nonlinear solver context
17909b94acceSBarry Smith 
17919b94acceSBarry Smith    Output Parameter:
17929a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
17939a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
17949b94acceSBarry Smith 
17959b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
17969b94acceSBarry Smith @*/
17974b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
17989b94acceSBarry Smith {
179906a9b0adSLois Curfman McInnes   int ierr;
1800*3a40ed3dSBarry Smith 
1801*3a40ed3dSBarry Smith   PetscFunctionBegin;
180237fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
18034b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
180437fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
1805*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
18069b94acceSBarry Smith }
18079b94acceSBarry Smith 
18085615d1e5SSatish Balay #undef __FUNC__
1809d4bb536fSBarry Smith #define __FUNC__ "SNESPrintTypes_Private"
18109b94acceSBarry Smith /*
18114b0e389bSBarry Smith    SNESPrintTypes_Private - Prints the SNES methods available from the
18129b94acceSBarry Smith    options database.
18139b94acceSBarry Smith 
18149b94acceSBarry Smith    Input Parameters:
18150752156aSBarry Smith .  comm   - communicator (usually PETSC_COMM_WORLD)
18169b94acceSBarry Smith .  prefix - prefix (usually "-")
18174b0e389bSBarry Smith .  name   - the options database name (by default "snes_type")
18189b94acceSBarry Smith */
181933455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
18209b94acceSBarry Smith {
18219b94acceSBarry Smith   FuncList *entry;
1822*3a40ed3dSBarry Smith 
1823*3a40ed3dSBarry Smith   PetscFunctionBegin;
182437fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
182537fcc0dbSBarry Smith   entry = __SNESList->head;
182677c4ece6SBarry Smith   PetscPrintf(comm," %s%s (one of)",prefix,name);
18279b94acceSBarry Smith   while (entry) {
182877c4ece6SBarry Smith     PetscPrintf(comm," %s",entry->name);
18299b94acceSBarry Smith     entry = entry->next;
18309b94acceSBarry Smith   }
183177c4ece6SBarry Smith   PetscPrintf(comm,"\n");
1832*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
18339b94acceSBarry Smith }
18349b94acceSBarry Smith 
18355615d1e5SSatish Balay #undef __FUNC__
1836d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution"
18379b94acceSBarry Smith /*@C
18389b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
18399b94acceSBarry Smith    stored.
18409b94acceSBarry Smith 
18419b94acceSBarry Smith    Input Parameter:
18429b94acceSBarry Smith .  snes - the SNES context
18439b94acceSBarry Smith 
18449b94acceSBarry Smith    Output Parameter:
18459b94acceSBarry Smith .  x - the solution
18469b94acceSBarry Smith 
18479b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
18489b94acceSBarry Smith 
1849a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
18509b94acceSBarry Smith @*/
18519b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
18529b94acceSBarry Smith {
1853*3a40ed3dSBarry Smith   PetscFunctionBegin;
185477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18559b94acceSBarry Smith   *x = snes->vec_sol_always;
1856*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
18579b94acceSBarry Smith }
18589b94acceSBarry Smith 
18595615d1e5SSatish Balay #undef __FUNC__
1860d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate"
18619b94acceSBarry Smith /*@C
18629b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
18639b94acceSBarry Smith    stored.
18649b94acceSBarry Smith 
18659b94acceSBarry Smith    Input Parameter:
18669b94acceSBarry Smith .  snes - the SNES context
18679b94acceSBarry Smith 
18689b94acceSBarry Smith    Output Parameter:
18699b94acceSBarry Smith .  x - the solution update
18709b94acceSBarry Smith 
18719b94acceSBarry Smith    Notes:
18729b94acceSBarry Smith    This vector is implementation dependent.
18739b94acceSBarry Smith 
18749b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
18759b94acceSBarry Smith 
18769b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
18779b94acceSBarry Smith @*/
18789b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
18799b94acceSBarry Smith {
1880*3a40ed3dSBarry Smith   PetscFunctionBegin;
188177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18829b94acceSBarry Smith   *x = snes->vec_sol_update_always;
1883*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
18849b94acceSBarry Smith }
18859b94acceSBarry Smith 
18865615d1e5SSatish Balay #undef __FUNC__
1887d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction"
18889b94acceSBarry Smith /*@C
18893638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
18909b94acceSBarry Smith 
18919b94acceSBarry Smith    Input Parameter:
18929b94acceSBarry Smith .  snes - the SNES context
18939b94acceSBarry Smith 
18949b94acceSBarry Smith    Output Parameter:
18953638b69dSLois Curfman McInnes .  r - the function
18969b94acceSBarry Smith 
18979b94acceSBarry Smith    Notes:
18989b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
18999b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
19009b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
19019b94acceSBarry Smith 
1902a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
19039b94acceSBarry Smith 
190431615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
190531615d04SLois Curfman McInnes           SNESGetGradient()
19069b94acceSBarry Smith @*/
19079b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
19089b94acceSBarry Smith {
1909*3a40ed3dSBarry Smith   PetscFunctionBegin;
191077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1911*3a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
19129b94acceSBarry Smith   *r = snes->vec_func_always;
1913*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
19149b94acceSBarry Smith }
19159b94acceSBarry Smith 
19165615d1e5SSatish Balay #undef __FUNC__
1917d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient"
19189b94acceSBarry Smith /*@C
19193638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
19209b94acceSBarry Smith 
19219b94acceSBarry Smith    Input Parameter:
19229b94acceSBarry Smith .  snes - the SNES context
19239b94acceSBarry Smith 
19249b94acceSBarry Smith    Output Parameter:
19259b94acceSBarry Smith .  r - the gradient
19269b94acceSBarry Smith 
19279b94acceSBarry Smith    Notes:
19289b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
19299b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
19309b94acceSBarry Smith    SNESGetFunction().
19319b94acceSBarry Smith 
19329b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
19339b94acceSBarry Smith 
193431615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
19359b94acceSBarry Smith @*/
19369b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
19379b94acceSBarry Smith {
1938*3a40ed3dSBarry Smith   PetscFunctionBegin;
193977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1940*3a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1941*3a40ed3dSBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1942*3a40ed3dSBarry Smith   }
19439b94acceSBarry Smith   *r = snes->vec_func_always;
1944*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
19459b94acceSBarry Smith }
19469b94acceSBarry Smith 
19475615d1e5SSatish Balay #undef __FUNC__
1948d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction"
19499b94acceSBarry Smith /*@
19509b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
19519b94acceSBarry Smith    unconstrained minimization problems.
19529b94acceSBarry Smith 
19539b94acceSBarry Smith    Input Parameter:
19549b94acceSBarry Smith .  snes - the SNES context
19559b94acceSBarry Smith 
19569b94acceSBarry Smith    Output Parameter:
19579b94acceSBarry Smith .  r - the function
19589b94acceSBarry Smith 
19599b94acceSBarry Smith    Notes:
19609b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
19619b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
19629b94acceSBarry Smith    SNESGetFunction().
19639b94acceSBarry Smith 
19649b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
19659b94acceSBarry Smith 
196631615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
19679b94acceSBarry Smith @*/
19689b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
19699b94acceSBarry Smith {
1970*3a40ed3dSBarry Smith   PetscFunctionBegin;
197177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
197274679c65SBarry Smith   PetscValidScalarPointer(r);
1973*3a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1974*3a40ed3dSBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1975*3a40ed3dSBarry Smith   }
19769b94acceSBarry Smith   *r = snes->fc;
1977*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
19789b94acceSBarry Smith }
19799b94acceSBarry Smith 
19805615d1e5SSatish Balay #undef __FUNC__
1981d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix"
19823c7409f5SSatish Balay /*@C
19833c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1984d850072dSLois Curfman McInnes    SNES options in the database.
19853c7409f5SSatish Balay 
19863c7409f5SSatish Balay    Input Parameter:
19873c7409f5SSatish Balay .  snes - the SNES context
19883c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
19893c7409f5SSatish Balay 
1990d850072dSLois Curfman McInnes    Notes:
1991a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1992a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1993a83b1b31SSatish Balay    hyphen.
1994d850072dSLois Curfman McInnes 
19953c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1996a86d99e1SLois Curfman McInnes 
1997a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
19983c7409f5SSatish Balay @*/
19993c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
20003c7409f5SSatish Balay {
20013c7409f5SSatish Balay   int ierr;
20023c7409f5SSatish Balay 
2003*3a40ed3dSBarry Smith   PetscFunctionBegin;
200477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2005639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20063c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2007*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
20083c7409f5SSatish Balay }
20093c7409f5SSatish Balay 
20105615d1e5SSatish Balay #undef __FUNC__
2011d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix"
20123c7409f5SSatish Balay /*@C
2013f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2014d850072dSLois Curfman McInnes    SNES options in the database.
20153c7409f5SSatish Balay 
20163c7409f5SSatish Balay    Input Parameter:
20173c7409f5SSatish Balay .  snes - the SNES context
20183c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
20193c7409f5SSatish Balay 
2020d850072dSLois Curfman McInnes    Notes:
2021a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2022a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
2023a83b1b31SSatish Balay    hyphen.
2024d850072dSLois Curfman McInnes 
20253c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2026a86d99e1SLois Curfman McInnes 
2027a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20283c7409f5SSatish Balay @*/
20293c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
20303c7409f5SSatish Balay {
20313c7409f5SSatish Balay   int ierr;
20323c7409f5SSatish Balay 
2033*3a40ed3dSBarry Smith   PetscFunctionBegin;
203477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2035639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20363c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2037*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
20383c7409f5SSatish Balay }
20393c7409f5SSatish Balay 
20405615d1e5SSatish Balay #undef __FUNC__
2041d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix"
2042c83590e2SSatish Balay /*@
20433c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20443c7409f5SSatish Balay    SNES options in the database.
20453c7409f5SSatish Balay 
20463c7409f5SSatish Balay    Input Parameter:
20473c7409f5SSatish Balay .  snes - the SNES context
20483c7409f5SSatish Balay 
20493c7409f5SSatish Balay    Output Parameter:
20503c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20513c7409f5SSatish Balay 
20523c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2053a86d99e1SLois Curfman McInnes 
2054a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20553c7409f5SSatish Balay @*/
20563c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
20573c7409f5SSatish Balay {
20583c7409f5SSatish Balay   int ierr;
20593c7409f5SSatish Balay 
2060*3a40ed3dSBarry Smith   PetscFunctionBegin;
206177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2062639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
2063*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
20643c7409f5SSatish Balay }
20653c7409f5SSatish Balay 
20663c7409f5SSatish Balay 
20679b94acceSBarry Smith 
20689b94acceSBarry Smith 
20699b94acceSBarry Smith 
2070*3a40ed3dSBarry Smith 
2071