xref: /petsc/src/snes/interface/snes.c (revision 63c41f6a1560bbb6cf7ee09697a660f5641fb9ab)
19b94acceSBarry Smith #ifndef lint
2*63c41f6aSSatish Balay static char vcid[] = "$Id: snes.c,v 1.101 1996/12/16 20:39:19 balay Exp balay $";
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 
14a847f771SSatish Balay #undef __FUNCTION__
15a847f771SSatish Balay #define __FUNCTION__ "SNESView"
169b94acceSBarry Smith /*@
179b94acceSBarry Smith    SNESView - Prints the SNES data structure.
189b94acceSBarry Smith 
199b94acceSBarry Smith    Input Parameters:
209b94acceSBarry Smith .  SNES - the SNES context
219b94acceSBarry Smith .  viewer - visualization context
229b94acceSBarry Smith 
239b94acceSBarry Smith    Options Database Key:
249b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
259b94acceSBarry Smith 
269b94acceSBarry Smith    Notes:
279b94acceSBarry Smith    The available visualization contexts include
286d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
296d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
309b94acceSBarry Smith $       output where only the first processor opens
319b94acceSBarry Smith $       the file.  All other processors send their
329b94acceSBarry Smith $       data to the first processor to print.
339b94acceSBarry Smith 
349b94acceSBarry Smith    The user can open alternative vistualization contexts with
35dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
369b94acceSBarry Smith 
379b94acceSBarry Smith .keywords: SNES, view
389b94acceSBarry Smith 
39dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
409b94acceSBarry Smith @*/
419b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
429b94acceSBarry Smith {
439b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
449b94acceSBarry Smith   FILE                *fd;
459b94acceSBarry Smith   int                 ierr;
469b94acceSBarry Smith   SLES                sles;
479b94acceSBarry Smith   char                *method;
4819bcc07fSBarry Smith   ViewerType          vtype;
499b94acceSBarry Smith 
5074679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5174679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5274679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5374679c65SBarry Smith 
5419bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5519bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5690ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5777c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
584b0e389bSBarry Smith     SNESGetType(snes,PETSC_NULL,&method);
5977c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
609b94acceSBarry Smith     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
6177c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
629b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
639b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6477c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
659b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
669b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
679b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
6877c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
699b94acceSBarry Smith     if (snes->ksp_ewconv) {
709b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
719b94acceSBarry Smith       if (kctx) {
7277c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
739b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
749b94acceSBarry Smith         kctx->version);
7577c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
769b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
779b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
7877c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
799b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
809b94acceSBarry Smith       }
819b94acceSBarry Smith     }
82c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
83c30f285eSLois Curfman McInnes     SNESGetType(snes,PETSC_NULL,&method);
84c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
8519bcc07fSBarry Smith   }
869b94acceSBarry Smith   SNESGetSLES(snes,&sles);
879b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
889b94acceSBarry Smith   return 0;
899b94acceSBarry Smith }
909b94acceSBarry Smith 
91639f9d9dSBarry Smith /*
92639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
93639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
94639f9d9dSBarry Smith */
95639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
96639f9d9dSBarry Smith static int numberofsetfromoptions;
97639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
98639f9d9dSBarry Smith 
99a847f771SSatish Balay #undef __FUNCTION__
100a847f771SSatish Balay #define __FUNCTION__ "SNESAddOptionsChecker"
101639f9d9dSBarry Smith /*@
102639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
103639f9d9dSBarry Smith 
104639f9d9dSBarry Smith   Input Parameter:
105639f9d9dSBarry Smith .   snescheck - function that checks for options
106639f9d9dSBarry Smith 
107639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
108639f9d9dSBarry Smith @*/
109639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
110639f9d9dSBarry Smith {
111639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
112*63c41f6aSSatish Balay     SETERRQ(1,"Too many options checkers, only 5 allowed");
113639f9d9dSBarry Smith   }
114639f9d9dSBarry Smith 
115639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
116639f9d9dSBarry Smith   return 0;
117639f9d9dSBarry Smith }
118639f9d9dSBarry Smith 
119a847f771SSatish Balay #undef __FUNCTION__
120a847f771SSatish Balay #define __FUNCTION__ "SNESSetFromOptions"
1219b94acceSBarry Smith /*@
122682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1239b94acceSBarry Smith 
1249b94acceSBarry Smith    Input Parameter:
1259b94acceSBarry Smith .  snes - the SNES context
1269b94acceSBarry Smith 
1279b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1289b94acceSBarry Smith 
129a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1309b94acceSBarry Smith @*/
1319b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1329b94acceSBarry Smith {
1334b0e389bSBarry Smith   SNESType method;
1349b94acceSBarry Smith   double   tmp;
1359b94acceSBarry Smith   SLES     sles;
13681f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1373c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1389b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1399b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1409b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1419b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1429b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1439b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1449b94acceSBarry Smith 
14581f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
14681f57ec7SBarry Smith 
14781f57ec7SBarry Smith 
14877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
149*63c41f6aSSatish Balay   if (snes->setup_called) SETERRQ(1,"Must call prior to SNESSetUp");
150052efed2SBarry Smith   ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr);
151052efed2SBarry Smith   if (flg) {
152052efed2SBarry Smith     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
1539b94acceSBarry Smith   }
1545334005bSBarry Smith   else if (!snes->set_method_called) {
1555334005bSBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
15640191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
1575334005bSBarry Smith     }
1585334005bSBarry Smith     else {
15940191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
1605334005bSBarry Smith     }
1615334005bSBarry Smith   }
1625334005bSBarry Smith 
1633c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1643c7409f5SSatish Balay   if (flg) { SNESPrintHelp(snes); }
1653c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
166d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);}
1673c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
168d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);}
1693c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
170d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); }
1713c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1723c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
173d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
174d7a720efSLois Curfman McInnes   if (flg) { SNESSetTrustRegionTolerance(snes,tmp); }
175d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
176d7a720efSLois Curfman McInnes   if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); }
1773c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1783c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
179b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
180b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
181b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
182b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
183b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
184b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
185b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
1869b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
1879b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
1885bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
1895bbad29bSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);}
1903c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
191639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
1923c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
193d31a3109SLois Curfman McInnes   if (flg) {
194639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
195d31a3109SLois Curfman McInnes   }
19681f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
1973c7409f5SSatish Balay   if (flg) {
19817699dbbSLois Curfman McInnes     int    rank = 0;
199d7e8b826SBarry Smith     DrawLG lg;
20017699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
20117699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
20217699dbbSLois Curfman McInnes     if (!rank) {
20381f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2049b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
205f8c826e1SBarry Smith       snes->xmonitor = lg;
2069b94acceSBarry Smith       PLogObjectParent(snes,lg);
2079b94acceSBarry Smith     }
2089b94acceSBarry Smith   }
2093c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2103c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2119b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2129b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
21394a424c1SBarry Smith     PLogInfo(snes,
21431615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
21531615d04SLois Curfman McInnes   }
21631615d04SLois Curfman McInnes   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
21731615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
21831615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
21994a424c1SBarry Smith     PLogInfo(snes,
22031615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2219b94acceSBarry Smith   }
222639f9d9dSBarry Smith 
223639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
224639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
225639f9d9dSBarry Smith   }
226639f9d9dSBarry Smith 
2279b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2289b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
2299b94acceSBarry Smith   if (!snes->setfromoptions) return 0;
2309b94acceSBarry Smith   return (*snes->setfromoptions)(snes);
2319b94acceSBarry Smith }
2329b94acceSBarry Smith 
233a847f771SSatish Balay #undef __FUNCTION__
234a847f771SSatish Balay #define __FUNCTION__ "SNESPrintHelp"
2359b94acceSBarry Smith /*@
2369b94acceSBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2379b94acceSBarry Smith 
2389b94acceSBarry Smith    Input Parameter:
2399b94acceSBarry Smith .  snes - the SNES context
2409b94acceSBarry Smith 
241a703fe33SLois Curfman McInnes    Options Database Keys:
242a703fe33SLois Curfman McInnes $  -help, -h
243a703fe33SLois Curfman McInnes 
2449b94acceSBarry Smith .keywords: SNES, nonlinear, help
2459b94acceSBarry Smith 
246682d7d0cSBarry Smith .seealso: SNESSetFromOptions()
2479b94acceSBarry Smith @*/
2489b94acceSBarry Smith int SNESPrintHelp(SNES snes)
2499b94acceSBarry Smith {
2506daaf66cSBarry Smith   char                p[64];
2516daaf66cSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2526daaf66cSBarry Smith 
25377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2546daaf66cSBarry Smith 
2556daaf66cSBarry Smith   PetscStrcpy(p,"-");
2566daaf66cSBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2576daaf66cSBarry Smith 
2586daaf66cSBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2596daaf66cSBarry Smith 
260d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
26133455573SSatish Balay   SNESPrintTypes_Private(snes->comm,p,"snes_type");
26277c4ece6SBarry Smith   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
263b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
264b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
265b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
266b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
267b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
268b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
2695bbad29bSBarry Smith   PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n");
2705bbad29bSBarry Smith   PetscPrintf(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
271d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor\n",p);
272d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
273b18e04deSLois Curfman McInnes   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
27477c4ece6SBarry Smith     PetscPrintf(snes->comm,
275d31a3109SLois Curfman McInnes      " Options for solving systems of nonlinear equations only:\n");
27677c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
27777c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
278ef1dfb11SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
2791650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p);
2801650c7b0SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p);
28177c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
28277c4ece6SBarry Smith     PetscPrintf(snes->comm,
283b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
28477c4ece6SBarry Smith     PetscPrintf(snes->comm,
285b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
28677c4ece6SBarry Smith     PetscPrintf(snes->comm,
287b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
28877c4ece6SBarry Smith     PetscPrintf(snes->comm,
289b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
29077c4ece6SBarry Smith     PetscPrintf(snes->comm,
291b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
29277c4ece6SBarry Smith     PetscPrintf(snes->comm,
293b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
29477c4ece6SBarry Smith     PetscPrintf(snes->comm,
295b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
296b18e04deSLois Curfman McInnes   }
297b18e04deSLois Curfman McInnes   else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
29877c4ece6SBarry Smith     PetscPrintf(snes->comm,
299d31a3109SLois Curfman McInnes      " Options for solving unconstrained minimization problems only:\n");
300b18e04deSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
30177c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
30277c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
303d31a3109SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p);
304d31a3109SLois Curfman McInnes      PetscPrintf(snes->comm,"   %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p);
305b18e04deSLois Curfman McInnes   }
3064537a946SLois Curfman McInnes   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
30777c4ece6SBarry Smith   PetscPrintf(snes->comm,"a particular method\n");
3086daaf66cSBarry Smith   if (snes->printhelp) (*snes->printhelp)(snes,p);
3099b94acceSBarry Smith   return 0;
3109b94acceSBarry Smith }
311a847f771SSatish Balay 
312a847f771SSatish Balay #undef __FUNCTION__
313a847f771SSatish Balay #define __FUNCTION__ "SNESSetApplicationContext"
3149b94acceSBarry Smith /*@
3159b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3169b94acceSBarry Smith    the nonlinear solvers.
3179b94acceSBarry Smith 
3189b94acceSBarry Smith    Input Parameters:
3199b94acceSBarry Smith .  snes - the SNES context
3209b94acceSBarry Smith .  usrP - optional user context
3219b94acceSBarry Smith 
3229b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3239b94acceSBarry Smith 
3249b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3259b94acceSBarry Smith @*/
3269b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3279b94acceSBarry Smith {
32877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3299b94acceSBarry Smith   snes->user		= usrP;
3309b94acceSBarry Smith   return 0;
3319b94acceSBarry Smith }
33274679c65SBarry Smith 
333a847f771SSatish Balay #undef __FUNCTION__
334a847f771SSatish Balay #define __FUNCTION__ "SNESGetApplicationContext"
3359b94acceSBarry Smith /*@C
3369b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3379b94acceSBarry Smith    nonlinear solvers.
3389b94acceSBarry Smith 
3399b94acceSBarry Smith    Input Parameter:
3409b94acceSBarry Smith .  snes - SNES context
3419b94acceSBarry Smith 
3429b94acceSBarry Smith    Output Parameter:
3439b94acceSBarry Smith .  usrP - user context
3449b94acceSBarry Smith 
3459b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3469b94acceSBarry Smith 
3479b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3489b94acceSBarry Smith @*/
3499b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3509b94acceSBarry Smith {
35177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3529b94acceSBarry Smith   *usrP = snes->user;
3539b94acceSBarry Smith   return 0;
3549b94acceSBarry Smith }
35574679c65SBarry Smith 
356a847f771SSatish Balay #undef __FUNCTION__
357a847f771SSatish Balay #define __FUNCTION__ "SNESGetIterationNumber"
3589b94acceSBarry Smith /*@
3599b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3609b94acceSBarry Smith    nonlinear solver.
3619b94acceSBarry Smith 
3629b94acceSBarry Smith    Input Parameter:
3639b94acceSBarry Smith .  snes - SNES context
3649b94acceSBarry Smith 
3659b94acceSBarry Smith    Output Parameter:
3669b94acceSBarry Smith .  iter - iteration number
3679b94acceSBarry Smith 
3689b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3699b94acceSBarry Smith @*/
3709b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3719b94acceSBarry Smith {
37277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
37374679c65SBarry Smith   PetscValidIntPointer(iter);
3749b94acceSBarry Smith   *iter = snes->iter;
3759b94acceSBarry Smith   return 0;
3769b94acceSBarry Smith }
37774679c65SBarry Smith 
378a847f771SSatish Balay #undef __FUNCTION__
379a847f771SSatish Balay #define __FUNCTION__ "SNESGetFunctionNorm"
3809b94acceSBarry Smith /*@
3819b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3829b94acceSBarry Smith    with SNESSSetFunction().
3839b94acceSBarry Smith 
3849b94acceSBarry Smith    Input Parameter:
3859b94acceSBarry Smith .  snes - SNES context
3869b94acceSBarry Smith 
3879b94acceSBarry Smith    Output Parameter:
3889b94acceSBarry Smith .  fnorm - 2-norm of function
3899b94acceSBarry Smith 
3909b94acceSBarry Smith    Note:
3919b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
392a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
393a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3949b94acceSBarry Smith 
3959b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
396a86d99e1SLois Curfman McInnes 
397a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3989b94acceSBarry Smith @*/
3999b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
4009b94acceSBarry Smith {
40177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
40274679c65SBarry Smith   PetscValidScalarPointer(fnorm);
40374679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
404*63c41f6aSSatish Balay     SETERRQ(1,"For SNES_NONLINEAR_EQUATIONS only");
40574679c65SBarry Smith   }
4069b94acceSBarry Smith   *fnorm = snes->norm;
4079b94acceSBarry Smith   return 0;
4089b94acceSBarry Smith }
40974679c65SBarry Smith 
410a847f771SSatish Balay #undef __FUNCTION__
411a847f771SSatish Balay #define __FUNCTION__ "SNESGetGradientNorm"
4129b94acceSBarry Smith /*@
4139b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
4149b94acceSBarry Smith    with SNESSSetGradient().
4159b94acceSBarry Smith 
4169b94acceSBarry Smith    Input Parameter:
4179b94acceSBarry Smith .  snes - SNES context
4189b94acceSBarry Smith 
4199b94acceSBarry Smith    Output Parameter:
4209b94acceSBarry Smith .  fnorm - 2-norm of gradient
4219b94acceSBarry Smith 
4229b94acceSBarry Smith    Note:
4239b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
424a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
425a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4269b94acceSBarry Smith 
4279b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
428a86d99e1SLois Curfman McInnes 
429a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4309b94acceSBarry Smith @*/
4319b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
4329b94acceSBarry Smith {
43377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43474679c65SBarry Smith   PetscValidScalarPointer(gnorm);
43574679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
436*63c41f6aSSatish Balay     SETERRQ(1,"For SNES_UNCONSTRAINED_MINIMIZATION only");
43774679c65SBarry Smith   }
4389b94acceSBarry Smith   *gnorm = snes->norm;
4399b94acceSBarry Smith   return 0;
4409b94acceSBarry Smith }
44174679c65SBarry Smith 
442a847f771SSatish Balay #undef __FUNCTION__
443a847f771SSatish Balay #define __FUNCTION__ "SNESGetNumberUnsuccessfulSteps"
4449b94acceSBarry Smith /*@
4459b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4469b94acceSBarry Smith    attempted by the nonlinear solver.
4479b94acceSBarry Smith 
4489b94acceSBarry Smith    Input Parameter:
4499b94acceSBarry Smith .  snes - SNES context
4509b94acceSBarry Smith 
4519b94acceSBarry Smith    Output Parameter:
4529b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4539b94acceSBarry Smith 
4549b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4559b94acceSBarry Smith @*/
4569b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4579b94acceSBarry Smith {
45877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
45974679c65SBarry Smith   PetscValidIntPointer(nfails);
4609b94acceSBarry Smith   *nfails = snes->nfailures;
4619b94acceSBarry Smith   return 0;
4629b94acceSBarry Smith }
463a847f771SSatish Balay 
464a847f771SSatish Balay #undef __FUNCTION__
465a847f771SSatish Balay #define __FUNCTION__ "SNESGetSLES"
4669b94acceSBarry Smith /*@C
4679b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4689b94acceSBarry Smith 
4699b94acceSBarry Smith    Input Parameter:
4709b94acceSBarry Smith .  snes - the SNES context
4719b94acceSBarry Smith 
4729b94acceSBarry Smith    Output Parameter:
4739b94acceSBarry Smith .  sles - the SLES context
4749b94acceSBarry Smith 
4759b94acceSBarry Smith    Notes:
4769b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4779b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4789b94acceSBarry Smith    KSP and PC contexts as well.
4799b94acceSBarry Smith 
4809b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4819b94acceSBarry Smith 
4829b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4839b94acceSBarry Smith @*/
4849b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4859b94acceSBarry Smith {
48677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4879b94acceSBarry Smith   *sles = snes->sles;
4889b94acceSBarry Smith   return 0;
4899b94acceSBarry Smith }
4909b94acceSBarry Smith /* -----------------------------------------------------------*/
491a847f771SSatish Balay #undef __FUNCTION__
492a847f771SSatish Balay #define __FUNCTION__ "SNESCreate"
4939b94acceSBarry Smith /*@C
4949b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
4959b94acceSBarry Smith 
4969b94acceSBarry Smith    Input Parameter:
4979b94acceSBarry Smith .  comm - MPI communicator
4989b94acceSBarry Smith .  type - type of method, one of
4999b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
5009b94acceSBarry Smith $      (for systems of nonlinear equations)
5019b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
5029b94acceSBarry Smith $      (for unconstrained minimization)
5039b94acceSBarry Smith 
5049b94acceSBarry Smith    Output Parameter:
5059b94acceSBarry Smith .  outsnes - the new SNES context
5069b94acceSBarry Smith 
507c1f60f51SBarry Smith    Options Database Key:
50819bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
50919bd6747SLois Curfman McInnes $              and no preconditioning matrix
51019bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
51119bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
51219bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
51319bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
514c1f60f51SBarry Smith 
5159b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5169b94acceSBarry Smith 
51763a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5189b94acceSBarry Smith @*/
5194b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5209b94acceSBarry Smith {
5219b94acceSBarry Smith   int                 ierr;
5229b94acceSBarry Smith   SNES                snes;
5239b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
52437fcc0dbSBarry Smith 
5259b94acceSBarry Smith   *outsnes = 0;
5262a0acf01SLois Curfman McInnes   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS)
527*63c41f6aSSatish Balay     SETERRQ(1,"incorrect method type");
5280452661fSBarry Smith   PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm);
5299b94acceSBarry Smith   PLogObjectCreate(snes);
5309b94acceSBarry Smith   snes->max_its           = 50;
5319b94acceSBarry Smith   snes->max_funcs	  = 1000;
5329b94acceSBarry Smith   snes->norm		  = 0.0;
5335a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
534b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
535b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5369b94acceSBarry Smith     snes->atol		  = 1.e-10;
5375a2d0531SBarry Smith   }
538b4874afaSBarry Smith   else {
539b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
540b4874afaSBarry Smith     snes->ttol            = 0.0;
541b4874afaSBarry Smith     snes->atol		  = 1.e-50;
542b4874afaSBarry Smith   }
5439b94acceSBarry Smith   snes->xtol		  = 1.e-8;
544b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5459b94acceSBarry Smith   snes->nfuncs            = 0;
5469b94acceSBarry Smith   snes->nfailures         = 0;
547639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5489b94acceSBarry Smith   snes->data              = 0;
5499b94acceSBarry Smith   snes->view              = 0;
5509b94acceSBarry Smith   snes->computeumfunction = 0;
5519b94acceSBarry Smith   snes->umfunP            = 0;
5529b94acceSBarry Smith   snes->fc                = 0;
5539b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5549b94acceSBarry Smith   snes->fmin              = -1.e30;
5559b94acceSBarry Smith   snes->method_class      = type;
5569b94acceSBarry Smith   snes->set_method_called = 0;
557a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
5589b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5597d1a2b2bSBarry Smith   snes->type              = -1;
5606f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5616f24a144SLois Curfman McInnes   snes->vwork             = 0;
5626f24a144SLois Curfman McInnes   snes->nwork             = 0;
5639b94acceSBarry Smith 
5649b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5650452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
5669b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5679b94acceSBarry Smith   kctx->version     = 2;
5689b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5699b94acceSBarry Smith                              this was too large for some test cases */
5709b94acceSBarry Smith   kctx->rtol_last   = 0;
5719b94acceSBarry Smith   kctx->rtol_max    = .9;
5729b94acceSBarry Smith   kctx->gamma       = 1.0;
5739b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5749b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5759b94acceSBarry Smith   kctx->threshold   = .1;
5769b94acceSBarry Smith   kctx->lresid_last = 0;
5779b94acceSBarry Smith   kctx->norm_last   = 0;
5789b94acceSBarry Smith 
5799b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5809b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5815334005bSBarry Smith 
5829b94acceSBarry Smith   *outsnes = snes;
5839b94acceSBarry Smith   return 0;
5849b94acceSBarry Smith }
5859b94acceSBarry Smith 
5869b94acceSBarry Smith /* --------------------------------------------------------------- */
587a847f771SSatish Balay #undef __FUNCTION__
588a847f771SSatish Balay #define __FUNCTION__ "SNESSetFunction"
5899b94acceSBarry Smith /*@C
5909b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
5919b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
5929b94acceSBarry Smith    equations.
5939b94acceSBarry Smith 
5949b94acceSBarry Smith    Input Parameters:
5959b94acceSBarry Smith .  snes - the SNES context
5969b94acceSBarry Smith .  func - function evaluation routine
5979b94acceSBarry Smith .  r - vector to store function value
5982cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
5992cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6009b94acceSBarry Smith 
6019b94acceSBarry Smith    Calling sequence of func:
602313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
6039b94acceSBarry Smith 
6049b94acceSBarry Smith .  x - input vector
605313e4042SLois Curfman McInnes .  f - function vector
6062cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
6079b94acceSBarry Smith 
6089b94acceSBarry Smith    Notes:
6099b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6109b94acceSBarry Smith $      f'(x) x = -f(x),
6119b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
6129b94acceSBarry Smith 
6139b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6149b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6159b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6169b94acceSBarry Smith 
6179b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6189b94acceSBarry Smith 
619a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6209b94acceSBarry Smith @*/
6215334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6229b94acceSBarry Smith {
62377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
6249b94acceSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
625*63c41f6aSSatish Balay     "For SNES_NONLINEAR_EQUATIONS only");
6269b94acceSBarry Smith   snes->computefunction     = func;
6279b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6289b94acceSBarry Smith   snes->funP                = ctx;
6299b94acceSBarry Smith   return 0;
6309b94acceSBarry Smith }
6319b94acceSBarry Smith 
632a847f771SSatish Balay #undef __FUNCTION__
633a847f771SSatish Balay #define __FUNCTION__ "SNESComputeFunction"
6349b94acceSBarry Smith /*@
6359b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6369b94acceSBarry Smith    SNESSetFunction().
6379b94acceSBarry Smith 
6389b94acceSBarry Smith    Input Parameters:
6399b94acceSBarry Smith .  snes - the SNES context
6409b94acceSBarry Smith .  x - input vector
6419b94acceSBarry Smith 
6429b94acceSBarry Smith    Output Parameter:
6433638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6449b94acceSBarry Smith 
6451bffabb2SLois Curfman McInnes    Notes:
6469b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6479b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6489b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6499b94acceSBarry Smith 
6509b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6519b94acceSBarry Smith 
652a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6539b94acceSBarry Smith @*/
6549b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6559b94acceSBarry Smith {
6569b94acceSBarry Smith   int    ierr;
6579b94acceSBarry Smith 
65874679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
659*63c41f6aSSatish Balay     SETERRQ(1," For SNES_NONLINEAR_EQUATIONS only");
6609b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
6619b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
6629b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6639b94acceSBarry Smith   return 0;
6649b94acceSBarry Smith }
6659b94acceSBarry Smith 
666a847f771SSatish Balay #undef __FUNCTION__
667a847f771SSatish Balay #define __FUNCTION__ "SNESSetMinimizationFunction"
6689b94acceSBarry Smith /*@C
6699b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6709b94acceSBarry Smith    unconstrained minimization.
6719b94acceSBarry Smith 
6729b94acceSBarry Smith    Input Parameters:
6739b94acceSBarry Smith .  snes - the SNES context
6749b94acceSBarry Smith .  func - function evaluation routine
675044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
676044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6779b94acceSBarry Smith 
6789b94acceSBarry Smith    Calling sequence of func:
6799b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
6809b94acceSBarry Smith 
6819b94acceSBarry Smith .  x - input vector
6829b94acceSBarry Smith .  f - function
683044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
6849b94acceSBarry Smith 
6859b94acceSBarry Smith    Notes:
6869b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
6879b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
6889b94acceSBarry Smith    SNESSetFunction().
6899b94acceSBarry Smith 
6909b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
6919b94acceSBarry Smith 
692a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
693a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
6949b94acceSBarry Smith @*/
6959b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
6969b94acceSBarry Smith                       void *ctx)
6979b94acceSBarry Smith {
69877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
6999b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
700*63c41f6aSSatish Balay     "Only for SNES_UNCONSTRAINED_MINIMIZATION");
7019b94acceSBarry Smith   snes->computeumfunction   = func;
7029b94acceSBarry Smith   snes->umfunP              = ctx;
7039b94acceSBarry Smith   return 0;
7049b94acceSBarry Smith }
7059b94acceSBarry Smith 
706a847f771SSatish Balay #undef __FUNCTION__
707a847f771SSatish Balay #define __FUNCTION__ "SNESComputeMinimizationFunction"
7089b94acceSBarry Smith /*@
7099b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7109b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7119b94acceSBarry Smith 
7129b94acceSBarry Smith    Input Parameters:
7139b94acceSBarry Smith .  snes - the SNES context
7149b94acceSBarry Smith .  x - input vector
7159b94acceSBarry Smith 
7169b94acceSBarry Smith    Output Parameter:
7179b94acceSBarry Smith .  y - function value
7189b94acceSBarry Smith 
7199b94acceSBarry Smith    Notes:
7209b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7219b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7229b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
723a86d99e1SLois Curfman McInnes 
724a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
725a86d99e1SLois Curfman McInnes 
726a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
727a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7289b94acceSBarry Smith @*/
7299b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7309b94acceSBarry Smith {
7319b94acceSBarry Smith   int    ierr;
7329b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
733*63c41f6aSSatish Balay     "Only for SNES_UNCONSTRAINED_MINIMIZATION");
7349b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
7359b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
7369b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7379b94acceSBarry Smith   return 0;
7389b94acceSBarry Smith }
7399b94acceSBarry Smith 
740a847f771SSatish Balay #undef __FUNCTION__
741a847f771SSatish Balay #define __FUNCTION__ "SNESSetGradient"
7429b94acceSBarry Smith /*@C
7439b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7449b94acceSBarry Smith    vector for use by the SNES routines.
7459b94acceSBarry Smith 
7469b94acceSBarry Smith    Input Parameters:
7479b94acceSBarry Smith .  snes - the SNES context
7489b94acceSBarry Smith .  func - function evaluation routine
749044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
750044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7519b94acceSBarry Smith .  r - vector to store gradient value
7529b94acceSBarry Smith 
7539b94acceSBarry Smith    Calling sequence of func:
7549b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7559b94acceSBarry Smith 
7569b94acceSBarry Smith .  x - input vector
7579b94acceSBarry Smith .  g - gradient vector
758044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7599b94acceSBarry Smith 
7609b94acceSBarry Smith    Notes:
7619b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7629b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7639b94acceSBarry Smith    SNESSetFunction().
7649b94acceSBarry Smith 
7659b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7669b94acceSBarry Smith 
767a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
768a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
7699b94acceSBarry Smith @*/
77074679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
7719b94acceSBarry Smith {
77277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
7739b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
774*63c41f6aSSatish Balay     "For SNES_UNCONSTRAINED_MINIMIZATION only");
7759b94acceSBarry Smith   snes->computefunction     = func;
7769b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7779b94acceSBarry Smith   snes->funP                = ctx;
7789b94acceSBarry Smith   return 0;
7799b94acceSBarry Smith }
7809b94acceSBarry Smith 
781a847f771SSatish Balay #undef __FUNCTION__
782a847f771SSatish Balay #define __FUNCTION__ "SNESComputeGradient"
7839b94acceSBarry Smith /*@
784a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
785a86d99e1SLois Curfman McInnes    SNESSetGradient().
7869b94acceSBarry Smith 
7879b94acceSBarry Smith    Input Parameters:
7889b94acceSBarry Smith .  snes - the SNES context
7899b94acceSBarry Smith .  x - input vector
7909b94acceSBarry Smith 
7919b94acceSBarry Smith    Output Parameter:
7929b94acceSBarry Smith .  y - gradient vector
7939b94acceSBarry Smith 
7949b94acceSBarry Smith    Notes:
7959b94acceSBarry Smith    SNESComputeGradient() is valid only for
7969b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7979b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
798a86d99e1SLois Curfman McInnes 
799a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
800a86d99e1SLois Curfman McInnes 
801a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
802a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8039b94acceSBarry Smith @*/
8049b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8059b94acceSBarry Smith {
8069b94acceSBarry Smith   int    ierr;
80774679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
808*63c41f6aSSatish Balay     SETERRQ(1,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8099b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
8109b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
8119b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8129b94acceSBarry Smith   return 0;
8139b94acceSBarry Smith }
8149b94acceSBarry Smith 
815a847f771SSatish Balay #undef __FUNCTION__
816a847f771SSatish Balay #define __FUNCTION__ "SNESComputeJacobian"
81762fef451SLois Curfman McInnes /*@
81862fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
81962fef451SLois Curfman McInnes    set with SNESSetJacobian().
82062fef451SLois Curfman McInnes 
82162fef451SLois Curfman McInnes    Input Parameters:
82262fef451SLois Curfman McInnes .  snes - the SNES context
82362fef451SLois Curfman McInnes .  x - input vector
82462fef451SLois Curfman McInnes 
82562fef451SLois Curfman McInnes    Output Parameters:
82662fef451SLois Curfman McInnes .  A - Jacobian matrix
82762fef451SLois Curfman McInnes .  B - optional preconditioning matrix
82862fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
82962fef451SLois Curfman McInnes 
83062fef451SLois Curfman McInnes    Notes:
83162fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
83262fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
83362fef451SLois Curfman McInnes 
834dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
835dc5a77f8SLois Curfman McInnes    flag parameter.
83662fef451SLois Curfman McInnes 
83762fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
83862fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
83962fef451SLois Curfman McInnes    methods is SNESComputeJacobian().
84062fef451SLois Curfman McInnes 
84162fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
84262fef451SLois Curfman McInnes 
84362fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
84462fef451SLois Curfman McInnes @*/
8459b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8469b94acceSBarry Smith {
8479b94acceSBarry Smith   int    ierr;
84874679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
849*63c41f6aSSatish Balay     SETERRQ(1,"For SNES_NONLINEAR_EQUATIONS only");
8509b94acceSBarry Smith   if (!snes->computejacobian) return 0;
8519b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
852c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
8539b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
8549b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8556d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
85677c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
85777c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8589b94acceSBarry Smith   return 0;
8599b94acceSBarry Smith }
8609b94acceSBarry Smith 
861a847f771SSatish Balay #undef __FUNCTION__
862a847f771SSatish Balay #define __FUNCTION__ "SNESComputeHessian"
86362fef451SLois Curfman McInnes /*@
86462fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
86562fef451SLois Curfman McInnes    set with SNESSetHessian().
86662fef451SLois Curfman McInnes 
86762fef451SLois Curfman McInnes    Input Parameters:
86862fef451SLois Curfman McInnes .  snes - the SNES context
86962fef451SLois Curfman McInnes .  x - input vector
87062fef451SLois Curfman McInnes 
87162fef451SLois Curfman McInnes    Output Parameters:
87262fef451SLois Curfman McInnes .  A - Hessian matrix
87362fef451SLois Curfman McInnes .  B - optional preconditioning matrix
87462fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
87562fef451SLois Curfman McInnes 
87662fef451SLois Curfman McInnes    Notes:
87762fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
87862fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
87962fef451SLois Curfman McInnes 
880dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
881dc5a77f8SLois Curfman McInnes    flag parameter.
88262fef451SLois Curfman McInnes 
88362fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
88462fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
88562fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
88662fef451SLois Curfman McInnes 
88762fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
88862fef451SLois Curfman McInnes 
889a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
890a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
89162fef451SLois Curfman McInnes @*/
89262fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
8939b94acceSBarry Smith {
8949b94acceSBarry Smith   int    ierr;
89574679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
896*63c41f6aSSatish Balay     SETERRQ(1,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8979b94acceSBarry Smith   if (!snes->computejacobian) return 0;
89862fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
8990f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
90062fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
90162fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9020f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
90377c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
90477c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9059b94acceSBarry Smith   return 0;
9069b94acceSBarry Smith }
9079b94acceSBarry Smith 
908a847f771SSatish Balay #undef __FUNCTION__
909a847f771SSatish Balay #define __FUNCTION__ "SNESSetJacobian"
9109b94acceSBarry Smith /*@C
9119b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
912044dda88SLois Curfman McInnes    location to store the matrix.
9139b94acceSBarry Smith 
9149b94acceSBarry Smith    Input Parameters:
9159b94acceSBarry Smith .  snes - the SNES context
9169b94acceSBarry Smith .  A - Jacobian matrix
9179b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9189b94acceSBarry Smith .  func - Jacobian evaluation routine
9192cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9202cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9219b94acceSBarry Smith 
9229b94acceSBarry Smith    Calling sequence of func:
923313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9249b94acceSBarry Smith 
9259b94acceSBarry Smith .  x - input vector
9269b94acceSBarry Smith .  A - Jacobian matrix
9279b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
928ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
929ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9302cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
9319b94acceSBarry Smith 
9329b94acceSBarry Smith    Notes:
933dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9342cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
935ac21db08SLois Curfman McInnes 
936ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9379b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9389b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9399b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9409b94acceSBarry Smith    throughout the global iterations.
9419b94acceSBarry Smith 
9429b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9439b94acceSBarry Smith 
944ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
9459b94acceSBarry Smith @*/
9469b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9479b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9489b94acceSBarry Smith {
94977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
95074679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
951*63c41f6aSSatish Balay     SETERRQ(1,"For SNES_NONLINEAR_EQUATIONS only");
9529b94acceSBarry Smith   snes->computejacobian = func;
9539b94acceSBarry Smith   snes->jacP            = ctx;
9549b94acceSBarry Smith   snes->jacobian        = A;
9559b94acceSBarry Smith   snes->jacobian_pre    = B;
9569b94acceSBarry Smith   return 0;
9579b94acceSBarry Smith }
95862fef451SLois Curfman McInnes 
959a847f771SSatish Balay #undef __FUNCTION__
960a847f771SSatish Balay #define __FUNCTION__ "SNESGetJacobian"
961b4fd4287SBarry Smith /*@
962b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
963b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
964b4fd4287SBarry Smith 
965b4fd4287SBarry Smith    Input Parameter:
966b4fd4287SBarry Smith .  snes - the nonlinear solver context
967b4fd4287SBarry Smith 
968b4fd4287SBarry Smith    Output Parameters:
969b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
970b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
971b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
972b4fd4287SBarry Smith 
973b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
974b4fd4287SBarry Smith @*/
975b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
976b4fd4287SBarry Smith {
97774679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
978*63c41f6aSSatish Balay     SETERRQ(1,"For SNES_NONLINEAR_EQUATIONS only");
979b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
980b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
981b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
982b4fd4287SBarry Smith   return 0;
983b4fd4287SBarry Smith }
984b4fd4287SBarry Smith 
985a847f771SSatish Balay #undef __FUNCTION__
986a847f771SSatish Balay #define __FUNCTION__ "SNESSetHessian"
9879b94acceSBarry Smith /*@C
9889b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
989044dda88SLois Curfman McInnes    location to store the matrix.
9909b94acceSBarry Smith 
9919b94acceSBarry Smith    Input Parameters:
9929b94acceSBarry Smith .  snes - the SNES context
9939b94acceSBarry Smith .  A - Hessian matrix
9949b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
9959b94acceSBarry Smith .  func - Jacobian evaluation routine
996313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
997313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
9989b94acceSBarry Smith 
9999b94acceSBarry Smith    Calling sequence of func:
1000313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10019b94acceSBarry Smith 
10029b94acceSBarry Smith .  x - input vector
10039b94acceSBarry Smith .  A - Hessian matrix
10049b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1005ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1006ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10072cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10089b94acceSBarry Smith 
10099b94acceSBarry Smith    Notes:
1010dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10112cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1012ac21db08SLois Curfman McInnes 
10139b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10149b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10159b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10169b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10179b94acceSBarry Smith    throughout the global iterations.
10189b94acceSBarry Smith 
10199b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10209b94acceSBarry Smith 
1021ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10229b94acceSBarry Smith @*/
10239b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10249b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10259b94acceSBarry Smith {
102677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
102774679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
1028*63c41f6aSSatish Balay     SETERRQ(1,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10299b94acceSBarry Smith   snes->computejacobian = func;
10309b94acceSBarry Smith   snes->jacP            = ctx;
10319b94acceSBarry Smith   snes->jacobian        = A;
10329b94acceSBarry Smith   snes->jacobian_pre    = B;
10339b94acceSBarry Smith   return 0;
10349b94acceSBarry Smith }
10359b94acceSBarry Smith 
1036a847f771SSatish Balay #undef __FUNCTION__
1037a847f771SSatish Balay #define __FUNCTION__ "SNESGetHessian"
103862fef451SLois Curfman McInnes /*@
103962fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
104062fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
104162fef451SLois Curfman McInnes 
104262fef451SLois Curfman McInnes    Input Parameter:
104362fef451SLois Curfman McInnes .  snes - the nonlinear solver context
104462fef451SLois Curfman McInnes 
104562fef451SLois Curfman McInnes    Output Parameters:
104662fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
104762fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
104862fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
104962fef451SLois Curfman McInnes 
105062fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
105162fef451SLois Curfman McInnes @*/
105262fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
105362fef451SLois Curfman McInnes {
105474679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
1055*63c41f6aSSatish Balay     SETERRQ(1,"For SNES_UNCONSTRAINED_MINIMIZATION only");
105662fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
105762fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
105862fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
105962fef451SLois Curfman McInnes   return 0;
106062fef451SLois Curfman McInnes }
106162fef451SLois Curfman McInnes 
10629b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10639b94acceSBarry Smith 
1064a847f771SSatish Balay #undef __FUNCTION__
1065a847f771SSatish Balay #define __FUNCTION__ "SNESSetUp"
10669b94acceSBarry Smith /*@
10679b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1068272ac6f2SLois Curfman McInnes    of a nonlinear solver.
10699b94acceSBarry Smith 
10709b94acceSBarry Smith    Input Parameter:
10719b94acceSBarry Smith .  snes - the SNES context
10728ddd3da0SLois Curfman McInnes .  x - the solution vector
10739b94acceSBarry Smith 
1074272ac6f2SLois Curfman McInnes    Notes:
1075272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1076272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1077272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1078272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1079272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1080272ac6f2SLois Curfman McInnes 
10819b94acceSBarry Smith .keywords: SNES, nonlinear, setup
10829b94acceSBarry Smith 
10839b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
10849b94acceSBarry Smith @*/
10858ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
10869b94acceSBarry Smith {
1087272ac6f2SLois Curfman McInnes   int ierr, flg;
108877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
108977c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
10908ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
10919b94acceSBarry Smith 
1092c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1093c1f60f51SBarry Smith   /*
1094c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1095dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1096c1f60f51SBarry Smith   */
1097c1f60f51SBarry Smith   if (flg) {
1098c1f60f51SBarry Smith     Mat J;
1099c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1100c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1101c1f60f51SBarry Smith     snes->mfshell = J;
1102c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1103c1f60f51SBarry Smith       snes->jacobian = J;
1104c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1105c1f60f51SBarry Smith     }
1106c1f60f51SBarry Smith     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1107c1f60f51SBarry Smith       snes->jacobian = J;
1108c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1109*63c41f6aSSatish Balay     } else SETERRQ(1,"Method class doesn't support matrix-free operator option");
1110c1f60f51SBarry Smith   }
1111272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1112c1f60f51SBarry Smith   /*
1113dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1114c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1115c1f60f51SBarry Smith    */
111631615d04SLois Curfman McInnes   if (flg) {
1117272ac6f2SLois Curfman McInnes     Mat J;
1118272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1119272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1120272ac6f2SLois Curfman McInnes     snes->mfshell = J;
112183e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
112283e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
112394a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
112483e56afdSLois Curfman McInnes     }
112583e56afdSLois Curfman McInnes     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
112683e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
112794a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1128*63c41f6aSSatish Balay     } else SETERRQ(1,"Method class doesn't support matrix-free option");
1129272ac6f2SLois Curfman McInnes   }
11309b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1131*63c41f6aSSatish Balay     if (!snes->vec_func) SETERRQ(1,"Must call SNESSetFunction() first");
1132*63c41f6aSSatish Balay     if (!snes->computefunction) SETERRQ(1,"Must call SNESSetFunction() first");
1133*63c41f6aSSatish Balay     if (!snes->jacobian) SETERRQ(1,"Must call SNESSetJacobian() first");
1134*63c41f6aSSatish Balay     if (snes->vec_func == snes->vec_sol) SETERRQ(1,"Solution vector cannot be function vector");
1135a703fe33SLois Curfman McInnes 
1136a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
113740191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1138a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1139a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1140a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1141a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1142a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1143a703fe33SLois Curfman McInnes     }
11449b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1145*63c41f6aSSatish Balay     if (!snes->vec_func) SETERRQ(1,"Must call SNESSetGradient() first");
1146*63c41f6aSSatish Balay     if (!snes->computefunction) SETERRQ(1,"Must call SNESSetGradient() first");
114737fcc0dbSBarry Smith     if (!snes->computeumfunction)
1148*63c41f6aSSatish Balay       SETERRQ(1,"Must call SNESSetMinimizationFunction() first");
1149*63c41f6aSSatish Balay     if (!snes->jacobian) SETERRQ(1,"Must call SNESSetHessian() first");
1150*63c41f6aSSatish Balay   } else SETERRQ(1,"Unknown method class");
1151a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1152a703fe33SLois Curfman McInnes   snes->setup_called = 1;
1153a703fe33SLois Curfman McInnes   return 0;
11549b94acceSBarry Smith }
11559b94acceSBarry Smith 
1156a847f771SSatish Balay #undef __FUNCTION__
1157a847f771SSatish Balay #define __FUNCTION__ "SNESDestroy"
11589b94acceSBarry Smith /*@C
11599b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11609b94acceSBarry Smith    with SNESCreate().
11619b94acceSBarry Smith 
11629b94acceSBarry Smith    Input Parameter:
11639b94acceSBarry Smith .  snes - the SNES context
11649b94acceSBarry Smith 
11659b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
11669b94acceSBarry Smith 
116763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
11689b94acceSBarry Smith @*/
11699b94acceSBarry Smith int SNESDestroy(SNES snes)
11709b94acceSBarry Smith {
11719b94acceSBarry Smith   int ierr;
117277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
11739b94acceSBarry Smith   ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);
11740452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
11759b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
11769b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
11773e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
11786f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
11799b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
11800452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
11819b94acceSBarry Smith   return 0;
11829b94acceSBarry Smith }
11839b94acceSBarry Smith 
11849b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11859b94acceSBarry Smith 
1186a847f771SSatish Balay #undef __FUNCTION__
1187a847f771SSatish Balay #define __FUNCTION__ "SNESSetTolerances"
11889b94acceSBarry Smith /*@
1189d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11909b94acceSBarry Smith 
11919b94acceSBarry Smith    Input Parameters:
11929b94acceSBarry Smith .  snes - the SNES context
119333174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
119433174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
119533174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
119633174efeSLois Curfman McInnes            of the change in the solution between steps
119733174efeSLois Curfman McInnes .  maxit - maximum number of iterations
119833174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
11999b94acceSBarry Smith 
120033174efeSLois Curfman McInnes    Options Database Keys:
120133174efeSLois Curfman McInnes $    -snes_atol <atol>
120233174efeSLois Curfman McInnes $    -snes_rtol <rtol>
120333174efeSLois Curfman McInnes $    -snes_stol <stol>
120433174efeSLois Curfman McInnes $    -snes_max_it <maxit>
120533174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
12069b94acceSBarry Smith 
1207d7a720efSLois Curfman McInnes    Notes:
12089b94acceSBarry Smith    The default maximum number of iterations is 50.
12099b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
12109b94acceSBarry Smith 
121133174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
12129b94acceSBarry Smith 
1213d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
12149b94acceSBarry Smith @*/
1215d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
12169b94acceSBarry Smith {
121777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1218d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1219d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1220d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1221d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1222d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
12239b94acceSBarry Smith   return 0;
12249b94acceSBarry Smith }
12259b94acceSBarry Smith 
1226a847f771SSatish Balay #undef __FUNCTION__
1227a847f771SSatish Balay #define __FUNCTION__ "SNESGetTolerances"
12289b94acceSBarry Smith /*@
122933174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
123033174efeSLois Curfman McInnes 
123133174efeSLois Curfman McInnes    Input Parameters:
123233174efeSLois Curfman McInnes .  snes - the SNES context
123333174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
123433174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
123533174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
123633174efeSLois Curfman McInnes            of the change in the solution between steps
123733174efeSLois Curfman McInnes .  maxit - maximum number of iterations
123833174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
123933174efeSLois Curfman McInnes 
124033174efeSLois Curfman McInnes    Notes:
124133174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
124233174efeSLois Curfman McInnes 
124333174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
124433174efeSLois Curfman McInnes 
124533174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
124633174efeSLois Curfman McInnes @*/
124733174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
124833174efeSLois Curfman McInnes {
124933174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
125033174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
125133174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
125233174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
125333174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
125433174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
125533174efeSLois Curfman McInnes   return 0;
125633174efeSLois Curfman McInnes }
125733174efeSLois Curfman McInnes 
1258a847f771SSatish Balay #undef __FUNCTION__
1259a847f771SSatish Balay #define __FUNCTION__ "SNESSetTrustRegionTolerance"
126033174efeSLois Curfman McInnes /*@
12619b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12629b94acceSBarry Smith 
12639b94acceSBarry Smith    Input Parameters:
12649b94acceSBarry Smith .  snes - the SNES context
12659b94acceSBarry Smith .  tol - tolerance
12669b94acceSBarry Smith 
12679b94acceSBarry Smith    Options Database Key:
1268d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
12699b94acceSBarry Smith 
12709b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12719b94acceSBarry Smith 
1272d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
12739b94acceSBarry Smith @*/
12749b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
12759b94acceSBarry Smith {
127677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12779b94acceSBarry Smith   snes->deltatol = tol;
12789b94acceSBarry Smith   return 0;
12799b94acceSBarry Smith }
12809b94acceSBarry Smith 
1281a847f771SSatish Balay #undef __FUNCTION__
1282a847f771SSatish Balay #define __FUNCTION__ "SNESSetMinimizationFunctionTolerance"
12839b94acceSBarry Smith /*@
128477c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
12859b94acceSBarry Smith    for unconstrained minimization solvers.
12869b94acceSBarry Smith 
12879b94acceSBarry Smith    Input Parameters:
12889b94acceSBarry Smith .  snes - the SNES context
12899b94acceSBarry Smith .  ftol - minimum function tolerance
12909b94acceSBarry Smith 
12919b94acceSBarry Smith    Options Database Key:
1292d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
12939b94acceSBarry Smith 
12949b94acceSBarry Smith    Note:
129577c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
12969b94acceSBarry Smith    methods only.
12979b94acceSBarry Smith 
12989b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
12999b94acceSBarry Smith 
1300d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
13019b94acceSBarry Smith @*/
130277c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
13039b94acceSBarry Smith {
130477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13059b94acceSBarry Smith   snes->fmin = ftol;
13069b94acceSBarry Smith   return 0;
13079b94acceSBarry Smith }
13089b94acceSBarry Smith 
13099b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
13109b94acceSBarry Smith 
1311a847f771SSatish Balay #undef __FUNCTION__
1312a847f771SSatish Balay #define __FUNCTION__ "SNESSetMonitor"
13139b94acceSBarry Smith /*@C
13149b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
13159b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
13169b94acceSBarry Smith    progress.
13179b94acceSBarry Smith 
13189b94acceSBarry Smith    Input Parameters:
13199b94acceSBarry Smith .  snes - the SNES context
13209b94acceSBarry Smith .  func - monitoring routine
1321044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1322044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
13239b94acceSBarry Smith 
13249b94acceSBarry Smith    Calling sequence of func:
1325682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
13269b94acceSBarry Smith 
13279b94acceSBarry Smith $    snes - the SNES context
13289b94acceSBarry Smith $    its - iteration number
1329044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
13309b94acceSBarry Smith $
13319b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13329b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
13339b94acceSBarry Smith $
13349b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13359b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
13369b94acceSBarry Smith 
13379665c990SLois Curfman McInnes    Options Database Keys:
13389665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
13399665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
13409665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
13419665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
13429665c990SLois Curfman McInnes $                           been hardwired into a code by
13439665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
13449665c990SLois Curfman McInnes $                           does not cancel those set via
13459665c990SLois Curfman McInnes $                           the options database.
13469665c990SLois Curfman McInnes 
13479665c990SLois Curfman McInnes 
1348639f9d9dSBarry Smith    Notes:
13496bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13506bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13516bc08f3fSLois Curfman McInnes    order in which they were set.
1352639f9d9dSBarry Smith 
13539b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13549b94acceSBarry Smith 
13559b94acceSBarry Smith .seealso: SNESDefaultMonitor()
13569b94acceSBarry Smith @*/
135774679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
13589b94acceSBarry Smith {
1359639f9d9dSBarry Smith   if (!func) {
1360639f9d9dSBarry Smith     snes->numbermonitors = 0;
1361639f9d9dSBarry Smith     return 0;
1362639f9d9dSBarry Smith   }
1363639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1364*63c41f6aSSatish Balay     SETERRQ(1,"Too many monitors set");
1365639f9d9dSBarry Smith   }
1366639f9d9dSBarry Smith 
1367639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1368639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
13699b94acceSBarry Smith   return 0;
13709b94acceSBarry Smith }
13719b94acceSBarry Smith 
1372a847f771SSatish Balay #undef __FUNCTION__
1373a847f771SSatish Balay #define __FUNCTION__ "SNESSetConvergenceTest"
13749b94acceSBarry Smith /*@C
13759b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
13769b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
13779b94acceSBarry Smith 
13789b94acceSBarry Smith    Input Parameters:
13799b94acceSBarry Smith .  snes - the SNES context
13809b94acceSBarry Smith .  func - routine to test for convergence
1381044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1382044dda88SLois Curfman McInnes           (may be PETSC_NULL)
13839b94acceSBarry Smith 
13849b94acceSBarry Smith    Calling sequence of func:
13859b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
13869b94acceSBarry Smith              double f,void *cctx)
13879b94acceSBarry Smith 
13889b94acceSBarry Smith $    snes - the SNES context
1389044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
13909b94acceSBarry Smith $    xnorm - 2-norm of current iterate
13919b94acceSBarry Smith $
13929b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13939b94acceSBarry Smith $    gnorm - 2-norm of current step
13949b94acceSBarry Smith $    f - 2-norm of function
13959b94acceSBarry Smith $
13969b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13979b94acceSBarry Smith $    gnorm - 2-norm of current gradient
13989b94acceSBarry Smith $    f - function value
13999b94acceSBarry Smith 
14009b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14019b94acceSBarry Smith 
140240191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
140340191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
14049b94acceSBarry Smith @*/
140574679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
14069b94acceSBarry Smith {
14079b94acceSBarry Smith   (snes)->converged = func;
14089b94acceSBarry Smith   (snes)->cnvP      = cctx;
14099b94acceSBarry Smith   return 0;
14109b94acceSBarry Smith }
14119b94acceSBarry Smith 
1412a847f771SSatish Balay #undef __FUNCTION__
1413a847f771SSatish Balay #define __FUNCTION__ "SNESScaleStep_Private"
14149b94acceSBarry Smith /*
14159b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
14169b94acceSBarry Smith    positive parameter delta.
14179b94acceSBarry Smith 
14189b94acceSBarry Smith     Input Parameters:
14199b94acceSBarry Smith .   snes - the SNES context
14209b94acceSBarry Smith .   y - approximate solution of linear system
14219b94acceSBarry Smith .   fnorm - 2-norm of current function
14229b94acceSBarry Smith .   delta - trust region size
14239b94acceSBarry Smith 
14249b94acceSBarry Smith     Output Parameters:
14259b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
14269b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
14279b94acceSBarry Smith     region, and exceeds zero otherwise.
14289b94acceSBarry Smith .   ynorm - 2-norm of the step
14299b94acceSBarry Smith 
14309b94acceSBarry Smith     Note:
143140191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
14329b94acceSBarry Smith     is set to be the maximum allowable step size.
14339b94acceSBarry Smith 
14349b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
14359b94acceSBarry Smith */
14369b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
14379b94acceSBarry Smith                           double *gpnorm,double *ynorm)
14389b94acceSBarry Smith {
14399b94acceSBarry Smith   double norm;
14409b94acceSBarry Smith   Scalar cnorm;
1441cddf8d76SBarry Smith   VecNorm(y,NORM_2, &norm );
14429b94acceSBarry Smith   if (norm > *delta) {
14439b94acceSBarry Smith      norm = *delta/norm;
14449b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
14459b94acceSBarry Smith      cnorm = norm;
14469b94acceSBarry Smith      VecScale( &cnorm, y );
14479b94acceSBarry Smith      *ynorm = *delta;
14489b94acceSBarry Smith   } else {
14499b94acceSBarry Smith      *gpnorm = 0.0;
14509b94acceSBarry Smith      *ynorm = norm;
14519b94acceSBarry Smith   }
14529b94acceSBarry Smith   return 0;
14539b94acceSBarry Smith }
14549b94acceSBarry Smith 
1455a847f771SSatish Balay #undef __FUNCTION__
1456a847f771SSatish Balay #define __FUNCTION__ "SNESSolve"
14579b94acceSBarry Smith /*@
14589b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
145963a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
14609b94acceSBarry Smith 
14619b94acceSBarry Smith    Input Parameter:
14629b94acceSBarry Smith .  snes - the SNES context
14638ddd3da0SLois Curfman McInnes .  x - the solution vector
14649b94acceSBarry Smith 
14659b94acceSBarry Smith    Output Parameter:
14669b94acceSBarry Smith    its - number of iterations until termination
14679b94acceSBarry Smith 
14688ddd3da0SLois Curfman McInnes    Note:
14698ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
14708ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
14718ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
14728ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
14738ddd3da0SLois Curfman McInnes 
14749b94acceSBarry Smith .keywords: SNES, nonlinear, solve
14759b94acceSBarry Smith 
147663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
14779b94acceSBarry Smith @*/
14788ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
14799b94acceSBarry Smith {
14803c7409f5SSatish Balay   int ierr, flg;
1481052efed2SBarry Smith 
148277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
148374679c65SBarry Smith   PetscValidIntPointer(its);
1484c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1485c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
14869b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
14879b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
14889b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
14893c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
14906d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
14919b94acceSBarry Smith   return 0;
14929b94acceSBarry Smith }
14939b94acceSBarry Smith 
14949b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
149537fcc0dbSBarry Smith static NRList *__SNESList = 0;
14969b94acceSBarry Smith 
1497a847f771SSatish Balay #undef __FUNCTION__
1498a847f771SSatish Balay #define __FUNCTION__ "SNESSetType"
14999b94acceSBarry Smith /*@
15004b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
15019b94acceSBarry Smith 
15029b94acceSBarry Smith    Input Parameters:
15039b94acceSBarry Smith .  snes - the SNES context
15049b94acceSBarry Smith .  method - a known method
15059b94acceSBarry Smith 
1506ae12b187SLois Curfman McInnes   Options Database Command:
1507ae12b187SLois Curfman McInnes $ -snes_type  <method>
1508ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1509ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1510ae12b187SLois Curfman McInnes 
15119b94acceSBarry Smith    Notes:
15129b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
15139b94acceSBarry Smith $  Systems of nonlinear equations:
151440191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
151540191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
15169b94acceSBarry Smith $  Unconstrained minimization:
151740191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
151840191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
15199b94acceSBarry Smith 
1520ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1521ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1522ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1523ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1524ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1525ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1526ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1527ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1528ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1529ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1530a703fe33SLois Curfman McInnes 
1531f536c699SSatish Balay .keywords: SNES, set, method
15329b94acceSBarry Smith @*/
15334b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
15349b94acceSBarry Smith {
15359b94acceSBarry Smith   int (*r)(SNES);
153677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15377d1a2b2bSBarry Smith   if (snes->type == (int) method) return 0;
15387d1a2b2bSBarry Smith 
15399b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
154037fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
1541*63c41f6aSSatish Balay   if (!__SNESList) {SETERRQ(1,"Could not get methods");}
154237fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1543*63c41f6aSSatish Balay   if (!r) {SETERRQ(1,"Unknown method");}
15440452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
15459b94acceSBarry Smith   snes->set_method_called = 1;
15469b94acceSBarry Smith   return (*r)(snes);
15479b94acceSBarry Smith }
15489b94acceSBarry Smith 
15499b94acceSBarry Smith /* --------------------------------------------------------------------- */
1550a847f771SSatish Balay #undef __FUNCTION__
1551a847f771SSatish Balay #define __FUNCTION__ "SNESRegister"
15529b94acceSBarry Smith /*@C
15539b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
15544b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
15559b94acceSBarry Smith 
15569b94acceSBarry Smith    Input Parameters:
155740191667SLois Curfman McInnes .  name - for instance SNES_EQ_LS, SNES_EQ_TR, ...
155840191667SLois Curfman McInnes .  sname - corresponding string for name
15599b94acceSBarry Smith .  create - routine to create method context
15609b94acceSBarry Smith 
15619b94acceSBarry Smith .keywords: SNES, nonlinear, register
15629b94acceSBarry Smith 
15639b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
15649b94acceSBarry Smith @*/
15659b94acceSBarry Smith int SNESRegister(int name, char *sname, int (*create)(SNES))
15669b94acceSBarry Smith {
15679b94acceSBarry Smith   int ierr;
156837fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
156937fcc0dbSBarry Smith   NRRegister( __SNESList, name, sname, (int (*)(void*))create );
15709b94acceSBarry Smith   return 0;
15719b94acceSBarry Smith }
1572a847f771SSatish Balay 
15739b94acceSBarry Smith /* --------------------------------------------------------------------- */
1574a847f771SSatish Balay #undef __FUNCTION__
1575a847f771SSatish Balay #define __FUNCTION__ "SNESRegisterDestroy"
15769b94acceSBarry Smith /*@C
15779b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
15789b94acceSBarry Smith    registered by SNESRegister().
15799b94acceSBarry Smith 
15809b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
15819b94acceSBarry Smith 
15829b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
15839b94acceSBarry Smith @*/
15849b94acceSBarry Smith int SNESRegisterDestroy()
15859b94acceSBarry Smith {
158637fcc0dbSBarry Smith   if (__SNESList) {
158737fcc0dbSBarry Smith     NRDestroy( __SNESList );
158837fcc0dbSBarry Smith     __SNESList = 0;
15899b94acceSBarry Smith   }
15909b94acceSBarry Smith   return 0;
15919b94acceSBarry Smith }
15929b94acceSBarry Smith 
1593a847f771SSatish Balay #undef __FUNCTION__
1594a847f771SSatish Balay #define __FUNCTION__ "SNESGetTypeFromOptions_Private"
15959b94acceSBarry Smith /*
15964b0e389bSBarry Smith    SNESGetTypeFromOptions_Private - Sets the selected method from the
15979b94acceSBarry Smith    options database.
15989b94acceSBarry Smith 
15999b94acceSBarry Smith    Input Parameter:
16009b94acceSBarry Smith .  ctx - the SNES context
16019b94acceSBarry Smith 
16029b94acceSBarry Smith    Output Parameter:
16039b94acceSBarry Smith .  method -  solver method
16049b94acceSBarry Smith 
16059b94acceSBarry Smith    Returns:
16069b94acceSBarry Smith    Returns 1 if the method is found; 0 otherwise.
16079b94acceSBarry Smith 
16089b94acceSBarry Smith    Options Database Key:
16094b0e389bSBarry Smith $  -snes_type  method
16109b94acceSBarry Smith */
1611052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
16129b94acceSBarry Smith {
1613052efed2SBarry Smith   int  ierr;
16149b94acceSBarry Smith   char sbuf[50];
16155baf8537SBarry Smith 
1616052efed2SBarry Smith   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1617052efed2SBarry Smith   if (*flg) {
1618052efed2SBarry Smith     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
161937fcc0dbSBarry Smith     *method = (SNESType)NRFindID( __SNESList, sbuf );
16209b94acceSBarry Smith   }
16219b94acceSBarry Smith   return 0;
16229b94acceSBarry Smith }
16239b94acceSBarry Smith 
1624a847f771SSatish Balay #undef __FUNCTION__
1625a847f771SSatish Balay #define __FUNCTION__ "SNESGetType"
16269b94acceSBarry Smith /*@C
16279a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
16289b94acceSBarry Smith 
16299b94acceSBarry Smith    Input Parameter:
16304b0e389bSBarry Smith .  snes - nonlinear solver context
16319b94acceSBarry Smith 
16329b94acceSBarry Smith    Output Parameter:
16339a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
16349a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
16359b94acceSBarry Smith 
16369b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
16379b94acceSBarry Smith @*/
16384b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
16399b94acceSBarry Smith {
164006a9b0adSLois Curfman McInnes   int ierr;
164137fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
16424b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
164337fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
16449b94acceSBarry Smith   return 0;
16459b94acceSBarry Smith }
16469b94acceSBarry Smith 
16479b94acceSBarry Smith #include <stdio.h>
1648a847f771SSatish Balay #undef __FUNCTION__
1649a847f771SSatish Balay #define __FUNCTION__ "SNESPrintTypes_Private"
16509b94acceSBarry Smith /*
16514b0e389bSBarry Smith    SNESPrintTypes_Private - Prints the SNES methods available from the
16529b94acceSBarry Smith    options database.
16539b94acceSBarry Smith 
16549b94acceSBarry Smith    Input Parameters:
165533455573SSatish Balay .  comm   - communicator (usually MPI_COMM_WORLD)
16569b94acceSBarry Smith .  prefix - prefix (usually "-")
16574b0e389bSBarry Smith .  name   - the options database name (by default "snes_type")
16589b94acceSBarry Smith */
165933455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
16609b94acceSBarry Smith {
16619b94acceSBarry Smith   FuncList *entry;
166237fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
166337fcc0dbSBarry Smith   entry = __SNESList->head;
166477c4ece6SBarry Smith   PetscPrintf(comm," %s%s (one of)",prefix,name);
16659b94acceSBarry Smith   while (entry) {
166677c4ece6SBarry Smith     PetscPrintf(comm," %s",entry->name);
16679b94acceSBarry Smith     entry = entry->next;
16689b94acceSBarry Smith   }
166977c4ece6SBarry Smith   PetscPrintf(comm,"\n");
16709b94acceSBarry Smith   return 0;
16719b94acceSBarry Smith }
16729b94acceSBarry Smith 
1673a847f771SSatish Balay #undef __FUNCTION__
1674a847f771SSatish Balay #define __FUNCTION__ "SNESGetSolution"
16759b94acceSBarry Smith /*@C
16769b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
16779b94acceSBarry Smith    stored.
16789b94acceSBarry Smith 
16799b94acceSBarry Smith    Input Parameter:
16809b94acceSBarry Smith .  snes - the SNES context
16819b94acceSBarry Smith 
16829b94acceSBarry Smith    Output Parameter:
16839b94acceSBarry Smith .  x - the solution
16849b94acceSBarry Smith 
16859b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
16869b94acceSBarry Smith 
1687a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
16889b94acceSBarry Smith @*/
16899b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
16909b94acceSBarry Smith {
169177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16929b94acceSBarry Smith   *x = snes->vec_sol_always;
16939b94acceSBarry Smith   return 0;
16949b94acceSBarry Smith }
16959b94acceSBarry Smith 
1696a847f771SSatish Balay #undef __FUNCTION__
1697a847f771SSatish Balay #define __FUNCTION__ "SNESGetSolutionUpdate"
16989b94acceSBarry Smith /*@C
16999b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
17009b94acceSBarry Smith    stored.
17019b94acceSBarry Smith 
17029b94acceSBarry Smith    Input Parameter:
17039b94acceSBarry Smith .  snes - the SNES context
17049b94acceSBarry Smith 
17059b94acceSBarry Smith    Output Parameter:
17069b94acceSBarry Smith .  x - the solution update
17079b94acceSBarry Smith 
17089b94acceSBarry Smith    Notes:
17099b94acceSBarry Smith    This vector is implementation dependent.
17109b94acceSBarry Smith 
17119b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
17129b94acceSBarry Smith 
17139b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
17149b94acceSBarry Smith @*/
17159b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
17169b94acceSBarry Smith {
171777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17189b94acceSBarry Smith   *x = snes->vec_sol_update_always;
17199b94acceSBarry Smith   return 0;
17209b94acceSBarry Smith }
17219b94acceSBarry Smith 
1722a847f771SSatish Balay #undef __FUNCTION__
1723a847f771SSatish Balay #define __FUNCTION__ "SNESGetFunction"
17249b94acceSBarry Smith /*@C
17253638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
17269b94acceSBarry Smith 
17279b94acceSBarry Smith    Input Parameter:
17289b94acceSBarry Smith .  snes - the SNES context
17299b94acceSBarry Smith 
17309b94acceSBarry Smith    Output Parameter:
17313638b69dSLois Curfman McInnes .  r - the function
17329b94acceSBarry Smith 
17339b94acceSBarry Smith    Notes:
17349b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
17359b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
17369b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
17379b94acceSBarry Smith 
1738a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
17399b94acceSBarry Smith 
174031615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
174131615d04SLois Curfman McInnes           SNESGetGradient()
17429b94acceSBarry Smith @*/
17439b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
17449b94acceSBarry Smith {
174577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17469b94acceSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
1747*63c41f6aSSatish Balay     For SNES_NONLINEAR_EQUATIONS only");
17489b94acceSBarry Smith   *r = snes->vec_func_always;
17499b94acceSBarry Smith   return 0;
17509b94acceSBarry Smith }
17519b94acceSBarry Smith 
1752a847f771SSatish Balay #undef __FUNCTION__
1753a847f771SSatish Balay #define __FUNCTION__ "SNESGetGradient"
17549b94acceSBarry Smith /*@C
17553638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
17569b94acceSBarry Smith 
17579b94acceSBarry Smith    Input Parameter:
17589b94acceSBarry Smith .  snes - the SNES context
17599b94acceSBarry Smith 
17609b94acceSBarry Smith    Output Parameter:
17619b94acceSBarry Smith .  r - the gradient
17629b94acceSBarry Smith 
17639b94acceSBarry Smith    Notes:
17649b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
17659b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
17669b94acceSBarry Smith    SNESGetFunction().
17679b94acceSBarry Smith 
17689b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
17699b94acceSBarry Smith 
177031615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
17719b94acceSBarry Smith @*/
17729b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
17739b94acceSBarry Smith {
177477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17759b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
1776*63c41f6aSSatish Balay     "For SNES_UNCONSTRAINED_MINIMIZATION only");
17779b94acceSBarry Smith   *r = snes->vec_func_always;
17789b94acceSBarry Smith   return 0;
17799b94acceSBarry Smith }
17809b94acceSBarry Smith 
1781a847f771SSatish Balay #undef __FUNCTION__
1782a847f771SSatish Balay #define __FUNCTION__ "SNESGetMinimizationFunction"
17839b94acceSBarry Smith /*@
17849b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
17859b94acceSBarry Smith    unconstrained minimization problems.
17869b94acceSBarry Smith 
17879b94acceSBarry Smith    Input Parameter:
17889b94acceSBarry Smith .  snes - the SNES context
17899b94acceSBarry Smith 
17909b94acceSBarry Smith    Output Parameter:
17919b94acceSBarry Smith .  r - the function
17929b94acceSBarry Smith 
17939b94acceSBarry Smith    Notes:
17949b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
17959b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
17969b94acceSBarry Smith    SNESGetFunction().
17979b94acceSBarry Smith 
17989b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
17999b94acceSBarry Smith 
180031615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
18019b94acceSBarry Smith @*/
18029b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
18039b94acceSBarry Smith {
180477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
180574679c65SBarry Smith   PetscValidScalarPointer(r);
18069b94acceSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,
1807*63c41f6aSSatish Balay     "For SNES_UNCONSTRAINED_MINIMIZATION only");
18089b94acceSBarry Smith   *r = snes->fc;
18099b94acceSBarry Smith   return 0;
18109b94acceSBarry Smith }
18119b94acceSBarry Smith 
1812a847f771SSatish Balay #undef __FUNCTION__
1813a847f771SSatish Balay #define __FUNCTION__ "SNESSetOptionsPrefix"
18143c7409f5SSatish Balay /*@C
18153c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1816639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1817639f9d9dSBarry Smith    the prefix name.
18183c7409f5SSatish Balay 
18193c7409f5SSatish Balay    Input Parameter:
18203c7409f5SSatish Balay .  snes - the SNES context
18213c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18223c7409f5SSatish Balay 
18233c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1824a86d99e1SLois Curfman McInnes 
1825a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
18263c7409f5SSatish Balay @*/
18273c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
18283c7409f5SSatish Balay {
18293c7409f5SSatish Balay   int ierr;
18303c7409f5SSatish Balay 
183177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1832639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18333c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
18343c7409f5SSatish Balay   return 0;
18353c7409f5SSatish Balay }
18363c7409f5SSatish Balay 
1837a847f771SSatish Balay #undef __FUNCTION__
1838a847f771SSatish Balay #define __FUNCTION__ "SNESAppendOptionsPrefix"
18393c7409f5SSatish Balay /*@C
1840f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1841639f9d9dSBarry Smith    SNES options in the database. You must NOT include the - at the beginning of
1842639f9d9dSBarry Smith    the prefix name.
18433c7409f5SSatish Balay 
18443c7409f5SSatish Balay    Input Parameter:
18453c7409f5SSatish Balay .  snes - the SNES context
18463c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18473c7409f5SSatish Balay 
18483c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1849a86d99e1SLois Curfman McInnes 
1850a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
18513c7409f5SSatish Balay @*/
18523c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
18533c7409f5SSatish Balay {
18543c7409f5SSatish Balay   int ierr;
18553c7409f5SSatish Balay 
185677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1857639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18583c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
18593c7409f5SSatish Balay   return 0;
18603c7409f5SSatish Balay }
18613c7409f5SSatish Balay 
1862a847f771SSatish Balay #undef __FUNCTION__
1863a847f771SSatish Balay #define __FUNCTION__ "SNESGetOptionsPrefix"
1864c83590e2SSatish Balay /*@
18653c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
18663c7409f5SSatish Balay    SNES options in the database.
18673c7409f5SSatish Balay 
18683c7409f5SSatish Balay    Input Parameter:
18693c7409f5SSatish Balay .  snes - the SNES context
18703c7409f5SSatish Balay 
18713c7409f5SSatish Balay    Output Parameter:
18723c7409f5SSatish Balay .  prefix - pointer to the prefix string used
18733c7409f5SSatish Balay 
18743c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1875a86d99e1SLois Curfman McInnes 
1876a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
18773c7409f5SSatish Balay @*/
18783c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
18793c7409f5SSatish Balay {
18803c7409f5SSatish Balay   int ierr;
18813c7409f5SSatish Balay 
188277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1883639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18843c7409f5SSatish Balay   return 0;
18853c7409f5SSatish Balay }
18863c7409f5SSatish Balay 
18873c7409f5SSatish Balay 
18889b94acceSBarry Smith 
18899b94acceSBarry Smith 
18909b94acceSBarry Smith 
1891