xref: /petsc/src/snes/interface/snes.c (revision ca161407697b5d73508d887e7f1927b4808f5f9d)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*ca161407SBarry Smith static char vcid[] = "$Id: snes.c,v 1.132 1997/10/28 14:24:42 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 
1084cb2905SBarry Smith int SNESRegisterAllCalled = 0;
11*ca161407SBarry Smith static NRList *__SNESList = 0;
1284cb2905SBarry Smith 
135615d1e5SSatish Balay #undef __FUNC__
14d4bb536fSBarry Smith #define __FUNC__ "SNESView"
159b94acceSBarry Smith /*@
169b94acceSBarry Smith    SNESView - Prints the SNES data structure.
179b94acceSBarry Smith 
189b94acceSBarry Smith    Input Parameters:
199b94acceSBarry Smith .  SNES - the SNES context
209b94acceSBarry Smith .  viewer - visualization context
219b94acceSBarry Smith 
229b94acceSBarry Smith    Options Database Key:
239b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
249b94acceSBarry Smith 
259b94acceSBarry Smith    Notes:
269b94acceSBarry Smith    The available visualization contexts include
276d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
286d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
299b94acceSBarry Smith $       output where only the first processor opens
309b94acceSBarry Smith $       the file.  All other processors send their
319b94acceSBarry Smith $       data to the first processor to print.
329b94acceSBarry Smith 
339b94acceSBarry Smith    The user can open alternative vistualization contexts with
34dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
359b94acceSBarry Smith 
369b94acceSBarry Smith .keywords: SNES, view
379b94acceSBarry Smith 
38dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
399b94acceSBarry Smith @*/
409b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
419b94acceSBarry Smith {
429b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
439b94acceSBarry Smith   FILE                *fd;
449b94acceSBarry Smith   int                 ierr;
459b94acceSBarry Smith   SLES                sles;
469b94acceSBarry Smith   char                *method;
4719bcc07fSBarry Smith   ViewerType          vtype;
489b94acceSBarry Smith 
493a40ed3dSBarry Smith   PetscFunctionBegin;
5074679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5174679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5274679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5374679c65SBarry Smith 
5419bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5519bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5690ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5777c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
584b0e389bSBarry Smith     SNESGetType(snes,PETSC_NULL,&method);
5977c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
609b94acceSBarry Smith     if (snes->view) (*snes->view)((PetscObject)snes,viewer);
6177c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
629b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
639b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6477c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
659b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
669b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
677a00f4a9SLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
687a00f4a9SLois Curfman McInnes     "  total number of linear solver iterations=%d\n",snes->linear_its);
69ae3c334cSLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
7068632166SLois Curfman McInnes      "  total number of function evaluations=%d\n",snes->nfuncs);
719b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
7277c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
739b94acceSBarry Smith     if (snes->ksp_ewconv) {
749b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
759b94acceSBarry Smith       if (kctx) {
7677c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
779b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
789b94acceSBarry Smith         kctx->version);
7977c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
809b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
819b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
8277c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
839b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
849b94acceSBarry Smith       }
859b94acceSBarry Smith     }
86c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
87c30f285eSLois Curfman McInnes     SNESGetType(snes,PETSC_NULL,&method);
88c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
8919bcc07fSBarry Smith   }
909b94acceSBarry Smith   SNESGetSLES(snes,&sles);
919b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
923a40ed3dSBarry Smith   PetscFunctionReturn(0);
939b94acceSBarry Smith }
949b94acceSBarry Smith 
95639f9d9dSBarry Smith /*
96639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
97639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
98639f9d9dSBarry Smith */
99639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
100639f9d9dSBarry Smith static int numberofsetfromoptions;
101639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
102639f9d9dSBarry Smith 
1035615d1e5SSatish Balay #undef __FUNC__
104d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker"
105639f9d9dSBarry Smith /*@
106639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
107639f9d9dSBarry Smith 
108639f9d9dSBarry Smith     Input Parameter:
109639f9d9dSBarry Smith .   snescheck - function that checks for options
110639f9d9dSBarry Smith 
111639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
112639f9d9dSBarry Smith @*/
113639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
114639f9d9dSBarry Smith {
1153a40ed3dSBarry Smith   PetscFunctionBegin;
116639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
117e3372554SBarry Smith     SETERRQ(1,0,"Too many options checkers, only 5 allowed");
118639f9d9dSBarry Smith   }
119639f9d9dSBarry Smith 
120639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
1213a40ed3dSBarry Smith   PetscFunctionReturn(0);
122639f9d9dSBarry Smith }
123639f9d9dSBarry Smith 
1245615d1e5SSatish Balay #undef __FUNC__
1255615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1269b94acceSBarry Smith /*@
127682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1289b94acceSBarry Smith 
1299b94acceSBarry Smith    Input Parameter:
1309b94acceSBarry Smith .  snes - the SNES context
1319b94acceSBarry Smith 
13211ca99fdSLois Curfman McInnes    Notes:
13311ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
13411ca99fdSLois Curfman McInnes    the users manual.
13583e2fdc7SBarry Smith 
1369b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1379b94acceSBarry Smith 
138a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1399b94acceSBarry Smith @*/
1409b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1419b94acceSBarry Smith {
1424b0e389bSBarry Smith   SNESType method;
1439b94acceSBarry Smith   double   tmp;
1449b94acceSBarry Smith   SLES     sles;
14581f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1463c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1479b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1489b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1499b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1509b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1519b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1529b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1539b94acceSBarry Smith 
1543a40ed3dSBarry Smith   PetscFunctionBegin;
15581f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
15681f57ec7SBarry Smith 
15777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
158e3372554SBarry Smith   if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp");
159*ca161407SBarry Smith 
160*ca161407SBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll();CHKERRQ(ierr);}
161*ca161407SBarry Smith   ierr = NRGetTypeFromOptions(snes->prefix,"-snes_type",__SNESList,&method,&flg);CHKERRQ(ierr);
162052efed2SBarry Smith   if (flg) {
163052efed2SBarry Smith     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
164*ca161407SBarry Smith   } else if (!snes->set_method_called) {
1655334005bSBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
16640191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
167d64ed03dSBarry Smith     } else {
16840191667SLois Curfman McInnes       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
1695334005bSBarry Smith     }
1705334005bSBarry Smith   }
1715334005bSBarry Smith 
1723c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
173d64ed03dSBarry Smith   if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);}
1743c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
175d64ed03dSBarry Smith   if (flg) {
176d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
177d64ed03dSBarry Smith   }
1783c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
179d64ed03dSBarry Smith   if (flg) {
180d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
181d64ed03dSBarry Smith   }
1823c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
183d64ed03dSBarry Smith   if (flg) {
184d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
185d64ed03dSBarry Smith   }
1863c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1873c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
188d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
189d64ed03dSBarry Smith   if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); }
190d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
191d64ed03dSBarry Smith   if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);  CHKERRQ(ierr);}
1923c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1933c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
194b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
195b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
196b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
197b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
198b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
199b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
200b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
20193c39befSBarry Smith 
20293c39befSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr);
20393c39befSBarry Smith   if (flg) {snes->converged = 0;}
20493c39befSBarry Smith 
2059b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
2069b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
2075bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
2085bbad29bSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);}
2093c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
210639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
2113c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
212d31a3109SLois Curfman McInnes   if (flg) {
213639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
214d31a3109SLois Curfman McInnes   }
21581f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
2163c7409f5SSatish Balay   if (flg) {
21717699dbbSLois Curfman McInnes     int    rank = 0;
218d7e8b826SBarry Smith     DrawLG lg;
21917699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
22017699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
22117699dbbSLois Curfman McInnes     if (!rank) {
22281f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2239b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
224f8c826e1SBarry Smith       snes->xmonitor = lg;
2259b94acceSBarry Smith       PLogObjectParent(snes,lg);
2269b94acceSBarry Smith     }
2279b94acceSBarry Smith   }
2283c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2293c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2309b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2319b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
23294a424c1SBarry Smith     PLogInfo(snes,
23331615d04SLois Curfman McInnes       "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
23431615d04SLois Curfman McInnes   }
23531615d04SLois Curfman McInnes   else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
23631615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
23731615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
238d64ed03dSBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2399b94acceSBarry Smith   }
240639f9d9dSBarry Smith 
241639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
242639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
243639f9d9dSBarry Smith   }
244639f9d9dSBarry Smith 
2459b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2469b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
2473a40ed3dSBarry Smith   if (snes->setfromoptions) {
2483a40ed3dSBarry Smith     ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
2493a40ed3dSBarry Smith   }
2503a40ed3dSBarry Smith   PetscFunctionReturn(0);
2519b94acceSBarry Smith }
2529b94acceSBarry Smith 
253a847f771SSatish Balay 
2545615d1e5SSatish Balay #undef __FUNC__
255d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext"
2569b94acceSBarry Smith /*@
2579b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2589b94acceSBarry Smith    the nonlinear solvers.
2599b94acceSBarry Smith 
2609b94acceSBarry Smith    Input Parameters:
2619b94acceSBarry Smith .  snes - the SNES context
2629b94acceSBarry Smith .  usrP - optional user context
2639b94acceSBarry Smith 
2649b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2659b94acceSBarry Smith 
2669b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2679b94acceSBarry Smith @*/
2689b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
2699b94acceSBarry Smith {
2703a40ed3dSBarry Smith   PetscFunctionBegin;
27177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2729b94acceSBarry Smith   snes->user		= usrP;
2733a40ed3dSBarry Smith   PetscFunctionReturn(0);
2749b94acceSBarry Smith }
27574679c65SBarry Smith 
2765615d1e5SSatish Balay #undef __FUNC__
277d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext"
2789b94acceSBarry Smith /*@C
2799b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
2809b94acceSBarry Smith    nonlinear solvers.
2819b94acceSBarry Smith 
2829b94acceSBarry Smith    Input Parameter:
2839b94acceSBarry Smith .  snes - SNES context
2849b94acceSBarry Smith 
2859b94acceSBarry Smith    Output Parameter:
2869b94acceSBarry Smith .  usrP - user context
2879b94acceSBarry Smith 
2889b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
2899b94acceSBarry Smith 
2909b94acceSBarry Smith .seealso: SNESSetApplicationContext()
2919b94acceSBarry Smith @*/
2929b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
2939b94acceSBarry Smith {
2943a40ed3dSBarry Smith   PetscFunctionBegin;
29577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2969b94acceSBarry Smith   *usrP = snes->user;
2973a40ed3dSBarry Smith   PetscFunctionReturn(0);
2989b94acceSBarry Smith }
29974679c65SBarry Smith 
3005615d1e5SSatish Balay #undef __FUNC__
301d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber"
3029b94acceSBarry Smith /*@
3039b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3049b94acceSBarry Smith    nonlinear solver.
3059b94acceSBarry Smith 
3069b94acceSBarry Smith    Input Parameter:
3079b94acceSBarry Smith .  snes - SNES context
3089b94acceSBarry Smith 
3099b94acceSBarry Smith    Output Parameter:
3109b94acceSBarry Smith .  iter - iteration number
3119b94acceSBarry Smith 
3129b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3139b94acceSBarry Smith @*/
3149b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3159b94acceSBarry Smith {
3163a40ed3dSBarry Smith   PetscFunctionBegin;
31777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
31874679c65SBarry Smith   PetscValidIntPointer(iter);
3199b94acceSBarry Smith   *iter = snes->iter;
3203a40ed3dSBarry Smith   PetscFunctionReturn(0);
3219b94acceSBarry Smith }
32274679c65SBarry Smith 
3235615d1e5SSatish Balay #undef __FUNC__
3245615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
3259b94acceSBarry Smith /*@
3269b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3279b94acceSBarry Smith    with SNESSSetFunction().
3289b94acceSBarry Smith 
3299b94acceSBarry Smith    Input Parameter:
3309b94acceSBarry Smith .  snes - SNES context
3319b94acceSBarry Smith 
3329b94acceSBarry Smith    Output Parameter:
3339b94acceSBarry Smith .  fnorm - 2-norm of function
3349b94acceSBarry Smith 
3359b94acceSBarry Smith    Note:
3369b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
337a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
338a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3399b94acceSBarry Smith 
3409b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
341a86d99e1SLois Curfman McInnes 
342a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3439b94acceSBarry Smith @*/
3449b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3459b94acceSBarry Smith {
3463a40ed3dSBarry Smith   PetscFunctionBegin;
34777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
34874679c65SBarry Smith   PetscValidScalarPointer(fnorm);
34974679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
350d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only");
35174679c65SBarry Smith   }
3529b94acceSBarry Smith   *fnorm = snes->norm;
3533a40ed3dSBarry Smith   PetscFunctionReturn(0);
3549b94acceSBarry Smith }
35574679c65SBarry Smith 
3565615d1e5SSatish Balay #undef __FUNC__
3575615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
3589b94acceSBarry Smith /*@
3599b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
3609b94acceSBarry Smith    with SNESSSetGradient().
3619b94acceSBarry Smith 
3629b94acceSBarry Smith    Input Parameter:
3639b94acceSBarry Smith .  snes - SNES context
3649b94acceSBarry Smith 
3659b94acceSBarry Smith    Output Parameter:
3669b94acceSBarry Smith .  fnorm - 2-norm of gradient
3679b94acceSBarry Smith 
3689b94acceSBarry Smith    Note:
3699b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
370a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
371a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
3729b94acceSBarry Smith 
3739b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
374a86d99e1SLois Curfman McInnes 
375a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
3769b94acceSBarry Smith @*/
3779b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
3789b94acceSBarry Smith {
3793a40ed3dSBarry Smith   PetscFunctionBegin;
38077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
38174679c65SBarry Smith   PetscValidScalarPointer(gnorm);
38274679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
383d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
38474679c65SBarry Smith   }
3859b94acceSBarry Smith   *gnorm = snes->norm;
3863a40ed3dSBarry Smith   PetscFunctionReturn(0);
3879b94acceSBarry Smith }
38874679c65SBarry Smith 
3895615d1e5SSatish Balay #undef __FUNC__
390d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps"
3919b94acceSBarry Smith /*@
3929b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
3939b94acceSBarry Smith    attempted by the nonlinear solver.
3949b94acceSBarry Smith 
3959b94acceSBarry Smith    Input Parameter:
3969b94acceSBarry Smith .  snes - SNES context
3979b94acceSBarry Smith 
3989b94acceSBarry Smith    Output Parameter:
3999b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4009b94acceSBarry Smith 
401c96a6f78SLois Curfman McInnes    Notes:
402c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
403c96a6f78SLois Curfman McInnes 
4049b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4059b94acceSBarry Smith @*/
4069b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4079b94acceSBarry Smith {
4083a40ed3dSBarry Smith   PetscFunctionBegin;
40977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
41074679c65SBarry Smith   PetscValidIntPointer(nfails);
4119b94acceSBarry Smith   *nfails = snes->nfailures;
4123a40ed3dSBarry Smith   PetscFunctionReturn(0);
4139b94acceSBarry Smith }
414a847f771SSatish Balay 
4155615d1e5SSatish Balay #undef __FUNC__
416d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations"
417c96a6f78SLois Curfman McInnes /*@
418c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
419c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
420c96a6f78SLois Curfman McInnes 
421c96a6f78SLois Curfman McInnes    Input Parameter:
422c96a6f78SLois Curfman McInnes .  snes - SNES context
423c96a6f78SLois Curfman McInnes 
424c96a6f78SLois Curfman McInnes    Output Parameter:
425c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
426c96a6f78SLois Curfman McInnes 
427c96a6f78SLois Curfman McInnes    Notes:
428c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
429c96a6f78SLois Curfman McInnes 
430c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
431c96a6f78SLois Curfman McInnes @*/
43286bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
433c96a6f78SLois Curfman McInnes {
4343a40ed3dSBarry Smith   PetscFunctionBegin;
435c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
436c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
437c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
4383a40ed3dSBarry Smith   PetscFunctionReturn(0);
439c96a6f78SLois Curfman McInnes }
440c96a6f78SLois Curfman McInnes 
441c96a6f78SLois Curfman McInnes #undef __FUNC__
442d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES"
4439b94acceSBarry Smith /*@C
4449b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4459b94acceSBarry Smith 
4469b94acceSBarry Smith    Input Parameter:
4479b94acceSBarry Smith .  snes - the SNES context
4489b94acceSBarry Smith 
4499b94acceSBarry Smith    Output Parameter:
4509b94acceSBarry Smith .  sles - the SLES context
4519b94acceSBarry Smith 
4529b94acceSBarry Smith    Notes:
4539b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4549b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4559b94acceSBarry Smith    KSP and PC contexts as well.
4569b94acceSBarry Smith 
4579b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4589b94acceSBarry Smith 
4599b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4609b94acceSBarry Smith @*/
4619b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4629b94acceSBarry Smith {
4633a40ed3dSBarry Smith   PetscFunctionBegin;
46477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4659b94acceSBarry Smith   *sles = snes->sles;
4663a40ed3dSBarry Smith   PetscFunctionReturn(0);
4679b94acceSBarry Smith }
4689b94acceSBarry Smith /* -----------------------------------------------------------*/
4695615d1e5SSatish Balay #undef __FUNC__
4705615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
4719b94acceSBarry Smith /*@C
4729b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
4739b94acceSBarry Smith 
4749b94acceSBarry Smith    Input Parameter:
4759b94acceSBarry Smith .  comm - MPI communicator
4769b94acceSBarry Smith .  type - type of method, one of
4779b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
4789b94acceSBarry Smith $      (for systems of nonlinear equations)
4799b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
4809b94acceSBarry Smith $      (for unconstrained minimization)
4819b94acceSBarry Smith 
4829b94acceSBarry Smith    Output Parameter:
4839b94acceSBarry Smith .  outsnes - the new SNES context
4849b94acceSBarry Smith 
485c1f60f51SBarry Smith    Options Database Key:
48619bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
48719bd6747SLois Curfman McInnes $              and no preconditioning matrix
48819bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
48919bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
49019bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
49119bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
492c1f60f51SBarry Smith 
4939b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
4949b94acceSBarry Smith 
49563a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
4969b94acceSBarry Smith @*/
4974b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
4989b94acceSBarry Smith {
4999b94acceSBarry Smith   int                 ierr;
5009b94acceSBarry Smith   SNES                snes;
5019b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
50237fcc0dbSBarry Smith 
5033a40ed3dSBarry Smith   PetscFunctionBegin;
5049b94acceSBarry Smith   *outsnes = 0;
505d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
506d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
507d64ed03dSBarry Smith   }
508d4bb536fSBarry Smith   PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm,SNESDestroy,SNESView);
5099b94acceSBarry Smith   PLogObjectCreate(snes);
5109b94acceSBarry Smith   snes->max_its           = 50;
5119750a799SBarry Smith   snes->max_funcs	  = 10000;
5129b94acceSBarry Smith   snes->norm		  = 0.0;
5135a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
514b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
515b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5169b94acceSBarry Smith     snes->atol		  = 1.e-10;
517d64ed03dSBarry Smith   } else {
518b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
519b4874afaSBarry Smith     snes->ttol            = 0.0;
520b4874afaSBarry Smith     snes->atol		  = 1.e-50;
521b4874afaSBarry Smith   }
5229b94acceSBarry Smith   snes->xtol		  = 1.e-8;
523b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5249b94acceSBarry Smith   snes->nfuncs            = 0;
5259b94acceSBarry Smith   snes->nfailures         = 0;
5267a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
527639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5289b94acceSBarry Smith   snes->data              = 0;
5299b94acceSBarry Smith   snes->view              = 0;
5309b94acceSBarry Smith   snes->computeumfunction = 0;
5319b94acceSBarry Smith   snes->umfunP            = 0;
5329b94acceSBarry Smith   snes->fc                = 0;
5339b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5349b94acceSBarry Smith   snes->fmin              = -1.e30;
5359b94acceSBarry Smith   snes->method_class      = type;
5369b94acceSBarry Smith   snes->set_method_called = 0;
537a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
5389b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5397d1a2b2bSBarry Smith   snes->type              = -1;
5406f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5416f24a144SLois Curfman McInnes   snes->vwork             = 0;
5426f24a144SLois Curfman McInnes   snes->nwork             = 0;
543c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
544c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
545c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
5469b94acceSBarry Smith 
5479b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5480452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
549eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
5509b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5519b94acceSBarry Smith   kctx->version     = 2;
5529b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5539b94acceSBarry Smith                              this was too large for some test cases */
5549b94acceSBarry Smith   kctx->rtol_last   = 0;
5559b94acceSBarry Smith   kctx->rtol_max    = .9;
5569b94acceSBarry Smith   kctx->gamma       = 1.0;
5579b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5589b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5599b94acceSBarry Smith   kctx->threshold   = .1;
5609b94acceSBarry Smith   kctx->lresid_last = 0;
5619b94acceSBarry Smith   kctx->norm_last   = 0;
5629b94acceSBarry Smith 
5639b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5649b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5655334005bSBarry Smith 
5669b94acceSBarry Smith   *outsnes = snes;
5673a40ed3dSBarry Smith   PetscFunctionReturn(0);
5689b94acceSBarry Smith }
5699b94acceSBarry Smith 
5709b94acceSBarry Smith /* --------------------------------------------------------------- */
5715615d1e5SSatish Balay #undef __FUNC__
572d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction"
5739b94acceSBarry Smith /*@C
5749b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
5759b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
5769b94acceSBarry Smith    equations.
5779b94acceSBarry Smith 
5789b94acceSBarry Smith    Input Parameters:
5799b94acceSBarry Smith .  snes - the SNES context
5809b94acceSBarry Smith .  func - function evaluation routine
5819b94acceSBarry Smith .  r - vector to store function value
5822cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
5832cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
5849b94acceSBarry Smith 
5859b94acceSBarry Smith    Calling sequence of func:
586313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
5879b94acceSBarry Smith 
5889b94acceSBarry Smith .  x - input vector
589313e4042SLois Curfman McInnes .  f - function vector
5902cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
5919b94acceSBarry Smith 
5929b94acceSBarry Smith    Notes:
5939b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
5949b94acceSBarry Smith $      f'(x) x = -f(x),
5959b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
5969b94acceSBarry Smith 
5979b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
5989b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
5999b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6009b94acceSBarry Smith 
6019b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6029b94acceSBarry Smith 
603a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6049b94acceSBarry Smith @*/
6055334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6069b94acceSBarry Smith {
6073a40ed3dSBarry Smith   PetscFunctionBegin;
60877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
609e3372554SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
6109b94acceSBarry Smith   snes->computefunction     = func;
6119b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6129b94acceSBarry Smith   snes->funP                = ctx;
6133a40ed3dSBarry Smith   PetscFunctionReturn(0);
6149b94acceSBarry Smith }
6159b94acceSBarry Smith 
6165615d1e5SSatish Balay #undef __FUNC__
6175615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
6189b94acceSBarry Smith /*@
6199b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6209b94acceSBarry Smith    SNESSetFunction().
6219b94acceSBarry Smith 
6229b94acceSBarry Smith    Input Parameters:
6239b94acceSBarry Smith .  snes - the SNES context
6249b94acceSBarry Smith .  x - input vector
6259b94acceSBarry Smith 
6269b94acceSBarry Smith    Output Parameter:
6273638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6289b94acceSBarry Smith 
6291bffabb2SLois Curfman McInnes    Notes:
6309b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6319b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6329b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6339b94acceSBarry Smith 
6349b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6359b94acceSBarry Smith 
636a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6379b94acceSBarry Smith @*/
6389b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6399b94acceSBarry Smith {
6409b94acceSBarry Smith   int    ierr;
6419b94acceSBarry Smith 
6423a40ed3dSBarry Smith   PetscFunctionBegin;
643d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
644e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
645d4bb536fSBarry Smith   }
6469b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
647d64ed03dSBarry Smith   PetscStackPush("SNES user function");
6489b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
649d64ed03dSBarry Smith   PetscStackPop;
650ae3c334cSLois Curfman McInnes   snes->nfuncs++;
6519b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6523a40ed3dSBarry Smith   PetscFunctionReturn(0);
6539b94acceSBarry Smith }
6549b94acceSBarry Smith 
6555615d1e5SSatish Balay #undef __FUNC__
656d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction"
6579b94acceSBarry Smith /*@C
6589b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6599b94acceSBarry Smith    unconstrained minimization.
6609b94acceSBarry Smith 
6619b94acceSBarry Smith    Input Parameters:
6629b94acceSBarry Smith .  snes - the SNES context
6639b94acceSBarry Smith .  func - function evaluation routine
664044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
665044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6669b94acceSBarry Smith 
6679b94acceSBarry Smith    Calling sequence of func:
6689b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
6699b94acceSBarry Smith 
6709b94acceSBarry Smith .  x - input vector
6719b94acceSBarry Smith .  f - function
672044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
6739b94acceSBarry Smith 
6749b94acceSBarry Smith    Notes:
6759b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
6769b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
6779b94acceSBarry Smith    SNESSetFunction().
6789b94acceSBarry Smith 
6799b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
6809b94acceSBarry Smith 
681a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
682a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
6839b94acceSBarry Smith @*/
6849b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
6859b94acceSBarry Smith                       void *ctx)
6869b94acceSBarry Smith {
6873a40ed3dSBarry Smith   PetscFunctionBegin;
68877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
689e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
6909b94acceSBarry Smith   snes->computeumfunction   = func;
6919b94acceSBarry Smith   snes->umfunP              = ctx;
6923a40ed3dSBarry Smith   PetscFunctionReturn(0);
6939b94acceSBarry Smith }
6949b94acceSBarry Smith 
6955615d1e5SSatish Balay #undef __FUNC__
6965615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
6979b94acceSBarry Smith /*@
6989b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
6999b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7009b94acceSBarry Smith 
7019b94acceSBarry Smith    Input Parameters:
7029b94acceSBarry Smith .  snes - the SNES context
7039b94acceSBarry Smith .  x - input vector
7049b94acceSBarry Smith 
7059b94acceSBarry Smith    Output Parameter:
7069b94acceSBarry Smith .  y - function value
7079b94acceSBarry Smith 
7089b94acceSBarry Smith    Notes:
7099b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7109b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7119b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
712a86d99e1SLois Curfman McInnes 
713a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
714a86d99e1SLois Curfman McInnes 
715a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
716a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7179b94acceSBarry Smith @*/
7189b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7199b94acceSBarry Smith {
7209b94acceSBarry Smith   int    ierr;
7213a40ed3dSBarry Smith 
7223a40ed3dSBarry Smith   PetscFunctionBegin;
723e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
7249b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
725d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
7269b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
727d64ed03dSBarry Smith   PetscStackPop;
728ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7299b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7303a40ed3dSBarry Smith   PetscFunctionReturn(0);
7319b94acceSBarry Smith }
7329b94acceSBarry Smith 
7335615d1e5SSatish Balay #undef __FUNC__
734d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient"
7359b94acceSBarry Smith /*@C
7369b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7379b94acceSBarry Smith    vector for use by the SNES routines.
7389b94acceSBarry Smith 
7399b94acceSBarry Smith    Input Parameters:
7409b94acceSBarry Smith .  snes - the SNES context
7419b94acceSBarry Smith .  func - function evaluation routine
742044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
743044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7449b94acceSBarry Smith .  r - vector to store gradient value
7459b94acceSBarry Smith 
7469b94acceSBarry Smith    Calling sequence of func:
7479b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7489b94acceSBarry Smith 
7499b94acceSBarry Smith .  x - input vector
7509b94acceSBarry Smith .  g - gradient vector
751044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7529b94acceSBarry Smith 
7539b94acceSBarry Smith    Notes:
7549b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7559b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7569b94acceSBarry Smith    SNESSetFunction().
7579b94acceSBarry Smith 
7589b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7599b94acceSBarry Smith 
760a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
761a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
7629b94acceSBarry Smith @*/
76374679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
7649b94acceSBarry Smith {
7653a40ed3dSBarry Smith   PetscFunctionBegin;
76677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
767e3372554SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
7689b94acceSBarry Smith   snes->computefunction     = func;
7699b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7709b94acceSBarry Smith   snes->funP                = ctx;
7713a40ed3dSBarry Smith   PetscFunctionReturn(0);
7729b94acceSBarry Smith }
7739b94acceSBarry Smith 
7745615d1e5SSatish Balay #undef __FUNC__
7755615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
7769b94acceSBarry Smith /*@
777a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
778a86d99e1SLois Curfman McInnes    SNESSetGradient().
7799b94acceSBarry Smith 
7809b94acceSBarry Smith    Input Parameters:
7819b94acceSBarry Smith .  snes - the SNES context
7829b94acceSBarry Smith .  x - input vector
7839b94acceSBarry Smith 
7849b94acceSBarry Smith    Output Parameter:
7859b94acceSBarry Smith .  y - gradient vector
7869b94acceSBarry Smith 
7879b94acceSBarry Smith    Notes:
7889b94acceSBarry Smith    SNESComputeGradient() is valid only for
7899b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7909b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
791a86d99e1SLois Curfman McInnes 
792a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
793a86d99e1SLois Curfman McInnes 
794a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
795a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
7969b94acceSBarry Smith @*/
7979b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
7989b94acceSBarry Smith {
7999b94acceSBarry Smith   int    ierr;
8003a40ed3dSBarry Smith 
8013a40ed3dSBarry Smith   PetscFunctionBegin;
8023a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
803e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8043a40ed3dSBarry Smith   }
8059b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
806d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
8079b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
808d64ed03dSBarry Smith   PetscStackPop;
8099b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8103a40ed3dSBarry Smith   PetscFunctionReturn(0);
8119b94acceSBarry Smith }
8129b94acceSBarry Smith 
8135615d1e5SSatish Balay #undef __FUNC__
8145615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
81562fef451SLois Curfman McInnes /*@
81662fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
81762fef451SLois Curfman McInnes    set with SNESSetJacobian().
81862fef451SLois Curfman McInnes 
81962fef451SLois Curfman McInnes    Input Parameters:
82062fef451SLois Curfman McInnes .  snes - the SNES context
82162fef451SLois Curfman McInnes .  x - input vector
82262fef451SLois Curfman McInnes 
82362fef451SLois Curfman McInnes    Output Parameters:
82462fef451SLois Curfman McInnes .  A - Jacobian matrix
82562fef451SLois Curfman McInnes .  B - optional preconditioning matrix
82662fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
82762fef451SLois Curfman McInnes 
82862fef451SLois Curfman McInnes    Notes:
82962fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
83062fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
83162fef451SLois Curfman McInnes 
832dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
833dc5a77f8SLois Curfman McInnes    flag parameter.
83462fef451SLois Curfman McInnes 
83562fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
83662fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
837005c665bSBarry Smith    methods is SNESComputeHessian().
83862fef451SLois Curfman McInnes 
83962fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
84062fef451SLois Curfman McInnes 
84162fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
84262fef451SLois Curfman McInnes @*/
8439b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8449b94acceSBarry Smith {
8459b94acceSBarry Smith   int    ierr;
8463a40ed3dSBarry Smith 
8473a40ed3dSBarry Smith   PetscFunctionBegin;
8483a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
849e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
8503a40ed3dSBarry Smith   }
8513a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
8529b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
853c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
854d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
8559b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
856d64ed03dSBarry Smith   PetscStackPop;
8579b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8586d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
85977c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
86077c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8613a40ed3dSBarry Smith   PetscFunctionReturn(0);
8629b94acceSBarry Smith }
8639b94acceSBarry Smith 
8645615d1e5SSatish Balay #undef __FUNC__
8655615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
86662fef451SLois Curfman McInnes /*@
86762fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
86862fef451SLois Curfman McInnes    set with SNESSetHessian().
86962fef451SLois Curfman McInnes 
87062fef451SLois Curfman McInnes    Input Parameters:
87162fef451SLois Curfman McInnes .  snes - the SNES context
87262fef451SLois Curfman McInnes .  x - input vector
87362fef451SLois Curfman McInnes 
87462fef451SLois Curfman McInnes    Output Parameters:
87562fef451SLois Curfman McInnes .  A - Hessian matrix
87662fef451SLois Curfman McInnes .  B - optional preconditioning matrix
87762fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
87862fef451SLois Curfman McInnes 
87962fef451SLois Curfman McInnes    Notes:
88062fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
88162fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
88262fef451SLois Curfman McInnes 
883dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
884dc5a77f8SLois Curfman McInnes    flag parameter.
88562fef451SLois Curfman McInnes 
88662fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
88762fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
88862fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
88962fef451SLois Curfman McInnes 
89062fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
89162fef451SLois Curfman McInnes 
892a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
893a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
89462fef451SLois Curfman McInnes @*/
89562fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
8969b94acceSBarry Smith {
8979b94acceSBarry Smith   int    ierr;
8983a40ed3dSBarry Smith 
8993a40ed3dSBarry Smith   PetscFunctionBegin;
9003a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
901e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9023a40ed3dSBarry Smith   }
9033a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
90462fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9050f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
906d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
90762fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
908d64ed03dSBarry Smith   PetscStackPop;
90962fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9100f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
91177c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
91277c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9133a40ed3dSBarry Smith   PetscFunctionReturn(0);
9149b94acceSBarry Smith }
9159b94acceSBarry Smith 
9165615d1e5SSatish Balay #undef __FUNC__
917d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian"
9189b94acceSBarry Smith /*@C
9199b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
920044dda88SLois Curfman McInnes    location to store the matrix.
9219b94acceSBarry Smith 
9229b94acceSBarry Smith    Input Parameters:
9239b94acceSBarry Smith .  snes - the SNES context
9249b94acceSBarry Smith .  A - Jacobian matrix
9259b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9269b94acceSBarry Smith .  func - Jacobian evaluation routine
9272cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9282cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9299b94acceSBarry Smith 
9309b94acceSBarry Smith    Calling sequence of func:
931313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9329b94acceSBarry Smith 
9339b94acceSBarry Smith .  x - input vector
9349b94acceSBarry Smith .  A - Jacobian matrix
9359b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
936ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
937ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9382cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
9399b94acceSBarry Smith 
9409b94acceSBarry Smith    Notes:
941dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9422cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
943ac21db08SLois Curfman McInnes 
944ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9459b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9469b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9479b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9489b94acceSBarry Smith    throughout the global iterations.
9499b94acceSBarry Smith 
9509b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9519b94acceSBarry Smith 
952ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
9539b94acceSBarry Smith @*/
9549b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9559b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9569b94acceSBarry Smith {
9573a40ed3dSBarry Smith   PetscFunctionBegin;
95877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
959d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
9609b94acceSBarry Smith   snes->computejacobian = func;
9619b94acceSBarry Smith   snes->jacP            = ctx;
9629b94acceSBarry Smith   snes->jacobian        = A;
9639b94acceSBarry Smith   snes->jacobian_pre    = B;
9643a40ed3dSBarry Smith   PetscFunctionReturn(0);
9659b94acceSBarry Smith }
96662fef451SLois Curfman McInnes 
9675615d1e5SSatish Balay #undef __FUNC__
968d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian"
969b4fd4287SBarry Smith /*@
970b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
971b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
972b4fd4287SBarry Smith 
973b4fd4287SBarry Smith    Input Parameter:
974b4fd4287SBarry Smith .  snes - the nonlinear solver context
975b4fd4287SBarry Smith 
976b4fd4287SBarry Smith    Output Parameters:
977b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
978b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
979b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
980b4fd4287SBarry Smith 
981b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
982b4fd4287SBarry Smith @*/
983b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
984b4fd4287SBarry Smith {
9853a40ed3dSBarry Smith   PetscFunctionBegin;
9863a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
987e3372554SBarry Smith     SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
9883a40ed3dSBarry Smith   }
989b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
990b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
991b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
9923a40ed3dSBarry Smith   PetscFunctionReturn(0);
993b4fd4287SBarry Smith }
994b4fd4287SBarry Smith 
9955615d1e5SSatish Balay #undef __FUNC__
996d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian"
9979b94acceSBarry Smith /*@C
9989b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
999044dda88SLois Curfman McInnes    location to store the matrix.
10009b94acceSBarry Smith 
10019b94acceSBarry Smith    Input Parameters:
10029b94acceSBarry Smith .  snes - the SNES context
10039b94acceSBarry Smith .  A - Hessian matrix
10049b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
10059b94acceSBarry Smith .  func - Jacobian evaluation routine
1006313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1007313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10089b94acceSBarry Smith 
10099b94acceSBarry Smith    Calling sequence of func:
1010313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10119b94acceSBarry Smith 
10129b94acceSBarry Smith .  x - input vector
10139b94acceSBarry Smith .  A - Hessian matrix
10149b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1015ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1016ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10172cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10189b94acceSBarry Smith 
10199b94acceSBarry Smith    Notes:
1020dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10212cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1022ac21db08SLois Curfman McInnes 
10239b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10249b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10259b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10269b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10279b94acceSBarry Smith    throughout the global iterations.
10289b94acceSBarry Smith 
10299b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10309b94acceSBarry Smith 
1031ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10329b94acceSBarry Smith @*/
10339b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10349b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10359b94acceSBarry Smith {
10363a40ed3dSBarry Smith   PetscFunctionBegin;
103777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1038d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1039e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1040d4bb536fSBarry Smith   }
10419b94acceSBarry Smith   snes->computejacobian = func;
10429b94acceSBarry Smith   snes->jacP            = ctx;
10439b94acceSBarry Smith   snes->jacobian        = A;
10449b94acceSBarry Smith   snes->jacobian_pre    = B;
10453a40ed3dSBarry Smith   PetscFunctionReturn(0);
10469b94acceSBarry Smith }
10479b94acceSBarry Smith 
10485615d1e5SSatish Balay #undef __FUNC__
1049d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian"
105062fef451SLois Curfman McInnes /*@
105162fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
105262fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
105362fef451SLois Curfman McInnes 
105462fef451SLois Curfman McInnes    Input Parameter:
105562fef451SLois Curfman McInnes .  snes - the nonlinear solver context
105662fef451SLois Curfman McInnes 
105762fef451SLois Curfman McInnes    Output Parameters:
105862fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
105962fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
106062fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
106162fef451SLois Curfman McInnes 
106262fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
106362fef451SLois Curfman McInnes @*/
106462fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
106562fef451SLois Curfman McInnes {
10663a40ed3dSBarry Smith   PetscFunctionBegin;
10673a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1068e3372554SBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10693a40ed3dSBarry Smith   }
107062fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
107162fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
107262fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
10733a40ed3dSBarry Smith   PetscFunctionReturn(0);
107462fef451SLois Curfman McInnes }
107562fef451SLois Curfman McInnes 
10769b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10779b94acceSBarry Smith 
10785615d1e5SSatish Balay #undef __FUNC__
10795615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
10809b94acceSBarry Smith /*@
10819b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1082272ac6f2SLois Curfman McInnes    of a nonlinear solver.
10839b94acceSBarry Smith 
10849b94acceSBarry Smith    Input Parameter:
10859b94acceSBarry Smith .  snes - the SNES context
10868ddd3da0SLois Curfman McInnes .  x - the solution vector
10879b94acceSBarry Smith 
1088272ac6f2SLois Curfman McInnes    Notes:
1089272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1090272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1091272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1092272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1093272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1094272ac6f2SLois Curfman McInnes 
10959b94acceSBarry Smith .keywords: SNES, nonlinear, setup
10969b94acceSBarry Smith 
10979b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
10989b94acceSBarry Smith @*/
10998ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
11009b94acceSBarry Smith {
1101272ac6f2SLois Curfman McInnes   int ierr, flg;
11023a40ed3dSBarry Smith 
11033a40ed3dSBarry Smith   PetscFunctionBegin;
110477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
110577c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
11068ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
11079b94acceSBarry Smith 
1108c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1109c1f60f51SBarry Smith   /*
1110c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1111dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1112c1f60f51SBarry Smith   */
1113c1f60f51SBarry Smith   if (flg) {
1114c1f60f51SBarry Smith     Mat J;
1115c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1116c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1117c1f60f51SBarry Smith     snes->mfshell = J;
1118c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1119c1f60f51SBarry Smith       snes->jacobian = J;
1120c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1121d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1122c1f60f51SBarry Smith       snes->jacobian = J;
1123c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1124d4bb536fSBarry Smith     } else {
1125d4bb536fSBarry Smith       SETERRQ(1,0,"Method class doesn't support matrix-free operator option");
1126d4bb536fSBarry Smith     }
1127c1f60f51SBarry Smith   }
1128272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1129c1f60f51SBarry Smith   /*
1130dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1131c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1132c1f60f51SBarry Smith    */
113331615d04SLois Curfman McInnes   if (flg) {
1134272ac6f2SLois Curfman McInnes     Mat J;
1135272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1136272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1137272ac6f2SLois Curfman McInnes     snes->mfshell = J;
113883e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
113983e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
114094a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1141d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
114283e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
114394a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1144d4bb536fSBarry Smith     } else {
1145d4bb536fSBarry Smith       SETERRQ(1,0,"Method class doesn't support matrix-free option");
1146d4bb536fSBarry Smith     }
1147272ac6f2SLois Curfman McInnes   }
11489b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1149e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first");
1150e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first");
1151e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first");
1152e3372554SBarry Smith     if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector");
1153a703fe33SLois Curfman McInnes 
1154a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
115540191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1156a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1157a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1158a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1159a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1160a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1161a703fe33SLois Curfman McInnes     }
11629b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1163e3372554SBarry Smith     if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first");
1164e3372554SBarry Smith     if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first");
1165d4bb536fSBarry Smith     if (!snes->computeumfunction) SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first");
1166e3372554SBarry Smith     if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first");
1167d4bb536fSBarry Smith   } else {
1168d4bb536fSBarry Smith     SETERRQ(1,0,"Unknown method class");
1169d4bb536fSBarry Smith   }
1170a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1171a703fe33SLois Curfman McInnes   snes->setup_called = 1;
11723a40ed3dSBarry Smith   PetscFunctionReturn(0);
11739b94acceSBarry Smith }
11749b94acceSBarry Smith 
11755615d1e5SSatish Balay #undef __FUNC__
1176d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy"
11779b94acceSBarry Smith /*@C
11789b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11799b94acceSBarry Smith    with SNESCreate().
11809b94acceSBarry Smith 
11819b94acceSBarry Smith    Input Parameter:
11829b94acceSBarry Smith .  snes - the SNES context
11839b94acceSBarry Smith 
11849b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
11859b94acceSBarry Smith 
118663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
11879b94acceSBarry Smith @*/
11889b94acceSBarry Smith int SNESDestroy(SNES snes)
11899b94acceSBarry Smith {
11909b94acceSBarry Smith   int ierr;
11913a40ed3dSBarry Smith 
11923a40ed3dSBarry Smith   PetscFunctionBegin;
119377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
11943a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1195d4bb536fSBarry Smith 
11969750a799SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);}
11970452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
11989b94acceSBarry Smith   if (snes->mfshell) MatDestroy(snes->mfshell);
11999b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
12003e2e452bSBarry Smith   if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor);
12016f24a144SLois Curfman McInnes   if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork);
12029b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
12030452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
12043a40ed3dSBarry Smith   PetscFunctionReturn(0);
12059b94acceSBarry Smith }
12069b94acceSBarry Smith 
12079b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12089b94acceSBarry Smith 
12095615d1e5SSatish Balay #undef __FUNC__
12105615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
12119b94acceSBarry Smith /*@
1212d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12139b94acceSBarry Smith 
12149b94acceSBarry Smith    Input Parameters:
12159b94acceSBarry Smith .  snes - the SNES context
121633174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
121733174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
121833174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
121933174efeSLois Curfman McInnes            of the change in the solution between steps
122033174efeSLois Curfman McInnes .  maxit - maximum number of iterations
122133174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
12229b94acceSBarry Smith 
122333174efeSLois Curfman McInnes    Options Database Keys:
122433174efeSLois Curfman McInnes $    -snes_atol <atol>
122533174efeSLois Curfman McInnes $    -snes_rtol <rtol>
122633174efeSLois Curfman McInnes $    -snes_stol <stol>
122733174efeSLois Curfman McInnes $    -snes_max_it <maxit>
122833174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
12299b94acceSBarry Smith 
1230d7a720efSLois Curfman McInnes    Notes:
12319b94acceSBarry Smith    The default maximum number of iterations is 50.
12329b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
12339b94acceSBarry Smith 
123433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
12359b94acceSBarry Smith 
1236d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
12379b94acceSBarry Smith @*/
1238d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
12399b94acceSBarry Smith {
12403a40ed3dSBarry Smith   PetscFunctionBegin;
124177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1242d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1243d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1244d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1245d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1246d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
12473a40ed3dSBarry Smith   PetscFunctionReturn(0);
12489b94acceSBarry Smith }
12499b94acceSBarry Smith 
12505615d1e5SSatish Balay #undef __FUNC__
12515615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
12529b94acceSBarry Smith /*@
125333174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
125433174efeSLois Curfman McInnes 
125533174efeSLois Curfman McInnes    Input Parameters:
125633174efeSLois Curfman McInnes .  snes - the SNES context
125733174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
125833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
125933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
126033174efeSLois Curfman McInnes            of the change in the solution between steps
126133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
126233174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
126333174efeSLois Curfman McInnes 
126433174efeSLois Curfman McInnes    Notes:
126533174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
126633174efeSLois Curfman McInnes 
126733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
126833174efeSLois Curfman McInnes 
126933174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
127033174efeSLois Curfman McInnes @*/
127133174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
127233174efeSLois Curfman McInnes {
12733a40ed3dSBarry Smith   PetscFunctionBegin;
127433174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
127533174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
127633174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
127733174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
127833174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
127933174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
12803a40ed3dSBarry Smith   PetscFunctionReturn(0);
128133174efeSLois Curfman McInnes }
128233174efeSLois Curfman McInnes 
12835615d1e5SSatish Balay #undef __FUNC__
12845615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
128533174efeSLois Curfman McInnes /*@
12869b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12879b94acceSBarry Smith 
12889b94acceSBarry Smith    Input Parameters:
12899b94acceSBarry Smith .  snes - the SNES context
12909b94acceSBarry Smith .  tol - tolerance
12919b94acceSBarry Smith 
12929b94acceSBarry Smith    Options Database Key:
1293d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
12949b94acceSBarry Smith 
12959b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12969b94acceSBarry Smith 
1297d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
12989b94acceSBarry Smith @*/
12999b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
13009b94acceSBarry Smith {
13013a40ed3dSBarry Smith   PetscFunctionBegin;
130277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13039b94acceSBarry Smith   snes->deltatol = tol;
13043a40ed3dSBarry Smith   PetscFunctionReturn(0);
13059b94acceSBarry Smith }
13069b94acceSBarry Smith 
13075615d1e5SSatish Balay #undef __FUNC__
13085615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
13099b94acceSBarry Smith /*@
131077c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
13119b94acceSBarry Smith    for unconstrained minimization solvers.
13129b94acceSBarry Smith 
13139b94acceSBarry Smith    Input Parameters:
13149b94acceSBarry Smith .  snes - the SNES context
13159b94acceSBarry Smith .  ftol - minimum function tolerance
13169b94acceSBarry Smith 
13179b94acceSBarry Smith    Options Database Key:
1318d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
13199b94acceSBarry Smith 
13209b94acceSBarry Smith    Note:
132177c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
13229b94acceSBarry Smith    methods only.
13239b94acceSBarry Smith 
13249b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
13259b94acceSBarry Smith 
1326d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
13279b94acceSBarry Smith @*/
132877c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
13299b94acceSBarry Smith {
13303a40ed3dSBarry Smith   PetscFunctionBegin;
133177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13329b94acceSBarry Smith   snes->fmin = ftol;
13333a40ed3dSBarry Smith   PetscFunctionReturn(0);
13349b94acceSBarry Smith }
13359b94acceSBarry Smith 
13369b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
13379b94acceSBarry Smith 
13385615d1e5SSatish Balay #undef __FUNC__
1339d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor"
13409b94acceSBarry Smith /*@C
13419b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
13429b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
13439b94acceSBarry Smith    progress.
13449b94acceSBarry Smith 
13459b94acceSBarry Smith    Input Parameters:
13469b94acceSBarry Smith .  snes - the SNES context
13479b94acceSBarry Smith .  func - monitoring routine
1348044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1349044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
13509b94acceSBarry Smith 
13519b94acceSBarry Smith    Calling sequence of func:
1352682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
13539b94acceSBarry Smith 
13549b94acceSBarry Smith $    snes - the SNES context
13559b94acceSBarry Smith $    its - iteration number
1356044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
13579b94acceSBarry Smith $
13589b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13599b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
13609b94acceSBarry Smith $
13619b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13629b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
13639b94acceSBarry Smith 
13649665c990SLois Curfman McInnes    Options Database Keys:
13659665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
13669665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
13679665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
13689665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
13699665c990SLois Curfman McInnes $                           been hardwired into a code by
13709665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
13719665c990SLois Curfman McInnes $                           does not cancel those set via
13729665c990SLois Curfman McInnes $                           the options database.
13739665c990SLois Curfman McInnes 
13749665c990SLois Curfman McInnes 
1375639f9d9dSBarry Smith    Notes:
13766bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13776bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13786bc08f3fSLois Curfman McInnes    order in which they were set.
1379639f9d9dSBarry Smith 
13809b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13819b94acceSBarry Smith 
13829b94acceSBarry Smith .seealso: SNESDefaultMonitor()
13839b94acceSBarry Smith @*/
138474679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
13859b94acceSBarry Smith {
13863a40ed3dSBarry Smith   PetscFunctionBegin;
1387639f9d9dSBarry Smith   if (!func) {
1388639f9d9dSBarry Smith     snes->numbermonitors = 0;
13893a40ed3dSBarry Smith     PetscFunctionReturn(0);
1390639f9d9dSBarry Smith   }
1391639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1392e3372554SBarry Smith     SETERRQ(1,0,"Too many monitors set");
1393639f9d9dSBarry Smith   }
1394639f9d9dSBarry Smith 
1395639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1396639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
13973a40ed3dSBarry Smith   PetscFunctionReturn(0);
13989b94acceSBarry Smith }
13999b94acceSBarry Smith 
14005615d1e5SSatish Balay #undef __FUNC__
1401d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest"
14029b94acceSBarry Smith /*@C
14039b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
14049b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
14059b94acceSBarry Smith 
14069b94acceSBarry Smith    Input Parameters:
14079b94acceSBarry Smith .  snes - the SNES context
14089b94acceSBarry Smith .  func - routine to test for convergence
1409044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1410044dda88SLois Curfman McInnes           (may be PETSC_NULL)
14119b94acceSBarry Smith 
14129b94acceSBarry Smith    Calling sequence of func:
14139b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
14149b94acceSBarry Smith              double f,void *cctx)
14159b94acceSBarry Smith 
14169b94acceSBarry Smith $    snes - the SNES context
1417044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
14189b94acceSBarry Smith $    xnorm - 2-norm of current iterate
14199b94acceSBarry Smith $
14209b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14219b94acceSBarry Smith $    gnorm - 2-norm of current step
14229b94acceSBarry Smith $    f - 2-norm of function
14239b94acceSBarry Smith $
14249b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14259b94acceSBarry Smith $    gnorm - 2-norm of current gradient
14269b94acceSBarry Smith $    f - function value
14279b94acceSBarry Smith 
14289b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14299b94acceSBarry Smith 
143040191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
143140191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
14329b94acceSBarry Smith @*/
143374679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
14349b94acceSBarry Smith {
14353a40ed3dSBarry Smith   PetscFunctionBegin;
14369b94acceSBarry Smith   (snes)->converged = func;
14379b94acceSBarry Smith   (snes)->cnvP      = cctx;
14383a40ed3dSBarry Smith   PetscFunctionReturn(0);
14399b94acceSBarry Smith }
14409b94acceSBarry Smith 
14415615d1e5SSatish Balay #undef __FUNC__
1442d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory"
1443c9005455SLois Curfman McInnes /*@
1444c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1445c9005455SLois Curfman McInnes 
1446c9005455SLois Curfman McInnes    Input Parameters:
1447c9005455SLois Curfman McInnes .  snes - iterative context obtained from SNESCreate()
1448c9005455SLois Curfman McInnes .  a   - array to hold history
1449c9005455SLois Curfman McInnes .  na  - size of a
1450c9005455SLois Curfman McInnes 
1451c9005455SLois Curfman McInnes    Notes:
1452c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1453c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1454c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1455c9005455SLois Curfman McInnes    at each step.
1456c9005455SLois Curfman McInnes 
1457c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1458c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1459c9005455SLois Curfman McInnes    during the section of code that is being timed.
1460c9005455SLois Curfman McInnes 
1461c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1462c9005455SLois Curfman McInnes @*/
1463c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1464c9005455SLois Curfman McInnes {
14653a40ed3dSBarry Smith   PetscFunctionBegin;
1466c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1467c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1468c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1469c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
14703a40ed3dSBarry Smith   PetscFunctionReturn(0);
1471c9005455SLois Curfman McInnes }
1472c9005455SLois Curfman McInnes 
1473c9005455SLois Curfman McInnes #undef __FUNC__
14745615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
14759b94acceSBarry Smith /*
14769b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
14779b94acceSBarry Smith    positive parameter delta.
14789b94acceSBarry Smith 
14799b94acceSBarry Smith     Input Parameters:
14809b94acceSBarry Smith .   snes - the SNES context
14819b94acceSBarry Smith .   y - approximate solution of linear system
14829b94acceSBarry Smith .   fnorm - 2-norm of current function
14839b94acceSBarry Smith .   delta - trust region size
14849b94acceSBarry Smith 
14859b94acceSBarry Smith     Output Parameters:
14869b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
14879b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
14889b94acceSBarry Smith     region, and exceeds zero otherwise.
14899b94acceSBarry Smith .   ynorm - 2-norm of the step
14909b94acceSBarry Smith 
14919b94acceSBarry Smith     Note:
149240191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
14939b94acceSBarry Smith     is set to be the maximum allowable step size.
14949b94acceSBarry Smith 
14959b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
14969b94acceSBarry Smith */
14979b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
14989b94acceSBarry Smith                           double *gpnorm,double *ynorm)
14999b94acceSBarry Smith {
15009b94acceSBarry Smith   double norm;
15019b94acceSBarry Smith   Scalar cnorm;
15023a40ed3dSBarry Smith   int    ierr;
15033a40ed3dSBarry Smith 
15043a40ed3dSBarry Smith   PetscFunctionBegin;
15053a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
15069b94acceSBarry Smith   if (norm > *delta) {
15079b94acceSBarry Smith      norm = *delta/norm;
15089b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
15099b94acceSBarry Smith      cnorm = norm;
15109b94acceSBarry Smith      VecScale( &cnorm, y );
15119b94acceSBarry Smith      *ynorm = *delta;
15129b94acceSBarry Smith   } else {
15139b94acceSBarry Smith      *gpnorm = 0.0;
15149b94acceSBarry Smith      *ynorm = norm;
15159b94acceSBarry Smith   }
15163a40ed3dSBarry Smith   PetscFunctionReturn(0);
15179b94acceSBarry Smith }
15189b94acceSBarry Smith 
15195615d1e5SSatish Balay #undef __FUNC__
15205615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
15219b94acceSBarry Smith /*@
15229b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
152363a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
15249b94acceSBarry Smith 
15259b94acceSBarry Smith    Input Parameter:
15269b94acceSBarry Smith .  snes - the SNES context
15278ddd3da0SLois Curfman McInnes .  x - the solution vector
15289b94acceSBarry Smith 
15299b94acceSBarry Smith    Output Parameter:
15309b94acceSBarry Smith    its - number of iterations until termination
15319b94acceSBarry Smith 
15328ddd3da0SLois Curfman McInnes    Note:
15338ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
15348ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
15358ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
15368ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
15378ddd3da0SLois Curfman McInnes 
15389b94acceSBarry Smith .keywords: SNES, nonlinear, solve
15399b94acceSBarry Smith 
154063a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
15419b94acceSBarry Smith @*/
15428ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
15439b94acceSBarry Smith {
15443c7409f5SSatish Balay   int ierr, flg;
1545052efed2SBarry Smith 
15463a40ed3dSBarry Smith   PetscFunctionBegin;
154777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
154874679c65SBarry Smith   PetscValidIntPointer(its);
1549c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1550c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
15519b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1552c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
15539b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
15549b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
15553c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
15566d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
15573a40ed3dSBarry Smith   PetscFunctionReturn(0);
15589b94acceSBarry Smith }
15599b94acceSBarry Smith 
15609b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
15619b94acceSBarry Smith 
15625615d1e5SSatish Balay #undef __FUNC__
15635615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
15649b94acceSBarry Smith /*@
15654b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
15669b94acceSBarry Smith 
15679b94acceSBarry Smith    Input Parameters:
15689b94acceSBarry Smith .  snes - the SNES context
15699b94acceSBarry Smith .  method - a known method
15709b94acceSBarry Smith 
1571ae12b187SLois Curfman McInnes   Options Database Command:
1572ae12b187SLois Curfman McInnes $ -snes_type  <method>
1573ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1574ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1575ae12b187SLois Curfman McInnes 
15769b94acceSBarry Smith    Notes:
15779b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
15789b94acceSBarry Smith $  Systems of nonlinear equations:
157940191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
158040191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
15819b94acceSBarry Smith $  Unconstrained minimization:
158240191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
158340191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
15849b94acceSBarry Smith 
1585ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1586ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1587ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1588ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1589ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1590ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1591ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1592ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1593ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1594ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1595a703fe33SLois Curfman McInnes 
1596f536c699SSatish Balay .keywords: SNES, set, method
15979b94acceSBarry Smith @*/
15984b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
15999b94acceSBarry Smith {
160084cb2905SBarry Smith   int ierr;
16019b94acceSBarry Smith   int (*r)(SNES);
16023a40ed3dSBarry Smith 
16033a40ed3dSBarry Smith   PetscFunctionBegin;
160477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16053a40ed3dSBarry Smith   if (snes->type == (int) method) PetscFunctionReturn(0);
16067d1a2b2bSBarry Smith 
16079b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
160884cb2905SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
1609e3372554SBarry Smith   if (!__SNESList) {SETERRQ(1,0,"Could not get methods");}
161037fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1611e3372554SBarry Smith   if (!r) {SETERRQ(1,0,"Unknown method");}
16120452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
16139b94acceSBarry Smith   snes->set_method_called = 1;
16143a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
16153a40ed3dSBarry Smith   PetscFunctionReturn(0);
16169b94acceSBarry Smith }
16179b94acceSBarry Smith 
16189b94acceSBarry Smith /* --------------------------------------------------------------------- */
16195615d1e5SSatish Balay #undef __FUNC__
1620d4bb536fSBarry Smith #define __FUNC__ "SNESRegister"
16219b94acceSBarry Smith /*@C
16229b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
16234b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
16249b94acceSBarry Smith 
16259b94acceSBarry Smith    Input Parameters:
16262d872ea7SLois Curfman McInnes .  name - either a predefined name such as SNES_EQ_LS, or SNES_NEW
16272d872ea7SLois Curfman McInnes           to indicate a new user-defined solver
162840191667SLois Curfman McInnes .  sname - corresponding string for name
16299b94acceSBarry Smith .  create - routine to create method context
16309b94acceSBarry Smith 
163184cb2905SBarry Smith    Output Parameter:
163284cb2905SBarry Smith .  oname - type associated with this new solver
163384cb2905SBarry Smith 
16342d872ea7SLois Curfman McInnes    Notes:
16352d872ea7SLois Curfman McInnes    Multiple user-defined nonlinear solvers can be added by calling
16362d872ea7SLois Curfman McInnes    SNESRegister() with the input parameter "name" set to be SNES_NEW;
16372d872ea7SLois Curfman McInnes    each call will return a unique solver type in the output
16382d872ea7SLois Curfman McInnes    parameter "oname".
16392d872ea7SLois Curfman McInnes 
16409b94acceSBarry Smith .keywords: SNES, nonlinear, register
16419b94acceSBarry Smith 
16429b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
16439b94acceSBarry Smith @*/
164484cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES))
16459b94acceSBarry Smith {
16469b94acceSBarry Smith   int        ierr;
164784cb2905SBarry Smith   static int numberregistered = 0;
164884cb2905SBarry Smith 
16493a40ed3dSBarry Smith   PetscFunctionBegin;
1650d252947aSBarry Smith   if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++);
165184cb2905SBarry Smith 
165284cb2905SBarry Smith   if (oname) *oname = name;
165337fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
165484cb2905SBarry Smith   NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create );
16553a40ed3dSBarry Smith   PetscFunctionReturn(0);
16569b94acceSBarry Smith }
1657a847f771SSatish Balay 
16589b94acceSBarry Smith /* --------------------------------------------------------------------- */
16595615d1e5SSatish Balay #undef __FUNC__
1660d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy"
16619b94acceSBarry Smith /*@C
16629b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
16639b94acceSBarry Smith    registered by SNESRegister().
16649b94acceSBarry Smith 
16659b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
16669b94acceSBarry Smith 
16679b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
16689b94acceSBarry Smith @*/
16699b94acceSBarry Smith int SNESRegisterDestroy()
16709b94acceSBarry Smith {
16713a40ed3dSBarry Smith   PetscFunctionBegin;
167237fcc0dbSBarry Smith   if (__SNESList) {
167337fcc0dbSBarry Smith     NRDestroy( __SNESList );
167437fcc0dbSBarry Smith     __SNESList = 0;
16759b94acceSBarry Smith   }
167684cb2905SBarry Smith   SNESRegisterAllCalled = 0;
16773a40ed3dSBarry Smith   PetscFunctionReturn(0);
16789b94acceSBarry Smith }
16799b94acceSBarry Smith 
16805615d1e5SSatish Balay #undef __FUNC__
1681d4bb536fSBarry Smith #define __FUNC__ "SNESGetType"
16829b94acceSBarry Smith /*@C
16839a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
16849b94acceSBarry Smith 
16859b94acceSBarry Smith    Input Parameter:
16864b0e389bSBarry Smith .  snes - nonlinear solver context
16879b94acceSBarry Smith 
16889b94acceSBarry Smith    Output Parameter:
16899a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
16909a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
16919b94acceSBarry Smith 
16929b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
16939b94acceSBarry Smith @*/
16944b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
16959b94acceSBarry Smith {
169606a9b0adSLois Curfman McInnes   int ierr;
16973a40ed3dSBarry Smith 
16983a40ed3dSBarry Smith   PetscFunctionBegin;
169937fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
17004b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
170137fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
17023a40ed3dSBarry Smith   PetscFunctionReturn(0);
17039b94acceSBarry Smith }
17049b94acceSBarry Smith 
17055615d1e5SSatish Balay #undef __FUNC__
1706d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution"
17079b94acceSBarry Smith /*@C
17089b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
17099b94acceSBarry Smith    stored.
17109b94acceSBarry Smith 
17119b94acceSBarry Smith    Input Parameter:
17129b94acceSBarry Smith .  snes - the SNES context
17139b94acceSBarry Smith 
17149b94acceSBarry Smith    Output Parameter:
17159b94acceSBarry Smith .  x - the solution
17169b94acceSBarry Smith 
17179b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
17189b94acceSBarry Smith 
1719a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
17209b94acceSBarry Smith @*/
17219b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
17229b94acceSBarry Smith {
17233a40ed3dSBarry Smith   PetscFunctionBegin;
172477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17259b94acceSBarry Smith   *x = snes->vec_sol_always;
17263a40ed3dSBarry Smith   PetscFunctionReturn(0);
17279b94acceSBarry Smith }
17289b94acceSBarry Smith 
17295615d1e5SSatish Balay #undef __FUNC__
1730d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate"
17319b94acceSBarry Smith /*@C
17329b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
17339b94acceSBarry Smith    stored.
17349b94acceSBarry Smith 
17359b94acceSBarry Smith    Input Parameter:
17369b94acceSBarry Smith .  snes - the SNES context
17379b94acceSBarry Smith 
17389b94acceSBarry Smith    Output Parameter:
17399b94acceSBarry Smith .  x - the solution update
17409b94acceSBarry Smith 
17419b94acceSBarry Smith    Notes:
17429b94acceSBarry Smith    This vector is implementation dependent.
17439b94acceSBarry Smith 
17449b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
17459b94acceSBarry Smith 
17469b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
17479b94acceSBarry Smith @*/
17489b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
17499b94acceSBarry Smith {
17503a40ed3dSBarry Smith   PetscFunctionBegin;
175177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17529b94acceSBarry Smith   *x = snes->vec_sol_update_always;
17533a40ed3dSBarry Smith   PetscFunctionReturn(0);
17549b94acceSBarry Smith }
17559b94acceSBarry Smith 
17565615d1e5SSatish Balay #undef __FUNC__
1757d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction"
17589b94acceSBarry Smith /*@C
17593638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
17609b94acceSBarry Smith 
17619b94acceSBarry Smith    Input Parameter:
17629b94acceSBarry Smith .  snes - the SNES context
17639b94acceSBarry Smith 
17649b94acceSBarry Smith    Output Parameter:
17653638b69dSLois Curfman McInnes .  r - the function
17669b94acceSBarry Smith 
17679b94acceSBarry Smith    Notes:
17689b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
17699b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
17709b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
17719b94acceSBarry Smith 
1772a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
17739b94acceSBarry Smith 
177431615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
177531615d04SLois Curfman McInnes           SNESGetGradient()
17769b94acceSBarry Smith @*/
17779b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
17789b94acceSBarry Smith {
17793a40ed3dSBarry Smith   PetscFunctionBegin;
178077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17813a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
17829b94acceSBarry Smith   *r = snes->vec_func_always;
17833a40ed3dSBarry Smith   PetscFunctionReturn(0);
17849b94acceSBarry Smith }
17859b94acceSBarry Smith 
17865615d1e5SSatish Balay #undef __FUNC__
1787d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient"
17889b94acceSBarry Smith /*@C
17893638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
17909b94acceSBarry Smith 
17919b94acceSBarry Smith    Input Parameter:
17929b94acceSBarry Smith .  snes - the SNES context
17939b94acceSBarry Smith 
17949b94acceSBarry Smith    Output Parameter:
17959b94acceSBarry Smith .  r - the gradient
17969b94acceSBarry Smith 
17979b94acceSBarry Smith    Notes:
17989b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
17999b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18009b94acceSBarry Smith    SNESGetFunction().
18019b94acceSBarry Smith 
18029b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
18039b94acceSBarry Smith 
180431615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
18059b94acceSBarry Smith @*/
18069b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
18079b94acceSBarry Smith {
18083a40ed3dSBarry Smith   PetscFunctionBegin;
180977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18103a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
18113a40ed3dSBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
18123a40ed3dSBarry Smith   }
18139b94acceSBarry Smith   *r = snes->vec_func_always;
18143a40ed3dSBarry Smith   PetscFunctionReturn(0);
18159b94acceSBarry Smith }
18169b94acceSBarry Smith 
18175615d1e5SSatish Balay #undef __FUNC__
1818d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction"
18199b94acceSBarry Smith /*@
18209b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
18219b94acceSBarry Smith    unconstrained minimization problems.
18229b94acceSBarry Smith 
18239b94acceSBarry Smith    Input Parameter:
18249b94acceSBarry Smith .  snes - the SNES context
18259b94acceSBarry Smith 
18269b94acceSBarry Smith    Output Parameter:
18279b94acceSBarry Smith .  r - the function
18289b94acceSBarry Smith 
18299b94acceSBarry Smith    Notes:
18309b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
18319b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18329b94acceSBarry Smith    SNESGetFunction().
18339b94acceSBarry Smith 
18349b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
18359b94acceSBarry Smith 
183631615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
18379b94acceSBarry Smith @*/
18389b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
18399b94acceSBarry Smith {
18403a40ed3dSBarry Smith   PetscFunctionBegin;
184177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
184274679c65SBarry Smith   PetscValidScalarPointer(r);
18433a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
18443a40ed3dSBarry Smith     SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
18453a40ed3dSBarry Smith   }
18469b94acceSBarry Smith   *r = snes->fc;
18473a40ed3dSBarry Smith   PetscFunctionReturn(0);
18489b94acceSBarry Smith }
18499b94acceSBarry Smith 
18505615d1e5SSatish Balay #undef __FUNC__
1851d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix"
18523c7409f5SSatish Balay /*@C
18533c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1854d850072dSLois Curfman McInnes    SNES options in the database.
18553c7409f5SSatish Balay 
18563c7409f5SSatish Balay    Input Parameter:
18573c7409f5SSatish Balay .  snes - the SNES context
18583c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18593c7409f5SSatish Balay 
1860d850072dSLois Curfman McInnes    Notes:
1861a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1862a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1863a83b1b31SSatish Balay    hyphen.
1864d850072dSLois Curfman McInnes 
18653c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1866a86d99e1SLois Curfman McInnes 
1867a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
18683c7409f5SSatish Balay @*/
18693c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
18703c7409f5SSatish Balay {
18713c7409f5SSatish Balay   int ierr;
18723c7409f5SSatish Balay 
18733a40ed3dSBarry Smith   PetscFunctionBegin;
187477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1875639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18763c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
18773a40ed3dSBarry Smith   PetscFunctionReturn(0);
18783c7409f5SSatish Balay }
18793c7409f5SSatish Balay 
18805615d1e5SSatish Balay #undef __FUNC__
1881d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix"
18823c7409f5SSatish Balay /*@C
1883f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1884d850072dSLois Curfman McInnes    SNES options in the database.
18853c7409f5SSatish Balay 
18863c7409f5SSatish Balay    Input Parameter:
18873c7409f5SSatish Balay .  snes - the SNES context
18883c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18893c7409f5SSatish Balay 
1890d850072dSLois Curfman McInnes    Notes:
1891a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1892a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1893a83b1b31SSatish Balay    hyphen.
1894d850072dSLois Curfman McInnes 
18953c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1896a86d99e1SLois Curfman McInnes 
1897a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
18983c7409f5SSatish Balay @*/
18993c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
19003c7409f5SSatish Balay {
19013c7409f5SSatish Balay   int ierr;
19023c7409f5SSatish Balay 
19033a40ed3dSBarry Smith   PetscFunctionBegin;
190477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1905639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19063c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19073a40ed3dSBarry Smith   PetscFunctionReturn(0);
19083c7409f5SSatish Balay }
19093c7409f5SSatish Balay 
19105615d1e5SSatish Balay #undef __FUNC__
1911d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix"
1912c83590e2SSatish Balay /*@
19133c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
19143c7409f5SSatish Balay    SNES options in the database.
19153c7409f5SSatish Balay 
19163c7409f5SSatish Balay    Input Parameter:
19173c7409f5SSatish Balay .  snes - the SNES context
19183c7409f5SSatish Balay 
19193c7409f5SSatish Balay    Output Parameter:
19203c7409f5SSatish Balay .  prefix - pointer to the prefix string used
19213c7409f5SSatish Balay 
19223c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1923a86d99e1SLois Curfman McInnes 
1924a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
19253c7409f5SSatish Balay @*/
19263c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
19273c7409f5SSatish Balay {
19283c7409f5SSatish Balay   int ierr;
19293c7409f5SSatish Balay 
19303a40ed3dSBarry Smith   PetscFunctionBegin;
193177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1932639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19333a40ed3dSBarry Smith   PetscFunctionReturn(0);
19343c7409f5SSatish Balay }
19353c7409f5SSatish Balay 
1936*ca161407SBarry Smith #undef __FUNC__
1937*ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp"
1938*ca161407SBarry Smith /*@
1939*ca161407SBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
1940*ca161407SBarry Smith 
1941*ca161407SBarry Smith    Input Parameter:
1942*ca161407SBarry Smith .  snes - the SNES context
1943*ca161407SBarry Smith 
1944*ca161407SBarry Smith    Options Database Keys:
1945*ca161407SBarry Smith $  -help, -h
1946*ca161407SBarry Smith 
1947*ca161407SBarry Smith .keywords: SNES, nonlinear, help
1948*ca161407SBarry Smith 
1949*ca161407SBarry Smith .seealso: SNESSetFromOptions()
1950*ca161407SBarry Smith @*/
1951*ca161407SBarry Smith int SNESPrintHelp(SNES snes)
1952*ca161407SBarry Smith {
1953*ca161407SBarry Smith   char                p[64];
1954*ca161407SBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
1955*ca161407SBarry Smith   int                 ierr;
1956*ca161407SBarry Smith 
1957*ca161407SBarry Smith   PetscFunctionBegin;
1958*ca161407SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1959*ca161407SBarry Smith 
1960*ca161407SBarry Smith   PetscStrcpy(p,"-");
1961*ca161407SBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
1962*ca161407SBarry Smith 
1963*ca161407SBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
1964*ca161407SBarry Smith 
1965*ca161407SBarry Smith   PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n");
1966*ca161407SBarry Smith   ierr = NRPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",__SNESList);CHKERRQ(ierr);
1967*ca161407SBarry Smith   PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
1968*ca161407SBarry Smith   PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
1969*ca161407SBarry Smith   PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
1970*ca161407SBarry Smith   PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
1971*ca161407SBarry Smith   PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
1972*ca161407SBarry Smith   PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
1973*ca161407SBarry Smith   PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
1974*ca161407SBarry Smith   PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n");
1975*ca161407SBarry Smith   PetscPrintf(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
1976*ca161407SBarry Smith   PetscPrintf(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
1977*ca161407SBarry Smith     residual norm at each iteration.\n",p);
1978*ca161407SBarry Smith   PetscPrintf(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
1979*ca161407SBarry Smith     residual norm for small residual norms. This is useful to conceal\n\
1980*ca161407SBarry Smith     meaningless digits that may be different on different machines.\n",p);
1981*ca161407SBarry Smith   PetscPrintf(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
1982*ca161407SBarry Smith   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
1983*ca161407SBarry Smith     PetscPrintf(snes->comm,
1984*ca161407SBarry Smith      " Options for solving systems of nonlinear equations only:\n");
1985*ca161407SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
1986*ca161407SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
1987*ca161407SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
1988*ca161407SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
1989*ca161407SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
1990*ca161407SBarry Smith     PetscPrintf(snes->comm,
1991*ca161407SBarry Smith      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
1992*ca161407SBarry Smith     PetscPrintf(snes->comm,
1993*ca161407SBarry Smith      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
1994*ca161407SBarry Smith     PetscPrintf(snes->comm,
1995*ca161407SBarry Smith      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
1996*ca161407SBarry Smith     PetscPrintf(snes->comm,
1997*ca161407SBarry Smith      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
1998*ca161407SBarry Smith     PetscPrintf(snes->comm,
1999*ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
2000*ca161407SBarry Smith     PetscPrintf(snes->comm,
2001*ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
2002*ca161407SBarry Smith     PetscPrintf(snes->comm,
2003*ca161407SBarry Smith      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
2004*ca161407SBarry Smith   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
2005*ca161407SBarry Smith     PetscPrintf(snes->comm,
2006*ca161407SBarry Smith      " Options for solving unconstrained minimization problems only:\n");
2007*ca161407SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
2008*ca161407SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
2009*ca161407SBarry Smith     PetscPrintf(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
2010*ca161407SBarry Smith   }
2011*ca161407SBarry Smith   PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
2012*ca161407SBarry Smith   PetscPrintf(snes->comm,"a particular method\n");
2013*ca161407SBarry Smith   if (snes->printhelp) {
2014*ca161407SBarry Smith     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2015*ca161407SBarry Smith   }
2016*ca161407SBarry Smith   PetscFunctionReturn(0);
2017*ca161407SBarry Smith }
20183c7409f5SSatish Balay 
20199b94acceSBarry Smith 
20209b94acceSBarry Smith 
20219b94acceSBarry Smith 
20223a40ed3dSBarry Smith 
2023