xref: /petsc/src/snes/interface/snes.c (revision f09e8eb94a771781a812a8d81a9ca3d36ec35eba)
19b94acceSBarry Smith #ifndef lint
2*f09e8eb9SSatish Balay static char vcid[] = "$Id: snes.c,v 1.121 1997/04/12 19:22:21 balay Exp balay $";
39b94acceSBarry Smith #endif
49b94acceSBarry Smith 
570f55243SBarry Smith #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
6f5eb4b81SSatish Balay #include "src/sys/nreg.h"
79b94acceSBarry Smith #include "pinclude/pviewer.h"
89b94acceSBarry Smith #include <math.h>
99b94acceSBarry Smith 
10052efed2SBarry Smith extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*);
1133455573SSatish Balay extern int SNESPrintTypes_Private(MPI_Comm,char*,char*);
129b94acceSBarry Smith 
1384cb2905SBarry Smith int SNESRegisterAllCalled = 0;
1484cb2905SBarry Smith 
155615d1e5SSatish Balay #undef __FUNC__
165eea60f9SBarry Smith #define __FUNC__ "SNESView" /* ADIC Ignore */
179b94acceSBarry Smith /*@
189b94acceSBarry Smith    SNESView - Prints the SNES data structure.
199b94acceSBarry Smith 
209b94acceSBarry Smith    Input Parameters:
219b94acceSBarry Smith .  SNES - the SNES context
229b94acceSBarry Smith .  viewer - visualization context
239b94acceSBarry Smith 
249b94acceSBarry Smith    Options Database Key:
259b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
269b94acceSBarry Smith 
279b94acceSBarry Smith    Notes:
289b94acceSBarry Smith    The available visualization contexts include
296d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
306d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
319b94acceSBarry Smith $       output where only the first processor opens
329b94acceSBarry Smith $       the file.  All other processors send their
339b94acceSBarry Smith $       data to the first processor to print.
349b94acceSBarry Smith 
359b94acceSBarry Smith    The user can open alternative vistualization contexts with
36dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
379b94acceSBarry Smith 
389b94acceSBarry Smith .keywords: SNES, view
399b94acceSBarry Smith 
40dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
419b94acceSBarry Smith @*/
429b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
439b94acceSBarry Smith {
449b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
459b94acceSBarry Smith   FILE                *fd;
469b94acceSBarry Smith   int                 ierr;
479b94acceSBarry Smith   SLES                sles;
489b94acceSBarry Smith   char                *method;
4919bcc07fSBarry Smith   ViewerType          vtype;
509b94acceSBarry Smith 
5174679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5274679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5374679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5474679c65SBarry Smith 
5519bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5619bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5790ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5877c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
594b0e389bSBarry Smith     SNESGetType(snes,PETSC_NULL,&method);
6077c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
619b94acceSBarry Smith     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
6277c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
639b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
649b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6577c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
669b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
679b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
687a00f4a9SLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
697a00f4a9SLois Curfman McInnes     "  total number of linear solver iterations=%d\n",snes->linear_its);
70ae3c334cSLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
7168632166SLois Curfman McInnes      "  total number of function evaluations=%d\n",snes->nfuncs);
729b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
7377c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
749b94acceSBarry Smith     if (snes->ksp_ewconv) {
759b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
769b94acceSBarry Smith       if (kctx) {
7777c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
789b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
799b94acceSBarry Smith         kctx->version);
8077c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
819b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
829b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
8377c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
849b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
859b94acceSBarry Smith       }
869b94acceSBarry Smith     }
87c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
88c30f285eSLois Curfman McInnes     SNESGetType(snes,PETSC_NULL,&method);
89c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
9019bcc07fSBarry Smith   }
919b94acceSBarry Smith   SNESGetSLES(snes,&sles);
929b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
939b94acceSBarry Smith   return 0;
949b94acceSBarry Smith }
959b94acceSBarry Smith 
96639f9d9dSBarry Smith /*
97639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
98639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
99639f9d9dSBarry Smith */
100639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
101639f9d9dSBarry Smith static int numberofsetfromoptions;
102639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
103639f9d9dSBarry Smith 
1045615d1e5SSatish Balay #undef __FUNC__
1055eea60f9SBarry Smith #define __FUNC__ "SNESAddOptionsChecker" /* ADIC Ignore */
106639f9d9dSBarry Smith /*@
107639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
108639f9d9dSBarry Smith 
109639f9d9dSBarry Smith   Input Parameter:
110639f9d9dSBarry Smith .   snescheck - function that checks for options
111639f9d9dSBarry Smith 
112639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
113639f9d9dSBarry Smith @*/
114639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
115639f9d9dSBarry Smith {
116639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
117e3372554SBarry Smith     SETERRQ(1,0,"Too many options checkers, only 5 allowed");
118639f9d9dSBarry Smith   }
119639f9d9dSBarry Smith 
120639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
121639f9d9dSBarry Smith   return 0;
122639f9d9dSBarry Smith }
123639f9d9dSBarry Smith 
1245615d1e5SSatish Balay #undef __FUNC__
1255615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1269b94acceSBarry Smith /*@
127682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1289b94acceSBarry Smith 
1299b94acceSBarry Smith    Input Parameter:
1309b94acceSBarry Smith .  snes - the SNES context
1319b94acceSBarry Smith 
1329b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1339b94acceSBarry Smith 
134a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1359b94acceSBarry Smith @*/
1369b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1379b94acceSBarry Smith {
1384b0e389bSBarry Smith   SNESType method;
1399b94acceSBarry Smith   double   tmp;
1409b94acceSBarry Smith   SLES     sles;
14181f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1423c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1439b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1449b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1459b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1469b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1479b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1489b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1499b94acceSBarry Smith 
15081f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
15181f57ec7SBarry Smith 
15281f57ec7SBarry Smith 
15377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
154e3372554SBarry Smith   if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp");
155052efed2SBarry Smith   ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr);
156052efed2SBarry Smith   if (flg) {
157052efed2SBarry Smith     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
1589b94acceSBarry Smith   }
1595334005bSBarry Smith   else if (!snes->set_method_called) {
1605334005bSBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
16140191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
1625334005bSBarry Smith     }
1635334005bSBarry Smith     else {
16440191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
1655334005bSBarry Smith     }
1665334005bSBarry Smith   }
1675334005bSBarry Smith 
1683c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1693c7409f5SSatish Balay   if (flg) { SNESPrintHelp(snes); }
1703c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
171d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);}
1723c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
173d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);}
1743c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
175d7a720efSLois Curfman McInnes   if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); }
1763c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1773c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
178d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
179d7a720efSLois Curfman McInnes   if (flg) { SNESSetTrustRegionTolerance(snes,tmp); }
180d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
181d7a720efSLois Curfman McInnes   if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); }
1823c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1833c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
184b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
185b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
186b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
187b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
188b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
189b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
190b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
1919b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
1929b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
1935bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
1945bbad29bSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);}
1953c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
196639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
1973c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
198d31a3109SLois Curfman McInnes   if (flg) {
199639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
200d31a3109SLois Curfman McInnes   }
20181f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
2023c7409f5SSatish Balay   if (flg) {
20317699dbbSLois Curfman McInnes     int    rank = 0;
204d7e8b826SBarry Smith     DrawLG lg;
20517699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
20617699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
20717699dbbSLois Curfman McInnes     if (!rank) {
20881f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2099b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
210f8c826e1SBarry Smith       snes->xmonitor = lg;
2119b94acceSBarry Smith       PLogObjectParent(snes,lg);
2129b94acceSBarry Smith     }
2139b94acceSBarry Smith   }
2143c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2153c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2169b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2179b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
21894a424c1SBarry Smith     PLogInfo(snes,
21931615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
22031615d04SLois Curfman McInnes   }
22131615d04SLois Curfman McInnes   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
22231615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
22331615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
22494a424c1SBarry Smith     PLogInfo(snes,
22531615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2269b94acceSBarry Smith   }
227639f9d9dSBarry Smith 
228639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
229639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
230639f9d9dSBarry Smith   }
231639f9d9dSBarry Smith 
2329b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2339b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
2349b94acceSBarry Smith   if (!snes->setfromoptions) return 0;
2359b94acceSBarry Smith   return (*snes->setfromoptions)(snes);
2369b94acceSBarry Smith }
2379b94acceSBarry Smith 
2385615d1e5SSatish Balay #undef __FUNC__
2395eea60f9SBarry Smith #define __FUNC__ "SNESPrintHelp" /* ADIC Ignore */
2409b94acceSBarry Smith /*@
2419b94acceSBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2429b94acceSBarry Smith 
2439b94acceSBarry Smith    Input Parameter:
2449b94acceSBarry Smith .  snes - the SNES context
2459b94acceSBarry Smith 
246a703fe33SLois Curfman McInnes    Options Database Keys:
247a703fe33SLois Curfman McInnes $  -help, -h
248a703fe33SLois Curfman McInnes 
2499b94acceSBarry Smith .keywords: SNES, nonlinear, help
2509b94acceSBarry Smith 
251682d7d0cSBarry Smith .seealso: SNESSetFromOptions()
2529b94acceSBarry Smith @*/
2539b94acceSBarry Smith int SNESPrintHelp(SNES snes)
2549b94acceSBarry Smith {
2556daaf66cSBarry Smith   char                p[64];
2566daaf66cSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2576daaf66cSBarry Smith 
25877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2596daaf66cSBarry Smith 
2606daaf66cSBarry Smith   PetscStrcpy(p,"-");
2616daaf66cSBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2626daaf66cSBarry Smith 
2636daaf66cSBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2646daaf66cSBarry Smith 
265d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
26633455573SSatish Balay   SNESPrintTypes_Private(snes->comm,p,"snes_type");
26777c4ece6SBarry Smith   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
268b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
269b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
270b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
271b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
272b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
273b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
2745bbad29bSBarry Smith   PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n");
2755bbad29bSBarry Smith   PetscPrintf(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
276d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor\n",p);
277d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
278b18e04deSLois Curfman McInnes   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
27977c4ece6SBarry Smith     PetscPrintf(snes->comm,
280d31a3109SLois Curfman McInnes      " Options for solving systems of nonlinear equations only:\n");
28177c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
28277c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
283ef1dfb11SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
28477c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
28577c4ece6SBarry Smith     PetscPrintf(snes->comm,
286b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
28777c4ece6SBarry Smith     PetscPrintf(snes->comm,
288b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
28977c4ece6SBarry Smith     PetscPrintf(snes->comm,
290b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
29177c4ece6SBarry Smith     PetscPrintf(snes->comm,
292b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
29377c4ece6SBarry Smith     PetscPrintf(snes->comm,
294b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
29577c4ece6SBarry Smith     PetscPrintf(snes->comm,
296b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
29777c4ece6SBarry Smith     PetscPrintf(snes->comm,
298b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
299b18e04deSLois Curfman McInnes   }
300b18e04deSLois Curfman McInnes   else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
30177c4ece6SBarry Smith     PetscPrintf(snes->comm,
302d31a3109SLois Curfman McInnes      " Options for solving unconstrained minimization problems only:\n");
303b18e04deSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
30477c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
30577c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
306b18e04deSLois Curfman McInnes   }
3074537a946SLois Curfman McInnes   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
30877c4ece6SBarry Smith   PetscPrintf(snes->comm,"a particular method\n");
3096daaf66cSBarry Smith   if (snes->printhelp) (*snes->printhelp)(snes,p);
3109b94acceSBarry Smith   return 0;
3119b94acceSBarry Smith }
312a847f771SSatish Balay 
3135615d1e5SSatish Balay #undef __FUNC__
3145eea60f9SBarry Smith #define __FUNC__ "SNESSetApplicationContext" /* ADIC Ignore */
3159b94acceSBarry Smith /*@
3169b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3179b94acceSBarry Smith    the nonlinear solvers.
3189b94acceSBarry Smith 
3199b94acceSBarry Smith    Input Parameters:
3209b94acceSBarry Smith .  snes - the SNES context
3219b94acceSBarry Smith .  usrP - optional user context
3229b94acceSBarry Smith 
3239b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3249b94acceSBarry Smith 
3259b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3269b94acceSBarry Smith @*/
3279b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3289b94acceSBarry Smith {
32977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3309b94acceSBarry Smith   snes->user		= usrP;
3319b94acceSBarry Smith   return 0;
3329b94acceSBarry Smith }
33374679c65SBarry Smith 
3345615d1e5SSatish Balay #undef __FUNC__
3355eea60f9SBarry Smith #define __FUNC__ "SNESGetApplicationContext" /* ADIC Ignore */
3369b94acceSBarry Smith /*@C
3379b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3389b94acceSBarry Smith    nonlinear solvers.
3399b94acceSBarry Smith 
3409b94acceSBarry Smith    Input Parameter:
3419b94acceSBarry Smith .  snes - SNES context
3429b94acceSBarry Smith 
3439b94acceSBarry Smith    Output Parameter:
3449b94acceSBarry Smith .  usrP - user context
3459b94acceSBarry Smith 
3469b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3479b94acceSBarry Smith 
3489b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3499b94acceSBarry Smith @*/
3509b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3519b94acceSBarry Smith {
35277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3539b94acceSBarry Smith   *usrP = snes->user;
3549b94acceSBarry Smith   return 0;
3559b94acceSBarry Smith }
35674679c65SBarry Smith 
3575615d1e5SSatish Balay #undef __FUNC__
3585eea60f9SBarry Smith #define __FUNC__ "SNESGetIterationNumber" /* ADIC Ignore */
3599b94acceSBarry Smith /*@
3609b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3619b94acceSBarry Smith    nonlinear solver.
3629b94acceSBarry Smith 
3639b94acceSBarry Smith    Input Parameter:
3649b94acceSBarry Smith .  snes - SNES context
3659b94acceSBarry Smith 
3669b94acceSBarry Smith    Output Parameter:
3679b94acceSBarry Smith .  iter - iteration number
3689b94acceSBarry Smith 
3699b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3709b94acceSBarry Smith @*/
3719b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3729b94acceSBarry Smith {
37377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
37474679c65SBarry Smith   PetscValidIntPointer(iter);
3759b94acceSBarry Smith   *iter = snes->iter;
3769b94acceSBarry Smith   return 0;
3779b94acceSBarry Smith }
37874679c65SBarry Smith 
3795615d1e5SSatish Balay #undef __FUNC__
3805615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
3819b94acceSBarry Smith /*@
3829b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3839b94acceSBarry Smith    with SNESSSetFunction().
3849b94acceSBarry Smith 
3859b94acceSBarry Smith    Input Parameter:
3869b94acceSBarry Smith .  snes - SNES context
3879b94acceSBarry Smith 
3889b94acceSBarry Smith    Output Parameter:
3899b94acceSBarry Smith .  fnorm - 2-norm of function
3909b94acceSBarry Smith 
3919b94acceSBarry Smith    Note:
3929b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
393a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
394a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3959b94acceSBarry Smith 
3969b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
397a86d99e1SLois Curfman McInnes 
398a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3999b94acceSBarry Smith @*/
4009b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
4019b94acceSBarry Smith {
40277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
40374679c65SBarry Smith   PetscValidScalarPointer(fnorm);
40474679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
405d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only");
40674679c65SBarry Smith   }
4079b94acceSBarry Smith   *fnorm = snes->norm;
4089b94acceSBarry Smith   return 0;
4099b94acceSBarry Smith }
41074679c65SBarry Smith 
4115615d1e5SSatish Balay #undef __FUNC__
4125615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
4139b94acceSBarry Smith /*@
4149b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
4159b94acceSBarry Smith    with SNESSSetGradient().
4169b94acceSBarry Smith 
4179b94acceSBarry Smith    Input Parameter:
4189b94acceSBarry Smith .  snes - SNES context
4199b94acceSBarry Smith 
4209b94acceSBarry Smith    Output Parameter:
4219b94acceSBarry Smith .  fnorm - 2-norm of gradient
4229b94acceSBarry Smith 
4239b94acceSBarry Smith    Note:
4249b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
425a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
426a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4279b94acceSBarry Smith 
4289b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
429a86d99e1SLois Curfman McInnes 
430a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4319b94acceSBarry Smith @*/
4329b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
4339b94acceSBarry Smith {
43477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43574679c65SBarry Smith   PetscValidScalarPointer(gnorm);
43674679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
437d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
43874679c65SBarry Smith   }
4399b94acceSBarry Smith   *gnorm = snes->norm;
4409b94acceSBarry Smith   return 0;
4419b94acceSBarry Smith }
44274679c65SBarry Smith 
4435615d1e5SSatish Balay #undef __FUNC__
4445eea60f9SBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" /* ADIC Ignore */
4459b94acceSBarry Smith /*@
4469b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4479b94acceSBarry Smith    attempted by the nonlinear solver.
4489b94acceSBarry Smith 
4499b94acceSBarry Smith    Input Parameter:
4509b94acceSBarry Smith .  snes - SNES context
4519b94acceSBarry Smith 
4529b94acceSBarry Smith    Output Parameter:
4539b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4549b94acceSBarry Smith 
455c96a6f78SLois Curfman McInnes    Notes:
456c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
457c96a6f78SLois Curfman McInnes 
4589b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4599b94acceSBarry Smith @*/
4609b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4619b94acceSBarry Smith {
46277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
46374679c65SBarry Smith   PetscValidIntPointer(nfails);
4649b94acceSBarry Smith   *nfails = snes->nfailures;
4659b94acceSBarry Smith   return 0;
4669b94acceSBarry Smith }
467a847f771SSatish Balay 
4685615d1e5SSatish Balay #undef __FUNC__
4695eea60f9SBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" /* ADIC Ignore */
470c96a6f78SLois Curfman McInnes /*@
471c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
472c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
473c96a6f78SLois Curfman McInnes 
474c96a6f78SLois Curfman McInnes    Input Parameter:
475c96a6f78SLois Curfman McInnes .  snes - SNES context
476c96a6f78SLois Curfman McInnes 
477c96a6f78SLois Curfman McInnes    Output Parameter:
478c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
479c96a6f78SLois Curfman McInnes 
480c96a6f78SLois Curfman McInnes    Notes:
481c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
482c96a6f78SLois Curfman McInnes 
483c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
484c96a6f78SLois Curfman McInnes @*/
48586bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
486c96a6f78SLois Curfman McInnes {
487c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
488c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
489c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
490c96a6f78SLois Curfman McInnes   return 0;
491c96a6f78SLois Curfman McInnes }
492c96a6f78SLois Curfman McInnes 
493c96a6f78SLois Curfman McInnes #undef __FUNC__
4945eea60f9SBarry Smith #define __FUNC__ "SNESGetSLES" /* ADIC Ignore */
4959b94acceSBarry Smith /*@C
4969b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4979b94acceSBarry Smith 
4989b94acceSBarry Smith    Input Parameter:
4999b94acceSBarry Smith .  snes - the SNES context
5009b94acceSBarry Smith 
5019b94acceSBarry Smith    Output Parameter:
5029b94acceSBarry Smith .  sles - the SLES context
5039b94acceSBarry Smith 
5049b94acceSBarry Smith    Notes:
5059b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
5069b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5079b94acceSBarry Smith    KSP and PC contexts as well.
5089b94acceSBarry Smith 
5099b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
5109b94acceSBarry Smith 
5119b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
5129b94acceSBarry Smith @*/
5139b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
5149b94acceSBarry Smith {
51577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5169b94acceSBarry Smith   *sles = snes->sles;
5179b94acceSBarry Smith   return 0;
5189b94acceSBarry Smith }
5199b94acceSBarry Smith /* -----------------------------------------------------------*/
5205615d1e5SSatish Balay #undef __FUNC__
5215615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
5229b94acceSBarry Smith /*@C
5239b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5249b94acceSBarry Smith 
5259b94acceSBarry Smith    Input Parameter:
5269b94acceSBarry Smith .  comm - MPI communicator
5279b94acceSBarry Smith .  type - type of method, one of
5289b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
5299b94acceSBarry Smith $      (for systems of nonlinear equations)
5309b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
5319b94acceSBarry Smith $      (for unconstrained minimization)
5329b94acceSBarry Smith 
5339b94acceSBarry Smith    Output Parameter:
5349b94acceSBarry Smith .  outsnes - the new SNES context
5359b94acceSBarry Smith 
536c1f60f51SBarry Smith    Options Database Key:
53719bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
53819bd6747SLois Curfman McInnes $              and no preconditioning matrix
53919bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
54019bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
54119bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
54219bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
543c1f60f51SBarry Smith 
5449b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5459b94acceSBarry Smith 
54663a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5479b94acceSBarry Smith @*/
5484b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5499b94acceSBarry Smith {
5509b94acceSBarry Smith   int                 ierr;
5519b94acceSBarry Smith   SNES                snes;
5529b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
55337fcc0dbSBarry Smith 
5549b94acceSBarry Smith   *outsnes = 0;
5552a0acf01SLois Curfman McInnes   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS)
556d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
557*f09e8eb9SSatish Balay   PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm);
5589b94acceSBarry Smith   PLogObjectCreate(snes);
5599b94acceSBarry Smith   snes->max_its           = 50;
5609750a799SBarry Smith   snes->max_funcs	  = 10000;
5619b94acceSBarry Smith   snes->norm		  = 0.0;
5625a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
563b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
564b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5659b94acceSBarry Smith     snes->atol		  = 1.e-10;
5665a2d0531SBarry Smith   }
567b4874afaSBarry Smith   else {
568b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
569b4874afaSBarry Smith     snes->ttol            = 0.0;
570b4874afaSBarry Smith     snes->atol		  = 1.e-50;
571b4874afaSBarry Smith   }
5729b94acceSBarry Smith   snes->xtol		  = 1.e-8;
573b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5749b94acceSBarry Smith   snes->nfuncs            = 0;
5759b94acceSBarry Smith   snes->nfailures         = 0;
5767a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
577639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5789b94acceSBarry Smith   snes->data              = 0;
5799b94acceSBarry Smith   snes->view              = 0;
5809b94acceSBarry Smith   snes->computeumfunction = 0;
5819b94acceSBarry Smith   snes->umfunP            = 0;
5829b94acceSBarry Smith   snes->fc                = 0;
5839b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5849b94acceSBarry Smith   snes->fmin              = -1.e30;
5859b94acceSBarry Smith   snes->method_class      = type;
5869b94acceSBarry Smith   snes->set_method_called = 0;
587a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
5889b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5897d1a2b2bSBarry Smith   snes->type              = -1;
5906f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5916f24a144SLois Curfman McInnes   snes->vwork             = 0;
5926f24a144SLois Curfman McInnes   snes->nwork             = 0;
593c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
594c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
595c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
5969b94acceSBarry Smith 
5979b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5980452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
5999b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6009b94acceSBarry Smith   kctx->version     = 2;
6019b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6029b94acceSBarry Smith                              this was too large for some test cases */
6039b94acceSBarry Smith   kctx->rtol_last   = 0;
6049b94acceSBarry Smith   kctx->rtol_max    = .9;
6059b94acceSBarry Smith   kctx->gamma       = 1.0;
6069b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6079b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6089b94acceSBarry Smith   kctx->threshold   = .1;
6099b94acceSBarry Smith   kctx->lresid_last = 0;
6109b94acceSBarry Smith   kctx->norm_last   = 0;
6119b94acceSBarry Smith 
6129b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
6139b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
6145334005bSBarry Smith 
6159b94acceSBarry Smith   *outsnes = snes;
6169b94acceSBarry Smith   return 0;
6179b94acceSBarry Smith }
6189b94acceSBarry Smith 
6199b94acceSBarry Smith /* --------------------------------------------------------------- */
6205615d1e5SSatish Balay #undef __FUNC__
6215eea60f9SBarry Smith #define __FUNC__ "SNESSetFunction" /* ADIC Ignore */
6229b94acceSBarry Smith /*@C
6239b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6249b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6259b94acceSBarry Smith    equations.
6269b94acceSBarry Smith 
6279b94acceSBarry Smith    Input Parameters:
6289b94acceSBarry Smith .  snes - the SNES context
6299b94acceSBarry Smith .  func - function evaluation routine
6309b94acceSBarry Smith .  r - vector to store function value
6312cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
6322cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6339b94acceSBarry Smith 
6349b94acceSBarry Smith    Calling sequence of func:
635313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
6369b94acceSBarry Smith 
6379b94acceSBarry Smith .  x - input vector
638313e4042SLois Curfman McInnes .  f - function vector
6392cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
6409b94acceSBarry Smith 
6419b94acceSBarry Smith    Notes:
6429b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6439b94acceSBarry Smith $      f'(x) x = -f(x),
6449b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
6459b94acceSBarry Smith 
6469b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6479b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6489b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6499b94acceSBarry Smith 
6509b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6519b94acceSBarry Smith 
652a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6539b94acceSBarry Smith @*/
6545334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6559b94acceSBarry Smith {
65677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
657e3372554SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
6589b94acceSBarry Smith   snes->computefunction     = func;
6599b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6609b94acceSBarry Smith   snes->funP                = ctx;
6619b94acceSBarry Smith   return 0;
6629b94acceSBarry Smith }
6639b94acceSBarry Smith 
6645615d1e5SSatish Balay #undef __FUNC__
6655615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
6669b94acceSBarry Smith /*@
6679b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6689b94acceSBarry Smith    SNESSetFunction().
6699b94acceSBarry Smith 
6709b94acceSBarry Smith    Input Parameters:
6719b94acceSBarry Smith .  snes - the SNES context
6729b94acceSBarry Smith .  x - input vector
6739b94acceSBarry Smith 
6749b94acceSBarry Smith    Output Parameter:
6753638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6769b94acceSBarry Smith 
6771bffabb2SLois Curfman McInnes    Notes:
6789b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6799b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6809b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6819b94acceSBarry Smith 
6829b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6839b94acceSBarry Smith 
684a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6859b94acceSBarry Smith @*/
6869b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6879b94acceSBarry Smith {
6889b94acceSBarry Smith   int    ierr;
6899b94acceSBarry Smith 
69074679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
691e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
6929b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
6939b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
694ae3c334cSLois Curfman McInnes   snes->nfuncs++;
6959b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6969b94acceSBarry Smith   return 0;
6979b94acceSBarry Smith }
6989b94acceSBarry Smith 
6995615d1e5SSatish Balay #undef __FUNC__
7005eea60f9SBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" /* ADIC Ignore */
7019b94acceSBarry Smith /*@C
7029b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
7039b94acceSBarry Smith    unconstrained minimization.
7049b94acceSBarry Smith 
7059b94acceSBarry Smith    Input Parameters:
7069b94acceSBarry Smith .  snes - the SNES context
7079b94acceSBarry Smith .  func - function evaluation routine
708044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
709044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7109b94acceSBarry Smith 
7119b94acceSBarry Smith    Calling sequence of func:
7129b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
7139b94acceSBarry Smith 
7149b94acceSBarry Smith .  x - input vector
7159b94acceSBarry Smith .  f - function
716044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
7179b94acceSBarry Smith 
7189b94acceSBarry Smith    Notes:
7199b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7209b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7219b94acceSBarry Smith    SNESSetFunction().
7229b94acceSBarry Smith 
7239b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
7249b94acceSBarry Smith 
725a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
726a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
7279b94acceSBarry Smith @*/
7289b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
7299b94acceSBarry Smith                       void *ctx)
7309b94acceSBarry Smith {
73177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
732e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
7339b94acceSBarry Smith   snes->computeumfunction   = func;
7349b94acceSBarry Smith   snes->umfunP              = ctx;
7359b94acceSBarry Smith   return 0;
7369b94acceSBarry Smith }
7379b94acceSBarry Smith 
7385615d1e5SSatish Balay #undef __FUNC__
7395615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7409b94acceSBarry Smith /*@
7419b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7429b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7439b94acceSBarry Smith 
7449b94acceSBarry Smith    Input Parameters:
7459b94acceSBarry Smith .  snes - the SNES context
7469b94acceSBarry Smith .  x - input vector
7479b94acceSBarry Smith 
7489b94acceSBarry Smith    Output Parameter:
7499b94acceSBarry Smith .  y - function value
7509b94acceSBarry Smith 
7519b94acceSBarry Smith    Notes:
7529b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7539b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7549b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
755a86d99e1SLois Curfman McInnes 
756a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
757a86d99e1SLois Curfman McInnes 
758a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
759a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7609b94acceSBarry Smith @*/
7619b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7629b94acceSBarry Smith {
7639b94acceSBarry Smith   int    ierr;
764e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
7659b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
7669b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
767ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7689b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7699b94acceSBarry Smith   return 0;
7709b94acceSBarry Smith }
7719b94acceSBarry Smith 
7725615d1e5SSatish Balay #undef __FUNC__
7735eea60f9SBarry Smith #define __FUNC__ "SNESSetGradient" /* ADIC Ignore */
7749b94acceSBarry Smith /*@C
7759b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7769b94acceSBarry Smith    vector for use by the SNES routines.
7779b94acceSBarry Smith 
7789b94acceSBarry Smith    Input Parameters:
7799b94acceSBarry Smith .  snes - the SNES context
7809b94acceSBarry Smith .  func - function evaluation routine
781044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
782044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7839b94acceSBarry Smith .  r - vector to store gradient value
7849b94acceSBarry Smith 
7859b94acceSBarry Smith    Calling sequence of func:
7869b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7879b94acceSBarry Smith 
7889b94acceSBarry Smith .  x - input vector
7899b94acceSBarry Smith .  g - gradient vector
790044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7919b94acceSBarry Smith 
7929b94acceSBarry Smith    Notes:
7939b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7949b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7959b94acceSBarry Smith    SNESSetFunction().
7969b94acceSBarry Smith 
7979b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7989b94acceSBarry Smith 
799a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
800a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
8019b94acceSBarry Smith @*/
80274679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
8039b94acceSBarry Smith {
80477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
805e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8069b94acceSBarry Smith   snes->computefunction     = func;
8079b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
8089b94acceSBarry Smith   snes->funP                = ctx;
8099b94acceSBarry Smith   return 0;
8109b94acceSBarry Smith }
8119b94acceSBarry Smith 
8125615d1e5SSatish Balay #undef __FUNC__
8135615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
8149b94acceSBarry Smith /*@
815a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
816a86d99e1SLois Curfman McInnes    SNESSetGradient().
8179b94acceSBarry Smith 
8189b94acceSBarry Smith    Input Parameters:
8199b94acceSBarry Smith .  snes - the SNES context
8209b94acceSBarry Smith .  x - input vector
8219b94acceSBarry Smith 
8229b94acceSBarry Smith    Output Parameter:
8239b94acceSBarry Smith .  y - gradient vector
8249b94acceSBarry Smith 
8259b94acceSBarry Smith    Notes:
8269b94acceSBarry Smith    SNESComputeGradient() is valid only for
8279b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8289b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
829a86d99e1SLois Curfman McInnes 
830a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
831a86d99e1SLois Curfman McInnes 
832a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
833a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8349b94acceSBarry Smith @*/
8359b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8369b94acceSBarry Smith {
8379b94acceSBarry Smith   int    ierr;
83874679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
839e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8409b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
8419b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
8429b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8439b94acceSBarry Smith   return 0;
8449b94acceSBarry Smith }
8459b94acceSBarry Smith 
8465615d1e5SSatish Balay #undef __FUNC__
8475615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
84862fef451SLois Curfman McInnes /*@
84962fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
85062fef451SLois Curfman McInnes    set with SNESSetJacobian().
85162fef451SLois Curfman McInnes 
85262fef451SLois Curfman McInnes    Input Parameters:
85362fef451SLois Curfman McInnes .  snes - the SNES context
85462fef451SLois Curfman McInnes .  x - input vector
85562fef451SLois Curfman McInnes 
85662fef451SLois Curfman McInnes    Output Parameters:
85762fef451SLois Curfman McInnes .  A - Jacobian matrix
85862fef451SLois Curfman McInnes .  B - optional preconditioning matrix
85962fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
86062fef451SLois Curfman McInnes 
86162fef451SLois Curfman McInnes    Notes:
86262fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
86362fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
86462fef451SLois Curfman McInnes 
865dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
866dc5a77f8SLois Curfman McInnes    flag parameter.
86762fef451SLois Curfman McInnes 
86862fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
86962fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
87062fef451SLois Curfman McInnes    methods is SNESComputeJacobian().
87162fef451SLois Curfman McInnes 
87262fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
87362fef451SLois Curfman McInnes 
87462fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
87562fef451SLois Curfman McInnes @*/
8769b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8779b94acceSBarry Smith {
8789b94acceSBarry Smith   int    ierr;
87974679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
880e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
8819b94acceSBarry Smith   if (!snes->computejacobian) return 0;
8829b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
883c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
8849b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
8859b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8866d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
88777c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
88877c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8899b94acceSBarry Smith   return 0;
8909b94acceSBarry Smith }
8919b94acceSBarry Smith 
8925615d1e5SSatish Balay #undef __FUNC__
8935615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
89462fef451SLois Curfman McInnes /*@
89562fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
89662fef451SLois Curfman McInnes    set with SNESSetHessian().
89762fef451SLois Curfman McInnes 
89862fef451SLois Curfman McInnes    Input Parameters:
89962fef451SLois Curfman McInnes .  snes - the SNES context
90062fef451SLois Curfman McInnes .  x - input vector
90162fef451SLois Curfman McInnes 
90262fef451SLois Curfman McInnes    Output Parameters:
90362fef451SLois Curfman McInnes .  A - Hessian matrix
90462fef451SLois Curfman McInnes .  B - optional preconditioning matrix
90562fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
90662fef451SLois Curfman McInnes 
90762fef451SLois Curfman McInnes    Notes:
90862fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
90962fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
91062fef451SLois Curfman McInnes 
911dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
912dc5a77f8SLois Curfman McInnes    flag parameter.
91362fef451SLois Curfman McInnes 
91462fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
91562fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
91662fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
91762fef451SLois Curfman McInnes 
91862fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
91962fef451SLois Curfman McInnes 
920a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
921a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
92262fef451SLois Curfman McInnes @*/
92362fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
9249b94acceSBarry Smith {
9259b94acceSBarry Smith   int    ierr;
92674679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
927e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9289b94acceSBarry Smith   if (!snes->computejacobian) return 0;
92962fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9300f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
93162fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
93262fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9330f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
93477c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
93577c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9369b94acceSBarry Smith   return 0;
9379b94acceSBarry Smith }
9389b94acceSBarry Smith 
9395615d1e5SSatish Balay #undef __FUNC__
9405eea60f9SBarry Smith #define __FUNC__ "SNESSetJacobian" /* ADIC Ignore */
9419b94acceSBarry Smith /*@C
9429b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
943044dda88SLois Curfman McInnes    location to store the matrix.
9449b94acceSBarry Smith 
9459b94acceSBarry Smith    Input Parameters:
9469b94acceSBarry Smith .  snes - the SNES context
9479b94acceSBarry Smith .  A - Jacobian matrix
9489b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9499b94acceSBarry Smith .  func - Jacobian evaluation routine
9502cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9512cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9529b94acceSBarry Smith 
9539b94acceSBarry Smith    Calling sequence of func:
954313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9559b94acceSBarry Smith 
9569b94acceSBarry Smith .  x - input vector
9579b94acceSBarry Smith .  A - Jacobian matrix
9589b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
959ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
960ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9612cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
9629b94acceSBarry Smith 
9639b94acceSBarry Smith    Notes:
964dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9652cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
966ac21db08SLois Curfman McInnes 
967ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9689b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9699b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9709b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9719b94acceSBarry Smith    throughout the global iterations.
9729b94acceSBarry Smith 
9739b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9749b94acceSBarry Smith 
975ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
9769b94acceSBarry Smith @*/
9779b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9789b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9799b94acceSBarry Smith {
98077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
98174679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
982e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
9839b94acceSBarry Smith   snes->computejacobian = func;
9849b94acceSBarry Smith   snes->jacP            = ctx;
9859b94acceSBarry Smith   snes->jacobian        = A;
9869b94acceSBarry Smith   snes->jacobian_pre    = B;
9879b94acceSBarry Smith   return 0;
9889b94acceSBarry Smith }
98962fef451SLois Curfman McInnes 
9905615d1e5SSatish Balay #undef __FUNC__
9915eea60f9SBarry Smith #define __FUNC__ "SNESGetJacobian" /* ADIC Ignore */
992b4fd4287SBarry Smith /*@
993b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
994b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
995b4fd4287SBarry Smith 
996b4fd4287SBarry Smith    Input Parameter:
997b4fd4287SBarry Smith .  snes - the nonlinear solver context
998b4fd4287SBarry Smith 
999b4fd4287SBarry Smith    Output Parameters:
1000b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
1001b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
1002b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1003b4fd4287SBarry Smith 
1004b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1005b4fd4287SBarry Smith @*/
1006b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
1007b4fd4287SBarry Smith {
100874679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS)
1009e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
1010b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
1011b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
1012b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
1013b4fd4287SBarry Smith   return 0;
1014b4fd4287SBarry Smith }
1015b4fd4287SBarry Smith 
10165615d1e5SSatish Balay #undef __FUNC__
10175eea60f9SBarry Smith #define __FUNC__ "SNESSetHessian" /* ADIC Ignore */
10189b94acceSBarry Smith /*@C
10199b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1020044dda88SLois Curfman McInnes    location to store the matrix.
10219b94acceSBarry Smith 
10229b94acceSBarry Smith    Input Parameters:
10239b94acceSBarry Smith .  snes - the SNES context
10249b94acceSBarry Smith .  A - Hessian matrix
10259b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
10269b94acceSBarry Smith .  func - Jacobian evaluation routine
1027313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1028313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10299b94acceSBarry Smith 
10309b94acceSBarry Smith    Calling sequence of func:
1031313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10329b94acceSBarry Smith 
10339b94acceSBarry Smith .  x - input vector
10349b94acceSBarry Smith .  A - Hessian matrix
10359b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1036ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1037ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10382cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10399b94acceSBarry Smith 
10409b94acceSBarry Smith    Notes:
1041dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10422cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1043ac21db08SLois Curfman McInnes 
10449b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10459b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10469b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10479b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10489b94acceSBarry Smith    throughout the global iterations.
10499b94acceSBarry Smith 
10509b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10519b94acceSBarry Smith 
1052ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10539b94acceSBarry Smith @*/
10549b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10559b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10569b94acceSBarry Smith {
105777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
105874679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
1059e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10609b94acceSBarry Smith   snes->computejacobian = func;
10619b94acceSBarry Smith   snes->jacP            = ctx;
10629b94acceSBarry Smith   snes->jacobian        = A;
10639b94acceSBarry Smith   snes->jacobian_pre    = B;
10649b94acceSBarry Smith   return 0;
10659b94acceSBarry Smith }
10669b94acceSBarry Smith 
10675615d1e5SSatish Balay #undef __FUNC__
10685eea60f9SBarry Smith #define __FUNC__ "SNESGetHessian" /* ADIC Ignore */
106962fef451SLois Curfman McInnes /*@
107062fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
107162fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
107262fef451SLois Curfman McInnes 
107362fef451SLois Curfman McInnes    Input Parameter:
107462fef451SLois Curfman McInnes .  snes - the nonlinear solver context
107562fef451SLois Curfman McInnes 
107662fef451SLois Curfman McInnes    Output Parameters:
107762fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
107862fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
107962fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
108062fef451SLois Curfman McInnes 
108162fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
108262fef451SLois Curfman McInnes @*/
108362fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
108462fef451SLois Curfman McInnes {
108574679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION)
1086e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
108762fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
108862fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
108962fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
109062fef451SLois Curfman McInnes   return 0;
109162fef451SLois Curfman McInnes }
109262fef451SLois Curfman McInnes 
10939b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10949b94acceSBarry Smith 
10955615d1e5SSatish Balay #undef __FUNC__
10965615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
10979b94acceSBarry Smith /*@
10989b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1099272ac6f2SLois Curfman McInnes    of a nonlinear solver.
11009b94acceSBarry Smith 
11019b94acceSBarry Smith    Input Parameter:
11029b94acceSBarry Smith .  snes - the SNES context
11038ddd3da0SLois Curfman McInnes .  x - the solution vector
11049b94acceSBarry Smith 
1105272ac6f2SLois Curfman McInnes    Notes:
1106272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1107272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1108272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1109272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1110272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1111272ac6f2SLois Curfman McInnes 
11129b94acceSBarry Smith .keywords: SNES, nonlinear, setup
11139b94acceSBarry Smith 
11149b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
11159b94acceSBarry Smith @*/
11168ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
11179b94acceSBarry Smith {
1118272ac6f2SLois Curfman McInnes   int ierr, flg;
111977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
112077c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
11218ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
11229b94acceSBarry Smith 
1123c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1124c1f60f51SBarry Smith   /*
1125c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1126dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1127c1f60f51SBarry Smith   */
1128c1f60f51SBarry Smith   if (flg) {
1129c1f60f51SBarry Smith     Mat J;
1130c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1131c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1132c1f60f51SBarry Smith     snes->mfshell = J;
1133c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1134c1f60f51SBarry Smith       snes->jacobian = J;
1135c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1136c1f60f51SBarry Smith     }
1137c1f60f51SBarry Smith     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1138c1f60f51SBarry Smith       snes->jacobian = J;
1139c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1140e3372554SBarry Smith     } else SETERRQ(1,0,"Method class doesn't support matrix-free operator option");
1141c1f60f51SBarry Smith   }
1142272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1143c1f60f51SBarry Smith   /*
1144dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1145c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1146c1f60f51SBarry Smith    */
114731615d04SLois Curfman McInnes   if (flg) {
1148272ac6f2SLois Curfman McInnes     Mat J;
1149272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1150272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1151272ac6f2SLois Curfman McInnes     snes->mfshell = J;
115283e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
115383e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
115494a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
115583e56afdSLois Curfman McInnes     }
115683e56afdSLois Curfman McInnes     else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
115783e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
115894a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1159e3372554SBarry Smith     } else SETERRQ(1,0,"Method class doesn't support matrix-free option");
1160272ac6f2SLois Curfman McInnes   }
11619b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1162e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first");
1163e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first");
1164e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first");
1165e3372554SBarry Smith     if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector");
1166a703fe33SLois Curfman McInnes 
1167a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
116840191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1169a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1170a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1171a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1172a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1173a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1174a703fe33SLois Curfman McInnes     }
11759b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1176e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first");
1177e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first");
117837fcc0dbSBarry Smith     if (!snes->computeumfunction)
1179e3372554SBarry Smith       SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first");
1180e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first");
1181e3372554SBarry Smith   } else SETERRQ(1,0,"Unknown method class");
1182a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1183a703fe33SLois Curfman McInnes   snes->setup_called = 1;
1184a703fe33SLois Curfman McInnes   return 0;
11859b94acceSBarry Smith }
11869b94acceSBarry Smith 
11875615d1e5SSatish Balay #undef __FUNC__
11885eea60f9SBarry Smith #define __FUNC__ "SNESDestroy" /* ADIC Ignore */
11899b94acceSBarry Smith /*@C
11909b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11919b94acceSBarry Smith    with SNESCreate().
11929b94acceSBarry Smith 
11939b94acceSBarry Smith    Input Parameter:
11949b94acceSBarry Smith .  snes - the SNES context
11959b94acceSBarry Smith 
11969b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
11979b94acceSBarry Smith 
119863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
11999b94acceSBarry Smith @*/
12009b94acceSBarry Smith int SNESDestroy(SNES snes)
12019b94acceSBarry Smith {
12029b94acceSBarry Smith   int ierr;
120377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12049750a799SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);}
12050452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
12069b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
12079b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
12083e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
12096f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
12109b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
12110452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
12129b94acceSBarry Smith   return 0;
12139b94acceSBarry Smith }
12149b94acceSBarry Smith 
12159b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12169b94acceSBarry Smith 
12175615d1e5SSatish Balay #undef __FUNC__
12185615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
12199b94acceSBarry Smith /*@
1220d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12219b94acceSBarry Smith 
12229b94acceSBarry Smith    Input Parameters:
12239b94acceSBarry Smith .  snes - the SNES context
122433174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
122533174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
122633174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
122733174efeSLois Curfman McInnes            of the change in the solution between steps
122833174efeSLois Curfman McInnes .  maxit - maximum number of iterations
122933174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
12309b94acceSBarry Smith 
123133174efeSLois Curfman McInnes    Options Database Keys:
123233174efeSLois Curfman McInnes $    -snes_atol <atol>
123333174efeSLois Curfman McInnes $    -snes_rtol <rtol>
123433174efeSLois Curfman McInnes $    -snes_stol <stol>
123533174efeSLois Curfman McInnes $    -snes_max_it <maxit>
123633174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
12379b94acceSBarry Smith 
1238d7a720efSLois Curfman McInnes    Notes:
12399b94acceSBarry Smith    The default maximum number of iterations is 50.
12409b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
12419b94acceSBarry Smith 
124233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
12439b94acceSBarry Smith 
1244d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
12459b94acceSBarry Smith @*/
1246d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
12479b94acceSBarry Smith {
124877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1249d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1250d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1251d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1252d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1253d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
12549b94acceSBarry Smith   return 0;
12559b94acceSBarry Smith }
12569b94acceSBarry Smith 
12575615d1e5SSatish Balay #undef __FUNC__
12585615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
12599b94acceSBarry Smith /*@
126033174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
126133174efeSLois Curfman McInnes 
126233174efeSLois Curfman McInnes    Input Parameters:
126333174efeSLois Curfman McInnes .  snes - the SNES context
126433174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
126533174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
126633174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
126733174efeSLois Curfman McInnes            of the change in the solution between steps
126833174efeSLois Curfman McInnes .  maxit - maximum number of iterations
126933174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
127033174efeSLois Curfman McInnes 
127133174efeSLois Curfman McInnes    Notes:
127233174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
127333174efeSLois Curfman McInnes 
127433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
127533174efeSLois Curfman McInnes 
127633174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
127733174efeSLois Curfman McInnes @*/
127833174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
127933174efeSLois Curfman McInnes {
128033174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
128133174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
128233174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
128333174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
128433174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
128533174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
128633174efeSLois Curfman McInnes   return 0;
128733174efeSLois Curfman McInnes }
128833174efeSLois Curfman McInnes 
12895615d1e5SSatish Balay #undef __FUNC__
12905615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
129133174efeSLois Curfman McInnes /*@
12929b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12939b94acceSBarry Smith 
12949b94acceSBarry Smith    Input Parameters:
12959b94acceSBarry Smith .  snes - the SNES context
12969b94acceSBarry Smith .  tol - tolerance
12979b94acceSBarry Smith 
12989b94acceSBarry Smith    Options Database Key:
1299d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
13009b94acceSBarry Smith 
13019b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
13029b94acceSBarry Smith 
1303d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
13049b94acceSBarry Smith @*/
13059b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
13069b94acceSBarry Smith {
130777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13089b94acceSBarry Smith   snes->deltatol = tol;
13099b94acceSBarry Smith   return 0;
13109b94acceSBarry Smith }
13119b94acceSBarry Smith 
13125615d1e5SSatish Balay #undef __FUNC__
13135615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
13149b94acceSBarry Smith /*@
131577c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
13169b94acceSBarry Smith    for unconstrained minimization solvers.
13179b94acceSBarry Smith 
13189b94acceSBarry Smith    Input Parameters:
13199b94acceSBarry Smith .  snes - the SNES context
13209b94acceSBarry Smith .  ftol - minimum function tolerance
13219b94acceSBarry Smith 
13229b94acceSBarry Smith    Options Database Key:
1323d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
13249b94acceSBarry Smith 
13259b94acceSBarry Smith    Note:
132677c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
13279b94acceSBarry Smith    methods only.
13289b94acceSBarry Smith 
13299b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
13309b94acceSBarry Smith 
1331d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
13329b94acceSBarry Smith @*/
133377c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
13349b94acceSBarry Smith {
133577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13369b94acceSBarry Smith   snes->fmin = ftol;
13379b94acceSBarry Smith   return 0;
13389b94acceSBarry Smith }
13399b94acceSBarry Smith 
13409b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
13419b94acceSBarry Smith 
13425615d1e5SSatish Balay #undef __FUNC__
13435eea60f9SBarry Smith #define __FUNC__ "SNESSetMonitor" /* ADIC Ignore */
13449b94acceSBarry Smith /*@C
13459b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
13469b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
13479b94acceSBarry Smith    progress.
13489b94acceSBarry Smith 
13499b94acceSBarry Smith    Input Parameters:
13509b94acceSBarry Smith .  snes - the SNES context
13519b94acceSBarry Smith .  func - monitoring routine
1352044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1353044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
13549b94acceSBarry Smith 
13559b94acceSBarry Smith    Calling sequence of func:
1356682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
13579b94acceSBarry Smith 
13589b94acceSBarry Smith $    snes - the SNES context
13599b94acceSBarry Smith $    its - iteration number
1360044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
13619b94acceSBarry Smith $
13629b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13639b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
13649b94acceSBarry Smith $
13659b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13669b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
13679b94acceSBarry Smith 
13689665c990SLois Curfman McInnes    Options Database Keys:
13699665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
13709665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
13719665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
13729665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
13739665c990SLois Curfman McInnes $                           been hardwired into a code by
13749665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
13759665c990SLois Curfman McInnes $                           does not cancel those set via
13769665c990SLois Curfman McInnes $                           the options database.
13779665c990SLois Curfman McInnes 
13789665c990SLois Curfman McInnes 
1379639f9d9dSBarry Smith    Notes:
13806bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13816bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13826bc08f3fSLois Curfman McInnes    order in which they were set.
1383639f9d9dSBarry Smith 
13849b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13859b94acceSBarry Smith 
13869b94acceSBarry Smith .seealso: SNESDefaultMonitor()
13879b94acceSBarry Smith @*/
138874679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
13899b94acceSBarry Smith {
1390639f9d9dSBarry Smith   if (!func) {
1391639f9d9dSBarry Smith     snes->numbermonitors = 0;
1392639f9d9dSBarry Smith     return 0;
1393639f9d9dSBarry Smith   }
1394639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1395e3372554SBarry Smith     SETERRQ(1,0,"Too many monitors set");
1396639f9d9dSBarry Smith   }
1397639f9d9dSBarry Smith 
1398639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1399639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
14009b94acceSBarry Smith   return 0;
14019b94acceSBarry Smith }
14029b94acceSBarry Smith 
14035615d1e5SSatish Balay #undef __FUNC__
14045eea60f9SBarry Smith #define __FUNC__ "SNESSetConvergenceTest" /* ADIC Ignore */
14059b94acceSBarry Smith /*@C
14069b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
14079b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
14089b94acceSBarry Smith 
14099b94acceSBarry Smith    Input Parameters:
14109b94acceSBarry Smith .  snes - the SNES context
14119b94acceSBarry Smith .  func - routine to test for convergence
1412044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1413044dda88SLois Curfman McInnes           (may be PETSC_NULL)
14149b94acceSBarry Smith 
14159b94acceSBarry Smith    Calling sequence of func:
14169b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
14179b94acceSBarry Smith              double f,void *cctx)
14189b94acceSBarry Smith 
14199b94acceSBarry Smith $    snes - the SNES context
1420044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
14219b94acceSBarry Smith $    xnorm - 2-norm of current iterate
14229b94acceSBarry Smith $
14239b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14249b94acceSBarry Smith $    gnorm - 2-norm of current step
14259b94acceSBarry Smith $    f - 2-norm of function
14269b94acceSBarry Smith $
14279b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14289b94acceSBarry Smith $    gnorm - 2-norm of current gradient
14299b94acceSBarry Smith $    f - function value
14309b94acceSBarry Smith 
14319b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14329b94acceSBarry Smith 
143340191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
143440191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
14359b94acceSBarry Smith @*/
143674679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
14379b94acceSBarry Smith {
14389b94acceSBarry Smith   (snes)->converged = func;
14399b94acceSBarry Smith   (snes)->cnvP      = cctx;
14409b94acceSBarry Smith   return 0;
14419b94acceSBarry Smith }
14429b94acceSBarry Smith 
14435615d1e5SSatish Balay #undef __FUNC__
14445eea60f9SBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" /* ADIC Ignore */
1445c9005455SLois Curfman McInnes /*@
1446c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1447c9005455SLois Curfman McInnes 
1448c9005455SLois Curfman McInnes    Input Parameters:
1449c9005455SLois Curfman McInnes .  snes - iterative context obtained from SNESCreate()
1450c9005455SLois Curfman McInnes .  a   - array to hold history
1451c9005455SLois Curfman McInnes .  na  - size of a
1452c9005455SLois Curfman McInnes 
1453c9005455SLois Curfman McInnes    Notes:
1454c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1455c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1456c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1457c9005455SLois Curfman McInnes    at each step.
1458c9005455SLois Curfman McInnes 
1459c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1460c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1461c9005455SLois Curfman McInnes    during the section of code that is being timed.
1462c9005455SLois Curfman McInnes 
1463c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1464c9005455SLois Curfman McInnes @*/
1465c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1466c9005455SLois Curfman McInnes {
1467c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1468c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1469c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1470c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
1471c9005455SLois Curfman McInnes   return 0;
1472c9005455SLois Curfman McInnes }
1473c9005455SLois Curfman McInnes 
1474c9005455SLois Curfman McInnes #undef __FUNC__
14755615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
14769b94acceSBarry Smith /*
14779b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
14789b94acceSBarry Smith    positive parameter delta.
14799b94acceSBarry Smith 
14809b94acceSBarry Smith     Input Parameters:
14819b94acceSBarry Smith .   snes - the SNES context
14829b94acceSBarry Smith .   y - approximate solution of linear system
14839b94acceSBarry Smith .   fnorm - 2-norm of current function
14849b94acceSBarry Smith .   delta - trust region size
14859b94acceSBarry Smith 
14869b94acceSBarry Smith     Output Parameters:
14879b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
14889b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
14899b94acceSBarry Smith     region, and exceeds zero otherwise.
14909b94acceSBarry Smith .   ynorm - 2-norm of the step
14919b94acceSBarry Smith 
14929b94acceSBarry Smith     Note:
149340191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
14949b94acceSBarry Smith     is set to be the maximum allowable step size.
14959b94acceSBarry Smith 
14969b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
14979b94acceSBarry Smith */
14989b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
14999b94acceSBarry Smith                           double *gpnorm,double *ynorm)
15009b94acceSBarry Smith {
15019b94acceSBarry Smith   double norm;
15029b94acceSBarry Smith   Scalar cnorm;
1503cddf8d76SBarry Smith   VecNorm(y,NORM_2, &norm );
15049b94acceSBarry Smith   if (norm > *delta) {
15059b94acceSBarry Smith      norm = *delta/norm;
15069b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
15079b94acceSBarry Smith      cnorm = norm;
15089b94acceSBarry Smith      VecScale( &cnorm, y );
15099b94acceSBarry Smith      *ynorm = *delta;
15109b94acceSBarry Smith   } else {
15119b94acceSBarry Smith      *gpnorm = 0.0;
15129b94acceSBarry Smith      *ynorm = norm;
15139b94acceSBarry Smith   }
15149b94acceSBarry Smith   return 0;
15159b94acceSBarry Smith }
15169b94acceSBarry Smith 
15175615d1e5SSatish Balay #undef __FUNC__
15185615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
15199b94acceSBarry Smith /*@
15209b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
152163a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
15229b94acceSBarry Smith 
15239b94acceSBarry Smith    Input Parameter:
15249b94acceSBarry Smith .  snes - the SNES context
15258ddd3da0SLois Curfman McInnes .  x - the solution vector
15269b94acceSBarry Smith 
15279b94acceSBarry Smith    Output Parameter:
15289b94acceSBarry Smith    its - number of iterations until termination
15299b94acceSBarry Smith 
15308ddd3da0SLois Curfman McInnes    Note:
15318ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
15328ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
15338ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
15348ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
15358ddd3da0SLois Curfman McInnes 
15369b94acceSBarry Smith .keywords: SNES, nonlinear, solve
15379b94acceSBarry Smith 
153863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
15399b94acceSBarry Smith @*/
15408ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
15419b94acceSBarry Smith {
15423c7409f5SSatish Balay   int ierr, flg;
1543052efed2SBarry Smith 
154477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
154574679c65SBarry Smith   PetscValidIntPointer(its);
1546c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1547c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
15489b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1549c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
15509b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
15519b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
15523c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
15536d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
15549b94acceSBarry Smith   return 0;
15559b94acceSBarry Smith }
15569b94acceSBarry Smith 
15579b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
155837fcc0dbSBarry Smith static NRList *__SNESList = 0;
15599b94acceSBarry Smith 
15605615d1e5SSatish Balay #undef __FUNC__
15615615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
15629b94acceSBarry Smith /*@
15634b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
15649b94acceSBarry Smith 
15659b94acceSBarry Smith    Input Parameters:
15669b94acceSBarry Smith .  snes - the SNES context
15679b94acceSBarry Smith .  method - a known method
15689b94acceSBarry Smith 
1569ae12b187SLois Curfman McInnes   Options Database Command:
1570ae12b187SLois Curfman McInnes $ -snes_type  <method>
1571ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1572ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1573ae12b187SLois Curfman McInnes 
15749b94acceSBarry Smith    Notes:
15759b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
15769b94acceSBarry Smith $  Systems of nonlinear equations:
157740191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
157840191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
15799b94acceSBarry Smith $  Unconstrained minimization:
158040191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
158140191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
15829b94acceSBarry Smith 
1583ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1584ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1585ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1586ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1587ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1588ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1589ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1590ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1591ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1592ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1593a703fe33SLois Curfman McInnes 
1594f536c699SSatish Balay .keywords: SNES, set, method
15959b94acceSBarry Smith @*/
15964b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
15979b94acceSBarry Smith {
159884cb2905SBarry Smith   int ierr;
15999b94acceSBarry Smith   int (*r)(SNES);
160077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16017d1a2b2bSBarry Smith   if (snes->type == (int) method) return 0;
16027d1a2b2bSBarry Smith 
16039b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
160484cb2905SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
1605e3372554SBarry Smith   if (!__SNESList) {SETERRQ(1,0,"Could not get methods");}
160637fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1607e3372554SBarry Smith   if (!r) {SETERRQ(1,0,"Unknown method");}
16080452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
16099b94acceSBarry Smith   snes->set_method_called = 1;
16109b94acceSBarry Smith   return (*r)(snes);
16119b94acceSBarry Smith }
16129b94acceSBarry Smith 
16139b94acceSBarry Smith /* --------------------------------------------------------------------- */
16145615d1e5SSatish Balay #undef __FUNC__
16155eea60f9SBarry Smith #define __FUNC__ "SNESRegister" /* ADIC Ignore */
16169b94acceSBarry Smith /*@C
16179b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
16184b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
16199b94acceSBarry Smith 
16209b94acceSBarry Smith    Input Parameters:
16212d872ea7SLois Curfman McInnes .  name - either a predefined name such as SNES_EQ_LS, or SNES_NEW
16222d872ea7SLois Curfman McInnes           to indicate a new user-defined solver
162340191667SLois Curfman McInnes .  sname - corresponding string for name
16249b94acceSBarry Smith .  create - routine to create method context
16259b94acceSBarry Smith 
162684cb2905SBarry Smith    Output Parameter:
162784cb2905SBarry Smith .  oname - type associated with this new solver
162884cb2905SBarry Smith 
16292d872ea7SLois Curfman McInnes    Notes:
16302d872ea7SLois Curfman McInnes    Multiple user-defined nonlinear solvers can be added by calling
16312d872ea7SLois Curfman McInnes    SNESRegister() with the input parameter "name" set to be SNES_NEW;
16322d872ea7SLois Curfman McInnes    each call will return a unique solver type in the output
16332d872ea7SLois Curfman McInnes    parameter "oname".
16342d872ea7SLois Curfman McInnes 
16359b94acceSBarry Smith .keywords: SNES, nonlinear, register
16369b94acceSBarry Smith 
16379b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
16389b94acceSBarry Smith @*/
163984cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES))
16409b94acceSBarry Smith {
16419b94acceSBarry Smith   int ierr;
164284cb2905SBarry Smith   static int numberregistered = 0;
164384cb2905SBarry Smith 
1644d252947aSBarry Smith   if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++);
164584cb2905SBarry Smith 
164684cb2905SBarry Smith   if (oname) *oname = name;
164737fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
164884cb2905SBarry Smith   NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create );
16499b94acceSBarry Smith   return 0;
16509b94acceSBarry Smith }
1651a847f771SSatish Balay 
16529b94acceSBarry Smith /* --------------------------------------------------------------------- */
16535615d1e5SSatish Balay #undef __FUNC__
16545eea60f9SBarry Smith #define __FUNC__ "SNESRegisterDestroy" /* ADIC Ignore */
16559b94acceSBarry Smith /*@C
16569b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
16579b94acceSBarry Smith    registered by SNESRegister().
16589b94acceSBarry Smith 
16599b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
16609b94acceSBarry Smith 
16619b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
16629b94acceSBarry Smith @*/
16639b94acceSBarry Smith int SNESRegisterDestroy()
16649b94acceSBarry Smith {
166537fcc0dbSBarry Smith   if (__SNESList) {
166637fcc0dbSBarry Smith     NRDestroy( __SNESList );
166737fcc0dbSBarry Smith     __SNESList = 0;
16689b94acceSBarry Smith   }
166984cb2905SBarry Smith   SNESRegisterAllCalled = 0;
16709b94acceSBarry Smith   return 0;
16719b94acceSBarry Smith }
16729b94acceSBarry Smith 
16735615d1e5SSatish Balay #undef __FUNC__
16745eea60f9SBarry Smith #define __FUNC__ "SNESGetTypeFromOptions_Private" /* ADIC Ignore */
16759b94acceSBarry Smith /*
16764b0e389bSBarry Smith    SNESGetTypeFromOptions_Private - Sets the selected method from the
16779b94acceSBarry Smith    options database.
16789b94acceSBarry Smith 
16799b94acceSBarry Smith    Input Parameter:
16809b94acceSBarry Smith .  ctx - the SNES context
16819b94acceSBarry Smith 
16829b94acceSBarry Smith    Output Parameter:
16839b94acceSBarry Smith .  method -  solver method
16849b94acceSBarry Smith 
16859b94acceSBarry Smith    Returns:
16869b94acceSBarry Smith    Returns 1 if the method is found; 0 otherwise.
16879b94acceSBarry Smith 
16889b94acceSBarry Smith    Options Database Key:
16894b0e389bSBarry Smith $  -snes_type  method
16909b94acceSBarry Smith */
1691052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
16929b94acceSBarry Smith {
1693052efed2SBarry Smith   int  ierr;
16949b94acceSBarry Smith   char sbuf[50];
16955baf8537SBarry Smith 
1696052efed2SBarry Smith   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1697052efed2SBarry Smith   if (*flg) {
1698052efed2SBarry Smith     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
169937fcc0dbSBarry Smith     *method = (SNESType)NRFindID( __SNESList, sbuf );
17009b94acceSBarry Smith   }
17019b94acceSBarry Smith   return 0;
17029b94acceSBarry Smith }
17039b94acceSBarry Smith 
17045615d1e5SSatish Balay #undef __FUNC__
17055eea60f9SBarry Smith #define __FUNC__ "SNESGetType" /* ADIC Ignore */
17069b94acceSBarry Smith /*@C
17079a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
17089b94acceSBarry Smith 
17099b94acceSBarry Smith    Input Parameter:
17104b0e389bSBarry Smith .  snes - nonlinear solver context
17119b94acceSBarry Smith 
17129b94acceSBarry Smith    Output Parameter:
17139a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
17149a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
17159b94acceSBarry Smith 
17169b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
17179b94acceSBarry Smith @*/
17184b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
17199b94acceSBarry Smith {
172006a9b0adSLois Curfman McInnes   int ierr;
172137fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
17224b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
172337fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
17249b94acceSBarry Smith   return 0;
17259b94acceSBarry Smith }
17269b94acceSBarry Smith 
17279b94acceSBarry Smith #include <stdio.h>
17285615d1e5SSatish Balay #undef __FUNC__
17295eea60f9SBarry Smith #define __FUNC__ "SNESPrintTypes_Private" /* ADIC Ignore */
17309b94acceSBarry Smith /*
17314b0e389bSBarry Smith    SNESPrintTypes_Private - Prints the SNES methods available from the
17329b94acceSBarry Smith    options database.
17339b94acceSBarry Smith 
17349b94acceSBarry Smith    Input Parameters:
173533455573SSatish Balay .  comm   - communicator (usually MPI_COMM_WORLD)
17369b94acceSBarry Smith .  prefix - prefix (usually "-")
17374b0e389bSBarry Smith .  name   - the options database name (by default "snes_type")
17389b94acceSBarry Smith */
173933455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
17409b94acceSBarry Smith {
17419b94acceSBarry Smith   FuncList *entry;
174237fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
174337fcc0dbSBarry Smith   entry = __SNESList->head;
174477c4ece6SBarry Smith   PetscPrintf(comm," %s%s (one of)",prefix,name);
17459b94acceSBarry Smith   while (entry) {
174677c4ece6SBarry Smith     PetscPrintf(comm," %s",entry->name);
17479b94acceSBarry Smith     entry = entry->next;
17489b94acceSBarry Smith   }
174977c4ece6SBarry Smith   PetscPrintf(comm,"\n");
17509b94acceSBarry Smith   return 0;
17519b94acceSBarry Smith }
17529b94acceSBarry Smith 
17535615d1e5SSatish Balay #undef __FUNC__
17545eea60f9SBarry Smith #define __FUNC__ "SNESGetSolution" /* ADIC Ignore */
17559b94acceSBarry Smith /*@C
17569b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
17579b94acceSBarry Smith    stored.
17589b94acceSBarry Smith 
17599b94acceSBarry Smith    Input Parameter:
17609b94acceSBarry Smith .  snes - the SNES context
17619b94acceSBarry Smith 
17629b94acceSBarry Smith    Output Parameter:
17639b94acceSBarry Smith .  x - the solution
17649b94acceSBarry Smith 
17659b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
17669b94acceSBarry Smith 
1767a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
17689b94acceSBarry Smith @*/
17699b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
17709b94acceSBarry Smith {
177177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17729b94acceSBarry Smith   *x = snes->vec_sol_always;
17739b94acceSBarry Smith   return 0;
17749b94acceSBarry Smith }
17759b94acceSBarry Smith 
17765615d1e5SSatish Balay #undef __FUNC__
17775eea60f9SBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" /* ADIC Ignore */
17789b94acceSBarry Smith /*@C
17799b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
17809b94acceSBarry Smith    stored.
17819b94acceSBarry Smith 
17829b94acceSBarry Smith    Input Parameter:
17839b94acceSBarry Smith .  snes - the SNES context
17849b94acceSBarry Smith 
17859b94acceSBarry Smith    Output Parameter:
17869b94acceSBarry Smith .  x - the solution update
17879b94acceSBarry Smith 
17889b94acceSBarry Smith    Notes:
17899b94acceSBarry Smith    This vector is implementation dependent.
17909b94acceSBarry Smith 
17919b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
17929b94acceSBarry Smith 
17939b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
17949b94acceSBarry Smith @*/
17959b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
17969b94acceSBarry Smith {
179777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17989b94acceSBarry Smith   *x = snes->vec_sol_update_always;
17999b94acceSBarry Smith   return 0;
18009b94acceSBarry Smith }
18019b94acceSBarry Smith 
18025615d1e5SSatish Balay #undef __FUNC__
18035eea60f9SBarry Smith #define __FUNC__ "SNESGetFunction" /* ADIC Ignore */
18049b94acceSBarry Smith /*@C
18053638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
18069b94acceSBarry Smith 
18079b94acceSBarry Smith    Input Parameter:
18089b94acceSBarry Smith .  snes - the SNES context
18099b94acceSBarry Smith 
18109b94acceSBarry Smith    Output Parameter:
18113638b69dSLois Curfman McInnes .  r - the function
18129b94acceSBarry Smith 
18139b94acceSBarry Smith    Notes:
18149b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
18159b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
18169b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
18179b94acceSBarry Smith 
1818a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
18199b94acceSBarry Smith 
182031615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
182131615d04SLois Curfman McInnes           SNESGetGradient()
18229b94acceSBarry Smith @*/
18239b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
18249b94acceSBarry Smith {
182577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1826e3372554SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,
18277c792091SSatish Balay     "For SNES_NONLINEAR_EQUATIONS only");
18289b94acceSBarry Smith   *r = snes->vec_func_always;
18299b94acceSBarry Smith   return 0;
18309b94acceSBarry Smith }
18319b94acceSBarry Smith 
18325615d1e5SSatish Balay #undef __FUNC__
18335eea60f9SBarry Smith #define __FUNC__ "SNESGetGradient" /* ADIC Ignore */
18349b94acceSBarry Smith /*@C
18353638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
18369b94acceSBarry Smith 
18379b94acceSBarry Smith    Input Parameter:
18389b94acceSBarry Smith .  snes - the SNES context
18399b94acceSBarry Smith 
18409b94acceSBarry Smith    Output Parameter:
18419b94acceSBarry Smith .  r - the gradient
18429b94acceSBarry Smith 
18439b94acceSBarry Smith    Notes:
18449b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
18459b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18469b94acceSBarry Smith    SNESGetFunction().
18479b94acceSBarry Smith 
18489b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
18499b94acceSBarry Smith 
185031615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
18519b94acceSBarry Smith @*/
18529b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
18539b94acceSBarry Smith {
185477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1855e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,
185663c41f6aSSatish Balay     "For SNES_UNCONSTRAINED_MINIMIZATION only");
18579b94acceSBarry Smith   *r = snes->vec_func_always;
18589b94acceSBarry Smith   return 0;
18599b94acceSBarry Smith }
18609b94acceSBarry Smith 
18615615d1e5SSatish Balay #undef __FUNC__
18625eea60f9SBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" /* ADIC Ignore */
18639b94acceSBarry Smith /*@
18649b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
18659b94acceSBarry Smith    unconstrained minimization problems.
18669b94acceSBarry Smith 
18679b94acceSBarry Smith    Input Parameter:
18689b94acceSBarry Smith .  snes - the SNES context
18699b94acceSBarry Smith 
18709b94acceSBarry Smith    Output Parameter:
18719b94acceSBarry Smith .  r - the function
18729b94acceSBarry Smith 
18739b94acceSBarry Smith    Notes:
18749b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
18759b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18769b94acceSBarry Smith    SNESGetFunction().
18779b94acceSBarry Smith 
18789b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
18799b94acceSBarry Smith 
188031615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
18819b94acceSBarry Smith @*/
18829b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
18839b94acceSBarry Smith {
188477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
188574679c65SBarry Smith   PetscValidScalarPointer(r);
1886e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,
188763c41f6aSSatish Balay     "For SNES_UNCONSTRAINED_MINIMIZATION only");
18889b94acceSBarry Smith   *r = snes->fc;
18899b94acceSBarry Smith   return 0;
18909b94acceSBarry Smith }
18919b94acceSBarry Smith 
18925615d1e5SSatish Balay #undef __FUNC__
18935eea60f9SBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" /* ADIC Ignore */
18943c7409f5SSatish Balay /*@C
18953c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1896d850072dSLois Curfman McInnes    SNES options in the database.
18973c7409f5SSatish Balay 
18983c7409f5SSatish Balay    Input Parameter:
18993c7409f5SSatish Balay .  snes - the SNES context
19003c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
19013c7409f5SSatish Balay 
1902d850072dSLois Curfman McInnes    Notes:
1903a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1904a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1905a83b1b31SSatish Balay    hyphen.
1906d850072dSLois Curfman McInnes 
19073c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1908a86d99e1SLois Curfman McInnes 
1909a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
19103c7409f5SSatish Balay @*/
19113c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
19123c7409f5SSatish Balay {
19133c7409f5SSatish Balay   int ierr;
19143c7409f5SSatish Balay 
191577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1916639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19173c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19183c7409f5SSatish Balay   return 0;
19193c7409f5SSatish Balay }
19203c7409f5SSatish Balay 
19215615d1e5SSatish Balay #undef __FUNC__
19225eea60f9SBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" /* ADIC Ignore */
19233c7409f5SSatish Balay /*@C
1924f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1925d850072dSLois Curfman McInnes    SNES options in the database.
19263c7409f5SSatish Balay 
19273c7409f5SSatish Balay    Input Parameter:
19283c7409f5SSatish Balay .  snes - the SNES context
19293c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
19303c7409f5SSatish Balay 
1931d850072dSLois Curfman McInnes    Notes:
1932a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1933a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1934a83b1b31SSatish Balay    hyphen.
1935d850072dSLois Curfman McInnes 
19363c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1937a86d99e1SLois Curfman McInnes 
1938a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
19393c7409f5SSatish Balay @*/
19403c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
19413c7409f5SSatish Balay {
19423c7409f5SSatish Balay   int ierr;
19433c7409f5SSatish Balay 
194477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1945639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19463c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19473c7409f5SSatish Balay   return 0;
19483c7409f5SSatish Balay }
19493c7409f5SSatish Balay 
19505615d1e5SSatish Balay #undef __FUNC__
19515eea60f9SBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" /* ADIC Ignore */
1952c83590e2SSatish Balay /*@
19533c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
19543c7409f5SSatish Balay    SNES options in the database.
19553c7409f5SSatish Balay 
19563c7409f5SSatish Balay    Input Parameter:
19573c7409f5SSatish Balay .  snes - the SNES context
19583c7409f5SSatish Balay 
19593c7409f5SSatish Balay    Output Parameter:
19603c7409f5SSatish Balay .  prefix - pointer to the prefix string used
19613c7409f5SSatish Balay 
19623c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1963a86d99e1SLois Curfman McInnes 
1964a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
19653c7409f5SSatish Balay @*/
19663c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
19673c7409f5SSatish Balay {
19683c7409f5SSatish Balay   int ierr;
19693c7409f5SSatish Balay 
197077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1971639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19723c7409f5SSatish Balay   return 0;
19733c7409f5SSatish Balay }
19743c7409f5SSatish Balay 
19753c7409f5SSatish Balay 
19769b94acceSBarry Smith 
19779b94acceSBarry Smith 
19789b94acceSBarry Smith 
1979