xref: /petsc/src/snes/interface/snes.c (revision 639f9d9dbbc54d6ac4e42e98283c540b41bb2cee)
19b94acceSBarry Smith #ifndef lint
2*639f9d9dSBarry Smith static char vcid[] = "$Id: snes.c,v 1.94 1996/10/15 21:58:38 balay Exp bsmith $";
39b94acceSBarry Smith #endif
49b94acceSBarry Smith 
59b94acceSBarry Smith #include "draw.h"          /*I "draw.h"  I*/
670f55243SBarry Smith #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
7f5eb4b81SSatish Balay #include "src/sys/nreg.h"
89b94acceSBarry Smith #include "pinclude/pviewer.h"
99b94acceSBarry Smith #include <math.h>
109b94acceSBarry Smith 
11052efed2SBarry Smith extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*);
1233455573SSatish Balay extern int SNESPrintTypes_Private(MPI_Comm,char*,char*);
139b94acceSBarry Smith 
149b94acceSBarry Smith /*@
159b94acceSBarry Smith    SNESView - Prints the SNES data structure.
169b94acceSBarry Smith 
179b94acceSBarry Smith    Input Parameters:
189b94acceSBarry Smith .  SNES - the SNES context
199b94acceSBarry Smith .  viewer - visualization context
209b94acceSBarry Smith 
219b94acceSBarry Smith    Options Database Key:
229b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
239b94acceSBarry Smith 
249b94acceSBarry Smith    Notes:
259b94acceSBarry Smith    The available visualization contexts include
266d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
276d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
289b94acceSBarry Smith $       output where only the first processor opens
299b94acceSBarry Smith $       the file.  All other processors send their
309b94acceSBarry Smith $       data to the first processor to print.
319b94acceSBarry Smith 
329b94acceSBarry Smith    The user can open alternative vistualization contexts with
33dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
349b94acceSBarry Smith 
359b94acceSBarry Smith .keywords: SNES, view
369b94acceSBarry Smith 
37dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
389b94acceSBarry Smith @*/
399b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
409b94acceSBarry Smith {
419b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
429b94acceSBarry Smith   FILE                *fd;
439b94acceSBarry Smith   int                 ierr;
449b94acceSBarry Smith   SLES                sles;
459b94acceSBarry Smith   char                *method;
4619bcc07fSBarry Smith   ViewerType          vtype;
479b94acceSBarry Smith 
4874679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4974679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5074679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5174679c65SBarry Smith 
5219bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5319bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5490ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5577c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
564b0e389bSBarry Smith     SNESGetType(snes,PETSC_NULL,&method);
5777c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
589b94acceSBarry Smith     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
5977c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
609b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
619b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6277c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
639b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
649b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
659b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
6677c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
679b94acceSBarry Smith     if (snes->ksp_ewconv) {
689b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
699b94acceSBarry Smith       if (kctx) {
7077c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
719b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
729b94acceSBarry Smith         kctx->version);
7377c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
749b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
759b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
7677c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
779b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
789b94acceSBarry Smith       }
799b94acceSBarry Smith     }
80c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
81c30f285eSLois Curfman McInnes     SNESGetType(snes,PETSC_NULL,&method);
82c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
8319bcc07fSBarry Smith   }
849b94acceSBarry Smith   SNESGetSLES(snes,&sles);
859b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
869b94acceSBarry Smith   return 0;
879b94acceSBarry Smith }
889b94acceSBarry Smith 
89*639f9d9dSBarry Smith /*
90*639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
91*639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
92*639f9d9dSBarry Smith */
93*639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
94*639f9d9dSBarry Smith static int numberofsetfromoptions;
95*639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
96*639f9d9dSBarry Smith 
97*639f9d9dSBarry Smith /*@
98*639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
99*639f9d9dSBarry Smith 
100*639f9d9dSBarry Smith   Input Parameter:
101*639f9d9dSBarry Smith .   snescheck - function that checks for options
102*639f9d9dSBarry Smith 
103*639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
104*639f9d9dSBarry Smith @*/
105*639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
106*639f9d9dSBarry Smith {
107*639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
108*639f9d9dSBarry Smith     SETERRQ(1,"SNESAddOptionsChecker:Too many options checkers, only 5 allowed");
109*639f9d9dSBarry Smith   }
110*639f9d9dSBarry Smith 
111*639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
112*639f9d9dSBarry Smith   return 0;
113*639f9d9dSBarry Smith }
114*639f9d9dSBarry Smith 
1159b94acceSBarry Smith /*@
116682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1179b94acceSBarry Smith 
1189b94acceSBarry Smith    Input Parameter:
1199b94acceSBarry Smith .  snes - the SNES context
1209b94acceSBarry Smith 
1219b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1229b94acceSBarry Smith 
123a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1249b94acceSBarry Smith @*/
1259b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1269b94acceSBarry Smith {
1274b0e389bSBarry Smith   SNESType method;
1289b94acceSBarry Smith   double   tmp;
1299b94acceSBarry Smith   SLES     sles;
130*639f9d9dSBarry Smith   int      ierr, flg,i;
1313c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1329b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1339b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1349b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1359b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1369b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1379b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1389b94acceSBarry Smith 
13977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
14074679c65SBarry Smith   if (snes->setup_called) SETERRQ(1,"SNESSetFromOptions:Must call prior to SNESSetUp");
141052efed2SBarry Smith   ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr);
142052efed2SBarry Smith   if (flg) {
143052efed2SBarry Smith     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
1449b94acceSBarry Smith   }
1455334005bSBarry Smith   else if (!snes->set_method_called) {
1465334005bSBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
14740191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
1485334005bSBarry Smith     }
1495334005bSBarry Smith     else {
15040191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
1515334005bSBarry Smith     }
1525334005bSBarry Smith   }
1535334005bSBarry Smith 
1543c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1553c7409f5SSatish Balay   if (flg) { SNESPrintHelp(snes); }
1563c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
157d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);}
1583c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
159d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);}
1603c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
161d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); }
1623c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1633c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
164d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
165d7a720efSLois Curfman McInnes   if (flg) { SNESSetTrustRegionTolerance(snes,tmp); }
166d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
167d7a720efSLois Curfman McInnes   if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); }
1683c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1693c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
170b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
171b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
172b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
173b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
174b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
175b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
176b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
1779b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
1789b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
1793c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
180*639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
1813c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
182d31a3109SLois Curfman McInnes   if (flg) {
183*639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
184d31a3109SLois Curfman McInnes   }
1853c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_xmonitor",&flg); CHKERRQ(ierr);
1863c7409f5SSatish Balay   if (flg) {
18717699dbbSLois Curfman McInnes     int    rank = 0;
188d7e8b826SBarry Smith     DrawLG lg;
18917699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
19017699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
19117699dbbSLois Curfman McInnes     if (!rank) {
1929b94acceSBarry Smith       ierr = SNESLGMonitorCreate(0,0,0,0,300,300,&lg); CHKERRQ(ierr);
1939b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
194f8c826e1SBarry Smith       snes->xmonitor = lg;
1959b94acceSBarry Smith       PLogObjectParent(snes,lg);
1969b94acceSBarry Smith     }
1979b94acceSBarry Smith   }
1983c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
1993c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2009b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2019b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
20294a424c1SBarry Smith     PLogInfo(snes,
20331615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
20431615d04SLois Curfman McInnes   }
20531615d04SLois Curfman McInnes   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
20631615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
20731615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
20894a424c1SBarry Smith     PLogInfo(snes,
20931615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2109b94acceSBarry Smith   }
211*639f9d9dSBarry Smith 
212*639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
213*639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
214*639f9d9dSBarry Smith   }
215*639f9d9dSBarry Smith 
2169b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2179b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
2189b94acceSBarry Smith   if (!snes->setfromoptions) return 0;
2199b94acceSBarry Smith   return (*snes->setfromoptions)(snes);
2209b94acceSBarry Smith }
2219b94acceSBarry Smith 
2229b94acceSBarry Smith /*@
2239b94acceSBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2249b94acceSBarry Smith 
2259b94acceSBarry Smith    Input Parameter:
2269b94acceSBarry Smith .  snes - the SNES context
2279b94acceSBarry Smith 
228a703fe33SLois Curfman McInnes    Options Database Keys:
229a703fe33SLois Curfman McInnes $  -help, -h
230a703fe33SLois Curfman McInnes 
2319b94acceSBarry Smith .keywords: SNES, nonlinear, help
2329b94acceSBarry Smith 
233682d7d0cSBarry Smith .seealso: SNESSetFromOptions()
2349b94acceSBarry Smith @*/
2359b94acceSBarry Smith int SNESPrintHelp(SNES snes)
2369b94acceSBarry Smith {
2376daaf66cSBarry Smith   char                p[64];
2386daaf66cSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2396daaf66cSBarry Smith 
24077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2416daaf66cSBarry Smith 
2426daaf66cSBarry Smith   PetscStrcpy(p,"-");
2436daaf66cSBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2446daaf66cSBarry Smith 
2456daaf66cSBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2466daaf66cSBarry Smith 
247d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
24833455573SSatish Balay   SNESPrintTypes_Private(snes->comm,p,"snes_type");
24977c4ece6SBarry Smith   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
250b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
251b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
252b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
253b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
254b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
255b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
256d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm," SNES Monitoring Options: Choose 1 of the following\n");
257d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor\n",p);
258d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
259b18e04deSLois Curfman McInnes   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
26077c4ece6SBarry Smith     PetscPrintf(snes->comm,
261d31a3109SLois Curfman McInnes      " Options for solving systems of nonlinear equations only:\n");
26277c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
26377c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
264ef1dfb11SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
2651650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p);
2661650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p);
26777c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
26877c4ece6SBarry Smith     PetscPrintf(snes->comm,
269b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
27077c4ece6SBarry Smith     PetscPrintf(snes->comm,
271b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
27277c4ece6SBarry Smith     PetscPrintf(snes->comm,
273b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
27477c4ece6SBarry Smith     PetscPrintf(snes->comm,
275b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
27677c4ece6SBarry Smith     PetscPrintf(snes->comm,
277b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
27877c4ece6SBarry Smith     PetscPrintf(snes->comm,
279b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
28077c4ece6SBarry Smith     PetscPrintf(snes->comm,
281b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
282b18e04deSLois Curfman McInnes   }
283b18e04deSLois Curfman McInnes   else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
28477c4ece6SBarry Smith     PetscPrintf(snes->comm,
285d31a3109SLois Curfman McInnes      " Options for solving unconstrained minimization problems only:\n");
286b18e04deSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
28777c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
28877c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
289d31a3109SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p);
290d31a3109SLois Curfman McInnes      PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p);
291b18e04deSLois Curfman McInnes   }
2924537a946SLois Curfman McInnes   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
29377c4ece6SBarry Smith   PetscPrintf(snes->comm,"a particular method\n");
2946daaf66cSBarry Smith   if (snes->printhelp) (*snes->printhelp)(snes,p);
2959b94acceSBarry Smith   return 0;
2969b94acceSBarry Smith }
2979b94acceSBarry Smith /*@
2989b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2999b94acceSBarry Smith    the nonlinear solvers.
3009b94acceSBarry Smith 
3019b94acceSBarry Smith    Input Parameters:
3029b94acceSBarry Smith .  snes - the SNES context
3039b94acceSBarry Smith .  usrP - optional user context
3049b94acceSBarry Smith 
3059b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3069b94acceSBarry Smith 
3079b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3089b94acceSBarry Smith @*/
3099b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3109b94acceSBarry Smith {
31177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3129b94acceSBarry Smith   snes->user		= usrP;
3139b94acceSBarry Smith   return 0;
3149b94acceSBarry Smith }
31574679c65SBarry Smith 
3169b94acceSBarry Smith /*@C
3179b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3189b94acceSBarry Smith    nonlinear solvers.
3199b94acceSBarry Smith 
3209b94acceSBarry Smith    Input Parameter:
3219b94acceSBarry Smith .  snes - SNES context
3229b94acceSBarry Smith 
3239b94acceSBarry Smith    Output Parameter:
3249b94acceSBarry Smith .  usrP - user context
3259b94acceSBarry Smith 
3269b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3279b94acceSBarry Smith 
3289b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3299b94acceSBarry Smith @*/
3309b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3319b94acceSBarry Smith {
33277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3339b94acceSBarry Smith   *usrP = snes->user;
3349b94acceSBarry Smith   return 0;
3359b94acceSBarry Smith }
33674679c65SBarry Smith 
3379b94acceSBarry Smith /*@
3389b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3399b94acceSBarry Smith    nonlinear solver.
3409b94acceSBarry Smith 
3419b94acceSBarry Smith    Input Parameter:
3429b94acceSBarry Smith .  snes - SNES context
3439b94acceSBarry Smith 
3449b94acceSBarry Smith    Output Parameter:
3459b94acceSBarry Smith .  iter - iteration number
3469b94acceSBarry Smith 
3479b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3489b94acceSBarry Smith @*/
3499b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3509b94acceSBarry Smith {
35177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
35274679c65SBarry Smith   PetscValidIntPointer(iter);
3539b94acceSBarry Smith   *iter = snes->iter;
3549b94acceSBarry Smith   return 0;
3559b94acceSBarry Smith }
35674679c65SBarry Smith 
3579b94acceSBarry Smith /*@
3589b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3599b94acceSBarry Smith    with SNESSSetFunction().
3609b94acceSBarry Smith 
3619b94acceSBarry Smith    Input Parameter:
3629b94acceSBarry Smith .  snes - SNES context
3639b94acceSBarry Smith 
3649b94acceSBarry Smith    Output Parameter:
3659b94acceSBarry Smith .  fnorm - 2-norm of function
3669b94acceSBarry Smith 
3679b94acceSBarry Smith    Note:
3689b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
369a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
370a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3719b94acceSBarry Smith 
3729b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
373a86d99e1SLois Curfman McInnes 
374a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3759b94acceSBarry Smith @*/
3769b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3779b94acceSBarry Smith {
37877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
37974679c65SBarry Smith   PetscValidScalarPointer(fnorm);
38074679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
38174679c65SBarry Smith     SETERRQ(1,"SNESGetFunctionNorm:For SNES_NONLINEAR_EQUATIONS only");
38274679c65SBarry Smith   }
3839b94acceSBarry Smith   *fnorm = snes->norm;
3849b94acceSBarry Smith   return 0;
3859b94acceSBarry Smith }
38674679c65SBarry Smith 
3879b94acceSBarry Smith /*@
3889b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
3899b94acceSBarry Smith    with SNESSSetGradient().
3909b94acceSBarry Smith 
3919b94acceSBarry Smith    Input Parameter:
3929b94acceSBarry Smith .  snes - SNES context
3939b94acceSBarry Smith 
3949b94acceSBarry Smith    Output Parameter:
3959b94acceSBarry Smith .  fnorm - 2-norm of gradient
3969b94acceSBarry Smith 
3979b94acceSBarry Smith    Note:
3989b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
399a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
400a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4019b94acceSBarry Smith 
4029b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
403a86d99e1SLois Curfman McInnes 
404a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4059b94acceSBarry Smith @*/
4069b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
4079b94acceSBarry Smith {
40877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
40974679c65SBarry Smith   PetscValidScalarPointer(gnorm);
41074679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
41174679c65SBarry Smith     SETERRQ(1,"SNESGetGradientNorm:For SNES_UNCONSTRAINED_MINIMIZATION only");
41274679c65SBarry Smith   }
4139b94acceSBarry Smith   *gnorm = snes->norm;
4149b94acceSBarry Smith   return 0;
4159b94acceSBarry Smith }
41674679c65SBarry Smith 
4179b94acceSBarry Smith /*@
4189b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4199b94acceSBarry Smith    attempted by the nonlinear solver.
4209b94acceSBarry Smith 
4219b94acceSBarry Smith    Input Parameter:
4229b94acceSBarry Smith .  snes - SNES context
4239b94acceSBarry Smith 
4249b94acceSBarry Smith    Output Parameter:
4259b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4269b94acceSBarry Smith 
4279b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4289b94acceSBarry Smith @*/
4299b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4309b94acceSBarry Smith {
43177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43274679c65SBarry Smith   PetscValidIntPointer(nfails);
4339b94acceSBarry Smith   *nfails = snes->nfailures;
4349b94acceSBarry Smith   return 0;
4359b94acceSBarry Smith }
4369b94acceSBarry Smith /*@C
4379b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4389b94acceSBarry Smith 
4399b94acceSBarry Smith    Input Parameter:
4409b94acceSBarry Smith .  snes - the SNES context
4419b94acceSBarry Smith 
4429b94acceSBarry Smith    Output Parameter:
4439b94acceSBarry Smith .  sles - the SLES context
4449b94acceSBarry Smith 
4459b94acceSBarry Smith    Notes:
4469b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4479b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4489b94acceSBarry Smith    KSP and PC contexts as well.
4499b94acceSBarry Smith 
4509b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4519b94acceSBarry Smith 
4529b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4539b94acceSBarry Smith @*/
4549b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4559b94acceSBarry Smith {
45677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4579b94acceSBarry Smith   *sles = snes->sles;
4589b94acceSBarry Smith   return 0;
4599b94acceSBarry Smith }
4609b94acceSBarry Smith /* -----------------------------------------------------------*/
4619b94acceSBarry Smith /*@C
4629b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
4639b94acceSBarry Smith 
4649b94acceSBarry Smith    Input Parameter:
4659b94acceSBarry Smith .  comm - MPI communicator
4669b94acceSBarry Smith .  type - type of method, one of
4679b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
4689b94acceSBarry Smith $      (for systems of nonlinear equations)
4699b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
4709b94acceSBarry Smith $      (for unconstrained minimization)
4719b94acceSBarry Smith 
4729b94acceSBarry Smith    Output Parameter:
4739b94acceSBarry Smith .  outsnes - the new SNES context
4749b94acceSBarry Smith 
475c1f60f51SBarry Smith    Options Database Key:
47619bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
47719bd6747SLois Curfman McInnes $              and no preconditioning matrix
47819bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
47919bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
48019bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
48119bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
482c1f60f51SBarry Smith 
4839b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
4849b94acceSBarry Smith 
48563a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
4869b94acceSBarry Smith @*/
4874b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
4889b94acceSBarry Smith {
4899b94acceSBarry Smith   int                 ierr;
4909b94acceSBarry Smith   SNES                snes;
4919b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
49237fcc0dbSBarry Smith 
4939b94acceSBarry Smith   *outsnes = 0;
4942a0acf01SLois Curfman McInnes   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS)
4952a0acf01SLois Curfman McInnes     SETERRQ(1,"SNESCreate:incorrect method type");
4960452661fSBarry Smith   PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm);
4979b94acceSBarry Smith   PLogObjectCreate(snes);
4989b94acceSBarry Smith   snes->max_its           = 50;
4999b94acceSBarry Smith   snes->max_funcs	  = 1000;
5009b94acceSBarry Smith   snes->norm		  = 0.0;
5015a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
502b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
503b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5049b94acceSBarry Smith     snes->atol		  = 1.e-10;
5055a2d0531SBarry Smith   }
506b4874afaSBarry Smith   else {
507b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
508b4874afaSBarry Smith     snes->ttol            = 0.0;
509b4874afaSBarry Smith     snes->atol		  = 1.e-50;
510b4874afaSBarry Smith   }
5119b94acceSBarry Smith   snes->xtol		  = 1.e-8;
512b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5139b94acceSBarry Smith   snes->nfuncs            = 0;
5149b94acceSBarry Smith   snes->nfailures         = 0;
515*639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5169b94acceSBarry Smith   snes->data              = 0;
5179b94acceSBarry Smith   snes->view              = 0;
5189b94acceSBarry Smith   snes->computeumfunction = 0;
5199b94acceSBarry Smith   snes->umfunP            = 0;
5209b94acceSBarry Smith   snes->fc                = 0;
5219b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5229b94acceSBarry Smith   snes->fmin              = -1.e30;
5239b94acceSBarry Smith   snes->method_class      = type;
5249b94acceSBarry Smith   snes->set_method_called = 0;
525a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
5269b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5277d1a2b2bSBarry Smith   snes->type              = -1;
5286f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5296f24a144SLois Curfman McInnes   snes->vwork             = 0;
5306f24a144SLois Curfman McInnes   snes->nwork             = 0;
5319b94acceSBarry Smith 
5329b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5330452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
5349b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5359b94acceSBarry Smith   kctx->version     = 2;
5369b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5379b94acceSBarry Smith                              this was too large for some test cases */
5389b94acceSBarry Smith   kctx->rtol_last   = 0;
5399b94acceSBarry Smith   kctx->rtol_max    = .9;
5409b94acceSBarry Smith   kctx->gamma       = 1.0;
5419b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5429b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5439b94acceSBarry Smith   kctx->threshold   = .1;
5449b94acceSBarry Smith   kctx->lresid_last = 0;
5459b94acceSBarry Smith   kctx->norm_last   = 0;
5469b94acceSBarry Smith 
5479b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5489b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5495334005bSBarry Smith 
5509b94acceSBarry Smith   *outsnes = snes;
5519b94acceSBarry Smith   return 0;
5529b94acceSBarry Smith }
5539b94acceSBarry Smith 
5549b94acceSBarry Smith /* --------------------------------------------------------------- */
5559b94acceSBarry Smith /*@C
5569b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
5579b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
5589b94acceSBarry Smith    equations.
5599b94acceSBarry Smith 
5609b94acceSBarry Smith    Input Parameters:
5619b94acceSBarry Smith .  snes - the SNES context
5629b94acceSBarry Smith .  func - function evaluation routine
5639b94acceSBarry Smith .  r - vector to store function value
5642cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
5652cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
5669b94acceSBarry Smith 
5679b94acceSBarry Smith    Calling sequence of func:
568313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
5699b94acceSBarry Smith 
5709b94acceSBarry Smith .  x - input vector
571313e4042SLois Curfman McInnes .  f - function vector
5722cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
5739b94acceSBarry Smith 
5749b94acceSBarry Smith    Notes:
5759b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
5769b94acceSBarry Smith $      f'(x) x = -f(x),
5779b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
5789b94acceSBarry Smith 
5799b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
5809b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
5819b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
5829b94acceSBarry Smith 
5839b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
5849b94acceSBarry Smith 
585a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
5869b94acceSBarry Smith @*/
5875334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
5889b94acceSBarry Smith {
58977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5909b94acceSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
59148d91487SBarry Smith     "SNESSetFunction:For SNES_NONLINEAR_EQUATIONS only");
5929b94acceSBarry Smith   snes->computefunction     = func;
5939b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
5949b94acceSBarry Smith   snes->funP                = ctx;
5959b94acceSBarry Smith   return 0;
5969b94acceSBarry Smith }
5979b94acceSBarry Smith 
5989b94acceSBarry Smith /*@
5999b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6009b94acceSBarry Smith    SNESSetFunction().
6019b94acceSBarry Smith 
6029b94acceSBarry Smith    Input Parameters:
6039b94acceSBarry Smith .  snes - the SNES context
6049b94acceSBarry Smith .  x - input vector
6059b94acceSBarry Smith 
6069b94acceSBarry Smith    Output Parameter:
6073638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6089b94acceSBarry Smith 
6091bffabb2SLois Curfman McInnes    Notes:
6109b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6119b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6129b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6139b94acceSBarry Smith 
6149b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6159b94acceSBarry Smith 
616a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6179b94acceSBarry Smith @*/
6189b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6199b94acceSBarry Smith {
6209b94acceSBarry Smith   int    ierr;
6219b94acceSBarry Smith 
62274679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
62374679c65SBarry Smith     SETERRQ(1,"SNESComputeFunction: For SNES_NONLINEAR_EQUATIONS only");
6249b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
6259b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
6269b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6279b94acceSBarry Smith   return 0;
6289b94acceSBarry Smith }
6299b94acceSBarry Smith 
6309b94acceSBarry Smith /*@C
6319b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6329b94acceSBarry Smith    unconstrained minimization.
6339b94acceSBarry Smith 
6349b94acceSBarry Smith    Input Parameters:
6359b94acceSBarry Smith .  snes - the SNES context
6369b94acceSBarry Smith .  func - function evaluation routine
637044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
638044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6399b94acceSBarry Smith 
6409b94acceSBarry Smith    Calling sequence of func:
6419b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
6429b94acceSBarry Smith 
6439b94acceSBarry Smith .  x - input vector
6449b94acceSBarry Smith .  f - function
645044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
6469b94acceSBarry Smith 
6479b94acceSBarry Smith    Notes:
6489b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
6499b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
6509b94acceSBarry Smith    SNESSetFunction().
6519b94acceSBarry Smith 
6529b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
6539b94acceSBarry Smith 
654a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
655a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
6569b94acceSBarry Smith @*/
6579b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
6589b94acceSBarry Smith                       void *ctx)
6599b94acceSBarry Smith {
66077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
6619b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
66248d91487SBarry Smith     "SNESSetMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION");
6639b94acceSBarry Smith   snes->computeumfunction   = func;
6649b94acceSBarry Smith   snes->umfunP              = ctx;
6659b94acceSBarry Smith   return 0;
6669b94acceSBarry Smith }
6679b94acceSBarry Smith 
6689b94acceSBarry Smith /*@
6699b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
6709b94acceSBarry Smith    set with SNESSetMinimizationFunction().
6719b94acceSBarry Smith 
6729b94acceSBarry Smith    Input Parameters:
6739b94acceSBarry Smith .  snes - the SNES context
6749b94acceSBarry Smith .  x - input vector
6759b94acceSBarry Smith 
6769b94acceSBarry Smith    Output Parameter:
6779b94acceSBarry Smith .  y - function value
6789b94acceSBarry Smith 
6799b94acceSBarry Smith    Notes:
6809b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
6819b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
6829b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
683a86d99e1SLois Curfman McInnes 
684a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
685a86d99e1SLois Curfman McInnes 
686a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
687a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
6889b94acceSBarry Smith @*/
6899b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
6909b94acceSBarry Smith {
6919b94acceSBarry Smith   int    ierr;
6929b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
69348d91487SBarry Smith     "SNESComputeMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION");
6949b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
6959b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
6969b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
6979b94acceSBarry Smith   return 0;
6989b94acceSBarry Smith }
6999b94acceSBarry Smith 
7009b94acceSBarry Smith /*@C
7019b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7029b94acceSBarry Smith    vector for use by the SNES routines.
7039b94acceSBarry Smith 
7049b94acceSBarry Smith    Input Parameters:
7059b94acceSBarry Smith .  snes - the SNES context
7069b94acceSBarry Smith .  func - function evaluation routine
707044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
708044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7099b94acceSBarry Smith .  r - vector to store gradient value
7109b94acceSBarry Smith 
7119b94acceSBarry Smith    Calling sequence of func:
7129b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7139b94acceSBarry Smith 
7149b94acceSBarry Smith .  x - input vector
7159b94acceSBarry Smith .  g - gradient vector
716044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7179b94acceSBarry Smith 
7189b94acceSBarry Smith    Notes:
7199b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7209b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7219b94acceSBarry Smith    SNESSetFunction().
7229b94acceSBarry Smith 
7239b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7249b94acceSBarry Smith 
725a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
726a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
7279b94acceSBarry Smith @*/
72874679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
7299b94acceSBarry Smith {
73077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
7319b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
73248d91487SBarry Smith     "SNESSetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
7339b94acceSBarry Smith   snes->computefunction     = func;
7349b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7359b94acceSBarry Smith   snes->funP                = ctx;
7369b94acceSBarry Smith   return 0;
7379b94acceSBarry Smith }
7389b94acceSBarry Smith 
7399b94acceSBarry Smith /*@
740a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
741a86d99e1SLois Curfman McInnes    SNESSetGradient().
7429b94acceSBarry Smith 
7439b94acceSBarry Smith    Input Parameters:
7449b94acceSBarry Smith .  snes - the SNES context
7459b94acceSBarry Smith .  x - input vector
7469b94acceSBarry Smith 
7479b94acceSBarry Smith    Output Parameter:
7489b94acceSBarry Smith .  y - gradient vector
7499b94acceSBarry Smith 
7509b94acceSBarry Smith    Notes:
7519b94acceSBarry Smith    SNESComputeGradient() is valid only for
7529b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7539b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
754a86d99e1SLois Curfman McInnes 
755a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
756a86d99e1SLois Curfman McInnes 
757a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
758a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
7599b94acceSBarry Smith @*/
7609b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
7619b94acceSBarry Smith {
7629b94acceSBarry Smith   int    ierr;
76374679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
76474679c65SBarry Smith     SETERRQ(1,"SNESComputeGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
7659b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
7669b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
7679b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
7689b94acceSBarry Smith   return 0;
7699b94acceSBarry Smith }
7709b94acceSBarry Smith 
77162fef451SLois Curfman McInnes /*@
77262fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
77362fef451SLois Curfman McInnes    set with SNESSetJacobian().
77462fef451SLois Curfman McInnes 
77562fef451SLois Curfman McInnes    Input Parameters:
77662fef451SLois Curfman McInnes .  snes - the SNES context
77762fef451SLois Curfman McInnes .  x - input vector
77862fef451SLois Curfman McInnes 
77962fef451SLois Curfman McInnes    Output Parameters:
78062fef451SLois Curfman McInnes .  A - Jacobian matrix
78162fef451SLois Curfman McInnes .  B - optional preconditioning matrix
78262fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
78362fef451SLois Curfman McInnes 
78462fef451SLois Curfman McInnes    Notes:
78562fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
78662fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
78762fef451SLois Curfman McInnes 
788dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
789dc5a77f8SLois Curfman McInnes    flag parameter.
79062fef451SLois Curfman McInnes 
79162fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
79262fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
79362fef451SLois Curfman McInnes    methods is SNESComputeJacobian().
79462fef451SLois Curfman McInnes 
79562fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
79662fef451SLois Curfman McInnes 
79762fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
79862fef451SLois Curfman McInnes @*/
7999b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8009b94acceSBarry Smith {
8019b94acceSBarry Smith   int    ierr;
80274679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
80374679c65SBarry Smith     SETERRQ(1,"SNESComputeJacobian: For SNES_NONLINEAR_EQUATIONS only");
8049b94acceSBarry Smith   if (!snes->computejacobian) return 0;
8059b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
806c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
8079b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
8089b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8096d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
81077c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
81177c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8129b94acceSBarry Smith   return 0;
8139b94acceSBarry Smith }
8149b94acceSBarry Smith 
81562fef451SLois Curfman McInnes /*@
81662fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
81762fef451SLois Curfman McInnes    set with SNESSetHessian().
81862fef451SLois Curfman McInnes 
81962fef451SLois Curfman McInnes    Input Parameters:
82062fef451SLois Curfman McInnes .  snes - the SNES context
82162fef451SLois Curfman McInnes .  x - input vector
82262fef451SLois Curfman McInnes 
82362fef451SLois Curfman McInnes    Output Parameters:
82462fef451SLois Curfman McInnes .  A - Hessian matrix
82562fef451SLois Curfman McInnes .  B - optional preconditioning matrix
82662fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
82762fef451SLois Curfman McInnes 
82862fef451SLois Curfman McInnes    Notes:
82962fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
83062fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
83162fef451SLois Curfman McInnes 
832dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
833dc5a77f8SLois Curfman McInnes    flag parameter.
83462fef451SLois Curfman McInnes 
83562fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
83662fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
83762fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
83862fef451SLois Curfman McInnes 
83962fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
84062fef451SLois Curfman McInnes 
841a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
842a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
84362fef451SLois Curfman McInnes @*/
84462fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
8459b94acceSBarry Smith {
8469b94acceSBarry Smith   int    ierr;
84774679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
84874679c65SBarry Smith     SETERRQ(1,"SNESComputeHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
8499b94acceSBarry Smith   if (!snes->computejacobian) return 0;
85062fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
8510f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
85262fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
85362fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
8540f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
85577c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
85677c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8579b94acceSBarry Smith   return 0;
8589b94acceSBarry Smith }
8599b94acceSBarry Smith 
8609b94acceSBarry Smith /*@C
8619b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
862044dda88SLois Curfman McInnes    location to store the matrix.
8639b94acceSBarry Smith 
8649b94acceSBarry Smith    Input Parameters:
8659b94acceSBarry Smith .  snes - the SNES context
8669b94acceSBarry Smith .  A - Jacobian matrix
8679b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
8689b94acceSBarry Smith .  func - Jacobian evaluation routine
8692cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
8702cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
8719b94acceSBarry Smith 
8729b94acceSBarry Smith    Calling sequence of func:
873313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
8749b94acceSBarry Smith 
8759b94acceSBarry Smith .  x - input vector
8769b94acceSBarry Smith .  A - Jacobian matrix
8779b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
878ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
879ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
8802cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
8819b94acceSBarry Smith 
8829b94acceSBarry Smith    Notes:
883dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
8842cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
885ac21db08SLois Curfman McInnes 
886ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
8879b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
8889b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
8899b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
8909b94acceSBarry Smith    throughout the global iterations.
8919b94acceSBarry Smith 
8929b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
8939b94acceSBarry Smith 
894ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
8959b94acceSBarry Smith @*/
8969b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
8979b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
8989b94acceSBarry Smith {
89977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
90074679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
90174679c65SBarry Smith     SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only");
9029b94acceSBarry Smith   snes->computejacobian = func;
9039b94acceSBarry Smith   snes->jacP            = ctx;
9049b94acceSBarry Smith   snes->jacobian        = A;
9059b94acceSBarry Smith   snes->jacobian_pre    = B;
9069b94acceSBarry Smith   return 0;
9079b94acceSBarry Smith }
90862fef451SLois Curfman McInnes 
909b4fd4287SBarry Smith /*@
910b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
911b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
912b4fd4287SBarry Smith 
913b4fd4287SBarry Smith    Input Parameter:
914b4fd4287SBarry Smith .  snes - the nonlinear solver context
915b4fd4287SBarry Smith 
916b4fd4287SBarry Smith    Output Parameters:
917b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
918b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
919b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
920b4fd4287SBarry Smith 
921b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
922b4fd4287SBarry Smith @*/
923b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
924b4fd4287SBarry Smith {
92574679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
92674679c65SBarry Smith     SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only");
927b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
928b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
929b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
930b4fd4287SBarry Smith   return 0;
931b4fd4287SBarry Smith }
932b4fd4287SBarry Smith 
9339b94acceSBarry Smith /*@C
9349b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
935044dda88SLois Curfman McInnes    location to store the matrix.
9369b94acceSBarry Smith 
9379b94acceSBarry Smith    Input Parameters:
9389b94acceSBarry Smith .  snes - the SNES context
9399b94acceSBarry Smith .  A - Hessian matrix
9409b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
9419b94acceSBarry Smith .  func - Jacobian evaluation routine
942313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
943313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
9449b94acceSBarry Smith 
9459b94acceSBarry Smith    Calling sequence of func:
946313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9479b94acceSBarry Smith 
9489b94acceSBarry Smith .  x - input vector
9499b94acceSBarry Smith .  A - Hessian matrix
9509b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
951ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
952ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9532cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
9549b94acceSBarry Smith 
9559b94acceSBarry Smith    Notes:
956dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9572cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
958ac21db08SLois Curfman McInnes 
9599b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
9609b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
9619b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9629b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9639b94acceSBarry Smith    throughout the global iterations.
9649b94acceSBarry Smith 
9659b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
9669b94acceSBarry Smith 
967ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
9689b94acceSBarry Smith @*/
9699b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9709b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9719b94acceSBarry Smith {
97277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
97374679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
97474679c65SBarry Smith     SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
9759b94acceSBarry Smith   snes->computejacobian = func;
9769b94acceSBarry Smith   snes->jacP            = ctx;
9779b94acceSBarry Smith   snes->jacobian        = A;
9789b94acceSBarry Smith   snes->jacobian_pre    = B;
9799b94acceSBarry Smith   return 0;
9809b94acceSBarry Smith }
9819b94acceSBarry Smith 
98262fef451SLois Curfman McInnes /*@
98362fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
98462fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
98562fef451SLois Curfman McInnes 
98662fef451SLois Curfman McInnes    Input Parameter:
98762fef451SLois Curfman McInnes .  snes - the nonlinear solver context
98862fef451SLois Curfman McInnes 
98962fef451SLois Curfman McInnes    Output Parameters:
99062fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
99162fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
99262fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
99362fef451SLois Curfman McInnes 
99462fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
99562fef451SLois Curfman McInnes @*/
99662fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
99762fef451SLois Curfman McInnes {
99874679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
99974679c65SBarry Smith     SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
100062fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
100162fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
100262fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
100362fef451SLois Curfman McInnes   return 0;
100462fef451SLois Curfman McInnes }
100562fef451SLois Curfman McInnes 
10069b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10079b94acceSBarry Smith 
10089b94acceSBarry Smith /*@
10099b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1010272ac6f2SLois Curfman McInnes    of a nonlinear solver.
10119b94acceSBarry Smith 
10129b94acceSBarry Smith    Input Parameter:
10139b94acceSBarry Smith .  snes - the SNES context
10148ddd3da0SLois Curfman McInnes .  x - the solution vector
10159b94acceSBarry Smith 
1016272ac6f2SLois Curfman McInnes    Notes:
1017272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1018272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1019272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1020272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1021272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1022272ac6f2SLois Curfman McInnes 
10239b94acceSBarry Smith .keywords: SNES, nonlinear, setup
10249b94acceSBarry Smith 
10259b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
10269b94acceSBarry Smith @*/
10278ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
10289b94acceSBarry Smith {
1029272ac6f2SLois Curfman McInnes   int ierr, flg;
103077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
103177c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
10328ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
10339b94acceSBarry Smith 
1034c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1035c1f60f51SBarry Smith   /*
1036c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1037dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1038c1f60f51SBarry Smith   */
1039c1f60f51SBarry Smith   if (flg) {
1040c1f60f51SBarry Smith     Mat J;
1041c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1042c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1043c1f60f51SBarry Smith     snes->mfshell = J;
1044c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1045c1f60f51SBarry Smith       snes->jacobian = J;
1046c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1047c1f60f51SBarry Smith     }
1048c1f60f51SBarry Smith     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1049c1f60f51SBarry Smith       snes->jacobian = J;
1050c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1051c1f60f51SBarry Smith     } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free operator option");
1052c1f60f51SBarry Smith   }
1053272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1054c1f60f51SBarry Smith   /*
1055dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1056c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1057c1f60f51SBarry Smith    */
105831615d04SLois Curfman McInnes   if (flg) {
1059272ac6f2SLois Curfman McInnes     Mat J;
1060272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1061272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1062272ac6f2SLois Curfman McInnes     snes->mfshell = J;
106383e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
106483e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
106594a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
106683e56afdSLois Curfman McInnes     }
106783e56afdSLois Curfman McInnes     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
106883e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
106994a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
107083e56afdSLois Curfman McInnes     } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free option");
1071272ac6f2SLois Curfman McInnes   }
10729b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
107337fcc0dbSBarry Smith     if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first");
107437fcc0dbSBarry Smith     if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first");
107537fcc0dbSBarry Smith     if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetJacobian() first");
1076d804570eSBarry Smith     if (snes->vec_func == snes->vec_sol) SETERRQ(1,"SNESSetUp:Solution vector cannot be function vector");
1077a703fe33SLois Curfman McInnes 
1078a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
107940191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1080a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1081a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1082a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1083a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1084a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1085a703fe33SLois Curfman McInnes     }
10869b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
108737fcc0dbSBarry Smith     if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first");
108837fcc0dbSBarry Smith     if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first");
108937fcc0dbSBarry Smith     if (!snes->computeumfunction)
109037fcc0dbSBarry Smith       SETERRQ(1,"SNESSetUp:Must call SNESSetMinimizationFunction() first");
109137fcc0dbSBarry Smith     if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetHessian() first");
10929b94acceSBarry Smith   } else SETERRQ(1,"SNESSetUp:Unknown method class");
1093a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1094a703fe33SLois Curfman McInnes   snes->setup_called = 1;
1095a703fe33SLois Curfman McInnes   return 0;
10969b94acceSBarry Smith }
10979b94acceSBarry Smith 
10989b94acceSBarry Smith /*@C
10999b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11009b94acceSBarry Smith    with SNESCreate().
11019b94acceSBarry Smith 
11029b94acceSBarry Smith    Input Parameter:
11039b94acceSBarry Smith .  snes - the SNES context
11049b94acceSBarry Smith 
11059b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
11069b94acceSBarry Smith 
110763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
11089b94acceSBarry Smith @*/
11099b94acceSBarry Smith int SNESDestroy(SNES snes)
11109b94acceSBarry Smith {
11119b94acceSBarry Smith   int ierr;
111277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
11139b94acceSBarry Smith   ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);
11140452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
11159b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
11169b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
11173e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
11186f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
11199b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
11200452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
11219b94acceSBarry Smith   return 0;
11229b94acceSBarry Smith }
11239b94acceSBarry Smith 
11249b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11259b94acceSBarry Smith 
11269b94acceSBarry Smith /*@
1127d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11289b94acceSBarry Smith 
11299b94acceSBarry Smith    Input Parameters:
11309b94acceSBarry Smith .  snes - the SNES context
113133174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
113233174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
113333174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
113433174efeSLois Curfman McInnes            of the change in the solution between steps
113533174efeSLois Curfman McInnes .  maxit - maximum number of iterations
113633174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
11379b94acceSBarry Smith 
113833174efeSLois Curfman McInnes    Options Database Keys:
113933174efeSLois Curfman McInnes $    -snes_atol <atol>
114033174efeSLois Curfman McInnes $    -snes_rtol <rtol>
114133174efeSLois Curfman McInnes $    -snes_stol <stol>
114233174efeSLois Curfman McInnes $    -snes_max_it <maxit>
114333174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
11449b94acceSBarry Smith 
1145d7a720efSLois Curfman McInnes    Notes:
11469b94acceSBarry Smith    The default maximum number of iterations is 50.
11479b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
11489b94acceSBarry Smith 
114933174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
11509b94acceSBarry Smith 
1151d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
11529b94acceSBarry Smith @*/
1153d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
11549b94acceSBarry Smith {
115577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1156d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1157d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1158d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1159d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1160d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
11619b94acceSBarry Smith   return 0;
11629b94acceSBarry Smith }
11639b94acceSBarry Smith 
11649b94acceSBarry Smith /*@
116533174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
116633174efeSLois Curfman McInnes 
116733174efeSLois Curfman McInnes    Input Parameters:
116833174efeSLois Curfman McInnes .  snes - the SNES context
116933174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
117033174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
117133174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
117233174efeSLois Curfman McInnes            of the change in the solution between steps
117333174efeSLois Curfman McInnes .  maxit - maximum number of iterations
117433174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
117533174efeSLois Curfman McInnes 
117633174efeSLois Curfman McInnes    Notes:
117733174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
117833174efeSLois Curfman McInnes 
117933174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
118033174efeSLois Curfman McInnes 
118133174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
118233174efeSLois Curfman McInnes @*/
118333174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
118433174efeSLois Curfman McInnes {
118533174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
118633174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
118733174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
118833174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
118933174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
119033174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
119133174efeSLois Curfman McInnes   return 0;
119233174efeSLois Curfman McInnes }
119333174efeSLois Curfman McInnes 
119433174efeSLois Curfman McInnes /*@
11959b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
11969b94acceSBarry Smith 
11979b94acceSBarry Smith    Input Parameters:
11989b94acceSBarry Smith .  snes - the SNES context
11999b94acceSBarry Smith .  tol - tolerance
12009b94acceSBarry Smith 
12019b94acceSBarry Smith    Options Database Key:
1202d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
12039b94acceSBarry Smith 
12049b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12059b94acceSBarry Smith 
1206d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
12079b94acceSBarry Smith @*/
12089b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
12099b94acceSBarry Smith {
121077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12119b94acceSBarry Smith   snes->deltatol = tol;
12129b94acceSBarry Smith   return 0;
12139b94acceSBarry Smith }
12149b94acceSBarry Smith 
12159b94acceSBarry Smith /*@
121677c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
12179b94acceSBarry Smith    for unconstrained minimization solvers.
12189b94acceSBarry Smith 
12199b94acceSBarry Smith    Input Parameters:
12209b94acceSBarry Smith .  snes - the SNES context
12219b94acceSBarry Smith .  ftol - minimum function tolerance
12229b94acceSBarry Smith 
12239b94acceSBarry Smith    Options Database Key:
1224d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
12259b94acceSBarry Smith 
12269b94acceSBarry Smith    Note:
122777c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
12289b94acceSBarry Smith    methods only.
12299b94acceSBarry Smith 
12309b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
12319b94acceSBarry Smith 
1232d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
12339b94acceSBarry Smith @*/
123477c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
12359b94acceSBarry Smith {
123677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12379b94acceSBarry Smith   snes->fmin = ftol;
12389b94acceSBarry Smith   return 0;
12399b94acceSBarry Smith }
12409b94acceSBarry Smith 
12419b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
12429b94acceSBarry Smith 
12439b94acceSBarry Smith /*@C
12449b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
12459b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
12469b94acceSBarry Smith    progress.
12479b94acceSBarry Smith 
12489b94acceSBarry Smith    Input Parameters:
12499b94acceSBarry Smith .  snes - the SNES context
12509b94acceSBarry Smith .  func - monitoring routine
1251044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1252044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
12539b94acceSBarry Smith 
12549b94acceSBarry Smith    Calling sequence of func:
1255682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
12569b94acceSBarry Smith 
12579b94acceSBarry Smith $    snes - the SNES context
12589b94acceSBarry Smith $    its - iteration number
1259044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
12609b94acceSBarry Smith $
12619b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
12629b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
12639b94acceSBarry Smith $
12649b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
12659b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
12669b94acceSBarry Smith 
1267*639f9d9dSBarry Smith    Notes:
1268*639f9d9dSBarry Smith    Several different monitor routines may be set and all will be called.
1269*639f9d9dSBarry Smith 
12709b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
12719b94acceSBarry Smith 
12729b94acceSBarry Smith .seealso: SNESDefaultMonitor()
12739b94acceSBarry Smith @*/
127474679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
12759b94acceSBarry Smith {
1276*639f9d9dSBarry Smith   if (!func) {
1277*639f9d9dSBarry Smith     snes->numbermonitors = 0;
1278*639f9d9dSBarry Smith     return 0;
1279*639f9d9dSBarry Smith   }
1280*639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1281*639f9d9dSBarry Smith     SETERRQ(1,"SNESSetMonitor:Too many monitors set");
1282*639f9d9dSBarry Smith   }
1283*639f9d9dSBarry Smith 
1284*639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1285*639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
12869b94acceSBarry Smith   return 0;
12879b94acceSBarry Smith }
12889b94acceSBarry Smith 
12899b94acceSBarry Smith /*@C
12909b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
12919b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
12929b94acceSBarry Smith 
12939b94acceSBarry Smith    Input Parameters:
12949b94acceSBarry Smith .  snes - the SNES context
12959b94acceSBarry Smith .  func - routine to test for convergence
1296044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1297044dda88SLois Curfman McInnes           (may be PETSC_NULL)
12989b94acceSBarry Smith 
12999b94acceSBarry Smith    Calling sequence of func:
13009b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
13019b94acceSBarry Smith              double f,void *cctx)
13029b94acceSBarry Smith 
13039b94acceSBarry Smith $    snes - the SNES context
1304044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
13059b94acceSBarry Smith $    xnorm - 2-norm of current iterate
13069b94acceSBarry Smith $
13079b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13089b94acceSBarry Smith $    gnorm - 2-norm of current step
13099b94acceSBarry Smith $    f - 2-norm of function
13109b94acceSBarry Smith $
13119b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13129b94acceSBarry Smith $    gnorm - 2-norm of current gradient
13139b94acceSBarry Smith $    f - function value
13149b94acceSBarry Smith 
13159b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
13169b94acceSBarry Smith 
131740191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
131840191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
13199b94acceSBarry Smith @*/
132074679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
13219b94acceSBarry Smith {
13229b94acceSBarry Smith   (snes)->converged = func;
13239b94acceSBarry Smith   (snes)->cnvP      = cctx;
13249b94acceSBarry Smith   return 0;
13259b94acceSBarry Smith }
13269b94acceSBarry Smith 
13279b94acceSBarry Smith /*
13289b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
13299b94acceSBarry Smith    positive parameter delta.
13309b94acceSBarry Smith 
13319b94acceSBarry Smith     Input Parameters:
13329b94acceSBarry Smith .   snes - the SNES context
13339b94acceSBarry Smith .   y - approximate solution of linear system
13349b94acceSBarry Smith .   fnorm - 2-norm of current function
13359b94acceSBarry Smith .   delta - trust region size
13369b94acceSBarry Smith 
13379b94acceSBarry Smith     Output Parameters:
13389b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
13399b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
13409b94acceSBarry Smith     region, and exceeds zero otherwise.
13419b94acceSBarry Smith .   ynorm - 2-norm of the step
13429b94acceSBarry Smith 
13439b94acceSBarry Smith     Note:
134440191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
13459b94acceSBarry Smith     is set to be the maximum allowable step size.
13469b94acceSBarry Smith 
13479b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
13489b94acceSBarry Smith */
13499b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
13509b94acceSBarry Smith                           double *gpnorm,double *ynorm)
13519b94acceSBarry Smith {
13529b94acceSBarry Smith   double norm;
13539b94acceSBarry Smith   Scalar cnorm;
1354cddf8d76SBarry Smith   VecNorm(y,NORM_2, &norm );
13559b94acceSBarry Smith   if (norm > *delta) {
13569b94acceSBarry Smith      norm = *delta/norm;
13579b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
13589b94acceSBarry Smith      cnorm = norm;
13599b94acceSBarry Smith      VecScale( &cnorm, y );
13609b94acceSBarry Smith      *ynorm = *delta;
13619b94acceSBarry Smith   } else {
13629b94acceSBarry Smith      *gpnorm = 0.0;
13639b94acceSBarry Smith      *ynorm = norm;
13649b94acceSBarry Smith   }
13659b94acceSBarry Smith   return 0;
13669b94acceSBarry Smith }
13679b94acceSBarry Smith 
13689b94acceSBarry Smith /*@
13699b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
137063a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
13719b94acceSBarry Smith 
13729b94acceSBarry Smith    Input Parameter:
13739b94acceSBarry Smith .  snes - the SNES context
13748ddd3da0SLois Curfman McInnes .  x - the solution vector
13759b94acceSBarry Smith 
13769b94acceSBarry Smith    Output Parameter:
13779b94acceSBarry Smith    its - number of iterations until termination
13789b94acceSBarry Smith 
13798ddd3da0SLois Curfman McInnes    Note:
13808ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
13818ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
13828ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
13838ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
13848ddd3da0SLois Curfman McInnes 
13859b94acceSBarry Smith .keywords: SNES, nonlinear, solve
13869b94acceSBarry Smith 
138763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
13889b94acceSBarry Smith @*/
13898ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
13909b94acceSBarry Smith {
13913c7409f5SSatish Balay   int ierr, flg;
1392052efed2SBarry Smith 
139377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
139474679c65SBarry Smith   PetscValidIntPointer(its);
1395c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1396c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
13979b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
13989b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
13999b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
14003c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
14016d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
14029b94acceSBarry Smith   return 0;
14039b94acceSBarry Smith }
14049b94acceSBarry Smith 
14059b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
140637fcc0dbSBarry Smith static NRList *__SNESList = 0;
14079b94acceSBarry Smith 
14089b94acceSBarry Smith /*@
14094b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
14109b94acceSBarry Smith 
14119b94acceSBarry Smith    Input Parameters:
14129b94acceSBarry Smith .  snes - the SNES context
14139b94acceSBarry Smith .  method - a known method
14149b94acceSBarry Smith 
14159b94acceSBarry Smith    Notes:
14169b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
14179b94acceSBarry Smith $  Systems of nonlinear equations:
141840191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
141940191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
14209b94acceSBarry Smith $  Unconstrained minimization:
142140191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
142240191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
14239b94acceSBarry Smith 
14249b94acceSBarry Smith   Options Database Command:
14254b0e389bSBarry Smith $ -snes_type  <method>
14269b94acceSBarry Smith $    Use -help for a list of available methods
14279b94acceSBarry Smith $    (for instance, ls or tr)
1428a703fe33SLois Curfman McInnes 
1429f536c699SSatish Balay .keywords: SNES, set, method
14309b94acceSBarry Smith @*/
14314b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
14329b94acceSBarry Smith {
14339b94acceSBarry Smith   int (*r)(SNES);
143477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
14357d1a2b2bSBarry Smith   if (snes->type == (int) method) return 0;
14367d1a2b2bSBarry Smith 
14379b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
143837fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
143937fcc0dbSBarry Smith   if (!__SNESList) {SETERRQ(1,"SNESSetType:Could not get methods");}
144037fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
14414b0e389bSBarry Smith   if (!r) {SETERRQ(1,"SNESSetType:Unknown method");}
14420452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
14439b94acceSBarry Smith   snes->set_method_called = 1;
14449b94acceSBarry Smith   return (*r)(snes);
14459b94acceSBarry Smith }
14469b94acceSBarry Smith 
14479b94acceSBarry Smith /* --------------------------------------------------------------------- */
14489b94acceSBarry Smith /*@C
14499b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
14504b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
14519b94acceSBarry Smith 
14529b94acceSBarry Smith    Input Parameters:
145340191667SLois Curfman McInnes .  name - for instance SNES_EQ_LS, SNES_EQ_TR, ...
145440191667SLois Curfman McInnes .  sname - corresponding string for name
14559b94acceSBarry Smith .  create - routine to create method context
14569b94acceSBarry Smith 
14579b94acceSBarry Smith .keywords: SNES, nonlinear, register
14589b94acceSBarry Smith 
14599b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
14609b94acceSBarry Smith @*/
14619b94acceSBarry Smith int SNESRegister(int name, char *sname, int (*create)(SNES))
14629b94acceSBarry Smith {
14639b94acceSBarry Smith   int ierr;
146437fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
146537fcc0dbSBarry Smith   NRRegister( __SNESList, name, sname, (int (*)(void*))create );
14669b94acceSBarry Smith   return 0;
14679b94acceSBarry Smith }
14689b94acceSBarry Smith /* --------------------------------------------------------------------- */
14699b94acceSBarry Smith /*@C
14709b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
14719b94acceSBarry Smith    registered by SNESRegister().
14729b94acceSBarry Smith 
14739b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
14749b94acceSBarry Smith 
14759b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
14769b94acceSBarry Smith @*/
14779b94acceSBarry Smith int SNESRegisterDestroy()
14789b94acceSBarry Smith {
147937fcc0dbSBarry Smith   if (__SNESList) {
148037fcc0dbSBarry Smith     NRDestroy( __SNESList );
148137fcc0dbSBarry Smith     __SNESList = 0;
14829b94acceSBarry Smith   }
14839b94acceSBarry Smith   return 0;
14849b94acceSBarry Smith }
14859b94acceSBarry Smith 
14869b94acceSBarry Smith /*
14874b0e389bSBarry Smith    SNESGetTypeFromOptions_Private - Sets the selected method from the
14889b94acceSBarry Smith    options database.
14899b94acceSBarry Smith 
14909b94acceSBarry Smith    Input Parameter:
14919b94acceSBarry Smith .  ctx - the SNES context
14929b94acceSBarry Smith 
14939b94acceSBarry Smith    Output Parameter:
14949b94acceSBarry Smith .  method -  solver method
14959b94acceSBarry Smith 
14969b94acceSBarry Smith    Returns:
14979b94acceSBarry Smith    Returns 1 if the method is found; 0 otherwise.
14989b94acceSBarry Smith 
14999b94acceSBarry Smith    Options Database Key:
15004b0e389bSBarry Smith $  -snes_type  method
15019b94acceSBarry Smith */
1502052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
15039b94acceSBarry Smith {
1504052efed2SBarry Smith   int  ierr;
15059b94acceSBarry Smith   char sbuf[50];
15065baf8537SBarry Smith 
1507052efed2SBarry Smith   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1508052efed2SBarry Smith   if (*flg) {
1509052efed2SBarry Smith     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
151037fcc0dbSBarry Smith     *method = (SNESType)NRFindID( __SNESList, sbuf );
15119b94acceSBarry Smith   }
15129b94acceSBarry Smith   return 0;
15139b94acceSBarry Smith }
15149b94acceSBarry Smith 
15159b94acceSBarry Smith /*@C
15169a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
15179b94acceSBarry Smith 
15189b94acceSBarry Smith    Input Parameter:
15194b0e389bSBarry Smith .  snes - nonlinear solver context
15209b94acceSBarry Smith 
15219b94acceSBarry Smith    Output Parameter:
15229a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
15239a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
15249b94acceSBarry Smith 
15259b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
15269b94acceSBarry Smith @*/
15274b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
15289b94acceSBarry Smith {
152906a9b0adSLois Curfman McInnes   int ierr;
153037fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
15314b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
153237fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
15339b94acceSBarry Smith   return 0;
15349b94acceSBarry Smith }
15359b94acceSBarry Smith 
15369b94acceSBarry Smith #include <stdio.h>
15379b94acceSBarry Smith /*
15384b0e389bSBarry Smith    SNESPrintTypes_Private - Prints the SNES methods available from the
15399b94acceSBarry Smith    options database.
15409b94acceSBarry Smith 
15419b94acceSBarry Smith    Input Parameters:
154233455573SSatish Balay .  comm   - communicator (usually MPI_COMM_WORLD)
15439b94acceSBarry Smith .  prefix - prefix (usually "-")
15444b0e389bSBarry Smith .  name   - the options database name (by default "snes_type")
15459b94acceSBarry Smith */
154633455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
15479b94acceSBarry Smith {
15489b94acceSBarry Smith   FuncList *entry;
154937fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
155037fcc0dbSBarry Smith   entry = __SNESList->head;
155177c4ece6SBarry Smith   PetscPrintf(comm," %s%s (one of)",prefix,name);
15529b94acceSBarry Smith   while (entry) {
155377c4ece6SBarry Smith     PetscPrintf(comm," %s",entry->name);
15549b94acceSBarry Smith     entry = entry->next;
15559b94acceSBarry Smith   }
155677c4ece6SBarry Smith   PetscPrintf(comm,"\n");
15579b94acceSBarry Smith   return 0;
15589b94acceSBarry Smith }
15599b94acceSBarry Smith 
15609b94acceSBarry Smith /*@C
15619b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
15629b94acceSBarry Smith    stored.
15639b94acceSBarry Smith 
15649b94acceSBarry Smith    Input Parameter:
15659b94acceSBarry Smith .  snes - the SNES context
15669b94acceSBarry Smith 
15679b94acceSBarry Smith    Output Parameter:
15689b94acceSBarry Smith .  x - the solution
15699b94acceSBarry Smith 
15709b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
15719b94acceSBarry Smith 
1572a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
15739b94acceSBarry Smith @*/
15749b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
15759b94acceSBarry Smith {
157677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15779b94acceSBarry Smith   *x = snes->vec_sol_always;
15789b94acceSBarry Smith   return 0;
15799b94acceSBarry Smith }
15809b94acceSBarry Smith 
15819b94acceSBarry Smith /*@C
15829b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
15839b94acceSBarry Smith    stored.
15849b94acceSBarry Smith 
15859b94acceSBarry Smith    Input Parameter:
15869b94acceSBarry Smith .  snes - the SNES context
15879b94acceSBarry Smith 
15889b94acceSBarry Smith    Output Parameter:
15899b94acceSBarry Smith .  x - the solution update
15909b94acceSBarry Smith 
15919b94acceSBarry Smith    Notes:
15929b94acceSBarry Smith    This vector is implementation dependent.
15939b94acceSBarry Smith 
15949b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
15959b94acceSBarry Smith 
15969b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
15979b94acceSBarry Smith @*/
15989b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
15999b94acceSBarry Smith {
160077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16019b94acceSBarry Smith   *x = snes->vec_sol_update_always;
16029b94acceSBarry Smith   return 0;
16039b94acceSBarry Smith }
16049b94acceSBarry Smith 
16059b94acceSBarry Smith /*@C
16063638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
16079b94acceSBarry Smith 
16089b94acceSBarry Smith    Input Parameter:
16099b94acceSBarry Smith .  snes - the SNES context
16109b94acceSBarry Smith 
16119b94acceSBarry Smith    Output Parameter:
16123638b69dSLois Curfman McInnes .  r - the function
16139b94acceSBarry Smith 
16149b94acceSBarry Smith    Notes:
16159b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
16169b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
16179b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
16189b94acceSBarry Smith 
1619a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
16209b94acceSBarry Smith 
162131615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
162231615d04SLois Curfman McInnes           SNESGetGradient()
16239b94acceSBarry Smith @*/
16249b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
16259b94acceSBarry Smith {
162677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16279b94acceSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
162848d91487SBarry Smith     "SNESGetFunction:For SNES_NONLINEAR_EQUATIONS only");
16299b94acceSBarry Smith   *r = snes->vec_func_always;
16309b94acceSBarry Smith   return 0;
16319b94acceSBarry Smith }
16329b94acceSBarry Smith 
16339b94acceSBarry Smith /*@C
16343638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
16359b94acceSBarry Smith 
16369b94acceSBarry Smith    Input Parameter:
16379b94acceSBarry Smith .  snes - the SNES context
16389b94acceSBarry Smith 
16399b94acceSBarry Smith    Output Parameter:
16409b94acceSBarry Smith .  r - the gradient
16419b94acceSBarry Smith 
16429b94acceSBarry Smith    Notes:
16439b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
16449b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
16459b94acceSBarry Smith    SNESGetFunction().
16469b94acceSBarry Smith 
16479b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
16489b94acceSBarry Smith 
164931615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
16509b94acceSBarry Smith @*/
16519b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
16529b94acceSBarry Smith {
165377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16549b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
165548d91487SBarry Smith     "SNESGetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
16569b94acceSBarry Smith   *r = snes->vec_func_always;
16579b94acceSBarry Smith   return 0;
16589b94acceSBarry Smith }
16599b94acceSBarry Smith 
16609b94acceSBarry Smith /*@
16619b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
16629b94acceSBarry Smith    unconstrained minimization problems.
16639b94acceSBarry Smith 
16649b94acceSBarry Smith    Input Parameter:
16659b94acceSBarry Smith .  snes - the SNES context
16669b94acceSBarry Smith 
16679b94acceSBarry Smith    Output Parameter:
16689b94acceSBarry Smith .  r - the function
16699b94acceSBarry Smith 
16709b94acceSBarry Smith    Notes:
16719b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
16729b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
16739b94acceSBarry Smith    SNESGetFunction().
16749b94acceSBarry Smith 
16759b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
16769b94acceSBarry Smith 
167731615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
16789b94acceSBarry Smith @*/
16799b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
16809b94acceSBarry Smith {
168177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
168274679c65SBarry Smith   PetscValidScalarPointer(r);
16839b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
168448d91487SBarry Smith     "SNESGetMinimizationFunction:For SNES_UNCONSTRAINED_MINIMIZATION only");
16859b94acceSBarry Smith   *r = snes->fc;
16869b94acceSBarry Smith   return 0;
16879b94acceSBarry Smith }
16889b94acceSBarry Smith 
16893c7409f5SSatish Balay /*@C
16903c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1691*639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1692*639f9d9dSBarry Smith    the prefix name.
16933c7409f5SSatish Balay 
16943c7409f5SSatish Balay    Input Parameter:
16953c7409f5SSatish Balay .  snes - the SNES context
16963c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
16973c7409f5SSatish Balay 
16983c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1699a86d99e1SLois Curfman McInnes 
1700a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
17013c7409f5SSatish Balay @*/
17023c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
17033c7409f5SSatish Balay {
17043c7409f5SSatish Balay   int ierr;
17053c7409f5SSatish Balay 
170677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1707*639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
17083c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
17093c7409f5SSatish Balay   return 0;
17103c7409f5SSatish Balay }
17113c7409f5SSatish Balay 
17123c7409f5SSatish Balay /*@C
1713f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1714*639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1715*639f9d9dSBarry Smith    the prefix name.
17163c7409f5SSatish Balay 
17173c7409f5SSatish Balay    Input Parameter:
17183c7409f5SSatish Balay .  snes - the SNES context
17193c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
17203c7409f5SSatish Balay 
17213c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1722a86d99e1SLois Curfman McInnes 
1723a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
17243c7409f5SSatish Balay @*/
17253c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
17263c7409f5SSatish Balay {
17273c7409f5SSatish Balay   int ierr;
17283c7409f5SSatish Balay 
172977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1730*639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
17313c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
17323c7409f5SSatish Balay   return 0;
17333c7409f5SSatish Balay }
17343c7409f5SSatish Balay 
1735c83590e2SSatish Balay /*@
17363c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
17373c7409f5SSatish Balay    SNES options in the database.
17383c7409f5SSatish Balay 
17393c7409f5SSatish Balay    Input Parameter:
17403c7409f5SSatish Balay .  snes - the SNES context
17413c7409f5SSatish Balay 
17423c7409f5SSatish Balay    Output Parameter:
17433c7409f5SSatish Balay .  prefix - pointer to the prefix string used
17443c7409f5SSatish Balay 
17453c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1746a86d99e1SLois Curfman McInnes 
1747a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
17483c7409f5SSatish Balay @*/
17493c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
17503c7409f5SSatish Balay {
17513c7409f5SSatish Balay   int ierr;
17523c7409f5SSatish Balay 
175377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1754*639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
17553c7409f5SSatish Balay   return 0;
17563c7409f5SSatish Balay }
17573c7409f5SSatish Balay 
17583c7409f5SSatish Balay 
17599b94acceSBarry Smith 
17609b94acceSBarry Smith 
17619b94acceSBarry Smith 
1762