xref: /petsc/src/snes/interface/snes.c (revision d64ed03d54a2fae6c0c58880c522fae65f427572)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*d64ed03dSBarry Smith static char vcid[] = "$Id: snes.c,v 1.131 1997/10/19 03:29:25 bsmith Exp bsmith $";
39b94acceSBarry Smith #endif
49b94acceSBarry Smith 
570f55243SBarry Smith #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
6f5eb4b81SSatish Balay #include "src/sys/nreg.h"
79b94acceSBarry Smith #include "pinclude/pviewer.h"
89b94acceSBarry Smith #include <math.h>
99b94acceSBarry Smith 
10052efed2SBarry Smith extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*);
1133455573SSatish Balay extern int SNESPrintTypes_Private(MPI_Comm,char*,char*);
129b94acceSBarry Smith 
1384cb2905SBarry Smith int SNESRegisterAllCalled = 0;
1484cb2905SBarry Smith 
155615d1e5SSatish Balay #undef __FUNC__
16d4bb536fSBarry Smith #define __FUNC__ "SNESView"
179b94acceSBarry Smith /*@
189b94acceSBarry Smith    SNESView - Prints the SNES data structure.
199b94acceSBarry Smith 
209b94acceSBarry Smith    Input Parameters:
219b94acceSBarry Smith .  SNES - the SNES context
229b94acceSBarry Smith .  viewer - visualization context
239b94acceSBarry Smith 
249b94acceSBarry Smith    Options Database Key:
259b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
269b94acceSBarry Smith 
279b94acceSBarry Smith    Notes:
289b94acceSBarry Smith    The available visualization contexts include
296d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
306d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
319b94acceSBarry Smith $       output where only the first processor opens
329b94acceSBarry Smith $       the file.  All other processors send their
339b94acceSBarry Smith $       data to the first processor to print.
349b94acceSBarry Smith 
359b94acceSBarry Smith    The user can open alternative vistualization contexts with
36dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
379b94acceSBarry Smith 
389b94acceSBarry Smith .keywords: SNES, view
399b94acceSBarry Smith 
40dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
419b94acceSBarry Smith @*/
429b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
439b94acceSBarry Smith {
449b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
459b94acceSBarry Smith   FILE                *fd;
469b94acceSBarry Smith   int                 ierr;
479b94acceSBarry Smith   SLES                sles;
489b94acceSBarry Smith   char                *method;
4919bcc07fSBarry Smith   ViewerType          vtype;
509b94acceSBarry Smith 
513a40ed3dSBarry Smith   PetscFunctionBegin;
5274679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5374679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5474679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5574679c65SBarry Smith 
5619bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5719bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5890ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5977c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
604b0e389bSBarry Smith     SNESGetType(snes,PETSC_NULL,&method);
6177c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
629b94acceSBarry Smith     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
6377c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
649b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
659b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6677c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
679b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
689b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
697a00f4a9SLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
707a00f4a9SLois Curfman McInnes     "  total number of linear solver iterations=%d\n",snes->linear_its);
71ae3c334cSLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
7268632166SLois Curfman McInnes      "  total number of function evaluations=%d\n",snes->nfuncs);
739b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
7477c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
759b94acceSBarry Smith     if (snes->ksp_ewconv) {
769b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
779b94acceSBarry Smith       if (kctx) {
7877c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
799b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
809b94acceSBarry Smith         kctx->version);
8177c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
829b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
839b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
8477c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
859b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
869b94acceSBarry Smith       }
879b94acceSBarry Smith     }
88c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
89c30f285eSLois Curfman McInnes     SNESGetType(snes,PETSC_NULL,&method);
90c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
9119bcc07fSBarry Smith   }
929b94acceSBarry Smith   SNESGetSLES(snes,&sles);
939b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
943a40ed3dSBarry Smith   PetscFunctionReturn(0);
959b94acceSBarry Smith }
969b94acceSBarry Smith 
97639f9d9dSBarry Smith /*
98639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
99639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
100639f9d9dSBarry Smith */
101639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
102639f9d9dSBarry Smith static int numberofsetfromoptions;
103639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
104639f9d9dSBarry Smith 
1055615d1e5SSatish Balay #undef __FUNC__
106d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker"
107639f9d9dSBarry Smith /*@
108639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
109639f9d9dSBarry Smith 
110639f9d9dSBarry Smith     Input Parameter:
111639f9d9dSBarry Smith .   snescheck - function that checks for options
112639f9d9dSBarry Smith 
113639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
114639f9d9dSBarry Smith @*/
115639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
116639f9d9dSBarry Smith {
1173a40ed3dSBarry Smith   PetscFunctionBegin;
118639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
119e3372554SBarry Smith     SETERRQ(1,0,"Too many options checkers, only 5 allowed");
120639f9d9dSBarry Smith   }
121639f9d9dSBarry Smith 
122639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
1233a40ed3dSBarry Smith   PetscFunctionReturn(0);
124639f9d9dSBarry Smith }
125639f9d9dSBarry Smith 
1265615d1e5SSatish Balay #undef __FUNC__
1275615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1289b94acceSBarry Smith /*@
129682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1309b94acceSBarry Smith 
1319b94acceSBarry Smith    Input Parameter:
1329b94acceSBarry Smith .  snes - the SNES context
1339b94acceSBarry Smith 
13411ca99fdSLois Curfman McInnes    Notes:
13511ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
13611ca99fdSLois Curfman McInnes    the users manual.
13783e2fdc7SBarry Smith 
1389b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1399b94acceSBarry Smith 
140a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1419b94acceSBarry Smith @*/
1429b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1439b94acceSBarry Smith {
1444b0e389bSBarry Smith   SNESType method;
1459b94acceSBarry Smith   double   tmp;
1469b94acceSBarry Smith   SLES     sles;
14781f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1483c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1499b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1509b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1519b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1529b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1539b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1549b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1559b94acceSBarry Smith 
1563a40ed3dSBarry Smith   PetscFunctionBegin;
15781f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
15881f57ec7SBarry Smith 
15977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
160e3372554SBarry Smith   if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp");
161052efed2SBarry Smith   ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr);
162052efed2SBarry Smith   if (flg) {
163052efed2SBarry Smith     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
1649b94acceSBarry Smith   }
1655334005bSBarry Smith   else if (!snes->set_method_called) {
1665334005bSBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
16740191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
168*d64ed03dSBarry Smith     } else {
16940191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
1705334005bSBarry Smith     }
1715334005bSBarry Smith   }
1725334005bSBarry Smith 
1733c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
174*d64ed03dSBarry Smith   if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);}
1753c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
176*d64ed03dSBarry Smith   if (flg) {
177*d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
178*d64ed03dSBarry Smith   }
1793c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
180*d64ed03dSBarry Smith   if (flg) {
181*d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
182*d64ed03dSBarry Smith   }
1833c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
184*d64ed03dSBarry Smith   if (flg) {
185*d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
186*d64ed03dSBarry Smith   }
1873c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1883c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
189d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
190*d64ed03dSBarry Smith   if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); }
191d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
192*d64ed03dSBarry Smith   if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);  CHKERRQ(ierr);}
1933c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1943c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
195b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
196b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
197b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
198b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
199b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
200b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
201b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
20293c39befSBarry Smith 
20393c39befSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr);
20493c39befSBarry Smith   if (flg) {snes->converged = 0;}
20593c39befSBarry Smith 
2069b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
2079b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
2085bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
2095bbad29bSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);}
2103c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
211639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
2123c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
213d31a3109SLois Curfman McInnes   if (flg) {
214639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
215d31a3109SLois Curfman McInnes   }
21681f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
2173c7409f5SSatish Balay   if (flg) {
21817699dbbSLois Curfman McInnes     int    rank = 0;
219d7e8b826SBarry Smith     DrawLG lg;
22017699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
22117699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
22217699dbbSLois Curfman McInnes     if (!rank) {
22381f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2249b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
225f8c826e1SBarry Smith       snes->xmonitor = lg;
2269b94acceSBarry Smith       PLogObjectParent(snes,lg);
2279b94acceSBarry Smith     }
2289b94acceSBarry Smith   }
2293c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2303c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2319b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2329b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
23394a424c1SBarry Smith     PLogInfo(snes,
23431615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
23531615d04SLois Curfman McInnes   }
23631615d04SLois Curfman McInnes   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
23731615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
23831615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
239*d64ed03dSBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2409b94acceSBarry Smith   }
241639f9d9dSBarry Smith 
242639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
243639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
244639f9d9dSBarry Smith   }
245639f9d9dSBarry Smith 
2469b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2479b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
2483a40ed3dSBarry Smith   if (snes->setfromoptions) {
2493a40ed3dSBarry Smith     ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
2503a40ed3dSBarry Smith   }
2513a40ed3dSBarry Smith   PetscFunctionReturn(0);
2529b94acceSBarry Smith }
2539b94acceSBarry Smith 
2545615d1e5SSatish Balay #undef __FUNC__
255d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp"
2569b94acceSBarry Smith /*@
2579b94acceSBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2589b94acceSBarry Smith 
2599b94acceSBarry Smith    Input Parameter:
2609b94acceSBarry Smith .  snes - the SNES context
2619b94acceSBarry Smith 
262a703fe33SLois Curfman McInnes    Options Database Keys:
263a703fe33SLois Curfman McInnes $  -help, -h
264a703fe33SLois Curfman McInnes 
2659b94acceSBarry Smith .keywords: SNES, nonlinear, help
2669b94acceSBarry Smith 
267682d7d0cSBarry Smith .seealso: SNESSetFromOptions()
2689b94acceSBarry Smith @*/
2699b94acceSBarry Smith int SNESPrintHelp(SNES snes)
2709b94acceSBarry Smith {
2716daaf66cSBarry Smith   char                p[64];
2726daaf66cSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
273*d64ed03dSBarry Smith   int                 ierr;
2746daaf66cSBarry Smith 
2753a40ed3dSBarry Smith   PetscFunctionBegin;
27677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2776daaf66cSBarry Smith 
2786daaf66cSBarry Smith   PetscStrcpy(p,"-");
2796daaf66cSBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2806daaf66cSBarry Smith 
2816daaf66cSBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2826daaf66cSBarry Smith 
283d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
28433455573SSatish Balay   SNESPrintTypes_Private(snes->comm,p,"snes_type");
28577c4ece6SBarry Smith   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
286b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
287b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
288b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
289b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
290b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
291b18e04deSLois Curfman McInnes   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
2925bbad29bSBarry Smith   PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n");
2935bbad29bSBarry Smith   PetscPrintf(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
29483e2fdc7SBarry Smith   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
29583e2fdc7SBarry Smith     residual norm at each iteration.\n",p);
29683e2fdc7SBarry Smith   PetscPrintf(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
29783e2fdc7SBarry Smith     residual norm for small residual norms. This is useful to conceal\n\
29883e2fdc7SBarry Smith     meaningless digits that may be different on different machines.\n",p);
299d31a3109SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
300b18e04deSLois Curfman McInnes   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
30177c4ece6SBarry Smith     PetscPrintf(snes->comm,
302d31a3109SLois Curfman McInnes      " Options for solving systems of nonlinear equations only:\n");
30377c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
30477c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
305ef1dfb11SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
30693c39befSBarry Smith     PetscPrintf(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
30777c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
30877c4ece6SBarry Smith     PetscPrintf(snes->comm,
309b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
31077c4ece6SBarry Smith     PetscPrintf(snes->comm,
311b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
31277c4ece6SBarry Smith     PetscPrintf(snes->comm,
313b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
31477c4ece6SBarry Smith     PetscPrintf(snes->comm,
315b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
31677c4ece6SBarry Smith     PetscPrintf(snes->comm,
317b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
31877c4ece6SBarry Smith     PetscPrintf(snes->comm,
319b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
32077c4ece6SBarry Smith     PetscPrintf(snes->comm,
321b18e04deSLois Curfman McInnes      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
322*d64ed03dSBarry Smith   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
32377c4ece6SBarry Smith     PetscPrintf(snes->comm,
324d31a3109SLois Curfman McInnes      " Options for solving unconstrained minimization problems only:\n");
325b18e04deSLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
32677c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
32777c4ece6SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
328b18e04deSLois Curfman McInnes   }
3294537a946SLois Curfman McInnes   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
33077c4ece6SBarry Smith   PetscPrintf(snes->comm,"a particular method\n");
331*d64ed03dSBarry Smith   if (snes->printhelp) {
332*d64ed03dSBarry Smith     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
333*d64ed03dSBarry Smith   }
3343a40ed3dSBarry Smith   PetscFunctionReturn(0);
3359b94acceSBarry Smith }
336a847f771SSatish Balay 
3375615d1e5SSatish Balay #undef __FUNC__
338d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext"
3399b94acceSBarry Smith /*@
3409b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3419b94acceSBarry Smith    the nonlinear solvers.
3429b94acceSBarry Smith 
3439b94acceSBarry Smith    Input Parameters:
3449b94acceSBarry Smith .  snes - the SNES context
3459b94acceSBarry Smith .  usrP - optional user context
3469b94acceSBarry Smith 
3479b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3489b94acceSBarry Smith 
3499b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3509b94acceSBarry Smith @*/
3519b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3529b94acceSBarry Smith {
3533a40ed3dSBarry Smith   PetscFunctionBegin;
35477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3559b94acceSBarry Smith   snes->user		= usrP;
3563a40ed3dSBarry Smith   PetscFunctionReturn(0);
3579b94acceSBarry Smith }
35874679c65SBarry Smith 
3595615d1e5SSatish Balay #undef __FUNC__
360d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext"
3619b94acceSBarry Smith /*@C
3629b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3639b94acceSBarry Smith    nonlinear solvers.
3649b94acceSBarry Smith 
3659b94acceSBarry Smith    Input Parameter:
3669b94acceSBarry Smith .  snes - SNES context
3679b94acceSBarry Smith 
3689b94acceSBarry Smith    Output Parameter:
3699b94acceSBarry Smith .  usrP - user context
3709b94acceSBarry Smith 
3719b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3729b94acceSBarry Smith 
3739b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3749b94acceSBarry Smith @*/
3759b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3769b94acceSBarry Smith {
3773a40ed3dSBarry Smith   PetscFunctionBegin;
37877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3799b94acceSBarry Smith   *usrP = snes->user;
3803a40ed3dSBarry Smith   PetscFunctionReturn(0);
3819b94acceSBarry Smith }
38274679c65SBarry Smith 
3835615d1e5SSatish Balay #undef __FUNC__
384d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber"
3859b94acceSBarry Smith /*@
3869b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3879b94acceSBarry Smith    nonlinear solver.
3889b94acceSBarry Smith 
3899b94acceSBarry Smith    Input Parameter:
3909b94acceSBarry Smith .  snes - SNES context
3919b94acceSBarry Smith 
3929b94acceSBarry Smith    Output Parameter:
3939b94acceSBarry Smith .  iter - iteration number
3949b94acceSBarry Smith 
3959b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3969b94acceSBarry Smith @*/
3979b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3989b94acceSBarry Smith {
3993a40ed3dSBarry Smith   PetscFunctionBegin;
40077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
40174679c65SBarry Smith   PetscValidIntPointer(iter);
4029b94acceSBarry Smith   *iter = snes->iter;
4033a40ed3dSBarry Smith   PetscFunctionReturn(0);
4049b94acceSBarry Smith }
40574679c65SBarry Smith 
4065615d1e5SSatish Balay #undef __FUNC__
4075615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
4089b94acceSBarry Smith /*@
4099b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
4109b94acceSBarry Smith    with SNESSSetFunction().
4119b94acceSBarry Smith 
4129b94acceSBarry Smith    Input Parameter:
4139b94acceSBarry Smith .  snes - SNES context
4149b94acceSBarry Smith 
4159b94acceSBarry Smith    Output Parameter:
4169b94acceSBarry Smith .  fnorm - 2-norm of function
4179b94acceSBarry Smith 
4189b94acceSBarry Smith    Note:
4199b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
420a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
421a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
4229b94acceSBarry Smith 
4239b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
424a86d99e1SLois Curfman McInnes 
425a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
4269b94acceSBarry Smith @*/
4279b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
4289b94acceSBarry Smith {
4293a40ed3dSBarry Smith   PetscFunctionBegin;
43077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43174679c65SBarry Smith   PetscValidScalarPointer(fnorm);
43274679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
433d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only");
43474679c65SBarry Smith   }
4359b94acceSBarry Smith   *fnorm = snes->norm;
4363a40ed3dSBarry Smith   PetscFunctionReturn(0);
4379b94acceSBarry Smith }
43874679c65SBarry Smith 
4395615d1e5SSatish Balay #undef __FUNC__
4405615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
4419b94acceSBarry Smith /*@
4429b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
4439b94acceSBarry Smith    with SNESSSetGradient().
4449b94acceSBarry Smith 
4459b94acceSBarry Smith    Input Parameter:
4469b94acceSBarry Smith .  snes - SNES context
4479b94acceSBarry Smith 
4489b94acceSBarry Smith    Output Parameter:
4499b94acceSBarry Smith .  fnorm - 2-norm of gradient
4509b94acceSBarry Smith 
4519b94acceSBarry Smith    Note:
4529b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
453a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
454a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4559b94acceSBarry Smith 
4569b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
457a86d99e1SLois Curfman McInnes 
458a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4599b94acceSBarry Smith @*/
4609b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
4619b94acceSBarry Smith {
4623a40ed3dSBarry Smith   PetscFunctionBegin;
46377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
46474679c65SBarry Smith   PetscValidScalarPointer(gnorm);
46574679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
466d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
46774679c65SBarry Smith   }
4689b94acceSBarry Smith   *gnorm = snes->norm;
4693a40ed3dSBarry Smith   PetscFunctionReturn(0);
4709b94acceSBarry Smith }
47174679c65SBarry Smith 
4725615d1e5SSatish Balay #undef __FUNC__
473d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps"
4749b94acceSBarry Smith /*@
4759b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4769b94acceSBarry Smith    attempted by the nonlinear solver.
4779b94acceSBarry Smith 
4789b94acceSBarry Smith    Input Parameter:
4799b94acceSBarry Smith .  snes - SNES context
4809b94acceSBarry Smith 
4819b94acceSBarry Smith    Output Parameter:
4829b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4839b94acceSBarry Smith 
484c96a6f78SLois Curfman McInnes    Notes:
485c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
486c96a6f78SLois Curfman McInnes 
4879b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4889b94acceSBarry Smith @*/
4899b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4909b94acceSBarry Smith {
4913a40ed3dSBarry Smith   PetscFunctionBegin;
49277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
49374679c65SBarry Smith   PetscValidIntPointer(nfails);
4949b94acceSBarry Smith   *nfails = snes->nfailures;
4953a40ed3dSBarry Smith   PetscFunctionReturn(0);
4969b94acceSBarry Smith }
497a847f771SSatish Balay 
4985615d1e5SSatish Balay #undef __FUNC__
499d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations"
500c96a6f78SLois Curfman McInnes /*@
501c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
502c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
503c96a6f78SLois Curfman McInnes 
504c96a6f78SLois Curfman McInnes    Input Parameter:
505c96a6f78SLois Curfman McInnes .  snes - SNES context
506c96a6f78SLois Curfman McInnes 
507c96a6f78SLois Curfman McInnes    Output Parameter:
508c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
509c96a6f78SLois Curfman McInnes 
510c96a6f78SLois Curfman McInnes    Notes:
511c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
512c96a6f78SLois Curfman McInnes 
513c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
514c96a6f78SLois Curfman McInnes @*/
51586bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
516c96a6f78SLois Curfman McInnes {
5173a40ed3dSBarry Smith   PetscFunctionBegin;
518c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
519c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
520c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
5213a40ed3dSBarry Smith   PetscFunctionReturn(0);
522c96a6f78SLois Curfman McInnes }
523c96a6f78SLois Curfman McInnes 
524c96a6f78SLois Curfman McInnes #undef __FUNC__
525d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES"
5269b94acceSBarry Smith /*@C
5279b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
5289b94acceSBarry Smith 
5299b94acceSBarry Smith    Input Parameter:
5309b94acceSBarry Smith .  snes - the SNES context
5319b94acceSBarry Smith 
5329b94acceSBarry Smith    Output Parameter:
5339b94acceSBarry Smith .  sles - the SLES context
5349b94acceSBarry Smith 
5359b94acceSBarry Smith    Notes:
5369b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
5379b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5389b94acceSBarry Smith    KSP and PC contexts as well.
5399b94acceSBarry Smith 
5409b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
5419b94acceSBarry Smith 
5429b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
5439b94acceSBarry Smith @*/
5449b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
5459b94acceSBarry Smith {
5463a40ed3dSBarry Smith   PetscFunctionBegin;
54777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5489b94acceSBarry Smith   *sles = snes->sles;
5493a40ed3dSBarry Smith   PetscFunctionReturn(0);
5509b94acceSBarry Smith }
5519b94acceSBarry Smith /* -----------------------------------------------------------*/
5525615d1e5SSatish Balay #undef __FUNC__
5535615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
5549b94acceSBarry Smith /*@C
5559b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5569b94acceSBarry Smith 
5579b94acceSBarry Smith    Input Parameter:
5589b94acceSBarry Smith .  comm - MPI communicator
5599b94acceSBarry Smith .  type - type of method, one of
5609b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
5619b94acceSBarry Smith $      (for systems of nonlinear equations)
5629b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
5639b94acceSBarry Smith $      (for unconstrained minimization)
5649b94acceSBarry Smith 
5659b94acceSBarry Smith    Output Parameter:
5669b94acceSBarry Smith .  outsnes - the new SNES context
5679b94acceSBarry Smith 
568c1f60f51SBarry Smith    Options Database Key:
56919bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
57019bd6747SLois Curfman McInnes $              and no preconditioning matrix
57119bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
57219bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
57319bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
57419bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
575c1f60f51SBarry Smith 
5769b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5779b94acceSBarry Smith 
57863a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5799b94acceSBarry Smith @*/
5804b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5819b94acceSBarry Smith {
5829b94acceSBarry Smith   int                 ierr;
5839b94acceSBarry Smith   SNES                snes;
5849b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
58537fcc0dbSBarry Smith 
5863a40ed3dSBarry Smith   PetscFunctionBegin;
5879b94acceSBarry Smith   *outsnes = 0;
588*d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
589d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
590*d64ed03dSBarry Smith   }
591d4bb536fSBarry Smith   PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm,SNESDestroy,SNESView);
5929b94acceSBarry Smith   PLogObjectCreate(snes);
5939b94acceSBarry Smith   snes->max_its           = 50;
5949750a799SBarry Smith   snes->max_funcs	  = 10000;
5959b94acceSBarry Smith   snes->norm		  = 0.0;
5965a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
597b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
598b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5999b94acceSBarry Smith     snes->atol		  = 1.e-10;
600*d64ed03dSBarry Smith   } else {
601b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
602b4874afaSBarry Smith     snes->ttol            = 0.0;
603b4874afaSBarry Smith     snes->atol		  = 1.e-50;
604b4874afaSBarry Smith   }
6059b94acceSBarry Smith   snes->xtol		  = 1.e-8;
606b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
6079b94acceSBarry Smith   snes->nfuncs            = 0;
6089b94acceSBarry Smith   snes->nfailures         = 0;
6097a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
610639f9d9dSBarry Smith   snes->numbermonitors    = 0;
6119b94acceSBarry Smith   snes->data              = 0;
6129b94acceSBarry Smith   snes->view              = 0;
6139b94acceSBarry Smith   snes->computeumfunction = 0;
6149b94acceSBarry Smith   snes->umfunP            = 0;
6159b94acceSBarry Smith   snes->fc                = 0;
6169b94acceSBarry Smith   snes->deltatol          = 1.e-12;
6179b94acceSBarry Smith   snes->fmin              = -1.e30;
6189b94acceSBarry Smith   snes->method_class      = type;
6199b94acceSBarry Smith   snes->set_method_called = 0;
620a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
6219b94acceSBarry Smith   snes->ksp_ewconv        = 0;
6227d1a2b2bSBarry Smith   snes->type              = -1;
6236f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
6246f24a144SLois Curfman McInnes   snes->vwork             = 0;
6256f24a144SLois Curfman McInnes   snes->nwork             = 0;
626c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
627c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
628c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
6299b94acceSBarry Smith 
6309b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
6310452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
632eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
6339b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6349b94acceSBarry Smith   kctx->version     = 2;
6359b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6369b94acceSBarry Smith                              this was too large for some test cases */
6379b94acceSBarry Smith   kctx->rtol_last   = 0;
6389b94acceSBarry Smith   kctx->rtol_max    = .9;
6399b94acceSBarry Smith   kctx->gamma       = 1.0;
6409b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6419b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6429b94acceSBarry Smith   kctx->threshold   = .1;
6439b94acceSBarry Smith   kctx->lresid_last = 0;
6449b94acceSBarry Smith   kctx->norm_last   = 0;
6459b94acceSBarry Smith 
6469b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
6479b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
6485334005bSBarry Smith 
6499b94acceSBarry Smith   *outsnes = snes;
6503a40ed3dSBarry Smith   PetscFunctionReturn(0);
6519b94acceSBarry Smith }
6529b94acceSBarry Smith 
6539b94acceSBarry Smith /* --------------------------------------------------------------- */
6545615d1e5SSatish Balay #undef __FUNC__
655d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction"
6569b94acceSBarry Smith /*@C
6579b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6589b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6599b94acceSBarry Smith    equations.
6609b94acceSBarry Smith 
6619b94acceSBarry Smith    Input Parameters:
6629b94acceSBarry Smith .  snes - the SNES context
6639b94acceSBarry Smith .  func - function evaluation routine
6649b94acceSBarry Smith .  r - vector to store function value
6652cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
6662cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6679b94acceSBarry Smith 
6689b94acceSBarry Smith    Calling sequence of func:
669313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
6709b94acceSBarry Smith 
6719b94acceSBarry Smith .  x - input vector
672313e4042SLois Curfman McInnes .  f - function vector
6732cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
6749b94acceSBarry Smith 
6759b94acceSBarry Smith    Notes:
6769b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6779b94acceSBarry Smith $      f'(x) x = -f(x),
6789b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
6799b94acceSBarry Smith 
6809b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6819b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6829b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6839b94acceSBarry Smith 
6849b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6859b94acceSBarry Smith 
686a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6879b94acceSBarry Smith @*/
6885334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6899b94acceSBarry Smith {
6903a40ed3dSBarry Smith   PetscFunctionBegin;
69177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
692e3372554SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
6939b94acceSBarry Smith   snes->computefunction     = func;
6949b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6959b94acceSBarry Smith   snes->funP                = ctx;
6963a40ed3dSBarry Smith   PetscFunctionReturn(0);
6979b94acceSBarry Smith }
6989b94acceSBarry Smith 
6995615d1e5SSatish Balay #undef __FUNC__
7005615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
7019b94acceSBarry Smith /*@
7029b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
7039b94acceSBarry Smith    SNESSetFunction().
7049b94acceSBarry Smith 
7059b94acceSBarry Smith    Input Parameters:
7069b94acceSBarry Smith .  snes - the SNES context
7079b94acceSBarry Smith .  x - input vector
7089b94acceSBarry Smith 
7099b94acceSBarry Smith    Output Parameter:
7103638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
7119b94acceSBarry Smith 
7121bffabb2SLois Curfman McInnes    Notes:
7139b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
7149b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
7159b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
7169b94acceSBarry Smith 
7179b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7189b94acceSBarry Smith 
719a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7209b94acceSBarry Smith @*/
7219b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
7229b94acceSBarry Smith {
7239b94acceSBarry Smith   int    ierr;
7249b94acceSBarry Smith 
7253a40ed3dSBarry Smith   PetscFunctionBegin;
726d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
727e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
728d4bb536fSBarry Smith   }
7299b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
730*d64ed03dSBarry Smith   PetscStackPush("SNES user function");
7319b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
732*d64ed03dSBarry Smith   PetscStackPop;
733ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7349b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
7353a40ed3dSBarry Smith   PetscFunctionReturn(0);
7369b94acceSBarry Smith }
7379b94acceSBarry Smith 
7385615d1e5SSatish Balay #undef __FUNC__
739d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction"
7409b94acceSBarry Smith /*@C
7419b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
7429b94acceSBarry Smith    unconstrained minimization.
7439b94acceSBarry Smith 
7449b94acceSBarry Smith    Input Parameters:
7459b94acceSBarry Smith .  snes - the SNES context
7469b94acceSBarry Smith .  func - function evaluation routine
747044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
748044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7499b94acceSBarry Smith 
7509b94acceSBarry Smith    Calling sequence of func:
7519b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
7529b94acceSBarry Smith 
7539b94acceSBarry Smith .  x - input vector
7549b94acceSBarry Smith .  f - function
755044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
7569b94acceSBarry Smith 
7579b94acceSBarry Smith    Notes:
7589b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7599b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7609b94acceSBarry Smith    SNESSetFunction().
7619b94acceSBarry Smith 
7629b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
7639b94acceSBarry Smith 
764a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
765a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
7669b94acceSBarry Smith @*/
7679b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
7689b94acceSBarry Smith                       void *ctx)
7699b94acceSBarry Smith {
7703a40ed3dSBarry Smith   PetscFunctionBegin;
77177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
772e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
7739b94acceSBarry Smith   snes->computeumfunction   = func;
7749b94acceSBarry Smith   snes->umfunP              = ctx;
7753a40ed3dSBarry Smith   PetscFunctionReturn(0);
7769b94acceSBarry Smith }
7779b94acceSBarry Smith 
7785615d1e5SSatish Balay #undef __FUNC__
7795615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7809b94acceSBarry Smith /*@
7819b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7829b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7839b94acceSBarry Smith 
7849b94acceSBarry Smith    Input Parameters:
7859b94acceSBarry Smith .  snes - the SNES context
7869b94acceSBarry Smith .  x - input vector
7879b94acceSBarry Smith 
7889b94acceSBarry Smith    Output Parameter:
7899b94acceSBarry Smith .  y - function value
7909b94acceSBarry Smith 
7919b94acceSBarry Smith    Notes:
7929b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7939b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7949b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
795a86d99e1SLois Curfman McInnes 
796a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
797a86d99e1SLois Curfman McInnes 
798a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
799a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
8009b94acceSBarry Smith @*/
8019b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
8029b94acceSBarry Smith {
8039b94acceSBarry Smith   int    ierr;
8043a40ed3dSBarry Smith 
8053a40ed3dSBarry Smith   PetscFunctionBegin;
806e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
8079b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
808*d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
8099b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
810*d64ed03dSBarry Smith   PetscStackPop;
811ae3c334cSLois Curfman McInnes   snes->nfuncs++;
8129b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
8133a40ed3dSBarry Smith   PetscFunctionReturn(0);
8149b94acceSBarry Smith }
8159b94acceSBarry Smith 
8165615d1e5SSatish Balay #undef __FUNC__
817d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient"
8189b94acceSBarry Smith /*@C
8199b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
8209b94acceSBarry Smith    vector for use by the SNES routines.
8219b94acceSBarry Smith 
8229b94acceSBarry Smith    Input Parameters:
8239b94acceSBarry Smith .  snes - the SNES context
8249b94acceSBarry Smith .  func - function evaluation routine
825044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
826044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
8279b94acceSBarry Smith .  r - vector to store gradient value
8289b94acceSBarry Smith 
8299b94acceSBarry Smith    Calling sequence of func:
8309b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
8319b94acceSBarry Smith 
8329b94acceSBarry Smith .  x - input vector
8339b94acceSBarry Smith .  g - gradient vector
834044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
8359b94acceSBarry Smith 
8369b94acceSBarry Smith    Notes:
8379b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
8389b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
8399b94acceSBarry Smith    SNESSetFunction().
8409b94acceSBarry Smith 
8419b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
8429b94acceSBarry Smith 
843a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
844a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
8459b94acceSBarry Smith @*/
84674679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
8479b94acceSBarry Smith {
8483a40ed3dSBarry Smith   PetscFunctionBegin;
84977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
850e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8519b94acceSBarry Smith   snes->computefunction     = func;
8529b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
8539b94acceSBarry Smith   snes->funP                = ctx;
8543a40ed3dSBarry Smith   PetscFunctionReturn(0);
8559b94acceSBarry Smith }
8569b94acceSBarry Smith 
8575615d1e5SSatish Balay #undef __FUNC__
8585615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
8599b94acceSBarry Smith /*@
860a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
861a86d99e1SLois Curfman McInnes    SNESSetGradient().
8629b94acceSBarry Smith 
8639b94acceSBarry Smith    Input Parameters:
8649b94acceSBarry Smith .  snes - the SNES context
8659b94acceSBarry Smith .  x - input vector
8669b94acceSBarry Smith 
8679b94acceSBarry Smith    Output Parameter:
8689b94acceSBarry Smith .  y - gradient vector
8699b94acceSBarry Smith 
8709b94acceSBarry Smith    Notes:
8719b94acceSBarry Smith    SNESComputeGradient() is valid only for
8729b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8739b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
874a86d99e1SLois Curfman McInnes 
875a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
876a86d99e1SLois Curfman McInnes 
877a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
878a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8799b94acceSBarry Smith @*/
8809b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8819b94acceSBarry Smith {
8829b94acceSBarry Smith   int    ierr;
8833a40ed3dSBarry Smith 
8843a40ed3dSBarry Smith   PetscFunctionBegin;
8853a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
886e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8873a40ed3dSBarry Smith   }
8889b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
889*d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
8909b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
891*d64ed03dSBarry Smith   PetscStackPop;
8929b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8933a40ed3dSBarry Smith   PetscFunctionReturn(0);
8949b94acceSBarry Smith }
8959b94acceSBarry Smith 
8965615d1e5SSatish Balay #undef __FUNC__
8975615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
89862fef451SLois Curfman McInnes /*@
89962fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
90062fef451SLois Curfman McInnes    set with SNESSetJacobian().
90162fef451SLois Curfman McInnes 
90262fef451SLois Curfman McInnes    Input Parameters:
90362fef451SLois Curfman McInnes .  snes - the SNES context
90462fef451SLois Curfman McInnes .  x - input vector
90562fef451SLois Curfman McInnes 
90662fef451SLois Curfman McInnes    Output Parameters:
90762fef451SLois Curfman McInnes .  A - Jacobian matrix
90862fef451SLois Curfman McInnes .  B - optional preconditioning matrix
90962fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
91062fef451SLois Curfman McInnes 
91162fef451SLois Curfman McInnes    Notes:
91262fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
91362fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
91462fef451SLois Curfman McInnes 
915dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
916dc5a77f8SLois Curfman McInnes    flag parameter.
91762fef451SLois Curfman McInnes 
91862fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
91962fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
920005c665bSBarry Smith    methods is SNESComputeHessian().
92162fef451SLois Curfman McInnes 
92262fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
92362fef451SLois Curfman McInnes 
92462fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
92562fef451SLois Curfman McInnes @*/
9269b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
9279b94acceSBarry Smith {
9289b94acceSBarry Smith   int    ierr;
9293a40ed3dSBarry Smith 
9303a40ed3dSBarry Smith   PetscFunctionBegin;
9313a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
932e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
9333a40ed3dSBarry Smith   }
9343a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
9359b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
936c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
937*d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
9389b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
939*d64ed03dSBarry Smith   PetscStackPop;
9409b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
9416d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
94277c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
94377c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9443a40ed3dSBarry Smith   PetscFunctionReturn(0);
9459b94acceSBarry Smith }
9469b94acceSBarry Smith 
9475615d1e5SSatish Balay #undef __FUNC__
9485615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
94962fef451SLois Curfman McInnes /*@
95062fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
95162fef451SLois Curfman McInnes    set with SNESSetHessian().
95262fef451SLois Curfman McInnes 
95362fef451SLois Curfman McInnes    Input Parameters:
95462fef451SLois Curfman McInnes .  snes - the SNES context
95562fef451SLois Curfman McInnes .  x - input vector
95662fef451SLois Curfman McInnes 
95762fef451SLois Curfman McInnes    Output Parameters:
95862fef451SLois Curfman McInnes .  A - Hessian matrix
95962fef451SLois Curfman McInnes .  B - optional preconditioning matrix
96062fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
96162fef451SLois Curfman McInnes 
96262fef451SLois Curfman McInnes    Notes:
96362fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
96462fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
96562fef451SLois Curfman McInnes 
966dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
967dc5a77f8SLois Curfman McInnes    flag parameter.
96862fef451SLois Curfman McInnes 
96962fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
97062fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
97162fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
97262fef451SLois Curfman McInnes 
97362fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
97462fef451SLois Curfman McInnes 
975a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
976a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
97762fef451SLois Curfman McInnes @*/
97862fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
9799b94acceSBarry Smith {
9809b94acceSBarry Smith   int    ierr;
9813a40ed3dSBarry Smith 
9823a40ed3dSBarry Smith   PetscFunctionBegin;
9833a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
984e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9853a40ed3dSBarry Smith   }
9863a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
98762fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9880f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
989*d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
99062fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
991*d64ed03dSBarry Smith   PetscStackPop;
99262fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9930f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
99477c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
99577c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9963a40ed3dSBarry Smith   PetscFunctionReturn(0);
9979b94acceSBarry Smith }
9989b94acceSBarry Smith 
9995615d1e5SSatish Balay #undef __FUNC__
1000d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian"
10019b94acceSBarry Smith /*@C
10029b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1003044dda88SLois Curfman McInnes    location to store the matrix.
10049b94acceSBarry Smith 
10059b94acceSBarry Smith    Input Parameters:
10069b94acceSBarry Smith .  snes - the SNES context
10079b94acceSBarry Smith .  A - Jacobian matrix
10089b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
10099b94acceSBarry Smith .  func - Jacobian evaluation routine
10102cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
10112cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
10129b94acceSBarry Smith 
10139b94acceSBarry Smith    Calling sequence of func:
1014313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10159b94acceSBarry Smith 
10169b94acceSBarry Smith .  x - input vector
10179b94acceSBarry Smith .  A - Jacobian matrix
10189b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1019ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1020ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10212cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
10229b94acceSBarry Smith 
10239b94acceSBarry Smith    Notes:
1024dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10252cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1026ac21db08SLois Curfman McInnes 
1027ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
10289b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
10299b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10309b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10319b94acceSBarry Smith    throughout the global iterations.
10329b94acceSBarry Smith 
10339b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
10349b94acceSBarry Smith 
1035ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
10369b94acceSBarry Smith @*/
10379b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10389b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10399b94acceSBarry Smith {
10403a40ed3dSBarry Smith   PetscFunctionBegin;
104177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1042d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
10439b94acceSBarry Smith   snes->computejacobian = func;
10449b94acceSBarry Smith   snes->jacP            = ctx;
10459b94acceSBarry Smith   snes->jacobian        = A;
10469b94acceSBarry Smith   snes->jacobian_pre    = B;
10473a40ed3dSBarry Smith   PetscFunctionReturn(0);
10489b94acceSBarry Smith }
104962fef451SLois Curfman McInnes 
10505615d1e5SSatish Balay #undef __FUNC__
1051d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian"
1052b4fd4287SBarry Smith /*@
1053b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1054b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1055b4fd4287SBarry Smith 
1056b4fd4287SBarry Smith    Input Parameter:
1057b4fd4287SBarry Smith .  snes - the nonlinear solver context
1058b4fd4287SBarry Smith 
1059b4fd4287SBarry Smith    Output Parameters:
1060b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
1061b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
1062b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1063b4fd4287SBarry Smith 
1064b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1065b4fd4287SBarry Smith @*/
1066b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
1067b4fd4287SBarry Smith {
10683a40ed3dSBarry Smith   PetscFunctionBegin;
10693a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1070e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
10713a40ed3dSBarry Smith   }
1072b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
1073b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
1074b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
10753a40ed3dSBarry Smith   PetscFunctionReturn(0);
1076b4fd4287SBarry Smith }
1077b4fd4287SBarry Smith 
10785615d1e5SSatish Balay #undef __FUNC__
1079d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian"
10809b94acceSBarry Smith /*@C
10819b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1082044dda88SLois Curfman McInnes    location to store the matrix.
10839b94acceSBarry Smith 
10849b94acceSBarry Smith    Input Parameters:
10859b94acceSBarry Smith .  snes - the SNES context
10869b94acceSBarry Smith .  A - Hessian matrix
10879b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
10889b94acceSBarry Smith .  func - Jacobian evaluation routine
1089313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1090313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10919b94acceSBarry Smith 
10929b94acceSBarry Smith    Calling sequence of func:
1093313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10949b94acceSBarry Smith 
10959b94acceSBarry Smith .  x - input vector
10969b94acceSBarry Smith .  A - Hessian matrix
10979b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1098ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1099ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
11002cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
11019b94acceSBarry Smith 
11029b94acceSBarry Smith    Notes:
1103dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
11042cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1105ac21db08SLois Curfman McInnes 
11069b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
11079b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
11089b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
11099b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
11109b94acceSBarry Smith    throughout the global iterations.
11119b94acceSBarry Smith 
11129b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
11139b94acceSBarry Smith 
1114ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
11159b94acceSBarry Smith @*/
11169b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
11179b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
11189b94acceSBarry Smith {
11193a40ed3dSBarry Smith   PetscFunctionBegin;
112077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1121d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1122e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1123d4bb536fSBarry Smith   }
11249b94acceSBarry Smith   snes->computejacobian = func;
11259b94acceSBarry Smith   snes->jacP            = ctx;
11269b94acceSBarry Smith   snes->jacobian        = A;
11279b94acceSBarry Smith   snes->jacobian_pre    = B;
11283a40ed3dSBarry Smith   PetscFunctionReturn(0);
11299b94acceSBarry Smith }
11309b94acceSBarry Smith 
11315615d1e5SSatish Balay #undef __FUNC__
1132d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian"
113362fef451SLois Curfman McInnes /*@
113462fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
113562fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
113662fef451SLois Curfman McInnes 
113762fef451SLois Curfman McInnes    Input Parameter:
113862fef451SLois Curfman McInnes .  snes - the nonlinear solver context
113962fef451SLois Curfman McInnes 
114062fef451SLois Curfman McInnes    Output Parameters:
114162fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
114262fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
114362fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
114462fef451SLois Curfman McInnes 
114562fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
114662fef451SLois Curfman McInnes @*/
114762fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
114862fef451SLois Curfman McInnes {
11493a40ed3dSBarry Smith   PetscFunctionBegin;
11503a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1151e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
11523a40ed3dSBarry Smith   }
115362fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
115462fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
115562fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
11563a40ed3dSBarry Smith   PetscFunctionReturn(0);
115762fef451SLois Curfman McInnes }
115862fef451SLois Curfman McInnes 
11599b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
11609b94acceSBarry Smith 
11615615d1e5SSatish Balay #undef __FUNC__
11625615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
11639b94acceSBarry Smith /*@
11649b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1165272ac6f2SLois Curfman McInnes    of a nonlinear solver.
11669b94acceSBarry Smith 
11679b94acceSBarry Smith    Input Parameter:
11689b94acceSBarry Smith .  snes - the SNES context
11698ddd3da0SLois Curfman McInnes .  x - the solution vector
11709b94acceSBarry Smith 
1171272ac6f2SLois Curfman McInnes    Notes:
1172272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1173272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1174272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1175272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1176272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1177272ac6f2SLois Curfman McInnes 
11789b94acceSBarry Smith .keywords: SNES, nonlinear, setup
11799b94acceSBarry Smith 
11809b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
11819b94acceSBarry Smith @*/
11828ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
11839b94acceSBarry Smith {
1184272ac6f2SLois Curfman McInnes   int ierr, flg;
11853a40ed3dSBarry Smith 
11863a40ed3dSBarry Smith   PetscFunctionBegin;
118777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
118877c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
11898ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
11909b94acceSBarry Smith 
1191c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1192c1f60f51SBarry Smith   /*
1193c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1194dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1195c1f60f51SBarry Smith   */
1196c1f60f51SBarry Smith   if (flg) {
1197c1f60f51SBarry Smith     Mat J;
1198c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1199c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1200c1f60f51SBarry Smith     snes->mfshell = J;
1201c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1202c1f60f51SBarry Smith       snes->jacobian = J;
1203c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1204*d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1205c1f60f51SBarry Smith       snes->jacobian = J;
1206c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1207d4bb536fSBarry Smith     } else {
1208d4bb536fSBarry Smith       SETERRQ(1,0,"Method class doesn't support matrix-free operator option");
1209d4bb536fSBarry Smith     }
1210c1f60f51SBarry Smith   }
1211272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1212c1f60f51SBarry Smith   /*
1213dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1214c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1215c1f60f51SBarry Smith    */
121631615d04SLois Curfman McInnes   if (flg) {
1217272ac6f2SLois Curfman McInnes     Mat J;
1218272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1219272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1220272ac6f2SLois Curfman McInnes     snes->mfshell = J;
122183e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
122283e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
122394a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1224*d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
122583e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
122694a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1227d4bb536fSBarry Smith     } else {
1228d4bb536fSBarry Smith       SETERRQ(1,0,"Method class doesn't support matrix-free option");
1229d4bb536fSBarry Smith     }
1230272ac6f2SLois Curfman McInnes   }
12319b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1232e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first");
1233e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first");
1234e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first");
1235e3372554SBarry Smith     if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector");
1236a703fe33SLois Curfman McInnes 
1237a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
123840191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1239a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1240a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1241a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1242a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1243a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1244a703fe33SLois Curfman McInnes     }
12459b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1246e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first");
1247e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first");
1248d4bb536fSBarry Smith     if (!snes->computeumfunction) SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first");
1249e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first");
1250d4bb536fSBarry Smith   } else {
1251d4bb536fSBarry Smith     SETERRQ(1,0,"Unknown method class");
1252d4bb536fSBarry Smith   }
1253a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1254a703fe33SLois Curfman McInnes   snes->setup_called = 1;
12553a40ed3dSBarry Smith   PetscFunctionReturn(0);
12569b94acceSBarry Smith }
12579b94acceSBarry Smith 
12585615d1e5SSatish Balay #undef __FUNC__
1259d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy"
12609b94acceSBarry Smith /*@C
12619b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
12629b94acceSBarry Smith    with SNESCreate().
12639b94acceSBarry Smith 
12649b94acceSBarry Smith    Input Parameter:
12659b94acceSBarry Smith .  snes - the SNES context
12669b94acceSBarry Smith 
12679b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
12689b94acceSBarry Smith 
126963a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
12709b94acceSBarry Smith @*/
12719b94acceSBarry Smith int SNESDestroy(SNES snes)
12729b94acceSBarry Smith {
12739b94acceSBarry Smith   int ierr;
12743a40ed3dSBarry Smith 
12753a40ed3dSBarry Smith   PetscFunctionBegin;
127677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12773a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1278d4bb536fSBarry Smith 
12799750a799SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);}
12800452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
12819b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
12829b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
12833e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
12846f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
12859b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
12860452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
12873a40ed3dSBarry Smith   PetscFunctionReturn(0);
12889b94acceSBarry Smith }
12899b94acceSBarry Smith 
12909b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12919b94acceSBarry Smith 
12925615d1e5SSatish Balay #undef __FUNC__
12935615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
12949b94acceSBarry Smith /*@
1295d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12969b94acceSBarry Smith 
12979b94acceSBarry Smith    Input Parameters:
12989b94acceSBarry Smith .  snes - the SNES context
129933174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
130033174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
130133174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
130233174efeSLois Curfman McInnes            of the change in the solution between steps
130333174efeSLois Curfman McInnes .  maxit - maximum number of iterations
130433174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
13059b94acceSBarry Smith 
130633174efeSLois Curfman McInnes    Options Database Keys:
130733174efeSLois Curfman McInnes $    -snes_atol <atol>
130833174efeSLois Curfman McInnes $    -snes_rtol <rtol>
130933174efeSLois Curfman McInnes $    -snes_stol <stol>
131033174efeSLois Curfman McInnes $    -snes_max_it <maxit>
131133174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
13129b94acceSBarry Smith 
1313d7a720efSLois Curfman McInnes    Notes:
13149b94acceSBarry Smith    The default maximum number of iterations is 50.
13159b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
13169b94acceSBarry Smith 
131733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
13189b94acceSBarry Smith 
1319d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
13209b94acceSBarry Smith @*/
1321d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
13229b94acceSBarry Smith {
13233a40ed3dSBarry Smith   PetscFunctionBegin;
132477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1325d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1326d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1327d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1328d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1329d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
13303a40ed3dSBarry Smith   PetscFunctionReturn(0);
13319b94acceSBarry Smith }
13329b94acceSBarry Smith 
13335615d1e5SSatish Balay #undef __FUNC__
13345615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
13359b94acceSBarry Smith /*@
133633174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
133733174efeSLois Curfman McInnes 
133833174efeSLois Curfman McInnes    Input Parameters:
133933174efeSLois Curfman McInnes .  snes - the SNES context
134033174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
134133174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
134233174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
134333174efeSLois Curfman McInnes            of the change in the solution between steps
134433174efeSLois Curfman McInnes .  maxit - maximum number of iterations
134533174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
134633174efeSLois Curfman McInnes 
134733174efeSLois Curfman McInnes    Notes:
134833174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
134933174efeSLois Curfman McInnes 
135033174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
135133174efeSLois Curfman McInnes 
135233174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
135333174efeSLois Curfman McInnes @*/
135433174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
135533174efeSLois Curfman McInnes {
13563a40ed3dSBarry Smith   PetscFunctionBegin;
135733174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
135833174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
135933174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
136033174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
136133174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
136233174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
13633a40ed3dSBarry Smith   PetscFunctionReturn(0);
136433174efeSLois Curfman McInnes }
136533174efeSLois Curfman McInnes 
13665615d1e5SSatish Balay #undef __FUNC__
13675615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
136833174efeSLois Curfman McInnes /*@
13699b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
13709b94acceSBarry Smith 
13719b94acceSBarry Smith    Input Parameters:
13729b94acceSBarry Smith .  snes - the SNES context
13739b94acceSBarry Smith .  tol - tolerance
13749b94acceSBarry Smith 
13759b94acceSBarry Smith    Options Database Key:
1376d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
13779b94acceSBarry Smith 
13789b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
13799b94acceSBarry Smith 
1380d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
13819b94acceSBarry Smith @*/
13829b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
13839b94acceSBarry Smith {
13843a40ed3dSBarry Smith   PetscFunctionBegin;
138577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13869b94acceSBarry Smith   snes->deltatol = tol;
13873a40ed3dSBarry Smith   PetscFunctionReturn(0);
13889b94acceSBarry Smith }
13899b94acceSBarry Smith 
13905615d1e5SSatish Balay #undef __FUNC__
13915615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
13929b94acceSBarry Smith /*@
139377c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
13949b94acceSBarry Smith    for unconstrained minimization solvers.
13959b94acceSBarry Smith 
13969b94acceSBarry Smith    Input Parameters:
13979b94acceSBarry Smith .  snes - the SNES context
13989b94acceSBarry Smith .  ftol - minimum function tolerance
13999b94acceSBarry Smith 
14009b94acceSBarry Smith    Options Database Key:
1401d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
14029b94acceSBarry Smith 
14039b94acceSBarry Smith    Note:
140477c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
14059b94acceSBarry Smith    methods only.
14069b94acceSBarry Smith 
14079b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
14089b94acceSBarry Smith 
1409d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
14109b94acceSBarry Smith @*/
141177c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
14129b94acceSBarry Smith {
14133a40ed3dSBarry Smith   PetscFunctionBegin;
141477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
14159b94acceSBarry Smith   snes->fmin = ftol;
14163a40ed3dSBarry Smith   PetscFunctionReturn(0);
14179b94acceSBarry Smith }
14189b94acceSBarry Smith 
14199b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
14209b94acceSBarry Smith 
14215615d1e5SSatish Balay #undef __FUNC__
1422d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor"
14239b94acceSBarry Smith /*@C
14249b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
14259b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
14269b94acceSBarry Smith    progress.
14279b94acceSBarry Smith 
14289b94acceSBarry Smith    Input Parameters:
14299b94acceSBarry Smith .  snes - the SNES context
14309b94acceSBarry Smith .  func - monitoring routine
1431044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1432044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
14339b94acceSBarry Smith 
14349b94acceSBarry Smith    Calling sequence of func:
1435682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
14369b94acceSBarry Smith 
14379b94acceSBarry Smith $    snes - the SNES context
14389b94acceSBarry Smith $    its - iteration number
1439044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
14409b94acceSBarry Smith $
14419b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14429b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
14439b94acceSBarry Smith $
14449b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14459b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
14469b94acceSBarry Smith 
14479665c990SLois Curfman McInnes    Options Database Keys:
14489665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
14499665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
14509665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
14519665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
14529665c990SLois Curfman McInnes $                           been hardwired into a code by
14539665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
14549665c990SLois Curfman McInnes $                           does not cancel those set via
14559665c990SLois Curfman McInnes $                           the options database.
14569665c990SLois Curfman McInnes 
14579665c990SLois Curfman McInnes 
1458639f9d9dSBarry Smith    Notes:
14596bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
14606bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
14616bc08f3fSLois Curfman McInnes    order in which they were set.
1462639f9d9dSBarry Smith 
14639b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
14649b94acceSBarry Smith 
14659b94acceSBarry Smith .seealso: SNESDefaultMonitor()
14669b94acceSBarry Smith @*/
146774679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
14689b94acceSBarry Smith {
14693a40ed3dSBarry Smith   PetscFunctionBegin;
1470639f9d9dSBarry Smith   if (!func) {
1471639f9d9dSBarry Smith     snes->numbermonitors = 0;
14723a40ed3dSBarry Smith     PetscFunctionReturn(0);
1473639f9d9dSBarry Smith   }
1474639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1475e3372554SBarry Smith     SETERRQ(1,0,"Too many monitors set");
1476639f9d9dSBarry Smith   }
1477639f9d9dSBarry Smith 
1478639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1479639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
14803a40ed3dSBarry Smith   PetscFunctionReturn(0);
14819b94acceSBarry Smith }
14829b94acceSBarry Smith 
14835615d1e5SSatish Balay #undef __FUNC__
1484d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest"
14859b94acceSBarry Smith /*@C
14869b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
14879b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
14889b94acceSBarry Smith 
14899b94acceSBarry Smith    Input Parameters:
14909b94acceSBarry Smith .  snes - the SNES context
14919b94acceSBarry Smith .  func - routine to test for convergence
1492044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1493044dda88SLois Curfman McInnes           (may be PETSC_NULL)
14949b94acceSBarry Smith 
14959b94acceSBarry Smith    Calling sequence of func:
14969b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
14979b94acceSBarry Smith              double f,void *cctx)
14989b94acceSBarry Smith 
14999b94acceSBarry Smith $    snes - the SNES context
1500044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
15019b94acceSBarry Smith $    xnorm - 2-norm of current iterate
15029b94acceSBarry Smith $
15039b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
15049b94acceSBarry Smith $    gnorm - 2-norm of current step
15059b94acceSBarry Smith $    f - 2-norm of function
15069b94acceSBarry Smith $
15079b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
15089b94acceSBarry Smith $    gnorm - 2-norm of current gradient
15099b94acceSBarry Smith $    f - function value
15109b94acceSBarry Smith 
15119b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
15129b94acceSBarry Smith 
151340191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
151440191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
15159b94acceSBarry Smith @*/
151674679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
15179b94acceSBarry Smith {
15183a40ed3dSBarry Smith   PetscFunctionBegin;
15199b94acceSBarry Smith   (snes)->converged = func;
15209b94acceSBarry Smith   (snes)->cnvP      = cctx;
15213a40ed3dSBarry Smith   PetscFunctionReturn(0);
15229b94acceSBarry Smith }
15239b94acceSBarry Smith 
15245615d1e5SSatish Balay #undef __FUNC__
1525d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory"
1526c9005455SLois Curfman McInnes /*@
1527c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1528c9005455SLois Curfman McInnes 
1529c9005455SLois Curfman McInnes    Input Parameters:
1530c9005455SLois Curfman McInnes .  snes - iterative context obtained from SNESCreate()
1531c9005455SLois Curfman McInnes .  a   - array to hold history
1532c9005455SLois Curfman McInnes .  na  - size of a
1533c9005455SLois Curfman McInnes 
1534c9005455SLois Curfman McInnes    Notes:
1535c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1536c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1537c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1538c9005455SLois Curfman McInnes    at each step.
1539c9005455SLois Curfman McInnes 
1540c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1541c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1542c9005455SLois Curfman McInnes    during the section of code that is being timed.
1543c9005455SLois Curfman McInnes 
1544c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1545c9005455SLois Curfman McInnes @*/
1546c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1547c9005455SLois Curfman McInnes {
15483a40ed3dSBarry Smith   PetscFunctionBegin;
1549c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1550c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1551c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1552c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
15533a40ed3dSBarry Smith   PetscFunctionReturn(0);
1554c9005455SLois Curfman McInnes }
1555c9005455SLois Curfman McInnes 
1556c9005455SLois Curfman McInnes #undef __FUNC__
15575615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
15589b94acceSBarry Smith /*
15599b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
15609b94acceSBarry Smith    positive parameter delta.
15619b94acceSBarry Smith 
15629b94acceSBarry Smith     Input Parameters:
15639b94acceSBarry Smith .   snes - the SNES context
15649b94acceSBarry Smith .   y - approximate solution of linear system
15659b94acceSBarry Smith .   fnorm - 2-norm of current function
15669b94acceSBarry Smith .   delta - trust region size
15679b94acceSBarry Smith 
15689b94acceSBarry Smith     Output Parameters:
15699b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
15709b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
15719b94acceSBarry Smith     region, and exceeds zero otherwise.
15729b94acceSBarry Smith .   ynorm - 2-norm of the step
15739b94acceSBarry Smith 
15749b94acceSBarry Smith     Note:
157540191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
15769b94acceSBarry Smith     is set to be the maximum allowable step size.
15779b94acceSBarry Smith 
15789b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
15799b94acceSBarry Smith */
15809b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
15819b94acceSBarry Smith                           double *gpnorm,double *ynorm)
15829b94acceSBarry Smith {
15839b94acceSBarry Smith   double norm;
15849b94acceSBarry Smith   Scalar cnorm;
15853a40ed3dSBarry Smith   int    ierr;
15863a40ed3dSBarry Smith 
15873a40ed3dSBarry Smith   PetscFunctionBegin;
15883a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
15899b94acceSBarry Smith   if (norm > *delta) {
15909b94acceSBarry Smith      norm = *delta/norm;
15919b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
15929b94acceSBarry Smith      cnorm = norm;
15939b94acceSBarry Smith      VecScale( &cnorm, y );
15949b94acceSBarry Smith      *ynorm = *delta;
15959b94acceSBarry Smith   } else {
15969b94acceSBarry Smith      *gpnorm = 0.0;
15979b94acceSBarry Smith      *ynorm = norm;
15989b94acceSBarry Smith   }
15993a40ed3dSBarry Smith   PetscFunctionReturn(0);
16009b94acceSBarry Smith }
16019b94acceSBarry Smith 
16025615d1e5SSatish Balay #undef __FUNC__
16035615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
16049b94acceSBarry Smith /*@
16059b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
160663a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
16079b94acceSBarry Smith 
16089b94acceSBarry Smith    Input Parameter:
16099b94acceSBarry Smith .  snes - the SNES context
16108ddd3da0SLois Curfman McInnes .  x - the solution vector
16119b94acceSBarry Smith 
16129b94acceSBarry Smith    Output Parameter:
16139b94acceSBarry Smith    its - number of iterations until termination
16149b94acceSBarry Smith 
16158ddd3da0SLois Curfman McInnes    Note:
16168ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
16178ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
16188ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
16198ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
16208ddd3da0SLois Curfman McInnes 
16219b94acceSBarry Smith .keywords: SNES, nonlinear, solve
16229b94acceSBarry Smith 
162363a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
16249b94acceSBarry Smith @*/
16258ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
16269b94acceSBarry Smith {
16273c7409f5SSatish Balay   int ierr, flg;
1628052efed2SBarry Smith 
16293a40ed3dSBarry Smith   PetscFunctionBegin;
163077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
163174679c65SBarry Smith   PetscValidIntPointer(its);
1632c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1633c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
16349b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1635c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
16369b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
16379b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
16383c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
16396d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
16403a40ed3dSBarry Smith   PetscFunctionReturn(0);
16419b94acceSBarry Smith }
16429b94acceSBarry Smith 
16439b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
164437fcc0dbSBarry Smith static NRList *__SNESList = 0;
16459b94acceSBarry Smith 
16465615d1e5SSatish Balay #undef __FUNC__
16475615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
16489b94acceSBarry Smith /*@
16494b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
16509b94acceSBarry Smith 
16519b94acceSBarry Smith    Input Parameters:
16529b94acceSBarry Smith .  snes - the SNES context
16539b94acceSBarry Smith .  method - a known method
16549b94acceSBarry Smith 
1655ae12b187SLois Curfman McInnes   Options Database Command:
1656ae12b187SLois Curfman McInnes $ -snes_type  <method>
1657ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1658ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1659ae12b187SLois Curfman McInnes 
16609b94acceSBarry Smith    Notes:
16619b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
16629b94acceSBarry Smith $  Systems of nonlinear equations:
166340191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
166440191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
16659b94acceSBarry Smith $  Unconstrained minimization:
166640191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
166740191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
16689b94acceSBarry Smith 
1669ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1670ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1671ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1672ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1673ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1674ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1675ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1676ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1677ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1678ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1679a703fe33SLois Curfman McInnes 
1680f536c699SSatish Balay .keywords: SNES, set, method
16819b94acceSBarry Smith @*/
16824b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
16839b94acceSBarry Smith {
168484cb2905SBarry Smith   int ierr;
16859b94acceSBarry Smith   int (*r)(SNES);
16863a40ed3dSBarry Smith 
16873a40ed3dSBarry Smith   PetscFunctionBegin;
168877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16893a40ed3dSBarry Smith   if (snes->type == (int) method) PetscFunctionReturn(0);
16907d1a2b2bSBarry Smith 
16919b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
169284cb2905SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
1693e3372554SBarry Smith   if (!__SNESList) {SETERRQ(1,0,"Could not get methods");}
169437fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1695e3372554SBarry Smith   if (!r) {SETERRQ(1,0,"Unknown method");}
16960452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
16979b94acceSBarry Smith   snes->set_method_called = 1;
16983a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
16993a40ed3dSBarry Smith   PetscFunctionReturn(0);
17009b94acceSBarry Smith }
17019b94acceSBarry Smith 
17029b94acceSBarry Smith /* --------------------------------------------------------------------- */
17035615d1e5SSatish Balay #undef __FUNC__
1704d4bb536fSBarry Smith #define __FUNC__ "SNESRegister"
17059b94acceSBarry Smith /*@C
17069b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
17074b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
17089b94acceSBarry Smith 
17099b94acceSBarry Smith    Input Parameters:
17102d872ea7SLois Curfman McInnes .  name - either a predefined name such as SNES_EQ_LS, or SNES_NEW
17112d872ea7SLois Curfman McInnes           to indicate a new user-defined solver
171240191667SLois Curfman McInnes .  sname - corresponding string for name
17139b94acceSBarry Smith .  create - routine to create method context
17149b94acceSBarry Smith 
171584cb2905SBarry Smith    Output Parameter:
171684cb2905SBarry Smith .  oname - type associated with this new solver
171784cb2905SBarry Smith 
17182d872ea7SLois Curfman McInnes    Notes:
17192d872ea7SLois Curfman McInnes    Multiple user-defined nonlinear solvers can be added by calling
17202d872ea7SLois Curfman McInnes    SNESRegister() with the input parameter "name" set to be SNES_NEW;
17212d872ea7SLois Curfman McInnes    each call will return a unique solver type in the output
17222d872ea7SLois Curfman McInnes    parameter "oname".
17232d872ea7SLois Curfman McInnes 
17249b94acceSBarry Smith .keywords: SNES, nonlinear, register
17259b94acceSBarry Smith 
17269b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
17279b94acceSBarry Smith @*/
172884cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES))
17299b94acceSBarry Smith {
17309b94acceSBarry Smith   int        ierr;
173184cb2905SBarry Smith   static int numberregistered = 0;
173284cb2905SBarry Smith 
17333a40ed3dSBarry Smith   PetscFunctionBegin;
1734d252947aSBarry Smith   if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++);
173584cb2905SBarry Smith 
173684cb2905SBarry Smith   if (oname) *oname = name;
173737fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
173884cb2905SBarry Smith   NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create );
17393a40ed3dSBarry Smith   PetscFunctionReturn(0);
17409b94acceSBarry Smith }
1741a847f771SSatish Balay 
17429b94acceSBarry Smith /* --------------------------------------------------------------------- */
17435615d1e5SSatish Balay #undef __FUNC__
1744d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy"
17459b94acceSBarry Smith /*@C
17469b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
17479b94acceSBarry Smith    registered by SNESRegister().
17489b94acceSBarry Smith 
17499b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
17509b94acceSBarry Smith 
17519b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
17529b94acceSBarry Smith @*/
17539b94acceSBarry Smith int SNESRegisterDestroy()
17549b94acceSBarry Smith {
17553a40ed3dSBarry Smith   PetscFunctionBegin;
175637fcc0dbSBarry Smith   if (__SNESList) {
175737fcc0dbSBarry Smith     NRDestroy( __SNESList );
175837fcc0dbSBarry Smith     __SNESList = 0;
17599b94acceSBarry Smith   }
176084cb2905SBarry Smith   SNESRegisterAllCalled = 0;
17613a40ed3dSBarry Smith   PetscFunctionReturn(0);
17629b94acceSBarry Smith }
17639b94acceSBarry Smith 
17645615d1e5SSatish Balay #undef __FUNC__
1765d4bb536fSBarry Smith #define __FUNC__ "SNESGetTypeFromOptions_Private"
17669b94acceSBarry Smith /*
17674b0e389bSBarry Smith    SNESGetTypeFromOptions_Private - Sets the selected method from the
17689b94acceSBarry Smith    options database.
17699b94acceSBarry Smith 
17709b94acceSBarry Smith    Input Parameter:
17719b94acceSBarry Smith .  ctx - the SNES context
17729b94acceSBarry Smith 
17739b94acceSBarry Smith    Output Parameter:
17749b94acceSBarry Smith .  method -  solver method
1775*d64ed03dSBarry Smith .  flg    - indicate if method found
17769b94acceSBarry Smith 
17779b94acceSBarry Smith    Options Database Key:
17784b0e389bSBarry Smith $  -snes_type  method
17799b94acceSBarry Smith */
1780052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg)
17819b94acceSBarry Smith {
1782052efed2SBarry Smith   int  ierr;
17839b94acceSBarry Smith   char sbuf[50];
17845baf8537SBarry Smith 
17853a40ed3dSBarry Smith   PetscFunctionBegin;
1786052efed2SBarry Smith   ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr);
1787052efed2SBarry Smith   if (*flg) {
1788052efed2SBarry Smith     if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
178937fcc0dbSBarry Smith     *method = (SNESType)NRFindID( __SNESList, sbuf );
1790*d64ed03dSBarry Smith     if (*method == (SNESType) -1) SETERRQ(1,1,"Invalid SNES type");
17919b94acceSBarry Smith   }
17923a40ed3dSBarry Smith   PetscFunctionReturn(0);
17939b94acceSBarry Smith }
17949b94acceSBarry Smith 
17955615d1e5SSatish Balay #undef __FUNC__
1796d4bb536fSBarry Smith #define __FUNC__ "SNESGetType"
17979b94acceSBarry Smith /*@C
17989a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
17999b94acceSBarry Smith 
18009b94acceSBarry Smith    Input Parameter:
18014b0e389bSBarry Smith .  snes - nonlinear solver context
18029b94acceSBarry Smith 
18039b94acceSBarry Smith    Output Parameter:
18049a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
18059a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
18069b94acceSBarry Smith 
18079b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
18089b94acceSBarry Smith @*/
18094b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
18109b94acceSBarry Smith {
181106a9b0adSLois Curfman McInnes   int ierr;
18123a40ed3dSBarry Smith 
18133a40ed3dSBarry Smith   PetscFunctionBegin;
181437fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
18154b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
181637fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
18173a40ed3dSBarry Smith   PetscFunctionReturn(0);
18189b94acceSBarry Smith }
18199b94acceSBarry Smith 
18205615d1e5SSatish Balay #undef __FUNC__
1821d4bb536fSBarry Smith #define __FUNC__ "SNESPrintTypes_Private"
18229b94acceSBarry Smith /*
18234b0e389bSBarry Smith    SNESPrintTypes_Private - Prints the SNES methods available from the
18249b94acceSBarry Smith    options database.
18259b94acceSBarry Smith 
18269b94acceSBarry Smith    Input Parameters:
18270752156aSBarry Smith .  comm   - communicator (usually PETSC_COMM_WORLD)
18289b94acceSBarry Smith .  prefix - prefix (usually "-")
18294b0e389bSBarry Smith .  name   - the options database name (by default "snes_type")
18309b94acceSBarry Smith */
183133455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name)
18329b94acceSBarry Smith {
18339b94acceSBarry Smith   FuncList *entry;
18343a40ed3dSBarry Smith 
18353a40ed3dSBarry Smith   PetscFunctionBegin;
183637fcc0dbSBarry Smith   if (!__SNESList) {SNESRegisterAll();}
183737fcc0dbSBarry Smith   entry = __SNESList->head;
183877c4ece6SBarry Smith   PetscPrintf(comm," %s%s (one of)",prefix,name);
18399b94acceSBarry Smith   while (entry) {
184077c4ece6SBarry Smith     PetscPrintf(comm," %s",entry->name);
18419b94acceSBarry Smith     entry = entry->next;
18429b94acceSBarry Smith   }
184377c4ece6SBarry Smith   PetscPrintf(comm,"\n");
18443a40ed3dSBarry Smith   PetscFunctionReturn(0);
18459b94acceSBarry Smith }
18469b94acceSBarry Smith 
18475615d1e5SSatish Balay #undef __FUNC__
1848d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution"
18499b94acceSBarry Smith /*@C
18509b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
18519b94acceSBarry Smith    stored.
18529b94acceSBarry Smith 
18539b94acceSBarry Smith    Input Parameter:
18549b94acceSBarry Smith .  snes - the SNES context
18559b94acceSBarry Smith 
18569b94acceSBarry Smith    Output Parameter:
18579b94acceSBarry Smith .  x - the solution
18589b94acceSBarry Smith 
18599b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
18609b94acceSBarry Smith 
1861a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
18629b94acceSBarry Smith @*/
18639b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
18649b94acceSBarry Smith {
18653a40ed3dSBarry Smith   PetscFunctionBegin;
186677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18679b94acceSBarry Smith   *x = snes->vec_sol_always;
18683a40ed3dSBarry Smith   PetscFunctionReturn(0);
18699b94acceSBarry Smith }
18709b94acceSBarry Smith 
18715615d1e5SSatish Balay #undef __FUNC__
1872d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate"
18739b94acceSBarry Smith /*@C
18749b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
18759b94acceSBarry Smith    stored.
18769b94acceSBarry Smith 
18779b94acceSBarry Smith    Input Parameter:
18789b94acceSBarry Smith .  snes - the SNES context
18799b94acceSBarry Smith 
18809b94acceSBarry Smith    Output Parameter:
18819b94acceSBarry Smith .  x - the solution update
18829b94acceSBarry Smith 
18839b94acceSBarry Smith    Notes:
18849b94acceSBarry Smith    This vector is implementation dependent.
18859b94acceSBarry Smith 
18869b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
18879b94acceSBarry Smith 
18889b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
18899b94acceSBarry Smith @*/
18909b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
18919b94acceSBarry Smith {
18923a40ed3dSBarry Smith   PetscFunctionBegin;
189377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18949b94acceSBarry Smith   *x = snes->vec_sol_update_always;
18953a40ed3dSBarry Smith   PetscFunctionReturn(0);
18969b94acceSBarry Smith }
18979b94acceSBarry Smith 
18985615d1e5SSatish Balay #undef __FUNC__
1899d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction"
19009b94acceSBarry Smith /*@C
19013638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
19029b94acceSBarry Smith 
19039b94acceSBarry Smith    Input Parameter:
19049b94acceSBarry Smith .  snes - the SNES context
19059b94acceSBarry Smith 
19069b94acceSBarry Smith    Output Parameter:
19073638b69dSLois Curfman McInnes .  r - the function
19089b94acceSBarry Smith 
19099b94acceSBarry Smith    Notes:
19109b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
19119b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
19129b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
19139b94acceSBarry Smith 
1914a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
19159b94acceSBarry Smith 
191631615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
191731615d04SLois Curfman McInnes           SNESGetGradient()
19189b94acceSBarry Smith @*/
19199b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
19209b94acceSBarry Smith {
19213a40ed3dSBarry Smith   PetscFunctionBegin;
192277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
19233a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
19249b94acceSBarry Smith   *r = snes->vec_func_always;
19253a40ed3dSBarry Smith   PetscFunctionReturn(0);
19269b94acceSBarry Smith }
19279b94acceSBarry Smith 
19285615d1e5SSatish Balay #undef __FUNC__
1929d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient"
19309b94acceSBarry Smith /*@C
19313638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
19329b94acceSBarry Smith 
19339b94acceSBarry Smith    Input Parameter:
19349b94acceSBarry Smith .  snes - the SNES context
19359b94acceSBarry Smith 
19369b94acceSBarry Smith    Output Parameter:
19379b94acceSBarry Smith .  r - the gradient
19389b94acceSBarry Smith 
19399b94acceSBarry Smith    Notes:
19409b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
19419b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
19429b94acceSBarry Smith    SNESGetFunction().
19439b94acceSBarry Smith 
19449b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
19459b94acceSBarry Smith 
194631615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
19479b94acceSBarry Smith @*/
19489b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
19499b94acceSBarry Smith {
19503a40ed3dSBarry Smith   PetscFunctionBegin;
195177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
19523a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
19533a40ed3dSBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
19543a40ed3dSBarry Smith   }
19559b94acceSBarry Smith   *r = snes->vec_func_always;
19563a40ed3dSBarry Smith   PetscFunctionReturn(0);
19579b94acceSBarry Smith }
19589b94acceSBarry Smith 
19595615d1e5SSatish Balay #undef __FUNC__
1960d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction"
19619b94acceSBarry Smith /*@
19629b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
19639b94acceSBarry Smith    unconstrained minimization problems.
19649b94acceSBarry Smith 
19659b94acceSBarry Smith    Input Parameter:
19669b94acceSBarry Smith .  snes - the SNES context
19679b94acceSBarry Smith 
19689b94acceSBarry Smith    Output Parameter:
19699b94acceSBarry Smith .  r - the function
19709b94acceSBarry Smith 
19719b94acceSBarry Smith    Notes:
19729b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
19739b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
19749b94acceSBarry Smith    SNESGetFunction().
19759b94acceSBarry Smith 
19769b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
19779b94acceSBarry Smith 
197831615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
19799b94acceSBarry Smith @*/
19809b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
19819b94acceSBarry Smith {
19823a40ed3dSBarry Smith   PetscFunctionBegin;
198377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
198474679c65SBarry Smith   PetscValidScalarPointer(r);
19853a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
19863a40ed3dSBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
19873a40ed3dSBarry Smith   }
19889b94acceSBarry Smith   *r = snes->fc;
19893a40ed3dSBarry Smith   PetscFunctionReturn(0);
19909b94acceSBarry Smith }
19919b94acceSBarry Smith 
19925615d1e5SSatish Balay #undef __FUNC__
1993d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix"
19943c7409f5SSatish Balay /*@C
19953c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1996d850072dSLois Curfman McInnes    SNES options in the database.
19973c7409f5SSatish Balay 
19983c7409f5SSatish Balay    Input Parameter:
19993c7409f5SSatish Balay .  snes - the SNES context
20003c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
20013c7409f5SSatish Balay 
2002d850072dSLois Curfman McInnes    Notes:
2003a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2004a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
2005a83b1b31SSatish Balay    hyphen.
2006d850072dSLois Curfman McInnes 
20073c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2008a86d99e1SLois Curfman McInnes 
2009a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
20103c7409f5SSatish Balay @*/
20113c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
20123c7409f5SSatish Balay {
20133c7409f5SSatish Balay   int ierr;
20143c7409f5SSatish Balay 
20153a40ed3dSBarry Smith   PetscFunctionBegin;
201677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2017639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20183c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
20193a40ed3dSBarry Smith   PetscFunctionReturn(0);
20203c7409f5SSatish Balay }
20213c7409f5SSatish Balay 
20225615d1e5SSatish Balay #undef __FUNC__
2023d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix"
20243c7409f5SSatish Balay /*@C
2025f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2026d850072dSLois Curfman McInnes    SNES options in the database.
20273c7409f5SSatish Balay 
20283c7409f5SSatish Balay    Input Parameter:
20293c7409f5SSatish Balay .  snes - the SNES context
20303c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
20313c7409f5SSatish Balay 
2032d850072dSLois Curfman McInnes    Notes:
2033a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2034a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
2035a83b1b31SSatish Balay    hyphen.
2036d850072dSLois Curfman McInnes 
20373c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2038a86d99e1SLois Curfman McInnes 
2039a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20403c7409f5SSatish Balay @*/
20413c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
20423c7409f5SSatish Balay {
20433c7409f5SSatish Balay   int ierr;
20443c7409f5SSatish Balay 
20453a40ed3dSBarry Smith   PetscFunctionBegin;
204677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2047639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20483c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
20493a40ed3dSBarry Smith   PetscFunctionReturn(0);
20503c7409f5SSatish Balay }
20513c7409f5SSatish Balay 
20525615d1e5SSatish Balay #undef __FUNC__
2053d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix"
2054c83590e2SSatish Balay /*@
20553c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20563c7409f5SSatish Balay    SNES options in the database.
20573c7409f5SSatish Balay 
20583c7409f5SSatish Balay    Input Parameter:
20593c7409f5SSatish Balay .  snes - the SNES context
20603c7409f5SSatish Balay 
20613c7409f5SSatish Balay    Output Parameter:
20623c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20633c7409f5SSatish Balay 
20643c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2065a86d99e1SLois Curfman McInnes 
2066a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20673c7409f5SSatish Balay @*/
20683c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
20693c7409f5SSatish Balay {
20703c7409f5SSatish Balay   int ierr;
20713c7409f5SSatish Balay 
20723a40ed3dSBarry Smith   PetscFunctionBegin;
207377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2074639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20753a40ed3dSBarry Smith   PetscFunctionReturn(0);
20763c7409f5SSatish Balay }
20773c7409f5SSatish Balay 
20783c7409f5SSatish Balay 
20799b94acceSBarry Smith 
20809b94acceSBarry Smith 
20819b94acceSBarry Smith 
20823a40ed3dSBarry Smith 
2083