xref: /petsc/src/snes/interface/snes.c (revision 5bbad29ba5a66636b67f9a8bba6f70788beb01fe)
19b94acceSBarry Smith #ifndef lint
2*5bbad29bSBarry Smith static char vcid[] = "$Id: snes.c,v 1.97 1996/11/13 15:41:22 curfman 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 
89639f9d9dSBarry Smith /*
90639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
91639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
92639f9d9dSBarry Smith */
93639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
94639f9d9dSBarry Smith static int numberofsetfromoptions;
95639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
96639f9d9dSBarry Smith 
97639f9d9dSBarry Smith /*@
98639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
99639f9d9dSBarry Smith 
100639f9d9dSBarry Smith   Input Parameter:
101639f9d9dSBarry Smith .   snescheck - function that checks for options
102639f9d9dSBarry Smith 
103639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
104639f9d9dSBarry Smith @*/
105639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
106639f9d9dSBarry Smith {
107639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
108639f9d9dSBarry Smith     SETERRQ(1,"SNESAddOptionsChecker:Too many options checkers, only 5 allowed");
109639f9d9dSBarry Smith   }
110639f9d9dSBarry Smith 
111639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
112639f9d9dSBarry Smith   return 0;
113639f9d9dSBarry Smith }
114639f9d9dSBarry 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;
130639f9d9dSBarry 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);
179*5bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
180*5bbad29bSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);}
1813c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
182639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
1833c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
184d31a3109SLois Curfman McInnes   if (flg) {
185639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
186d31a3109SLois Curfman McInnes   }
1873c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_xmonitor",&flg); CHKERRQ(ierr);
1883c7409f5SSatish Balay   if (flg) {
18917699dbbSLois Curfman McInnes     int    rank = 0;
190d7e8b826SBarry Smith     DrawLG lg;
19117699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
19217699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
19317699dbbSLois Curfman McInnes     if (!rank) {
1949b94acceSBarry Smith       ierr = SNESLGMonitorCreate(0,0,0,0,300,300,&lg); CHKERRQ(ierr);
1959b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
196f8c826e1SBarry Smith       snes->xmonitor = lg;
1979b94acceSBarry Smith       PLogObjectParent(snes,lg);
1989b94acceSBarry Smith     }
1999b94acceSBarry Smith   }
2003c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2013c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2029b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2039b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
20494a424c1SBarry Smith     PLogInfo(snes,
20531615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
20631615d04SLois Curfman McInnes   }
20731615d04SLois Curfman McInnes   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
20831615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
20931615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
21094a424c1SBarry Smith     PLogInfo(snes,
21131615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2129b94acceSBarry Smith   }
213639f9d9dSBarry Smith 
214639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
215639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
216639f9d9dSBarry Smith   }
217639f9d9dSBarry Smith 
2189b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2199b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
2209b94acceSBarry Smith   if (!snes->setfromoptions) return 0;
2219b94acceSBarry Smith   return (*snes->setfromoptions)(snes);
2229b94acceSBarry Smith }
2239b94acceSBarry Smith 
2249b94acceSBarry Smith /*@
2259b94acceSBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2269b94acceSBarry Smith 
2279b94acceSBarry Smith    Input Parameter:
2289b94acceSBarry Smith .  snes - the SNES context
2299b94acceSBarry Smith 
230a703fe33SLois Curfman McInnes    Options Database Keys:
231a703fe33SLois Curfman McInnes $  -help, -h
232a703fe33SLois Curfman McInnes 
2339b94acceSBarry Smith .keywords: SNES, nonlinear, help
2349b94acceSBarry Smith 
235682d7d0cSBarry Smith .seealso: SNESSetFromOptions()
2369b94acceSBarry Smith @*/
2379b94acceSBarry Smith int SNESPrintHelp(SNES snes)
2389b94acceSBarry Smith {
2396daaf66cSBarry Smith   char                p[64];
2406daaf66cSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2416daaf66cSBarry Smith 
24277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2436daaf66cSBarry Smith 
2446daaf66cSBarry Smith   PetscStrcpy(p,"-");
2456daaf66cSBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2466daaf66cSBarry Smith 
2476daaf66cSBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2486daaf66cSBarry Smith 
249d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
25033455573SSatish Balay   SNESPrintTypes_Private(snes->comm,p,"snes_type");
25177c4ece6SBarry Smith   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
252b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
253b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
254b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
255b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
256b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
257b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
258*5bbad29bSBarry Smith   PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n");
259*5bbad29bSBarry Smith   PetscPrintf(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
260d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor\n",p);
261d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
262b18e04deSLois Curfman McInnes   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
26377c4ece6SBarry Smith     PetscPrintf(snes->comm,
264d31a3109SLois Curfman McInnes      " Options for solving systems of nonlinear equations only:\n");
26577c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
26677c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
267ef1dfb11SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
2681650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p);
2691650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p);
27077c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
27177c4ece6SBarry Smith     PetscPrintf(snes->comm,
272b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
27377c4ece6SBarry Smith     PetscPrintf(snes->comm,
274b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
27577c4ece6SBarry Smith     PetscPrintf(snes->comm,
276b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
27777c4ece6SBarry Smith     PetscPrintf(snes->comm,
278b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
27977c4ece6SBarry Smith     PetscPrintf(snes->comm,
280b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
28177c4ece6SBarry Smith     PetscPrintf(snes->comm,
282b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
28377c4ece6SBarry Smith     PetscPrintf(snes->comm,
284b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
285b18e04deSLois Curfman McInnes   }
286b18e04deSLois Curfman McInnes   else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
28777c4ece6SBarry Smith     PetscPrintf(snes->comm,
288d31a3109SLois Curfman McInnes      " Options for solving unconstrained minimization problems only:\n");
289b18e04deSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
29077c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
29177c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
292d31a3109SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p);
293d31a3109SLois Curfman McInnes      PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p);
294b18e04deSLois Curfman McInnes   }
2954537a946SLois Curfman McInnes   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
29677c4ece6SBarry Smith   PetscPrintf(snes->comm,"a particular method\n");
2976daaf66cSBarry Smith   if (snes->printhelp) (*snes->printhelp)(snes,p);
2989b94acceSBarry Smith   return 0;
2999b94acceSBarry Smith }
3009b94acceSBarry Smith /*@
3019b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3029b94acceSBarry Smith    the nonlinear solvers.
3039b94acceSBarry Smith 
3049b94acceSBarry Smith    Input Parameters:
3059b94acceSBarry Smith .  snes - the SNES context
3069b94acceSBarry Smith .  usrP - optional user context
3079b94acceSBarry Smith 
3089b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3099b94acceSBarry Smith 
3109b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3119b94acceSBarry Smith @*/
3129b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3139b94acceSBarry Smith {
31477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3159b94acceSBarry Smith   snes->user		= usrP;
3169b94acceSBarry Smith   return 0;
3179b94acceSBarry Smith }
31874679c65SBarry Smith 
3199b94acceSBarry Smith /*@C
3209b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3219b94acceSBarry Smith    nonlinear solvers.
3229b94acceSBarry Smith 
3239b94acceSBarry Smith    Input Parameter:
3249b94acceSBarry Smith .  snes - SNES context
3259b94acceSBarry Smith 
3269b94acceSBarry Smith    Output Parameter:
3279b94acceSBarry Smith .  usrP - user context
3289b94acceSBarry Smith 
3299b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3309b94acceSBarry Smith 
3319b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3329b94acceSBarry Smith @*/
3339b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3349b94acceSBarry Smith {
33577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3369b94acceSBarry Smith   *usrP = snes->user;
3379b94acceSBarry Smith   return 0;
3389b94acceSBarry Smith }
33974679c65SBarry Smith 
3409b94acceSBarry Smith /*@
3419b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3429b94acceSBarry Smith    nonlinear solver.
3439b94acceSBarry Smith 
3449b94acceSBarry Smith    Input Parameter:
3459b94acceSBarry Smith .  snes - SNES context
3469b94acceSBarry Smith 
3479b94acceSBarry Smith    Output Parameter:
3489b94acceSBarry Smith .  iter - iteration number
3499b94acceSBarry Smith 
3509b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3519b94acceSBarry Smith @*/
3529b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3539b94acceSBarry Smith {
35477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
35574679c65SBarry Smith   PetscValidIntPointer(iter);
3569b94acceSBarry Smith   *iter = snes->iter;
3579b94acceSBarry Smith   return 0;
3589b94acceSBarry Smith }
35974679c65SBarry Smith 
3609b94acceSBarry Smith /*@
3619b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3629b94acceSBarry Smith    with SNESSSetFunction().
3639b94acceSBarry Smith 
3649b94acceSBarry Smith    Input Parameter:
3659b94acceSBarry Smith .  snes - SNES context
3669b94acceSBarry Smith 
3679b94acceSBarry Smith    Output Parameter:
3689b94acceSBarry Smith .  fnorm - 2-norm of function
3699b94acceSBarry Smith 
3709b94acceSBarry Smith    Note:
3719b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
372a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
373a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3749b94acceSBarry Smith 
3759b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
376a86d99e1SLois Curfman McInnes 
377a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3789b94acceSBarry Smith @*/
3799b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3809b94acceSBarry Smith {
38177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
38274679c65SBarry Smith   PetscValidScalarPointer(fnorm);
38374679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
38474679c65SBarry Smith     SETERRQ(1,"SNESGetFunctionNorm:For SNES_NONLINEAR_EQUATIONS only");
38574679c65SBarry Smith   }
3869b94acceSBarry Smith   *fnorm = snes->norm;
3879b94acceSBarry Smith   return 0;
3889b94acceSBarry Smith }
38974679c65SBarry Smith 
3909b94acceSBarry Smith /*@
3919b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
3929b94acceSBarry Smith    with SNESSSetGradient().
3939b94acceSBarry Smith 
3949b94acceSBarry Smith    Input Parameter:
3959b94acceSBarry Smith .  snes - SNES context
3969b94acceSBarry Smith 
3979b94acceSBarry Smith    Output Parameter:
3989b94acceSBarry Smith .  fnorm - 2-norm of gradient
3999b94acceSBarry Smith 
4009b94acceSBarry Smith    Note:
4019b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
402a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
403a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4049b94acceSBarry Smith 
4059b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
406a86d99e1SLois Curfman McInnes 
407a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4089b94acceSBarry Smith @*/
4099b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
4109b94acceSBarry Smith {
41177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
41274679c65SBarry Smith   PetscValidScalarPointer(gnorm);
41374679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
41474679c65SBarry Smith     SETERRQ(1,"SNESGetGradientNorm:For SNES_UNCONSTRAINED_MINIMIZATION only");
41574679c65SBarry Smith   }
4169b94acceSBarry Smith   *gnorm = snes->norm;
4179b94acceSBarry Smith   return 0;
4189b94acceSBarry Smith }
41974679c65SBarry Smith 
4209b94acceSBarry Smith /*@
4219b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4229b94acceSBarry Smith    attempted by the nonlinear solver.
4239b94acceSBarry Smith 
4249b94acceSBarry Smith    Input Parameter:
4259b94acceSBarry Smith .  snes - SNES context
4269b94acceSBarry Smith 
4279b94acceSBarry Smith    Output Parameter:
4289b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4299b94acceSBarry Smith 
4309b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4319b94acceSBarry Smith @*/
4329b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4339b94acceSBarry Smith {
43477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43574679c65SBarry Smith   PetscValidIntPointer(nfails);
4369b94acceSBarry Smith   *nfails = snes->nfailures;
4379b94acceSBarry Smith   return 0;
4389b94acceSBarry Smith }
4399b94acceSBarry Smith /*@C
4409b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4419b94acceSBarry Smith 
4429b94acceSBarry Smith    Input Parameter:
4439b94acceSBarry Smith .  snes - the SNES context
4449b94acceSBarry Smith 
4459b94acceSBarry Smith    Output Parameter:
4469b94acceSBarry Smith .  sles - the SLES context
4479b94acceSBarry Smith 
4489b94acceSBarry Smith    Notes:
4499b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4509b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4519b94acceSBarry Smith    KSP and PC contexts as well.
4529b94acceSBarry Smith 
4539b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4549b94acceSBarry Smith 
4559b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4569b94acceSBarry Smith @*/
4579b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4589b94acceSBarry Smith {
45977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4609b94acceSBarry Smith   *sles = snes->sles;
4619b94acceSBarry Smith   return 0;
4629b94acceSBarry Smith }
4639b94acceSBarry Smith /* -----------------------------------------------------------*/
4649b94acceSBarry Smith /*@C
4659b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
4669b94acceSBarry Smith 
4679b94acceSBarry Smith    Input Parameter:
4689b94acceSBarry Smith .  comm - MPI communicator
4699b94acceSBarry Smith .  type - type of method, one of
4709b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
4719b94acceSBarry Smith $      (for systems of nonlinear equations)
4729b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
4739b94acceSBarry Smith $      (for unconstrained minimization)
4749b94acceSBarry Smith 
4759b94acceSBarry Smith    Output Parameter:
4769b94acceSBarry Smith .  outsnes - the new SNES context
4779b94acceSBarry Smith 
478c1f60f51SBarry Smith    Options Database Key:
47919bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
48019bd6747SLois Curfman McInnes $              and no preconditioning matrix
48119bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
48219bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
48319bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
48419bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
485c1f60f51SBarry Smith 
4869b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
4879b94acceSBarry Smith 
48863a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
4899b94acceSBarry Smith @*/
4904b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
4919b94acceSBarry Smith {
4929b94acceSBarry Smith   int                 ierr;
4939b94acceSBarry Smith   SNES                snes;
4949b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
49537fcc0dbSBarry Smith 
4969b94acceSBarry Smith   *outsnes = 0;
4972a0acf01SLois Curfman McInnes   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS)
4982a0acf01SLois Curfman McInnes     SETERRQ(1,"SNESCreate:incorrect method type");
4990452661fSBarry Smith   PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm);
5009b94acceSBarry Smith   PLogObjectCreate(snes);
5019b94acceSBarry Smith   snes->max_its           = 50;
5029b94acceSBarry Smith   snes->max_funcs	  = 1000;
5039b94acceSBarry Smith   snes->norm		  = 0.0;
5045a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
505b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
506b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5079b94acceSBarry Smith     snes->atol		  = 1.e-10;
5085a2d0531SBarry Smith   }
509b4874afaSBarry Smith   else {
510b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
511b4874afaSBarry Smith     snes->ttol            = 0.0;
512b4874afaSBarry Smith     snes->atol		  = 1.e-50;
513b4874afaSBarry Smith   }
5149b94acceSBarry Smith   snes->xtol		  = 1.e-8;
515b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5169b94acceSBarry Smith   snes->nfuncs            = 0;
5179b94acceSBarry Smith   snes->nfailures         = 0;
518639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5199b94acceSBarry Smith   snes->data              = 0;
5209b94acceSBarry Smith   snes->view              = 0;
5219b94acceSBarry Smith   snes->computeumfunction = 0;
5229b94acceSBarry Smith   snes->umfunP            = 0;
5239b94acceSBarry Smith   snes->fc                = 0;
5249b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5259b94acceSBarry Smith   snes->fmin              = -1.e30;
5269b94acceSBarry Smith   snes->method_class      = type;
5279b94acceSBarry Smith   snes->set_method_called = 0;
528a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
5299b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5307d1a2b2bSBarry Smith   snes->type              = -1;
5316f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5326f24a144SLois Curfman McInnes   snes->vwork             = 0;
5336f24a144SLois Curfman McInnes   snes->nwork             = 0;
5349b94acceSBarry Smith 
5359b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5360452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
5379b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5389b94acceSBarry Smith   kctx->version     = 2;
5399b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5409b94acceSBarry Smith                              this was too large for some test cases */
5419b94acceSBarry Smith   kctx->rtol_last   = 0;
5429b94acceSBarry Smith   kctx->rtol_max    = .9;
5439b94acceSBarry Smith   kctx->gamma       = 1.0;
5449b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5459b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5469b94acceSBarry Smith   kctx->threshold   = .1;
5479b94acceSBarry Smith   kctx->lresid_last = 0;
5489b94acceSBarry Smith   kctx->norm_last   = 0;
5499b94acceSBarry Smith 
5509b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5519b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5525334005bSBarry Smith 
5539b94acceSBarry Smith   *outsnes = snes;
5549b94acceSBarry Smith   return 0;
5559b94acceSBarry Smith }
5569b94acceSBarry Smith 
5579b94acceSBarry Smith /* --------------------------------------------------------------- */
5589b94acceSBarry Smith /*@C
5599b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
5609b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
5619b94acceSBarry Smith    equations.
5629b94acceSBarry Smith 
5639b94acceSBarry Smith    Input Parameters:
5649b94acceSBarry Smith .  snes - the SNES context
5659b94acceSBarry Smith .  func - function evaluation routine
5669b94acceSBarry Smith .  r - vector to store function value
5672cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
5682cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
5699b94acceSBarry Smith 
5709b94acceSBarry Smith    Calling sequence of func:
571313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
5729b94acceSBarry Smith 
5739b94acceSBarry Smith .  x - input vector
574313e4042SLois Curfman McInnes .  f - function vector
5752cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
5769b94acceSBarry Smith 
5779b94acceSBarry Smith    Notes:
5789b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
5799b94acceSBarry Smith $      f'(x) x = -f(x),
5809b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
5819b94acceSBarry Smith 
5829b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
5839b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
5849b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
5859b94acceSBarry Smith 
5869b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
5879b94acceSBarry Smith 
588a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
5899b94acceSBarry Smith @*/
5905334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
5919b94acceSBarry Smith {
59277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5939b94acceSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
59448d91487SBarry Smith     "SNESSetFunction:For SNES_NONLINEAR_EQUATIONS only");
5959b94acceSBarry Smith   snes->computefunction     = func;
5969b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
5979b94acceSBarry Smith   snes->funP                = ctx;
5989b94acceSBarry Smith   return 0;
5999b94acceSBarry Smith }
6009b94acceSBarry Smith 
6019b94acceSBarry Smith /*@
6029b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6039b94acceSBarry Smith    SNESSetFunction().
6049b94acceSBarry Smith 
6059b94acceSBarry Smith    Input Parameters:
6069b94acceSBarry Smith .  snes - the SNES context
6079b94acceSBarry Smith .  x - input vector
6089b94acceSBarry Smith 
6099b94acceSBarry Smith    Output Parameter:
6103638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6119b94acceSBarry Smith 
6121bffabb2SLois Curfman McInnes    Notes:
6139b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6149b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6159b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6169b94acceSBarry Smith 
6179b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6189b94acceSBarry Smith 
619a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6209b94acceSBarry Smith @*/
6219b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6229b94acceSBarry Smith {
6239b94acceSBarry Smith   int    ierr;
6249b94acceSBarry Smith 
62574679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
62674679c65SBarry Smith     SETERRQ(1,"SNESComputeFunction: For SNES_NONLINEAR_EQUATIONS only");
6279b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
6289b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
6299b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6309b94acceSBarry Smith   return 0;
6319b94acceSBarry Smith }
6329b94acceSBarry Smith 
6339b94acceSBarry Smith /*@C
6349b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6359b94acceSBarry Smith    unconstrained minimization.
6369b94acceSBarry Smith 
6379b94acceSBarry Smith    Input Parameters:
6389b94acceSBarry Smith .  snes - the SNES context
6399b94acceSBarry Smith .  func - function evaluation routine
640044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
641044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6429b94acceSBarry Smith 
6439b94acceSBarry Smith    Calling sequence of func:
6449b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
6459b94acceSBarry Smith 
6469b94acceSBarry Smith .  x - input vector
6479b94acceSBarry Smith .  f - function
648044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
6499b94acceSBarry Smith 
6509b94acceSBarry Smith    Notes:
6519b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
6529b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
6539b94acceSBarry Smith    SNESSetFunction().
6549b94acceSBarry Smith 
6559b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
6569b94acceSBarry Smith 
657a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
658a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
6599b94acceSBarry Smith @*/
6609b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
6619b94acceSBarry Smith                       void *ctx)
6629b94acceSBarry Smith {
66377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
6649b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
66548d91487SBarry Smith     "SNESSetMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION");
6669b94acceSBarry Smith   snes->computeumfunction   = func;
6679b94acceSBarry Smith   snes->umfunP              = ctx;
6689b94acceSBarry Smith   return 0;
6699b94acceSBarry Smith }
6709b94acceSBarry Smith 
6719b94acceSBarry Smith /*@
6729b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
6739b94acceSBarry Smith    set with SNESSetMinimizationFunction().
6749b94acceSBarry Smith 
6759b94acceSBarry Smith    Input Parameters:
6769b94acceSBarry Smith .  snes - the SNES context
6779b94acceSBarry Smith .  x - input vector
6789b94acceSBarry Smith 
6799b94acceSBarry Smith    Output Parameter:
6809b94acceSBarry Smith .  y - function value
6819b94acceSBarry Smith 
6829b94acceSBarry Smith    Notes:
6839b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
6849b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
6859b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
686a86d99e1SLois Curfman McInnes 
687a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
688a86d99e1SLois Curfman McInnes 
689a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
690a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
6919b94acceSBarry Smith @*/
6929b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
6939b94acceSBarry Smith {
6949b94acceSBarry Smith   int    ierr;
6959b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
69648d91487SBarry Smith     "SNESComputeMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION");
6979b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
6989b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
6999b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7009b94acceSBarry Smith   return 0;
7019b94acceSBarry Smith }
7029b94acceSBarry Smith 
7039b94acceSBarry Smith /*@C
7049b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7059b94acceSBarry Smith    vector for use by the SNES routines.
7069b94acceSBarry Smith 
7079b94acceSBarry Smith    Input Parameters:
7089b94acceSBarry Smith .  snes - the SNES context
7099b94acceSBarry Smith .  func - function evaluation routine
710044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
711044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7129b94acceSBarry Smith .  r - vector to store gradient value
7139b94acceSBarry Smith 
7149b94acceSBarry Smith    Calling sequence of func:
7159b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7169b94acceSBarry Smith 
7179b94acceSBarry Smith .  x - input vector
7189b94acceSBarry Smith .  g - gradient vector
719044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7209b94acceSBarry Smith 
7219b94acceSBarry Smith    Notes:
7229b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7239b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7249b94acceSBarry Smith    SNESSetFunction().
7259b94acceSBarry Smith 
7269b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7279b94acceSBarry Smith 
728a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
729a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
7309b94acceSBarry Smith @*/
73174679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
7329b94acceSBarry Smith {
73377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
7349b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
73548d91487SBarry Smith     "SNESSetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
7369b94acceSBarry Smith   snes->computefunction     = func;
7379b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7389b94acceSBarry Smith   snes->funP                = ctx;
7399b94acceSBarry Smith   return 0;
7409b94acceSBarry Smith }
7419b94acceSBarry Smith 
7429b94acceSBarry Smith /*@
743a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
744a86d99e1SLois Curfman McInnes    SNESSetGradient().
7459b94acceSBarry Smith 
7469b94acceSBarry Smith    Input Parameters:
7479b94acceSBarry Smith .  snes - the SNES context
7489b94acceSBarry Smith .  x - input vector
7499b94acceSBarry Smith 
7509b94acceSBarry Smith    Output Parameter:
7519b94acceSBarry Smith .  y - gradient vector
7529b94acceSBarry Smith 
7539b94acceSBarry Smith    Notes:
7549b94acceSBarry Smith    SNESComputeGradient() is valid only for
7559b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7569b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
757a86d99e1SLois Curfman McInnes 
758a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
759a86d99e1SLois Curfman McInnes 
760a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
761a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
7629b94acceSBarry Smith @*/
7639b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
7649b94acceSBarry Smith {
7659b94acceSBarry Smith   int    ierr;
76674679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
76774679c65SBarry Smith     SETERRQ(1,"SNESComputeGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
7689b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
7699b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
7709b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
7719b94acceSBarry Smith   return 0;
7729b94acceSBarry Smith }
7739b94acceSBarry Smith 
77462fef451SLois Curfman McInnes /*@
77562fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
77662fef451SLois Curfman McInnes    set with SNESSetJacobian().
77762fef451SLois Curfman McInnes 
77862fef451SLois Curfman McInnes    Input Parameters:
77962fef451SLois Curfman McInnes .  snes - the SNES context
78062fef451SLois Curfman McInnes .  x - input vector
78162fef451SLois Curfman McInnes 
78262fef451SLois Curfman McInnes    Output Parameters:
78362fef451SLois Curfman McInnes .  A - Jacobian matrix
78462fef451SLois Curfman McInnes .  B - optional preconditioning matrix
78562fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
78662fef451SLois Curfman McInnes 
78762fef451SLois Curfman McInnes    Notes:
78862fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
78962fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
79062fef451SLois Curfman McInnes 
791dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
792dc5a77f8SLois Curfman McInnes    flag parameter.
79362fef451SLois Curfman McInnes 
79462fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
79562fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
79662fef451SLois Curfman McInnes    methods is SNESComputeJacobian().
79762fef451SLois Curfman McInnes 
79862fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
79962fef451SLois Curfman McInnes 
80062fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
80162fef451SLois Curfman McInnes @*/
8029b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8039b94acceSBarry Smith {
8049b94acceSBarry Smith   int    ierr;
80574679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
80674679c65SBarry Smith     SETERRQ(1,"SNESComputeJacobian: For SNES_NONLINEAR_EQUATIONS only");
8079b94acceSBarry Smith   if (!snes->computejacobian) return 0;
8089b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
809c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
8109b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
8119b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8126d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
81377c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
81477c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8159b94acceSBarry Smith   return 0;
8169b94acceSBarry Smith }
8179b94acceSBarry Smith 
81862fef451SLois Curfman McInnes /*@
81962fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
82062fef451SLois Curfman McInnes    set with SNESSetHessian().
82162fef451SLois Curfman McInnes 
82262fef451SLois Curfman McInnes    Input Parameters:
82362fef451SLois Curfman McInnes .  snes - the SNES context
82462fef451SLois Curfman McInnes .  x - input vector
82562fef451SLois Curfman McInnes 
82662fef451SLois Curfman McInnes    Output Parameters:
82762fef451SLois Curfman McInnes .  A - Hessian matrix
82862fef451SLois Curfman McInnes .  B - optional preconditioning matrix
82962fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
83062fef451SLois Curfman McInnes 
83162fef451SLois Curfman McInnes    Notes:
83262fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
83362fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
83462fef451SLois Curfman McInnes 
835dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
836dc5a77f8SLois Curfman McInnes    flag parameter.
83762fef451SLois Curfman McInnes 
83862fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
83962fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
84062fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
84162fef451SLois Curfman McInnes 
84262fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
84362fef451SLois Curfman McInnes 
844a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
845a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
84662fef451SLois Curfman McInnes @*/
84762fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
8489b94acceSBarry Smith {
8499b94acceSBarry Smith   int    ierr;
85074679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
85174679c65SBarry Smith     SETERRQ(1,"SNESComputeHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
8529b94acceSBarry Smith   if (!snes->computejacobian) return 0;
85362fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
8540f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
85562fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
85662fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
8570f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
85877c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
85977c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8609b94acceSBarry Smith   return 0;
8619b94acceSBarry Smith }
8629b94acceSBarry Smith 
8639b94acceSBarry Smith /*@C
8649b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
865044dda88SLois Curfman McInnes    location to store the matrix.
8669b94acceSBarry Smith 
8679b94acceSBarry Smith    Input Parameters:
8689b94acceSBarry Smith .  snes - the SNES context
8699b94acceSBarry Smith .  A - Jacobian matrix
8709b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
8719b94acceSBarry Smith .  func - Jacobian evaluation routine
8722cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
8732cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
8749b94acceSBarry Smith 
8759b94acceSBarry Smith    Calling sequence of func:
876313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
8779b94acceSBarry Smith 
8789b94acceSBarry Smith .  x - input vector
8799b94acceSBarry Smith .  A - Jacobian matrix
8809b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
881ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
882ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
8832cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
8849b94acceSBarry Smith 
8859b94acceSBarry Smith    Notes:
886dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
8872cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
888ac21db08SLois Curfman McInnes 
889ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
8909b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
8919b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
8929b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
8939b94acceSBarry Smith    throughout the global iterations.
8949b94acceSBarry Smith 
8959b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
8969b94acceSBarry Smith 
897ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
8989b94acceSBarry Smith @*/
8999b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9009b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9019b94acceSBarry Smith {
90277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
90374679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
90474679c65SBarry Smith     SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only");
9059b94acceSBarry Smith   snes->computejacobian = func;
9069b94acceSBarry Smith   snes->jacP            = ctx;
9079b94acceSBarry Smith   snes->jacobian        = A;
9089b94acceSBarry Smith   snes->jacobian_pre    = B;
9099b94acceSBarry Smith   return 0;
9109b94acceSBarry Smith }
91162fef451SLois Curfman McInnes 
912b4fd4287SBarry Smith /*@
913b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
914b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
915b4fd4287SBarry Smith 
916b4fd4287SBarry Smith    Input Parameter:
917b4fd4287SBarry Smith .  snes - the nonlinear solver context
918b4fd4287SBarry Smith 
919b4fd4287SBarry Smith    Output Parameters:
920b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
921b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
922b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
923b4fd4287SBarry Smith 
924b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
925b4fd4287SBarry Smith @*/
926b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
927b4fd4287SBarry Smith {
92874679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
92974679c65SBarry Smith     SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only");
930b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
931b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
932b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
933b4fd4287SBarry Smith   return 0;
934b4fd4287SBarry Smith }
935b4fd4287SBarry Smith 
9369b94acceSBarry Smith /*@C
9379b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
938044dda88SLois Curfman McInnes    location to store the matrix.
9399b94acceSBarry Smith 
9409b94acceSBarry Smith    Input Parameters:
9419b94acceSBarry Smith .  snes - the SNES context
9429b94acceSBarry Smith .  A - Hessian matrix
9439b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
9449b94acceSBarry Smith .  func - Jacobian evaluation routine
945313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
946313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
9479b94acceSBarry Smith 
9489b94acceSBarry Smith    Calling sequence of func:
949313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9509b94acceSBarry Smith 
9519b94acceSBarry Smith .  x - input vector
9529b94acceSBarry Smith .  A - Hessian matrix
9539b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
954ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
955ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9562cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
9579b94acceSBarry Smith 
9589b94acceSBarry Smith    Notes:
959dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9602cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
961ac21db08SLois Curfman McInnes 
9629b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
9639b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
9649b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9659b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9669b94acceSBarry Smith    throughout the global iterations.
9679b94acceSBarry Smith 
9689b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
9699b94acceSBarry Smith 
970ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
9719b94acceSBarry Smith @*/
9729b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9739b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9749b94acceSBarry Smith {
97577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
97674679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
97774679c65SBarry Smith     SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
9789b94acceSBarry Smith   snes->computejacobian = func;
9799b94acceSBarry Smith   snes->jacP            = ctx;
9809b94acceSBarry Smith   snes->jacobian        = A;
9819b94acceSBarry Smith   snes->jacobian_pre    = B;
9829b94acceSBarry Smith   return 0;
9839b94acceSBarry Smith }
9849b94acceSBarry Smith 
98562fef451SLois Curfman McInnes /*@
98662fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
98762fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
98862fef451SLois Curfman McInnes 
98962fef451SLois Curfman McInnes    Input Parameter:
99062fef451SLois Curfman McInnes .  snes - the nonlinear solver context
99162fef451SLois Curfman McInnes 
99262fef451SLois Curfman McInnes    Output Parameters:
99362fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
99462fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
99562fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
99662fef451SLois Curfman McInnes 
99762fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
99862fef451SLois Curfman McInnes @*/
99962fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
100062fef451SLois Curfman McInnes {
100174679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
100274679c65SBarry Smith     SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only");
100362fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
100462fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
100562fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
100662fef451SLois Curfman McInnes   return 0;
100762fef451SLois Curfman McInnes }
100862fef451SLois Curfman McInnes 
10099b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10109b94acceSBarry Smith 
10119b94acceSBarry Smith /*@
10129b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1013272ac6f2SLois Curfman McInnes    of a nonlinear solver.
10149b94acceSBarry Smith 
10159b94acceSBarry Smith    Input Parameter:
10169b94acceSBarry Smith .  snes - the SNES context
10178ddd3da0SLois Curfman McInnes .  x - the solution vector
10189b94acceSBarry Smith 
1019272ac6f2SLois Curfman McInnes    Notes:
1020272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1021272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1022272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1023272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1024272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1025272ac6f2SLois Curfman McInnes 
10269b94acceSBarry Smith .keywords: SNES, nonlinear, setup
10279b94acceSBarry Smith 
10289b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
10299b94acceSBarry Smith @*/
10308ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
10319b94acceSBarry Smith {
1032272ac6f2SLois Curfman McInnes   int ierr, flg;
103377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
103477c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
10358ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
10369b94acceSBarry Smith 
1037c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1038c1f60f51SBarry Smith   /*
1039c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1040dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1041c1f60f51SBarry Smith   */
1042c1f60f51SBarry Smith   if (flg) {
1043c1f60f51SBarry Smith     Mat J;
1044c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1045c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1046c1f60f51SBarry Smith     snes->mfshell = J;
1047c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1048c1f60f51SBarry Smith       snes->jacobian = J;
1049c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1050c1f60f51SBarry Smith     }
1051c1f60f51SBarry Smith     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1052c1f60f51SBarry Smith       snes->jacobian = J;
1053c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1054c1f60f51SBarry Smith     } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free operator option");
1055c1f60f51SBarry Smith   }
1056272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1057c1f60f51SBarry Smith   /*
1058dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1059c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1060c1f60f51SBarry Smith    */
106131615d04SLois Curfman McInnes   if (flg) {
1062272ac6f2SLois Curfman McInnes     Mat J;
1063272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1064272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1065272ac6f2SLois Curfman McInnes     snes->mfshell = J;
106683e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
106783e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
106894a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
106983e56afdSLois Curfman McInnes     }
107083e56afdSLois Curfman McInnes     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
107183e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
107294a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
107383e56afdSLois Curfman McInnes     } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free option");
1074272ac6f2SLois Curfman McInnes   }
10759b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
107637fcc0dbSBarry Smith     if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first");
107737fcc0dbSBarry Smith     if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first");
107837fcc0dbSBarry Smith     if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetJacobian() first");
1079d804570eSBarry Smith     if (snes->vec_func == snes->vec_sol) SETERRQ(1,"SNESSetUp:Solution vector cannot be function vector");
1080a703fe33SLois Curfman McInnes 
1081a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
108240191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1083a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1084a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1085a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1086a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1087a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1088a703fe33SLois Curfman McInnes     }
10899b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
109037fcc0dbSBarry Smith     if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first");
109137fcc0dbSBarry Smith     if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first");
109237fcc0dbSBarry Smith     if (!snes->computeumfunction)
109337fcc0dbSBarry Smith       SETERRQ(1,"SNESSetUp:Must call SNESSetMinimizationFunction() first");
109437fcc0dbSBarry Smith     if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetHessian() first");
10959b94acceSBarry Smith   } else SETERRQ(1,"SNESSetUp:Unknown method class");
1096a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1097a703fe33SLois Curfman McInnes   snes->setup_called = 1;
1098a703fe33SLois Curfman McInnes   return 0;
10999b94acceSBarry Smith }
11009b94acceSBarry Smith 
11019b94acceSBarry Smith /*@C
11029b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11039b94acceSBarry Smith    with SNESCreate().
11049b94acceSBarry Smith 
11059b94acceSBarry Smith    Input Parameter:
11069b94acceSBarry Smith .  snes - the SNES context
11079b94acceSBarry Smith 
11089b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
11099b94acceSBarry Smith 
111063a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
11119b94acceSBarry Smith @*/
11129b94acceSBarry Smith int SNESDestroy(SNES snes)
11139b94acceSBarry Smith {
11149b94acceSBarry Smith   int ierr;
111577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
11169b94acceSBarry Smith   ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);
11170452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
11189b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
11199b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
11203e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
11216f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
11229b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
11230452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
11249b94acceSBarry Smith   return 0;
11259b94acceSBarry Smith }
11269b94acceSBarry Smith 
11279b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11289b94acceSBarry Smith 
11299b94acceSBarry Smith /*@
1130d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11319b94acceSBarry Smith 
11329b94acceSBarry Smith    Input Parameters:
11339b94acceSBarry Smith .  snes - the SNES context
113433174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
113533174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
113633174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
113733174efeSLois Curfman McInnes            of the change in the solution between steps
113833174efeSLois Curfman McInnes .  maxit - maximum number of iterations
113933174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
11409b94acceSBarry Smith 
114133174efeSLois Curfman McInnes    Options Database Keys:
114233174efeSLois Curfman McInnes $    -snes_atol <atol>
114333174efeSLois Curfman McInnes $    -snes_rtol <rtol>
114433174efeSLois Curfman McInnes $    -snes_stol <stol>
114533174efeSLois Curfman McInnes $    -snes_max_it <maxit>
114633174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
11479b94acceSBarry Smith 
1148d7a720efSLois Curfman McInnes    Notes:
11499b94acceSBarry Smith    The default maximum number of iterations is 50.
11509b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
11519b94acceSBarry Smith 
115233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
11539b94acceSBarry Smith 
1154d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
11559b94acceSBarry Smith @*/
1156d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
11579b94acceSBarry Smith {
115877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1159d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1160d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1161d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1162d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1163d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
11649b94acceSBarry Smith   return 0;
11659b94acceSBarry Smith }
11669b94acceSBarry Smith 
11679b94acceSBarry Smith /*@
116833174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
116933174efeSLois Curfman McInnes 
117033174efeSLois Curfman McInnes    Input Parameters:
117133174efeSLois Curfman McInnes .  snes - the SNES context
117233174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
117333174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
117433174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
117533174efeSLois Curfman McInnes            of the change in the solution between steps
117633174efeSLois Curfman McInnes .  maxit - maximum number of iterations
117733174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
117833174efeSLois Curfman McInnes 
117933174efeSLois Curfman McInnes    Notes:
118033174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
118133174efeSLois Curfman McInnes 
118233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
118333174efeSLois Curfman McInnes 
118433174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
118533174efeSLois Curfman McInnes @*/
118633174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
118733174efeSLois Curfman McInnes {
118833174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
118933174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
119033174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
119133174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
119233174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
119333174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
119433174efeSLois Curfman McInnes   return 0;
119533174efeSLois Curfman McInnes }
119633174efeSLois Curfman McInnes 
119733174efeSLois Curfman McInnes /*@
11989b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
11999b94acceSBarry Smith 
12009b94acceSBarry Smith    Input Parameters:
12019b94acceSBarry Smith .  snes - the SNES context
12029b94acceSBarry Smith .  tol - tolerance
12039b94acceSBarry Smith 
12049b94acceSBarry Smith    Options Database Key:
1205d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
12069b94acceSBarry Smith 
12079b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12089b94acceSBarry Smith 
1209d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
12109b94acceSBarry Smith @*/
12119b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
12129b94acceSBarry Smith {
121377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12149b94acceSBarry Smith   snes->deltatol = tol;
12159b94acceSBarry Smith   return 0;
12169b94acceSBarry Smith }
12179b94acceSBarry Smith 
12189b94acceSBarry Smith /*@
121977c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
12209b94acceSBarry Smith    for unconstrained minimization solvers.
12219b94acceSBarry Smith 
12229b94acceSBarry Smith    Input Parameters:
12239b94acceSBarry Smith .  snes - the SNES context
12249b94acceSBarry Smith .  ftol - minimum function tolerance
12259b94acceSBarry Smith 
12269b94acceSBarry Smith    Options Database Key:
1227d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
12289b94acceSBarry Smith 
12299b94acceSBarry Smith    Note:
123077c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
12319b94acceSBarry Smith    methods only.
12329b94acceSBarry Smith 
12339b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
12349b94acceSBarry Smith 
1235d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
12369b94acceSBarry Smith @*/
123777c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
12389b94acceSBarry Smith {
123977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12409b94acceSBarry Smith   snes->fmin = ftol;
12419b94acceSBarry Smith   return 0;
12429b94acceSBarry Smith }
12439b94acceSBarry Smith 
12449b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
12459b94acceSBarry Smith 
12469b94acceSBarry Smith /*@C
12479b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
12489b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
12499b94acceSBarry Smith    progress.
12509b94acceSBarry Smith 
12519b94acceSBarry Smith    Input Parameters:
12529b94acceSBarry Smith .  snes - the SNES context
12539b94acceSBarry Smith .  func - monitoring routine
1254044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1255044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
12569b94acceSBarry Smith 
12579b94acceSBarry Smith    Calling sequence of func:
1258682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
12599b94acceSBarry Smith 
12609b94acceSBarry Smith $    snes - the SNES context
12619b94acceSBarry Smith $    its - iteration number
1262044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
12639b94acceSBarry Smith $
12649b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
12659b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
12669b94acceSBarry Smith $
12679b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
12689b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
12699b94acceSBarry Smith 
1270639f9d9dSBarry Smith    Notes:
12716bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
12726bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
12736bc08f3fSLois Curfman McInnes    order in which they were set.
1274639f9d9dSBarry Smith 
12759b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
12769b94acceSBarry Smith 
12779b94acceSBarry Smith .seealso: SNESDefaultMonitor()
12789b94acceSBarry Smith @*/
127974679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
12809b94acceSBarry Smith {
1281639f9d9dSBarry Smith   if (!func) {
1282639f9d9dSBarry Smith     snes->numbermonitors = 0;
1283639f9d9dSBarry Smith     return 0;
1284639f9d9dSBarry Smith   }
1285639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1286639f9d9dSBarry Smith     SETERRQ(1,"SNESSetMonitor:Too many monitors set");
1287639f9d9dSBarry Smith   }
1288639f9d9dSBarry Smith 
1289639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1290639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
12919b94acceSBarry Smith   return 0;
12929b94acceSBarry Smith }
12939b94acceSBarry Smith 
12949b94acceSBarry Smith /*@C
12959b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
12969b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
12979b94acceSBarry Smith 
12989b94acceSBarry Smith    Input Parameters:
12999b94acceSBarry Smith .  snes - the SNES context
13009b94acceSBarry Smith .  func - routine to test for convergence
1301044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1302044dda88SLois Curfman McInnes           (may be PETSC_NULL)
13039b94acceSBarry Smith 
13049b94acceSBarry Smith    Calling sequence of func:
13059b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
13069b94acceSBarry Smith              double f,void *cctx)
13079b94acceSBarry Smith 
13089b94acceSBarry Smith $    snes - the SNES context
1309044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
13109b94acceSBarry Smith $    xnorm - 2-norm of current iterate
13119b94acceSBarry Smith $
13129b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13139b94acceSBarry Smith $    gnorm - 2-norm of current step
13149b94acceSBarry Smith $    f - 2-norm of function
13159b94acceSBarry Smith $
13169b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13179b94acceSBarry Smith $    gnorm - 2-norm of current gradient
13189b94acceSBarry Smith $    f - function value
13199b94acceSBarry Smith 
13209b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
13219b94acceSBarry Smith 
132240191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
132340191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
13249b94acceSBarry Smith @*/
132574679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
13269b94acceSBarry Smith {
13279b94acceSBarry Smith   (snes)->converged = func;
13289b94acceSBarry Smith   (snes)->cnvP      = cctx;
13299b94acceSBarry Smith   return 0;
13309b94acceSBarry Smith }
13319b94acceSBarry Smith 
13329b94acceSBarry Smith /*
13339b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
13349b94acceSBarry Smith    positive parameter delta.
13359b94acceSBarry Smith 
13369b94acceSBarry Smith     Input Parameters:
13379b94acceSBarry Smith .   snes - the SNES context
13389b94acceSBarry Smith .   y - approximate solution of linear system
13399b94acceSBarry Smith .   fnorm - 2-norm of current function
13409b94acceSBarry Smith .   delta - trust region size
13419b94acceSBarry Smith 
13429b94acceSBarry Smith     Output Parameters:
13439b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
13449b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
13459b94acceSBarry Smith     region, and exceeds zero otherwise.
13469b94acceSBarry Smith .   ynorm - 2-norm of the step
13479b94acceSBarry Smith 
13489b94acceSBarry Smith     Note:
134940191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
13509b94acceSBarry Smith     is set to be the maximum allowable step size.
13519b94acceSBarry Smith 
13529b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
13539b94acceSBarry Smith */
13549b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
13559b94acceSBarry Smith                           double *gpnorm,double *ynorm)
13569b94acceSBarry Smith {
13579b94acceSBarry Smith   double norm;
13589b94acceSBarry Smith   Scalar cnorm;
1359cddf8d76SBarry Smith   VecNorm(y,NORM_2, &norm );
13609b94acceSBarry Smith   if (norm > *delta) {
13619b94acceSBarry Smith      norm = *delta/norm;
13629b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
13639b94acceSBarry Smith      cnorm = norm;
13649b94acceSBarry Smith      VecScale( &cnorm, y );
13659b94acceSBarry Smith      *ynorm = *delta;
13669b94acceSBarry Smith   } else {
13679b94acceSBarry Smith      *gpnorm = 0.0;
13689b94acceSBarry Smith      *ynorm = norm;
13699b94acceSBarry Smith   }
13709b94acceSBarry Smith   return 0;
13719b94acceSBarry Smith }
13729b94acceSBarry Smith 
13739b94acceSBarry Smith /*@
13749b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
137563a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
13769b94acceSBarry Smith 
13779b94acceSBarry Smith    Input Parameter:
13789b94acceSBarry Smith .  snes - the SNES context
13798ddd3da0SLois Curfman McInnes .  x - the solution vector
13809b94acceSBarry Smith 
13819b94acceSBarry Smith    Output Parameter:
13829b94acceSBarry Smith    its - number of iterations until termination
13839b94acceSBarry Smith 
13848ddd3da0SLois Curfman McInnes    Note:
13858ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
13868ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
13878ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
13888ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
13898ddd3da0SLois Curfman McInnes 
13909b94acceSBarry Smith .keywords: SNES, nonlinear, solve
13919b94acceSBarry Smith 
139263a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
13939b94acceSBarry Smith @*/
13948ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
13959b94acceSBarry Smith {
13963c7409f5SSatish Balay   int ierr, flg;
1397052efed2SBarry Smith 
139877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
139974679c65SBarry Smith   PetscValidIntPointer(its);
1400c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1401c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
14029b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
14039b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
14049b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
14053c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
14066d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
14079b94acceSBarry Smith   return 0;
14089b94acceSBarry Smith }
14099b94acceSBarry Smith 
14109b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
141137fcc0dbSBarry Smith static NRList *__SNESList = 0;
14129b94acceSBarry Smith 
14139b94acceSBarry Smith /*@
14144b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
14159b94acceSBarry Smith 
14169b94acceSBarry Smith    Input Parameters:
14179b94acceSBarry Smith .  snes - the SNES context
14189b94acceSBarry Smith .  method - a known method
14199b94acceSBarry Smith 
1420ae12b187SLois Curfman McInnes   Options Database Command:
1421ae12b187SLois Curfman McInnes $ -snes_type  <method>
1422ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1423ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1424ae12b187SLois Curfman McInnes 
14259b94acceSBarry Smith    Notes:
14269b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
14279b94acceSBarry Smith $  Systems of nonlinear equations:
142840191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
142940191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
14309b94acceSBarry Smith $  Unconstrained minimization:
143140191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
143240191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
14339b94acceSBarry Smith 
1434ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1435ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1436ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1437ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1438ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1439ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1440ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1441ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1442ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1443ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1444a703fe33SLois Curfman McInnes 
1445f536c699SSatish Balay .keywords: SNES, set, method
14469b94acceSBarry Smith @*/
14474b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
14489b94acceSBarry Smith {
14499b94acceSBarry Smith   int (*r)(SNES);
145077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
14517d1a2b2bSBarry Smith   if (snes->type == (int) method) return 0;
14527d1a2b2bSBarry Smith 
14539b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
145437fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
145537fcc0dbSBarry Smith   if (!__SNESList) {SETERRQ(1,"SNESSetType:Could not get methods");}
145637fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
14574b0e389bSBarry Smith   if (!r) {SETERRQ(1,"SNESSetType:Unknown method");}
14580452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
14599b94acceSBarry Smith   snes->set_method_called = 1;
14609b94acceSBarry Smith   return (*r)(snes);
14619b94acceSBarry Smith }
14629b94acceSBarry Smith 
14639b94acceSBarry Smith /* --------------------------------------------------------------------- */
14649b94acceSBarry Smith /*@C
14659b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
14664b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
14679b94acceSBarry Smith 
14689b94acceSBarry Smith    Input Parameters:
146940191667SLois Curfman McInnes .  name - for instance SNES_EQ_LS, SNES_EQ_TR, ...
147040191667SLois Curfman McInnes .  sname - corresponding string for name
14719b94acceSBarry Smith .  create - routine to create method context
14729b94acceSBarry Smith 
14739b94acceSBarry Smith .keywords: SNES, nonlinear, register
14749b94acceSBarry Smith 
14759b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
14769b94acceSBarry Smith @*/
14779b94acceSBarry Smith int SNESRegister(int name, char *sname, int (*create)(SNES))
14789b94acceSBarry Smith {
14799b94acceSBarry Smith   int ierr;
148037fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
148137fcc0dbSBarry Smith   NRRegister( __SNESList, name, sname, (int (*)(void*))create );
14829b94acceSBarry Smith   return 0;
14839b94acceSBarry Smith }
14849b94acceSBarry Smith /* --------------------------------------------------------------------- */
14859b94acceSBarry Smith /*@C
14869b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
14879b94acceSBarry Smith    registered by SNESRegister().
14889b94acceSBarry Smith 
14899b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
14909b94acceSBarry Smith 
14919b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
14929b94acceSBarry Smith @*/
14939b94acceSBarry Smith int SNESRegisterDestroy()
14949b94acceSBarry Smith {
149537fcc0dbSBarry Smith   if (__SNESList) {
149637fcc0dbSBarry Smith     NRDestroy( __SNESList );
149737fcc0dbSBarry Smith     __SNESList = 0;
14989b94acceSBarry Smith   }
14999b94acceSBarry Smith   return 0;
15009b94acceSBarry Smith }
15019b94acceSBarry Smith 
15029b94acceSBarry Smith /*
15034b0e389bSBarry Smith    SNESGetTypeFromOptions_Private - Sets the selected method from the
15049b94acceSBarry Smith    options database.
15059b94acceSBarry Smith 
15069b94acceSBarry Smith    Input Parameter:
15079b94acceSBarry Smith .  ctx - the SNES context
15089b94acceSBarry Smith 
15099b94acceSBarry Smith    Output Parameter:
15109b94acceSBarry Smith .  method -  solver method
15119b94acceSBarry Smith 
15129b94acceSBarry Smith    Returns:
15139b94acceSBarry Smith    Returns 1 if the method is found; 0 otherwise.
15149b94acceSBarry Smith 
15159b94acceSBarry Smith    Options Database Key:
15164b0e389bSBarry Smith $  -snes_type  method
15179b94acceSBarry Smith */
1518052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
15199b94acceSBarry Smith {
1520052efed2SBarry Smith   int  ierr;
15219b94acceSBarry Smith   char sbuf[50];
15225baf8537SBarry Smith 
1523052efed2SBarry Smith   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1524052efed2SBarry Smith   if (*flg) {
1525052efed2SBarry Smith     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
152637fcc0dbSBarry Smith     *method = (SNESType)NRFindID( __SNESList, sbuf );
15279b94acceSBarry Smith   }
15289b94acceSBarry Smith   return 0;
15299b94acceSBarry Smith }
15309b94acceSBarry Smith 
15319b94acceSBarry Smith /*@C
15329a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
15339b94acceSBarry Smith 
15349b94acceSBarry Smith    Input Parameter:
15354b0e389bSBarry Smith .  snes - nonlinear solver context
15369b94acceSBarry Smith 
15379b94acceSBarry Smith    Output Parameter:
15389a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
15399a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
15409b94acceSBarry Smith 
15419b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
15429b94acceSBarry Smith @*/
15434b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
15449b94acceSBarry Smith {
154506a9b0adSLois Curfman McInnes   int ierr;
154637fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
15474b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
154837fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
15499b94acceSBarry Smith   return 0;
15509b94acceSBarry Smith }
15519b94acceSBarry Smith 
15529b94acceSBarry Smith #include <stdio.h>
15539b94acceSBarry Smith /*
15544b0e389bSBarry Smith    SNESPrintTypes_Private - Prints the SNES methods available from the
15559b94acceSBarry Smith    options database.
15569b94acceSBarry Smith 
15579b94acceSBarry Smith    Input Parameters:
155833455573SSatish Balay .  comm   - communicator (usually MPI_COMM_WORLD)
15599b94acceSBarry Smith .  prefix - prefix (usually "-")
15604b0e389bSBarry Smith .  name   - the options database name (by default "snes_type")
15619b94acceSBarry Smith */
156233455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
15639b94acceSBarry Smith {
15649b94acceSBarry Smith   FuncList *entry;
156537fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
156637fcc0dbSBarry Smith   entry = __SNESList->head;
156777c4ece6SBarry Smith   PetscPrintf(comm," %s%s (one of)",prefix,name);
15689b94acceSBarry Smith   while (entry) {
156977c4ece6SBarry Smith     PetscPrintf(comm," %s",entry->name);
15709b94acceSBarry Smith     entry = entry->next;
15719b94acceSBarry Smith   }
157277c4ece6SBarry Smith   PetscPrintf(comm,"\n");
15739b94acceSBarry Smith   return 0;
15749b94acceSBarry Smith }
15759b94acceSBarry Smith 
15769b94acceSBarry Smith /*@C
15779b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
15789b94acceSBarry Smith    stored.
15799b94acceSBarry Smith 
15809b94acceSBarry Smith    Input Parameter:
15819b94acceSBarry Smith .  snes - the SNES context
15829b94acceSBarry Smith 
15839b94acceSBarry Smith    Output Parameter:
15849b94acceSBarry Smith .  x - the solution
15859b94acceSBarry Smith 
15869b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
15879b94acceSBarry Smith 
1588a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
15899b94acceSBarry Smith @*/
15909b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
15919b94acceSBarry Smith {
159277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15939b94acceSBarry Smith   *x = snes->vec_sol_always;
15949b94acceSBarry Smith   return 0;
15959b94acceSBarry Smith }
15969b94acceSBarry Smith 
15979b94acceSBarry Smith /*@C
15989b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
15999b94acceSBarry Smith    stored.
16009b94acceSBarry Smith 
16019b94acceSBarry Smith    Input Parameter:
16029b94acceSBarry Smith .  snes - the SNES context
16039b94acceSBarry Smith 
16049b94acceSBarry Smith    Output Parameter:
16059b94acceSBarry Smith .  x - the solution update
16069b94acceSBarry Smith 
16079b94acceSBarry Smith    Notes:
16089b94acceSBarry Smith    This vector is implementation dependent.
16099b94acceSBarry Smith 
16109b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
16119b94acceSBarry Smith 
16129b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
16139b94acceSBarry Smith @*/
16149b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
16159b94acceSBarry Smith {
161677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16179b94acceSBarry Smith   *x = snes->vec_sol_update_always;
16189b94acceSBarry Smith   return 0;
16199b94acceSBarry Smith }
16209b94acceSBarry Smith 
16219b94acceSBarry Smith /*@C
16223638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
16239b94acceSBarry Smith 
16249b94acceSBarry Smith    Input Parameter:
16259b94acceSBarry Smith .  snes - the SNES context
16269b94acceSBarry Smith 
16279b94acceSBarry Smith    Output Parameter:
16283638b69dSLois Curfman McInnes .  r - the function
16299b94acceSBarry Smith 
16309b94acceSBarry Smith    Notes:
16319b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
16329b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
16339b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
16349b94acceSBarry Smith 
1635a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
16369b94acceSBarry Smith 
163731615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
163831615d04SLois Curfman McInnes           SNESGetGradient()
16399b94acceSBarry Smith @*/
16409b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
16419b94acceSBarry Smith {
164277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16439b94acceSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
164448d91487SBarry Smith     "SNESGetFunction:For SNES_NONLINEAR_EQUATIONS only");
16459b94acceSBarry Smith   *r = snes->vec_func_always;
16469b94acceSBarry Smith   return 0;
16479b94acceSBarry Smith }
16489b94acceSBarry Smith 
16499b94acceSBarry Smith /*@C
16503638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
16519b94acceSBarry Smith 
16529b94acceSBarry Smith    Input Parameter:
16539b94acceSBarry Smith .  snes - the SNES context
16549b94acceSBarry Smith 
16559b94acceSBarry Smith    Output Parameter:
16569b94acceSBarry Smith .  r - the gradient
16579b94acceSBarry Smith 
16589b94acceSBarry Smith    Notes:
16599b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
16609b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
16619b94acceSBarry Smith    SNESGetFunction().
16629b94acceSBarry Smith 
16639b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
16649b94acceSBarry Smith 
166531615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
16669b94acceSBarry Smith @*/
16679b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
16689b94acceSBarry Smith {
166977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16709b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
167148d91487SBarry Smith     "SNESGetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only");
16729b94acceSBarry Smith   *r = snes->vec_func_always;
16739b94acceSBarry Smith   return 0;
16749b94acceSBarry Smith }
16759b94acceSBarry Smith 
16769b94acceSBarry Smith /*@
16779b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
16789b94acceSBarry Smith    unconstrained minimization problems.
16799b94acceSBarry Smith 
16809b94acceSBarry Smith    Input Parameter:
16819b94acceSBarry Smith .  snes - the SNES context
16829b94acceSBarry Smith 
16839b94acceSBarry Smith    Output Parameter:
16849b94acceSBarry Smith .  r - the function
16859b94acceSBarry Smith 
16869b94acceSBarry Smith    Notes:
16879b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
16889b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
16899b94acceSBarry Smith    SNESGetFunction().
16909b94acceSBarry Smith 
16919b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
16929b94acceSBarry Smith 
169331615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
16949b94acceSBarry Smith @*/
16959b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
16969b94acceSBarry Smith {
169777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
169874679c65SBarry Smith   PetscValidScalarPointer(r);
16999b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
170048d91487SBarry Smith     "SNESGetMinimizationFunction:For SNES_UNCONSTRAINED_MINIMIZATION only");
17019b94acceSBarry Smith   *r = snes->fc;
17029b94acceSBarry Smith   return 0;
17039b94acceSBarry Smith }
17049b94acceSBarry Smith 
17053c7409f5SSatish Balay /*@C
17063c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1707639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1708639f9d9dSBarry Smith    the prefix name.
17093c7409f5SSatish Balay 
17103c7409f5SSatish Balay    Input Parameter:
17113c7409f5SSatish Balay .  snes - the SNES context
17123c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
17133c7409f5SSatish Balay 
17143c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1715a86d99e1SLois Curfman McInnes 
1716a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
17173c7409f5SSatish Balay @*/
17183c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
17193c7409f5SSatish Balay {
17203c7409f5SSatish Balay   int ierr;
17213c7409f5SSatish Balay 
172277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1723639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
17243c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
17253c7409f5SSatish Balay   return 0;
17263c7409f5SSatish Balay }
17273c7409f5SSatish Balay 
17283c7409f5SSatish Balay /*@C
1729f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1730639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1731639f9d9dSBarry Smith    the prefix name.
17323c7409f5SSatish Balay 
17333c7409f5SSatish Balay    Input Parameter:
17343c7409f5SSatish Balay .  snes - the SNES context
17353c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
17363c7409f5SSatish Balay 
17373c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1738a86d99e1SLois Curfman McInnes 
1739a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
17403c7409f5SSatish Balay @*/
17413c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
17423c7409f5SSatish Balay {
17433c7409f5SSatish Balay   int ierr;
17443c7409f5SSatish Balay 
174577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1746639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
17473c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
17483c7409f5SSatish Balay   return 0;
17493c7409f5SSatish Balay }
17503c7409f5SSatish Balay 
1751c83590e2SSatish Balay /*@
17523c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
17533c7409f5SSatish Balay    SNES options in the database.
17543c7409f5SSatish Balay 
17553c7409f5SSatish Balay    Input Parameter:
17563c7409f5SSatish Balay .  snes - the SNES context
17573c7409f5SSatish Balay 
17583c7409f5SSatish Balay    Output Parameter:
17593c7409f5SSatish Balay .  prefix - pointer to the prefix string used
17603c7409f5SSatish Balay 
17613c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1762a86d99e1SLois Curfman McInnes 
1763a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
17643c7409f5SSatish Balay @*/
17653c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
17663c7409f5SSatish Balay {
17673c7409f5SSatish Balay   int ierr;
17683c7409f5SSatish Balay 
176977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1770639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
17713c7409f5SSatish Balay   return 0;
17723c7409f5SSatish Balay }
17733c7409f5SSatish Balay 
17743c7409f5SSatish Balay 
17759b94acceSBarry Smith 
17769b94acceSBarry Smith 
17779b94acceSBarry Smith 
1778