xref: /petsc/src/snes/interface/snes.c (revision 76be9ce4a233aaa47cda2bc7f5a27cd7faabecaa) !
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*76be9ce4SBarry Smith static char vcid[] = "$Id: snes.c,v 1.136 1998/01/06 20:12:19 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;
11ca161407SBarry 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) {
117a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,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);
158a8c6a408SBarry Smith   if (snes->setup_called) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp");
159ca161407SBarry Smith 
160ca161407SBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll();CHKERRQ(ierr);}
161ca161407SBarry Smith   ierr = NRGetTypeFromOptions(snes->prefix,"-snes_type",__SNESList,&method,&flg);CHKERRQ(ierr);
162052efed2SBarry Smith   if (flg) {
163052efed2SBarry Smith     ierr = SNESSetType(snes,method); CHKERRQ(ierr);
164ca161407SBarry 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);
232981c4779SBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
233981c4779SBarry Smith   } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
23431615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
23531615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
236d64ed03dSBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2379b94acceSBarry Smith   }
238639f9d9dSBarry Smith 
239639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
240639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
241639f9d9dSBarry Smith   }
242639f9d9dSBarry Smith 
2439b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2449b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
2453a40ed3dSBarry Smith   if (snes->setfromoptions) {
2463a40ed3dSBarry Smith     ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
2473a40ed3dSBarry Smith   }
2483a40ed3dSBarry Smith   PetscFunctionReturn(0);
2499b94acceSBarry Smith }
2509b94acceSBarry Smith 
251a847f771SSatish Balay 
2525615d1e5SSatish Balay #undef __FUNC__
253d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext"
2549b94acceSBarry Smith /*@
2559b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2569b94acceSBarry Smith    the nonlinear solvers.
2579b94acceSBarry Smith 
2589b94acceSBarry Smith    Input Parameters:
2599b94acceSBarry Smith .  snes - the SNES context
2609b94acceSBarry Smith .  usrP - optional user context
2619b94acceSBarry Smith 
2629b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2639b94acceSBarry Smith 
2649b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2659b94acceSBarry Smith @*/
2669b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
2679b94acceSBarry Smith {
2683a40ed3dSBarry Smith   PetscFunctionBegin;
26977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2709b94acceSBarry Smith   snes->user		= usrP;
2713a40ed3dSBarry Smith   PetscFunctionReturn(0);
2729b94acceSBarry Smith }
27374679c65SBarry Smith 
2745615d1e5SSatish Balay #undef __FUNC__
275d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext"
2769b94acceSBarry Smith /*@C
2779b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
2789b94acceSBarry Smith    nonlinear solvers.
2799b94acceSBarry Smith 
2809b94acceSBarry Smith    Input Parameter:
2819b94acceSBarry Smith .  snes - SNES context
2829b94acceSBarry Smith 
2839b94acceSBarry Smith    Output Parameter:
2849b94acceSBarry Smith .  usrP - user context
2859b94acceSBarry Smith 
2869b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
2879b94acceSBarry Smith 
2889b94acceSBarry Smith .seealso: SNESSetApplicationContext()
2899b94acceSBarry Smith @*/
2909b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
2919b94acceSBarry Smith {
2923a40ed3dSBarry Smith   PetscFunctionBegin;
29377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2949b94acceSBarry Smith   *usrP = snes->user;
2953a40ed3dSBarry Smith   PetscFunctionReturn(0);
2969b94acceSBarry Smith }
29774679c65SBarry Smith 
2985615d1e5SSatish Balay #undef __FUNC__
299d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber"
3009b94acceSBarry Smith /*@
3019b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3029b94acceSBarry Smith    nonlinear solver.
3039b94acceSBarry Smith 
3049b94acceSBarry Smith    Input Parameter:
3059b94acceSBarry Smith .  snes - SNES context
3069b94acceSBarry Smith 
3079b94acceSBarry Smith    Output Parameter:
3089b94acceSBarry Smith .  iter - iteration number
3099b94acceSBarry Smith 
3109b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3119b94acceSBarry Smith @*/
3129b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3139b94acceSBarry Smith {
3143a40ed3dSBarry Smith   PetscFunctionBegin;
31577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
31674679c65SBarry Smith   PetscValidIntPointer(iter);
3179b94acceSBarry Smith   *iter = snes->iter;
3183a40ed3dSBarry Smith   PetscFunctionReturn(0);
3199b94acceSBarry Smith }
32074679c65SBarry Smith 
3215615d1e5SSatish Balay #undef __FUNC__
3225615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
3239b94acceSBarry Smith /*@
3249b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3259b94acceSBarry Smith    with SNESSSetFunction().
3269b94acceSBarry Smith 
3279b94acceSBarry Smith    Input Parameter:
3289b94acceSBarry Smith .  snes - SNES context
3299b94acceSBarry Smith 
3309b94acceSBarry Smith    Output Parameter:
3319b94acceSBarry Smith .  fnorm - 2-norm of function
3329b94acceSBarry Smith 
3339b94acceSBarry Smith    Note:
3349b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
335a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
336a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3379b94acceSBarry Smith 
3389b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
339a86d99e1SLois Curfman McInnes 
340a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3419b94acceSBarry Smith @*/
3429b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3439b94acceSBarry Smith {
3443a40ed3dSBarry Smith   PetscFunctionBegin;
34577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
34674679c65SBarry Smith   PetscValidScalarPointer(fnorm);
34774679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
348d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only");
34974679c65SBarry Smith   }
3509b94acceSBarry Smith   *fnorm = snes->norm;
3513a40ed3dSBarry Smith   PetscFunctionReturn(0);
3529b94acceSBarry Smith }
35374679c65SBarry Smith 
3545615d1e5SSatish Balay #undef __FUNC__
3555615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
3569b94acceSBarry Smith /*@
3579b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
3589b94acceSBarry Smith    with SNESSSetGradient().
3599b94acceSBarry Smith 
3609b94acceSBarry Smith    Input Parameter:
3619b94acceSBarry Smith .  snes - SNES context
3629b94acceSBarry Smith 
3639b94acceSBarry Smith    Output Parameter:
3649b94acceSBarry Smith .  fnorm - 2-norm of gradient
3659b94acceSBarry Smith 
3669b94acceSBarry Smith    Note:
3679b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
368a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
369a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
3709b94acceSBarry Smith 
3719b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
372a86d99e1SLois Curfman McInnes 
373a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
3749b94acceSBarry Smith @*/
3759b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
3769b94acceSBarry Smith {
3773a40ed3dSBarry Smith   PetscFunctionBegin;
37877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
37974679c65SBarry Smith   PetscValidScalarPointer(gnorm);
38074679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
381d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
38274679c65SBarry Smith   }
3839b94acceSBarry Smith   *gnorm = snes->norm;
3843a40ed3dSBarry Smith   PetscFunctionReturn(0);
3859b94acceSBarry Smith }
38674679c65SBarry Smith 
3875615d1e5SSatish Balay #undef __FUNC__
388d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps"
3899b94acceSBarry Smith /*@
3909b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
3919b94acceSBarry Smith    attempted by the nonlinear solver.
3929b94acceSBarry Smith 
3939b94acceSBarry Smith    Input Parameter:
3949b94acceSBarry Smith .  snes - SNES context
3959b94acceSBarry Smith 
3969b94acceSBarry Smith    Output Parameter:
3979b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
3989b94acceSBarry Smith 
399c96a6f78SLois Curfman McInnes    Notes:
400c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
401c96a6f78SLois Curfman McInnes 
4029b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4039b94acceSBarry Smith @*/
4049b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4059b94acceSBarry Smith {
4063a40ed3dSBarry Smith   PetscFunctionBegin;
40777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
40874679c65SBarry Smith   PetscValidIntPointer(nfails);
4099b94acceSBarry Smith   *nfails = snes->nfailures;
4103a40ed3dSBarry Smith   PetscFunctionReturn(0);
4119b94acceSBarry Smith }
412a847f771SSatish Balay 
4135615d1e5SSatish Balay #undef __FUNC__
414d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations"
415c96a6f78SLois Curfman McInnes /*@
416c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
417c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
418c96a6f78SLois Curfman McInnes 
419c96a6f78SLois Curfman McInnes    Input Parameter:
420c96a6f78SLois Curfman McInnes .  snes - SNES context
421c96a6f78SLois Curfman McInnes 
422c96a6f78SLois Curfman McInnes    Output Parameter:
423c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
424c96a6f78SLois Curfman McInnes 
425c96a6f78SLois Curfman McInnes    Notes:
426c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
427c96a6f78SLois Curfman McInnes 
428c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
429c96a6f78SLois Curfman McInnes @*/
43086bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
431c96a6f78SLois Curfman McInnes {
4323a40ed3dSBarry Smith   PetscFunctionBegin;
433c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
434c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
435c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
4363a40ed3dSBarry Smith   PetscFunctionReturn(0);
437c96a6f78SLois Curfman McInnes }
438c96a6f78SLois Curfman McInnes 
439c96a6f78SLois Curfman McInnes #undef __FUNC__
440d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES"
4419b94acceSBarry Smith /*@C
4429b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4439b94acceSBarry Smith 
4449b94acceSBarry Smith    Input Parameter:
4459b94acceSBarry Smith .  snes - the SNES context
4469b94acceSBarry Smith 
4479b94acceSBarry Smith    Output Parameter:
4489b94acceSBarry Smith .  sles - the SLES context
4499b94acceSBarry Smith 
4509b94acceSBarry Smith    Notes:
4519b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4529b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4539b94acceSBarry Smith    KSP and PC contexts as well.
4549b94acceSBarry Smith 
4559b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4569b94acceSBarry Smith 
4579b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4589b94acceSBarry Smith @*/
4599b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4609b94acceSBarry Smith {
4613a40ed3dSBarry Smith   PetscFunctionBegin;
46277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4639b94acceSBarry Smith   *sles = snes->sles;
4643a40ed3dSBarry Smith   PetscFunctionReturn(0);
4659b94acceSBarry Smith }
4669b94acceSBarry Smith /* -----------------------------------------------------------*/
4675615d1e5SSatish Balay #undef __FUNC__
4685615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
4699b94acceSBarry Smith /*@C
4709b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
4719b94acceSBarry Smith 
4729b94acceSBarry Smith    Input Parameter:
4739b94acceSBarry Smith .  comm - MPI communicator
4749b94acceSBarry Smith .  type - type of method, one of
4759b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
4769b94acceSBarry Smith $      (for systems of nonlinear equations)
4779b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
4789b94acceSBarry Smith $      (for unconstrained minimization)
4799b94acceSBarry Smith 
4809b94acceSBarry Smith    Output Parameter:
4819b94acceSBarry Smith .  outsnes - the new SNES context
4829b94acceSBarry Smith 
483c1f60f51SBarry Smith    Options Database Key:
48419bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
48519bd6747SLois Curfman McInnes $              and no preconditioning matrix
48619bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
48719bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
48819bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
48919bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
490c1f60f51SBarry Smith 
4919b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
4929b94acceSBarry Smith 
49363a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
4949b94acceSBarry Smith @*/
4954b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
4969b94acceSBarry Smith {
4979b94acceSBarry Smith   int                 ierr;
4989b94acceSBarry Smith   SNES                snes;
4999b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
50037fcc0dbSBarry Smith 
5013a40ed3dSBarry Smith   PetscFunctionBegin;
5029b94acceSBarry Smith   *outsnes = 0;
503d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
504d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
505d64ed03dSBarry Smith   }
506d4bb536fSBarry Smith   PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm,SNESDestroy,SNESView);
5079b94acceSBarry Smith   PLogObjectCreate(snes);
5089b94acceSBarry Smith   snes->max_its           = 50;
5099750a799SBarry Smith   snes->max_funcs	  = 10000;
5109b94acceSBarry Smith   snes->norm		  = 0.0;
5115a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
512b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
513b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5149b94acceSBarry Smith     snes->atol		  = 1.e-10;
515d64ed03dSBarry Smith   } else {
516b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
517b4874afaSBarry Smith     snes->ttol            = 0.0;
518b4874afaSBarry Smith     snes->atol		  = 1.e-50;
519b4874afaSBarry Smith   }
5209b94acceSBarry Smith   snes->xtol		  = 1.e-8;
521b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5229b94acceSBarry Smith   snes->nfuncs            = 0;
5239b94acceSBarry Smith   snes->nfailures         = 0;
5247a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
525639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5269b94acceSBarry Smith   snes->data              = 0;
5279b94acceSBarry Smith   snes->view              = 0;
5289b94acceSBarry Smith   snes->computeumfunction = 0;
5299b94acceSBarry Smith   snes->umfunP            = 0;
5309b94acceSBarry Smith   snes->fc                = 0;
5319b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5329b94acceSBarry Smith   snes->fmin              = -1.e30;
5339b94acceSBarry Smith   snes->method_class      = type;
5349b94acceSBarry Smith   snes->set_method_called = 0;
535a703fe33SLois Curfman McInnes   snes->setup_called      = 0;
5369b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5377d1a2b2bSBarry Smith   snes->type              = -1;
5386f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5396f24a144SLois Curfman McInnes   snes->vwork             = 0;
5406f24a144SLois Curfman McInnes   snes->nwork             = 0;
541c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
542c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
543c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
5449b94acceSBarry Smith 
5459b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5460452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
547eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
5489b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5499b94acceSBarry Smith   kctx->version     = 2;
5509b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5519b94acceSBarry Smith                              this was too large for some test cases */
5529b94acceSBarry Smith   kctx->rtol_last   = 0;
5539b94acceSBarry Smith   kctx->rtol_max    = .9;
5549b94acceSBarry Smith   kctx->gamma       = 1.0;
5559b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5569b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5579b94acceSBarry Smith   kctx->threshold   = .1;
5589b94acceSBarry Smith   kctx->lresid_last = 0;
5599b94acceSBarry Smith   kctx->norm_last   = 0;
5609b94acceSBarry Smith 
5619b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5629b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5635334005bSBarry Smith 
5649b94acceSBarry Smith   *outsnes = snes;
5653a40ed3dSBarry Smith   PetscFunctionReturn(0);
5669b94acceSBarry Smith }
5679b94acceSBarry Smith 
5689b94acceSBarry Smith /* --------------------------------------------------------------- */
5695615d1e5SSatish Balay #undef __FUNC__
570d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction"
5719b94acceSBarry Smith /*@C
5729b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
5739b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
5749b94acceSBarry Smith    equations.
5759b94acceSBarry Smith 
5769b94acceSBarry Smith    Input Parameters:
5779b94acceSBarry Smith .  snes - the SNES context
5789b94acceSBarry Smith .  func - function evaluation routine
5799b94acceSBarry Smith .  r - vector to store function value
5802cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
5812cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
5829b94acceSBarry Smith 
5839b94acceSBarry Smith    Calling sequence of func:
584313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
5859b94acceSBarry Smith 
5869b94acceSBarry Smith .  x - input vector
587313e4042SLois Curfman McInnes .  f - function vector
5882cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
5899b94acceSBarry Smith 
5909b94acceSBarry Smith    Notes:
5919b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
5929b94acceSBarry Smith $      f'(x) x = -f(x),
5939b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
5949b94acceSBarry Smith 
5959b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
5969b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
5979b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
5989b94acceSBarry Smith 
5999b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6009b94acceSBarry Smith 
601a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6029b94acceSBarry Smith @*/
6035334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6049b94acceSBarry Smith {
6053a40ed3dSBarry Smith   PetscFunctionBegin;
60677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
607a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
608a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
609a8c6a408SBarry Smith   }
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) {
644a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,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);
689a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
690a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
691a8c6a408SBarry Smith   }
6929b94acceSBarry Smith   snes->computeumfunction   = func;
6939b94acceSBarry Smith   snes->umfunP              = ctx;
6943a40ed3dSBarry Smith   PetscFunctionReturn(0);
6959b94acceSBarry Smith }
6969b94acceSBarry Smith 
6975615d1e5SSatish Balay #undef __FUNC__
6985615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
6999b94acceSBarry Smith /*@
7009b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7019b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7029b94acceSBarry Smith 
7039b94acceSBarry Smith    Input Parameters:
7049b94acceSBarry Smith .  snes - the SNES context
7059b94acceSBarry Smith .  x - input vector
7069b94acceSBarry Smith 
7079b94acceSBarry Smith    Output Parameter:
7089b94acceSBarry Smith .  y - function value
7099b94acceSBarry Smith 
7109b94acceSBarry Smith    Notes:
7119b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7129b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7139b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
714a86d99e1SLois Curfman McInnes 
715a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
716a86d99e1SLois Curfman McInnes 
717a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
718a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7199b94acceSBarry Smith @*/
7209b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7219b94acceSBarry Smith {
7229b94acceSBarry Smith   int    ierr;
7233a40ed3dSBarry Smith 
7243a40ed3dSBarry Smith   PetscFunctionBegin;
725a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
726a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
727a8c6a408SBarry Smith   }
7289b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
729d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
7309b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
731d64ed03dSBarry Smith   PetscStackPop;
732ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7339b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7343a40ed3dSBarry Smith   PetscFunctionReturn(0);
7359b94acceSBarry Smith }
7369b94acceSBarry Smith 
7375615d1e5SSatish Balay #undef __FUNC__
738d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient"
7399b94acceSBarry Smith /*@C
7409b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7419b94acceSBarry Smith    vector for use by the SNES routines.
7429b94acceSBarry Smith 
7439b94acceSBarry Smith    Input Parameters:
7449b94acceSBarry Smith .  snes - the SNES context
7459b94acceSBarry Smith .  func - function evaluation routine
746044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
747044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7489b94acceSBarry Smith .  r - vector to store gradient value
7499b94acceSBarry Smith 
7509b94acceSBarry Smith    Calling sequence of func:
7519b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7529b94acceSBarry Smith 
7539b94acceSBarry Smith .  x - input vector
7549b94acceSBarry Smith .  g - gradient vector
755044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7569b94acceSBarry Smith 
7579b94acceSBarry Smith    Notes:
7589b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7599b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7609b94acceSBarry Smith    SNESSetFunction().
7619b94acceSBarry Smith 
7629b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7639b94acceSBarry Smith 
764a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
765a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
7669b94acceSBarry Smith @*/
76774679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
7689b94acceSBarry Smith {
7693a40ed3dSBarry Smith   PetscFunctionBegin;
77077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
771a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
772a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
773a8c6a408SBarry Smith   }
7749b94acceSBarry Smith   snes->computefunction     = func;
7759b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7769b94acceSBarry Smith   snes->funP                = ctx;
7773a40ed3dSBarry Smith   PetscFunctionReturn(0);
7789b94acceSBarry Smith }
7799b94acceSBarry Smith 
7805615d1e5SSatish Balay #undef __FUNC__
7815615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
7829b94acceSBarry Smith /*@
783a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
784a86d99e1SLois Curfman McInnes    SNESSetGradient().
7859b94acceSBarry Smith 
7869b94acceSBarry Smith    Input Parameters:
7879b94acceSBarry Smith .  snes - the SNES context
7889b94acceSBarry Smith .  x - input vector
7899b94acceSBarry Smith 
7909b94acceSBarry Smith    Output Parameter:
7919b94acceSBarry Smith .  y - gradient vector
7929b94acceSBarry Smith 
7939b94acceSBarry Smith    Notes:
7949b94acceSBarry Smith    SNESComputeGradient() is valid only for
7959b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7969b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
797a86d99e1SLois Curfman McInnes 
798a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
799a86d99e1SLois Curfman McInnes 
800a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
801a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8029b94acceSBarry Smith @*/
8039b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8049b94acceSBarry Smith {
8059b94acceSBarry Smith   int    ierr;
8063a40ed3dSBarry Smith 
8073a40ed3dSBarry Smith   PetscFunctionBegin;
8083a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
809a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8103a40ed3dSBarry Smith   }
8119b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
812d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
8139b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
814d64ed03dSBarry Smith   PetscStackPop;
8159b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8163a40ed3dSBarry Smith   PetscFunctionReturn(0);
8179b94acceSBarry Smith }
8189b94acceSBarry Smith 
8195615d1e5SSatish Balay #undef __FUNC__
8205615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
82162fef451SLois Curfman McInnes /*@
82262fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
82362fef451SLois Curfman McInnes    set with SNESSetJacobian().
82462fef451SLois Curfman McInnes 
82562fef451SLois Curfman McInnes    Input Parameters:
82662fef451SLois Curfman McInnes .  snes - the SNES context
82762fef451SLois Curfman McInnes .  x - input vector
82862fef451SLois Curfman McInnes 
82962fef451SLois Curfman McInnes    Output Parameters:
83062fef451SLois Curfman McInnes .  A - Jacobian matrix
83162fef451SLois Curfman McInnes .  B - optional preconditioning matrix
83262fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
83362fef451SLois Curfman McInnes 
83462fef451SLois Curfman McInnes    Notes:
83562fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
83662fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
83762fef451SLois Curfman McInnes 
838dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
839dc5a77f8SLois Curfman McInnes    flag parameter.
84062fef451SLois Curfman McInnes 
84162fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
84262fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
843005c665bSBarry Smith    methods is SNESComputeHessian().
84462fef451SLois Curfman McInnes 
84562fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
84662fef451SLois Curfman McInnes 
84762fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
84862fef451SLois Curfman McInnes @*/
8499b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8509b94acceSBarry Smith {
8519b94acceSBarry Smith   int    ierr;
8523a40ed3dSBarry Smith 
8533a40ed3dSBarry Smith   PetscFunctionBegin;
8543a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
855a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
8563a40ed3dSBarry Smith   }
8573a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
8589b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
859c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
860d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
8619b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
862d64ed03dSBarry Smith   PetscStackPop;
8639b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8646d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
86577c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
86677c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8673a40ed3dSBarry Smith   PetscFunctionReturn(0);
8689b94acceSBarry Smith }
8699b94acceSBarry Smith 
8705615d1e5SSatish Balay #undef __FUNC__
8715615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
87262fef451SLois Curfman McInnes /*@
87362fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
87462fef451SLois Curfman McInnes    set with SNESSetHessian().
87562fef451SLois Curfman McInnes 
87662fef451SLois Curfman McInnes    Input Parameters:
87762fef451SLois Curfman McInnes .  snes - the SNES context
87862fef451SLois Curfman McInnes .  x - input vector
87962fef451SLois Curfman McInnes 
88062fef451SLois Curfman McInnes    Output Parameters:
88162fef451SLois Curfman McInnes .  A - Hessian matrix
88262fef451SLois Curfman McInnes .  B - optional preconditioning matrix
88362fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
88462fef451SLois Curfman McInnes 
88562fef451SLois Curfman McInnes    Notes:
88662fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
88762fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
88862fef451SLois Curfman McInnes 
889dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
890dc5a77f8SLois Curfman McInnes    flag parameter.
89162fef451SLois Curfman McInnes 
89262fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
89362fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
89462fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
89562fef451SLois Curfman McInnes 
89662fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
89762fef451SLois Curfman McInnes 
898a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
899a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
90062fef451SLois Curfman McInnes @*/
90162fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
9029b94acceSBarry Smith {
9039b94acceSBarry Smith   int    ierr;
9043a40ed3dSBarry Smith 
9053a40ed3dSBarry Smith   PetscFunctionBegin;
9063a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
907a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9083a40ed3dSBarry Smith   }
9093a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
91062fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9110f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
912d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
91362fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
914d64ed03dSBarry Smith   PetscStackPop;
91562fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9160f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
91777c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
91877c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9193a40ed3dSBarry Smith   PetscFunctionReturn(0);
9209b94acceSBarry Smith }
9219b94acceSBarry Smith 
9225615d1e5SSatish Balay #undef __FUNC__
923d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian"
9249b94acceSBarry Smith /*@C
9259b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
926044dda88SLois Curfman McInnes    location to store the matrix.
9279b94acceSBarry Smith 
9289b94acceSBarry Smith    Input Parameters:
9299b94acceSBarry Smith .  snes - the SNES context
9309b94acceSBarry Smith .  A - Jacobian matrix
9319b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9329b94acceSBarry Smith .  func - Jacobian evaluation routine
9332cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9342cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9359b94acceSBarry Smith 
9369b94acceSBarry Smith    Calling sequence of func:
937313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9389b94acceSBarry Smith 
9399b94acceSBarry Smith .  x - input vector
9409b94acceSBarry Smith .  A - Jacobian matrix
9419b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
942ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
943ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9442cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
9459b94acceSBarry Smith 
9469b94acceSBarry Smith    Notes:
947dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9482cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
949ac21db08SLois Curfman McInnes 
950ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9519b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9529b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9539b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9549b94acceSBarry Smith    throughout the global iterations.
9559b94acceSBarry Smith 
9569b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9579b94acceSBarry Smith 
958ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
9599b94acceSBarry Smith @*/
9609b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9619b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9629b94acceSBarry Smith {
9633a40ed3dSBarry Smith   PetscFunctionBegin;
96477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
965a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
966a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
967a8c6a408SBarry Smith   }
9689b94acceSBarry Smith   snes->computejacobian = func;
9699b94acceSBarry Smith   snes->jacP            = ctx;
9709b94acceSBarry Smith   snes->jacobian        = A;
9719b94acceSBarry Smith   snes->jacobian_pre    = B;
9723a40ed3dSBarry Smith   PetscFunctionReturn(0);
9739b94acceSBarry Smith }
97462fef451SLois Curfman McInnes 
9755615d1e5SSatish Balay #undef __FUNC__
976d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian"
977b4fd4287SBarry Smith /*@
978b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
979b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
980b4fd4287SBarry Smith 
981b4fd4287SBarry Smith    Input Parameter:
982b4fd4287SBarry Smith .  snes - the nonlinear solver context
983b4fd4287SBarry Smith 
984b4fd4287SBarry Smith    Output Parameters:
985b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
986b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
987b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
988b4fd4287SBarry Smith 
989b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
990b4fd4287SBarry Smith @*/
991b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
992b4fd4287SBarry Smith {
9933a40ed3dSBarry Smith   PetscFunctionBegin;
9943a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
995a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
9963a40ed3dSBarry Smith   }
997b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
998b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
999b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
10003a40ed3dSBarry Smith   PetscFunctionReturn(0);
1001b4fd4287SBarry Smith }
1002b4fd4287SBarry Smith 
10035615d1e5SSatish Balay #undef __FUNC__
1004d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian"
10059b94acceSBarry Smith /*@C
10069b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1007044dda88SLois Curfman McInnes    location to store the matrix.
10089b94acceSBarry Smith 
10099b94acceSBarry Smith    Input Parameters:
10109b94acceSBarry Smith .  snes - the SNES context
10119b94acceSBarry Smith .  A - Hessian matrix
10129b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
10139b94acceSBarry Smith .  func - Jacobian evaluation routine
1014313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1015313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10169b94acceSBarry Smith 
10179b94acceSBarry Smith    Calling sequence of func:
1018313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10199b94acceSBarry Smith 
10209b94acceSBarry Smith .  x - input vector
10219b94acceSBarry Smith .  A - Hessian matrix
10229b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1023ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1024ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10252cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10269b94acceSBarry Smith 
10279b94acceSBarry Smith    Notes:
1028dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10292cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1030ac21db08SLois Curfman McInnes 
10319b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10329b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10339b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10349b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10359b94acceSBarry Smith    throughout the global iterations.
10369b94acceSBarry Smith 
10379b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10389b94acceSBarry Smith 
1039ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10409b94acceSBarry Smith @*/
10419b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10429b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10439b94acceSBarry Smith {
10443a40ed3dSBarry Smith   PetscFunctionBegin;
104577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1046d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1047a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1048d4bb536fSBarry Smith   }
10499b94acceSBarry Smith   snes->computejacobian = func;
10509b94acceSBarry Smith   snes->jacP            = ctx;
10519b94acceSBarry Smith   snes->jacobian        = A;
10529b94acceSBarry Smith   snes->jacobian_pre    = B;
10533a40ed3dSBarry Smith   PetscFunctionReturn(0);
10549b94acceSBarry Smith }
10559b94acceSBarry Smith 
10565615d1e5SSatish Balay #undef __FUNC__
1057d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian"
105862fef451SLois Curfman McInnes /*@
105962fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
106062fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
106162fef451SLois Curfman McInnes 
106262fef451SLois Curfman McInnes    Input Parameter:
106362fef451SLois Curfman McInnes .  snes - the nonlinear solver context
106462fef451SLois Curfman McInnes 
106562fef451SLois Curfman McInnes    Output Parameters:
106662fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
106762fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
106862fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
106962fef451SLois Curfman McInnes 
107062fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
107162fef451SLois Curfman McInnes @*/
107262fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
107362fef451SLois Curfman McInnes {
10743a40ed3dSBarry Smith   PetscFunctionBegin;
10753a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1076a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10773a40ed3dSBarry Smith   }
107862fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
107962fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
108062fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
10813a40ed3dSBarry Smith   PetscFunctionReturn(0);
108262fef451SLois Curfman McInnes }
108362fef451SLois Curfman McInnes 
10849b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10859b94acceSBarry Smith 
10865615d1e5SSatish Balay #undef __FUNC__
10875615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
10889b94acceSBarry Smith /*@
10899b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1090272ac6f2SLois Curfman McInnes    of a nonlinear solver.
10919b94acceSBarry Smith 
10929b94acceSBarry Smith    Input Parameter:
10939b94acceSBarry Smith .  snes - the SNES context
10948ddd3da0SLois Curfman McInnes .  x - the solution vector
10959b94acceSBarry Smith 
1096272ac6f2SLois Curfman McInnes    Notes:
1097272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1098272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1099272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1100272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1101272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1102272ac6f2SLois Curfman McInnes 
11039b94acceSBarry Smith .keywords: SNES, nonlinear, setup
11049b94acceSBarry Smith 
11059b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
11069b94acceSBarry Smith @*/
11078ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
11089b94acceSBarry Smith {
1109272ac6f2SLois Curfman McInnes   int ierr, flg;
11103a40ed3dSBarry Smith 
11113a40ed3dSBarry Smith   PetscFunctionBegin;
111277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
111377c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
11148ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
11159b94acceSBarry Smith 
1116c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1117c1f60f51SBarry Smith   /*
1118c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1119dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1120c1f60f51SBarry Smith   */
1121c1f60f51SBarry Smith   if (flg) {
1122c1f60f51SBarry Smith     Mat J;
1123c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1124c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1125c1f60f51SBarry Smith     snes->mfshell = J;
1126c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1127c1f60f51SBarry Smith       snes->jacobian = J;
1128c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1129d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1130c1f60f51SBarry Smith       snes->jacobian = J;
1131c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1132d4bb536fSBarry Smith     } else {
1133a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option");
1134d4bb536fSBarry Smith     }
1135c1f60f51SBarry Smith   }
1136272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1137c1f60f51SBarry Smith   /*
1138dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1139c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1140c1f60f51SBarry Smith    */
114131615d04SLois Curfman McInnes   if (flg) {
1142272ac6f2SLois Curfman McInnes     Mat J;
1143272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1144272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1145272ac6f2SLois Curfman McInnes     snes->mfshell = J;
114683e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
114783e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
114894a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1149d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
115083e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
115194a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1152d4bb536fSBarry Smith     } else {
1153a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option");
1154d4bb536fSBarry Smith     }
1155272ac6f2SLois Curfman McInnes   }
11569b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1157a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1158a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1159a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first");
1160a8c6a408SBarry Smith     if (snes->vec_func == snes->vec_sol) {
1161a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector");
1162a8c6a408SBarry Smith     }
1163a703fe33SLois Curfman McInnes 
1164a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
116540191667SLois Curfman McInnes     if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) {
1166a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1167a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1168a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1169a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1170a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1171a703fe33SLois Curfman McInnes     }
11729b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1173a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1174a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1175a8c6a408SBarry Smith     if (!snes->computeumfunction) {
1176a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first");
1177a8c6a408SBarry Smith     }
1178a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first");
1179d4bb536fSBarry Smith   } else {
1180a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
1181d4bb536fSBarry Smith   }
1182a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
1183a703fe33SLois Curfman McInnes   snes->setup_called = 1;
11843a40ed3dSBarry Smith   PetscFunctionReturn(0);
11859b94acceSBarry Smith }
11869b94acceSBarry Smith 
11875615d1e5SSatish Balay #undef __FUNC__
1188d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy"
11899b94acceSBarry Smith /*@C
11909b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11919b94acceSBarry Smith    with SNESCreate().
11929b94acceSBarry Smith 
11939b94acceSBarry Smith    Input Parameter:
11949b94acceSBarry Smith .  snes - the SNES context
11959b94acceSBarry Smith 
11969b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
11979b94acceSBarry Smith 
119863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
11999b94acceSBarry Smith @*/
12009b94acceSBarry Smith int SNESDestroy(SNES snes)
12019b94acceSBarry Smith {
12029b94acceSBarry Smith   int ierr;
12033a40ed3dSBarry Smith 
12043a40ed3dSBarry Smith   PetscFunctionBegin;
120577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12063a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1207d4bb536fSBarry Smith 
12089750a799SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);}
12090452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
1210522c5e43SBarry Smith   if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);}
12119b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
1212522c5e43SBarry Smith   if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);}
1213522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
12149b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
12150452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
12163a40ed3dSBarry Smith   PetscFunctionReturn(0);
12179b94acceSBarry Smith }
12189b94acceSBarry Smith 
12199b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12209b94acceSBarry Smith 
12215615d1e5SSatish Balay #undef __FUNC__
12225615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
12239b94acceSBarry Smith /*@
1224d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12259b94acceSBarry Smith 
12269b94acceSBarry Smith    Input Parameters:
12279b94acceSBarry Smith .  snes - the SNES context
122833174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
122933174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
123033174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
123133174efeSLois Curfman McInnes            of the change in the solution between steps
123233174efeSLois Curfman McInnes .  maxit - maximum number of iterations
123333174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
12349b94acceSBarry Smith 
123533174efeSLois Curfman McInnes    Options Database Keys:
123633174efeSLois Curfman McInnes $    -snes_atol <atol>
123733174efeSLois Curfman McInnes $    -snes_rtol <rtol>
123833174efeSLois Curfman McInnes $    -snes_stol <stol>
123933174efeSLois Curfman McInnes $    -snes_max_it <maxit>
124033174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
12419b94acceSBarry Smith 
1242d7a720efSLois Curfman McInnes    Notes:
12439b94acceSBarry Smith    The default maximum number of iterations is 50.
12449b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
12459b94acceSBarry Smith 
124633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
12479b94acceSBarry Smith 
1248d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
12499b94acceSBarry Smith @*/
1250d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
12519b94acceSBarry Smith {
12523a40ed3dSBarry Smith   PetscFunctionBegin;
125377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1254d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1255d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1256d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1257d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1258d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
12593a40ed3dSBarry Smith   PetscFunctionReturn(0);
12609b94acceSBarry Smith }
12619b94acceSBarry Smith 
12625615d1e5SSatish Balay #undef __FUNC__
12635615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
12649b94acceSBarry Smith /*@
126533174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
126633174efeSLois Curfman McInnes 
126733174efeSLois Curfman McInnes    Input Parameters:
126833174efeSLois Curfman McInnes .  snes - the SNES context
126933174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
127033174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
127133174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
127233174efeSLois Curfman McInnes            of the change in the solution between steps
127333174efeSLois Curfman McInnes .  maxit - maximum number of iterations
127433174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
127533174efeSLois Curfman McInnes 
127633174efeSLois Curfman McInnes    Notes:
127733174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
127833174efeSLois Curfman McInnes 
127933174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
128033174efeSLois Curfman McInnes 
128133174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
128233174efeSLois Curfman McInnes @*/
128333174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
128433174efeSLois Curfman McInnes {
12853a40ed3dSBarry Smith   PetscFunctionBegin;
128633174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
128733174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
128833174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
128933174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
129033174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
129133174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
12923a40ed3dSBarry Smith   PetscFunctionReturn(0);
129333174efeSLois Curfman McInnes }
129433174efeSLois Curfman McInnes 
12955615d1e5SSatish Balay #undef __FUNC__
12965615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
129733174efeSLois Curfman McInnes /*@
12989b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12999b94acceSBarry Smith 
13009b94acceSBarry Smith    Input Parameters:
13019b94acceSBarry Smith .  snes - the SNES context
13029b94acceSBarry Smith .  tol - tolerance
13039b94acceSBarry Smith 
13049b94acceSBarry Smith    Options Database Key:
1305d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
13069b94acceSBarry Smith 
13079b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
13089b94acceSBarry Smith 
1309d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
13109b94acceSBarry Smith @*/
13119b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
13129b94acceSBarry Smith {
13133a40ed3dSBarry Smith   PetscFunctionBegin;
131477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13159b94acceSBarry Smith   snes->deltatol = tol;
13163a40ed3dSBarry Smith   PetscFunctionReturn(0);
13179b94acceSBarry Smith }
13189b94acceSBarry Smith 
13195615d1e5SSatish Balay #undef __FUNC__
13205615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
13219b94acceSBarry Smith /*@
132277c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
13239b94acceSBarry Smith    for unconstrained minimization solvers.
13249b94acceSBarry Smith 
13259b94acceSBarry Smith    Input Parameters:
13269b94acceSBarry Smith .  snes - the SNES context
13279b94acceSBarry Smith .  ftol - minimum function tolerance
13289b94acceSBarry Smith 
13299b94acceSBarry Smith    Options Database Key:
1330d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
13319b94acceSBarry Smith 
13329b94acceSBarry Smith    Note:
133377c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
13349b94acceSBarry Smith    methods only.
13359b94acceSBarry Smith 
13369b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
13379b94acceSBarry Smith 
1338d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
13399b94acceSBarry Smith @*/
134077c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
13419b94acceSBarry Smith {
13423a40ed3dSBarry Smith   PetscFunctionBegin;
134377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13449b94acceSBarry Smith   snes->fmin = ftol;
13453a40ed3dSBarry Smith   PetscFunctionReturn(0);
13469b94acceSBarry Smith }
13479b94acceSBarry Smith 
13489b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
13499b94acceSBarry Smith 
13505615d1e5SSatish Balay #undef __FUNC__
1351d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor"
13529b94acceSBarry Smith /*@C
13539b94acceSBarry Smith    SNESSetMonitor - Sets the function that is to be used at every
13549b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
13559b94acceSBarry Smith    progress.
13569b94acceSBarry Smith 
13579b94acceSBarry Smith    Input Parameters:
13589b94acceSBarry Smith .  snes - the SNES context
13599b94acceSBarry Smith .  func - monitoring routine
1360044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1361044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
13629b94acceSBarry Smith 
13639b94acceSBarry Smith    Calling sequence of func:
1364682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
13659b94acceSBarry Smith 
13669b94acceSBarry Smith $    snes - the SNES context
13679b94acceSBarry Smith $    its - iteration number
1368044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
13699b94acceSBarry Smith $
13709b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13719b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
13729b94acceSBarry Smith $
13739b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13749b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
13759b94acceSBarry Smith 
13769665c990SLois Curfman McInnes    Options Database Keys:
13779665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
13789665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
13799665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
13809665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
13819665c990SLois Curfman McInnes $                           been hardwired into a code by
13829665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
13839665c990SLois Curfman McInnes $                           does not cancel those set via
13849665c990SLois Curfman McInnes $                           the options database.
13859665c990SLois Curfman McInnes 
13869665c990SLois Curfman McInnes 
1387639f9d9dSBarry Smith    Notes:
13886bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13896bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13906bc08f3fSLois Curfman McInnes    order in which they were set.
1391639f9d9dSBarry Smith 
13929b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13939b94acceSBarry Smith 
13949b94acceSBarry Smith .seealso: SNESDefaultMonitor()
13959b94acceSBarry Smith @*/
139674679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
13979b94acceSBarry Smith {
13983a40ed3dSBarry Smith   PetscFunctionBegin;
1399639f9d9dSBarry Smith   if (!func) {
1400639f9d9dSBarry Smith     snes->numbermonitors = 0;
14013a40ed3dSBarry Smith     PetscFunctionReturn(0);
1402639f9d9dSBarry Smith   }
1403639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1404a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
1405639f9d9dSBarry Smith   }
1406639f9d9dSBarry Smith 
1407639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1408639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
14093a40ed3dSBarry Smith   PetscFunctionReturn(0);
14109b94acceSBarry Smith }
14119b94acceSBarry Smith 
14125615d1e5SSatish Balay #undef __FUNC__
1413d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest"
14149b94acceSBarry Smith /*@C
14159b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
14169b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
14179b94acceSBarry Smith 
14189b94acceSBarry Smith    Input Parameters:
14199b94acceSBarry Smith .  snes - the SNES context
14209b94acceSBarry Smith .  func - routine to test for convergence
1421044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1422044dda88SLois Curfman McInnes           (may be PETSC_NULL)
14239b94acceSBarry Smith 
14249b94acceSBarry Smith    Calling sequence of func:
14259b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
14269b94acceSBarry Smith              double f,void *cctx)
14279b94acceSBarry Smith 
14289b94acceSBarry Smith $    snes - the SNES context
1429044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
14309b94acceSBarry Smith $    xnorm - 2-norm of current iterate
14319b94acceSBarry Smith $
14329b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14339b94acceSBarry Smith $    gnorm - 2-norm of current step
14349b94acceSBarry Smith $    f - 2-norm of function
14359b94acceSBarry Smith $
14369b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14379b94acceSBarry Smith $    gnorm - 2-norm of current gradient
14389b94acceSBarry Smith $    f - function value
14399b94acceSBarry Smith 
14409b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14419b94acceSBarry Smith 
144240191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
144340191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
14449b94acceSBarry Smith @*/
144574679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
14469b94acceSBarry Smith {
14473a40ed3dSBarry Smith   PetscFunctionBegin;
14489b94acceSBarry Smith   (snes)->converged = func;
14499b94acceSBarry Smith   (snes)->cnvP      = cctx;
14503a40ed3dSBarry Smith   PetscFunctionReturn(0);
14519b94acceSBarry Smith }
14529b94acceSBarry Smith 
14535615d1e5SSatish Balay #undef __FUNC__
1454d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory"
1455c9005455SLois Curfman McInnes /*@
1456c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1457c9005455SLois Curfman McInnes 
1458c9005455SLois Curfman McInnes    Input Parameters:
1459c9005455SLois Curfman McInnes .  snes - iterative context obtained from SNESCreate()
1460c9005455SLois Curfman McInnes .  a   - array to hold history
1461c9005455SLois Curfman McInnes .  na  - size of a
1462c9005455SLois Curfman McInnes 
1463c9005455SLois Curfman McInnes    Notes:
1464c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1465c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1466c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1467c9005455SLois Curfman McInnes    at each step.
1468c9005455SLois Curfman McInnes 
1469c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1470c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1471c9005455SLois Curfman McInnes    during the section of code that is being timed.
1472c9005455SLois Curfman McInnes 
1473c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1474c9005455SLois Curfman McInnes @*/
1475c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1476c9005455SLois Curfman McInnes {
14773a40ed3dSBarry Smith   PetscFunctionBegin;
1478c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1479c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1480c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1481c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
14823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1483c9005455SLois Curfman McInnes }
1484c9005455SLois Curfman McInnes 
1485c9005455SLois Curfman McInnes #undef __FUNC__
14865615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
14879b94acceSBarry Smith /*
14889b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
14899b94acceSBarry Smith    positive parameter delta.
14909b94acceSBarry Smith 
14919b94acceSBarry Smith     Input Parameters:
14929b94acceSBarry Smith .   snes - the SNES context
14939b94acceSBarry Smith .   y - approximate solution of linear system
14949b94acceSBarry Smith .   fnorm - 2-norm of current function
14959b94acceSBarry Smith .   delta - trust region size
14969b94acceSBarry Smith 
14979b94acceSBarry Smith     Output Parameters:
14989b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
14999b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
15009b94acceSBarry Smith     region, and exceeds zero otherwise.
15019b94acceSBarry Smith .   ynorm - 2-norm of the step
15029b94acceSBarry Smith 
15039b94acceSBarry Smith     Note:
150440191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
15059b94acceSBarry Smith     is set to be the maximum allowable step size.
15069b94acceSBarry Smith 
15079b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
15089b94acceSBarry Smith */
15099b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
15109b94acceSBarry Smith                           double *gpnorm,double *ynorm)
15119b94acceSBarry Smith {
15129b94acceSBarry Smith   double norm;
15139b94acceSBarry Smith   Scalar cnorm;
15143a40ed3dSBarry Smith   int    ierr;
15153a40ed3dSBarry Smith 
15163a40ed3dSBarry Smith   PetscFunctionBegin;
15173a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
15189b94acceSBarry Smith   if (norm > *delta) {
15199b94acceSBarry Smith      norm = *delta/norm;
15209b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
15219b94acceSBarry Smith      cnorm = norm;
15229b94acceSBarry Smith      VecScale( &cnorm, y );
15239b94acceSBarry Smith      *ynorm = *delta;
15249b94acceSBarry Smith   } else {
15259b94acceSBarry Smith      *gpnorm = 0.0;
15269b94acceSBarry Smith      *ynorm = norm;
15279b94acceSBarry Smith   }
15283a40ed3dSBarry Smith   PetscFunctionReturn(0);
15299b94acceSBarry Smith }
15309b94acceSBarry Smith 
15315615d1e5SSatish Balay #undef __FUNC__
15325615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
15339b94acceSBarry Smith /*@
15349b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
153563a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
15369b94acceSBarry Smith 
15379b94acceSBarry Smith    Input Parameter:
15389b94acceSBarry Smith .  snes - the SNES context
15398ddd3da0SLois Curfman McInnes .  x - the solution vector
15409b94acceSBarry Smith 
15419b94acceSBarry Smith    Output Parameter:
15429b94acceSBarry Smith    its - number of iterations until termination
15439b94acceSBarry Smith 
15448ddd3da0SLois Curfman McInnes    Note:
15458ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
15468ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
15478ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
15488ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
15498ddd3da0SLois Curfman McInnes 
15509b94acceSBarry Smith .keywords: SNES, nonlinear, solve
15519b94acceSBarry Smith 
155263a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
15539b94acceSBarry Smith @*/
15548ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
15559b94acceSBarry Smith {
15563c7409f5SSatish Balay   int ierr, flg;
1557052efed2SBarry Smith 
15583a40ed3dSBarry Smith   PetscFunctionBegin;
155977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
156074679c65SBarry Smith   PetscValidIntPointer(its);
1561c4fc05e7SBarry Smith   if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1562c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
15639b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1564c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
15659b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
15669b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
15673c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
15686d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
15693a40ed3dSBarry Smith   PetscFunctionReturn(0);
15709b94acceSBarry Smith }
15719b94acceSBarry Smith 
15729b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
15739b94acceSBarry Smith 
15745615d1e5SSatish Balay #undef __FUNC__
15755615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
15769b94acceSBarry Smith /*@
15774b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
15789b94acceSBarry Smith 
15799b94acceSBarry Smith    Input Parameters:
15809b94acceSBarry Smith .  snes - the SNES context
15819b94acceSBarry Smith .  method - a known method
15829b94acceSBarry Smith 
1583ae12b187SLois Curfman McInnes   Options Database Command:
1584ae12b187SLois Curfman McInnes $ -snes_type  <method>
1585ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1586ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1587ae12b187SLois Curfman McInnes 
15889b94acceSBarry Smith    Notes:
15899b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
15909b94acceSBarry Smith $  Systems of nonlinear equations:
159140191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
159240191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
15939b94acceSBarry Smith $  Unconstrained minimization:
159440191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
159540191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
15969b94acceSBarry Smith 
1597ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1598ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1599ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1600ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1601ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1602ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1603ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1604ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1605ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1606ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1607a703fe33SLois Curfman McInnes 
1608f536c699SSatish Balay .keywords: SNES, set, method
16099b94acceSBarry Smith @*/
16104b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
16119b94acceSBarry Smith {
161284cb2905SBarry Smith   int ierr;
16139b94acceSBarry Smith   int (*r)(SNES);
16143a40ed3dSBarry Smith 
16153a40ed3dSBarry Smith   PetscFunctionBegin;
161677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16173a40ed3dSBarry Smith   if (snes->type == (int) method) PetscFunctionReturn(0);
16187d1a2b2bSBarry Smith 
16199b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
162084cb2905SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
162137fcc0dbSBarry Smith   r =  (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 );
1622a8c6a408SBarry Smith   if (!r) {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method");}
16230452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
16249b94acceSBarry Smith   snes->set_method_called = 1;
16253a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
16263a40ed3dSBarry Smith   PetscFunctionReturn(0);
16279b94acceSBarry Smith }
16289b94acceSBarry Smith 
16299b94acceSBarry Smith /* --------------------------------------------------------------------- */
16305615d1e5SSatish Balay #undef __FUNC__
1631d4bb536fSBarry Smith #define __FUNC__ "SNESRegister"
16329b94acceSBarry Smith /*@C
16339b94acceSBarry Smith    SNESRegister - Adds the method to the nonlinear solver package, given
16344b0e389bSBarry Smith    a function pointer and a nonlinear solver name of the type SNESType.
16359b94acceSBarry Smith 
16369b94acceSBarry Smith    Input Parameters:
16372d872ea7SLois Curfman McInnes .  name - either a predefined name such as SNES_EQ_LS, or SNES_NEW
16382d872ea7SLois Curfman McInnes           to indicate a new user-defined solver
163940191667SLois Curfman McInnes .  sname - corresponding string for name
16409b94acceSBarry Smith .  create - routine to create method context
16419b94acceSBarry Smith 
164284cb2905SBarry Smith    Output Parameter:
164384cb2905SBarry Smith .  oname - type associated with this new solver
164484cb2905SBarry Smith 
16452d872ea7SLois Curfman McInnes    Notes:
16462d872ea7SLois Curfman McInnes    Multiple user-defined nonlinear solvers can be added by calling
16472d872ea7SLois Curfman McInnes    SNESRegister() with the input parameter "name" set to be SNES_NEW;
16482d872ea7SLois Curfman McInnes    each call will return a unique solver type in the output
16492d872ea7SLois Curfman McInnes    parameter "oname".
16502d872ea7SLois Curfman McInnes 
16519b94acceSBarry Smith .keywords: SNES, nonlinear, register
16529b94acceSBarry Smith 
16539b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy()
16549b94acceSBarry Smith @*/
165584cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES))
16569b94acceSBarry Smith {
16579b94acceSBarry Smith   int        ierr;
165884cb2905SBarry Smith   static int numberregistered = 0;
165984cb2905SBarry Smith 
16603a40ed3dSBarry Smith   PetscFunctionBegin;
1661d252947aSBarry Smith   if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++);
166284cb2905SBarry Smith 
166384cb2905SBarry Smith   if (oname) *oname = name;
166437fcc0dbSBarry Smith   if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);}
166584cb2905SBarry Smith   NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create );
16663a40ed3dSBarry Smith   PetscFunctionReturn(0);
16679b94acceSBarry Smith }
1668a847f771SSatish Balay 
16699b94acceSBarry Smith /* --------------------------------------------------------------------- */
16705615d1e5SSatish Balay #undef __FUNC__
1671d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy"
16729b94acceSBarry Smith /*@C
16739b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
16749b94acceSBarry Smith    registered by SNESRegister().
16759b94acceSBarry Smith 
16769b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
16779b94acceSBarry Smith 
16789b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
16799b94acceSBarry Smith @*/
16809b94acceSBarry Smith int SNESRegisterDestroy()
16819b94acceSBarry Smith {
16823a40ed3dSBarry Smith   PetscFunctionBegin;
168337fcc0dbSBarry Smith   if (__SNESList) {
168437fcc0dbSBarry Smith     NRDestroy( __SNESList );
168537fcc0dbSBarry Smith     __SNESList = 0;
16869b94acceSBarry Smith   }
168784cb2905SBarry Smith   SNESRegisterAllCalled = 0;
16883a40ed3dSBarry Smith   PetscFunctionReturn(0);
16899b94acceSBarry Smith }
16909b94acceSBarry Smith 
16915615d1e5SSatish Balay #undef __FUNC__
1692d4bb536fSBarry Smith #define __FUNC__ "SNESGetType"
16939b94acceSBarry Smith /*@C
16949a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
16959b94acceSBarry Smith 
16969b94acceSBarry Smith    Input Parameter:
16974b0e389bSBarry Smith .  snes - nonlinear solver context
16989b94acceSBarry Smith 
16999b94acceSBarry Smith    Output Parameter:
17009a28b0a6SLois Curfman McInnes .  method - SNES method (or use PETSC_NULL)
17019a28b0a6SLois Curfman McInnes .  name - name of SNES method (or use PETSC_NULL)
17029b94acceSBarry Smith 
17039b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
17049b94acceSBarry Smith @*/
17054b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name)
17069b94acceSBarry Smith {
170706a9b0adSLois Curfman McInnes   int ierr;
17083a40ed3dSBarry Smith 
17093a40ed3dSBarry Smith   PetscFunctionBegin;
171037fcc0dbSBarry Smith   if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);}
17114b0e389bSBarry Smith   if (method) *method = (SNESType) snes->type;
171237fcc0dbSBarry Smith   if (name)  *name = NRFindName( __SNESList, (int) snes->type );
17133a40ed3dSBarry Smith   PetscFunctionReturn(0);
17149b94acceSBarry Smith }
17159b94acceSBarry Smith 
17165615d1e5SSatish Balay #undef __FUNC__
1717d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution"
17189b94acceSBarry Smith /*@C
17199b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
17209b94acceSBarry Smith    stored.
17219b94acceSBarry Smith 
17229b94acceSBarry Smith    Input Parameter:
17239b94acceSBarry Smith .  snes - the SNES context
17249b94acceSBarry Smith 
17259b94acceSBarry Smith    Output Parameter:
17269b94acceSBarry Smith .  x - the solution
17279b94acceSBarry Smith 
17289b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
17299b94acceSBarry Smith 
1730a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
17319b94acceSBarry Smith @*/
17329b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
17339b94acceSBarry Smith {
17343a40ed3dSBarry Smith   PetscFunctionBegin;
173577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17369b94acceSBarry Smith   *x = snes->vec_sol_always;
17373a40ed3dSBarry Smith   PetscFunctionReturn(0);
17389b94acceSBarry Smith }
17399b94acceSBarry Smith 
17405615d1e5SSatish Balay #undef __FUNC__
1741d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate"
17429b94acceSBarry Smith /*@C
17439b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
17449b94acceSBarry Smith    stored.
17459b94acceSBarry Smith 
17469b94acceSBarry Smith    Input Parameter:
17479b94acceSBarry Smith .  snes - the SNES context
17489b94acceSBarry Smith 
17499b94acceSBarry Smith    Output Parameter:
17509b94acceSBarry Smith .  x - the solution update
17519b94acceSBarry Smith 
17529b94acceSBarry Smith    Notes:
17539b94acceSBarry Smith    This vector is implementation dependent.
17549b94acceSBarry Smith 
17559b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
17569b94acceSBarry Smith 
17579b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
17589b94acceSBarry Smith @*/
17599b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
17609b94acceSBarry Smith {
17613a40ed3dSBarry Smith   PetscFunctionBegin;
176277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17639b94acceSBarry Smith   *x = snes->vec_sol_update_always;
17643a40ed3dSBarry Smith   PetscFunctionReturn(0);
17659b94acceSBarry Smith }
17669b94acceSBarry Smith 
17675615d1e5SSatish Balay #undef __FUNC__
1768d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction"
17699b94acceSBarry Smith /*@C
17703638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
17719b94acceSBarry Smith 
17729b94acceSBarry Smith    Input Parameter:
17739b94acceSBarry Smith .  snes - the SNES context
17749b94acceSBarry Smith 
17759b94acceSBarry Smith    Output Parameter:
17763638b69dSLois Curfman McInnes .  r - the function
17779b94acceSBarry Smith 
17789b94acceSBarry Smith    Notes:
17799b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
17809b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
17819b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
17829b94acceSBarry Smith 
1783a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
17849b94acceSBarry Smith 
178531615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
178631615d04SLois Curfman McInnes           SNESGetGradient()
17879b94acceSBarry Smith @*/
17889b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
17899b94acceSBarry Smith {
17903a40ed3dSBarry Smith   PetscFunctionBegin;
179177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1792a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1793a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1794a8c6a408SBarry Smith   }
17959b94acceSBarry Smith   *r = snes->vec_func_always;
17963a40ed3dSBarry Smith   PetscFunctionReturn(0);
17979b94acceSBarry Smith }
17989b94acceSBarry Smith 
17995615d1e5SSatish Balay #undef __FUNC__
1800d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient"
18019b94acceSBarry Smith /*@C
18023638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
18039b94acceSBarry Smith 
18049b94acceSBarry Smith    Input Parameter:
18059b94acceSBarry Smith .  snes - the SNES context
18069b94acceSBarry Smith 
18079b94acceSBarry Smith    Output Parameter:
18089b94acceSBarry Smith .  r - the gradient
18099b94acceSBarry Smith 
18109b94acceSBarry Smith    Notes:
18119b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
18129b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18139b94acceSBarry Smith    SNESGetFunction().
18149b94acceSBarry Smith 
18159b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
18169b94acceSBarry Smith 
181731615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
18189b94acceSBarry Smith @*/
18199b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
18209b94acceSBarry Smith {
18213a40ed3dSBarry Smith   PetscFunctionBegin;
182277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18233a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1824a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
18253a40ed3dSBarry Smith   }
18269b94acceSBarry Smith   *r = snes->vec_func_always;
18273a40ed3dSBarry Smith   PetscFunctionReturn(0);
18289b94acceSBarry Smith }
18299b94acceSBarry Smith 
18305615d1e5SSatish Balay #undef __FUNC__
1831d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction"
18329b94acceSBarry Smith /*@
18339b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
18349b94acceSBarry Smith    unconstrained minimization problems.
18359b94acceSBarry Smith 
18369b94acceSBarry Smith    Input Parameter:
18379b94acceSBarry Smith .  snes - the SNES context
18389b94acceSBarry Smith 
18399b94acceSBarry Smith    Output Parameter:
18409b94acceSBarry Smith .  r - the function
18419b94acceSBarry Smith 
18429b94acceSBarry Smith    Notes:
18439b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
18449b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18459b94acceSBarry Smith    SNESGetFunction().
18469b94acceSBarry Smith 
18479b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
18489b94acceSBarry Smith 
184931615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
18509b94acceSBarry Smith @*/
18519b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
18529b94acceSBarry Smith {
18533a40ed3dSBarry Smith   PetscFunctionBegin;
185477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
185574679c65SBarry Smith   PetscValidScalarPointer(r);
18563a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1857a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
18583a40ed3dSBarry Smith   }
18599b94acceSBarry Smith   *r = snes->fc;
18603a40ed3dSBarry Smith   PetscFunctionReturn(0);
18619b94acceSBarry Smith }
18629b94acceSBarry Smith 
18635615d1e5SSatish Balay #undef __FUNC__
1864d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix"
18653c7409f5SSatish Balay /*@C
18663c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1867d850072dSLois Curfman McInnes    SNES options in the database.
18683c7409f5SSatish Balay 
18693c7409f5SSatish Balay    Input Parameter:
18703c7409f5SSatish Balay .  snes - the SNES context
18713c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18723c7409f5SSatish Balay 
1873d850072dSLois Curfman McInnes    Notes:
1874a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1875a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1876a83b1b31SSatish Balay    hyphen.
1877d850072dSLois Curfman McInnes 
18783c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1879a86d99e1SLois Curfman McInnes 
1880a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
18813c7409f5SSatish Balay @*/
18823c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
18833c7409f5SSatish Balay {
18843c7409f5SSatish Balay   int ierr;
18853c7409f5SSatish Balay 
18863a40ed3dSBarry Smith   PetscFunctionBegin;
188777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1888639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18893c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
18903a40ed3dSBarry Smith   PetscFunctionReturn(0);
18913c7409f5SSatish Balay }
18923c7409f5SSatish Balay 
18935615d1e5SSatish Balay #undef __FUNC__
1894d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix"
18953c7409f5SSatish Balay /*@C
1896f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1897d850072dSLois Curfman McInnes    SNES options in the database.
18983c7409f5SSatish Balay 
18993c7409f5SSatish Balay    Input Parameter:
19003c7409f5SSatish Balay .  snes - the SNES context
19013c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
19023c7409f5SSatish Balay 
1903d850072dSLois Curfman McInnes    Notes:
1904a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1905a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1906a83b1b31SSatish Balay    hyphen.
1907d850072dSLois Curfman McInnes 
19083c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1909a86d99e1SLois Curfman McInnes 
1910a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
19113c7409f5SSatish Balay @*/
19123c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
19133c7409f5SSatish Balay {
19143c7409f5SSatish Balay   int ierr;
19153c7409f5SSatish Balay 
19163a40ed3dSBarry Smith   PetscFunctionBegin;
191777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1918639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19193c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19203a40ed3dSBarry Smith   PetscFunctionReturn(0);
19213c7409f5SSatish Balay }
19223c7409f5SSatish Balay 
19235615d1e5SSatish Balay #undef __FUNC__
1924d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix"
1925c83590e2SSatish Balay /*@
19263c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
19273c7409f5SSatish Balay    SNES options in the database.
19283c7409f5SSatish Balay 
19293c7409f5SSatish Balay    Input Parameter:
19303c7409f5SSatish Balay .  snes - the SNES context
19313c7409f5SSatish Balay 
19323c7409f5SSatish Balay    Output Parameter:
19333c7409f5SSatish Balay .  prefix - pointer to the prefix string used
19343c7409f5SSatish Balay 
19353c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1936a86d99e1SLois Curfman McInnes 
1937a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
19383c7409f5SSatish Balay @*/
19393c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
19403c7409f5SSatish Balay {
19413c7409f5SSatish Balay   int ierr;
19423c7409f5SSatish Balay 
19433a40ed3dSBarry Smith   PetscFunctionBegin;
194477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1945639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19463a40ed3dSBarry Smith   PetscFunctionReturn(0);
19473c7409f5SSatish Balay }
19483c7409f5SSatish Balay 
1949ca161407SBarry Smith #undef __FUNC__
1950ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp"
1951ca161407SBarry Smith /*@
1952ca161407SBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
1953ca161407SBarry Smith 
1954ca161407SBarry Smith    Input Parameter:
1955ca161407SBarry Smith .  snes - the SNES context
1956ca161407SBarry Smith 
1957ca161407SBarry Smith    Options Database Keys:
1958ca161407SBarry Smith $  -help, -h
1959ca161407SBarry Smith 
1960ca161407SBarry Smith .keywords: SNES, nonlinear, help
1961ca161407SBarry Smith 
1962ca161407SBarry Smith .seealso: SNESSetFromOptions()
1963ca161407SBarry Smith @*/
1964ca161407SBarry Smith int SNESPrintHelp(SNES snes)
1965ca161407SBarry Smith {
1966ca161407SBarry Smith   char                p[64];
1967ca161407SBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
1968ca161407SBarry Smith   int                 ierr;
1969ca161407SBarry Smith 
1970ca161407SBarry Smith   PetscFunctionBegin;
1971ca161407SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1972ca161407SBarry Smith 
1973ca161407SBarry Smith   PetscStrcpy(p,"-");
1974ca161407SBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
1975ca161407SBarry Smith 
1976ca161407SBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
1977ca161407SBarry Smith 
1978*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");
1979ca161407SBarry Smith   ierr = NRPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",__SNESList);CHKERRQ(ierr);
1980*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
1981*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
1982*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
1983*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
1984*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
1985*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
1986*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
1987*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");
1988*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
1989*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
1990ca161407SBarry Smith     residual norm at each iteration.\n",p);
1991*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
1992ca161407SBarry Smith     residual norm for small residual norms. This is useful to conceal\n\
1993ca161407SBarry Smith     meaningless digits that may be different on different machines.\n",p);
1994*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
1995ca161407SBarry Smith   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
1996*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
1997ca161407SBarry Smith      " Options for solving systems of nonlinear equations only:\n");
1998*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
1999*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
2000*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
2001*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
2002*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
2003*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2004ca161407SBarry Smith      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
2005*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2006ca161407SBarry Smith      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
2007*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2008ca161407SBarry Smith      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
2009*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2010ca161407SBarry Smith      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
2011*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2012ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
2013*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2014ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
2015*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2016ca161407SBarry Smith      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
2017ca161407SBarry Smith   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
2018*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2019ca161407SBarry Smith      " Options for solving unconstrained minimization problems only:\n");
2020*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
2021*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
2022*76be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
2023ca161407SBarry Smith   }
2024*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
2025*76be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"a particular method\n");
2026ca161407SBarry Smith   if (snes->printhelp) {
2027ca161407SBarry Smith     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2028ca161407SBarry Smith   }
2029ca161407SBarry Smith   PetscFunctionReturn(0);
2030ca161407SBarry Smith }
20313c7409f5SSatish Balay 
20329b94acceSBarry Smith 
20339b94acceSBarry Smith 
20349b94acceSBarry Smith 
20353a40ed3dSBarry Smith 
2036