xref: /petsc/src/snes/interface/snes.c (revision 5cd905551efc76d042fcd2cc746dd1c50567705b)
1*5cd90555SBarry Smith 
2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
3*5cd90555SBarry Smith static char vcid[] = "$Id: snes.c,v 1.142 1998/04/03 23:17:44 bsmith Exp bsmith $";
49b94acceSBarry Smith #endif
59b94acceSBarry Smith 
670f55243SBarry Smith #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
7f5eb4b81SSatish Balay #include "src/sys/nreg.h"
89b94acceSBarry Smith #include "pinclude/pviewer.h"
99b94acceSBarry Smith #include <math.h>
109b94acceSBarry Smith 
1184cb2905SBarry Smith int SNESRegisterAllCalled = 0;
1282bf6240SBarry Smith DLList SNESList = 0;
1382bf6240SBarry Smith 
1484cb2905SBarry Smith 
155615d1e5SSatish Balay #undef __FUNC__
16d4bb536fSBarry Smith #define __FUNC__ "SNESView"
179b94acceSBarry Smith /*@
189b94acceSBarry Smith    SNESView - Prints the SNES data structure.
199b94acceSBarry Smith 
209b94acceSBarry Smith    Input Parameters:
219b94acceSBarry Smith .  SNES - the SNES context
229b94acceSBarry Smith .  viewer - visualization context
239b94acceSBarry Smith 
249b94acceSBarry Smith    Options Database Key:
259b94acceSBarry Smith $  -snes_view : calls SNESView() at end of SNESSolve()
269b94acceSBarry Smith 
279b94acceSBarry Smith    Notes:
289b94acceSBarry Smith    The available visualization contexts include
296d4a8577SBarry Smith $     VIEWER_STDOUT_SELF - standard output (default)
306d4a8577SBarry Smith $     VIEWER_STDOUT_WORLD - synchronized standard
319b94acceSBarry Smith $       output where only the first processor opens
329b94acceSBarry Smith $       the file.  All other processors send their
339b94acceSBarry Smith $       data to the first processor to print.
349b94acceSBarry Smith 
359b94acceSBarry Smith    The user can open alternative vistualization contexts with
36dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
379b94acceSBarry Smith 
389b94acceSBarry Smith .keywords: SNES, view
399b94acceSBarry Smith 
40dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
419b94acceSBarry Smith @*/
429b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
439b94acceSBarry Smith {
449b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
459b94acceSBarry Smith   FILE                *fd;
469b94acceSBarry Smith   int                 ierr;
479b94acceSBarry Smith   SLES                sles;
489b94acceSBarry Smith   char                *method;
4919bcc07fSBarry Smith   ViewerType          vtype;
509b94acceSBarry Smith 
513a40ed3dSBarry Smith   PetscFunctionBegin;
5274679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5374679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5474679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5574679c65SBarry Smith 
5619bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5719bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5890ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
5977c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
6082bf6240SBarry Smith     SNESGetType(snes,&method);
6177c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
62e1311b90SBarry Smith     if (snes->view) (*snes->view)(snes,viewer);
6377c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
649b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
659b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6677c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
679b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
689b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
697a00f4a9SLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
707a00f4a9SLois Curfman McInnes     "  total number of linear solver iterations=%d\n",snes->linear_its);
71ae3c334cSLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
7268632166SLois Curfman McInnes      "  total number of function evaluations=%d\n",snes->nfuncs);
739b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
7477c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
759b94acceSBarry Smith     if (snes->ksp_ewconv) {
769b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
779b94acceSBarry Smith       if (kctx) {
7877c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
799b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
809b94acceSBarry Smith         kctx->version);
8177c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
829b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
839b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
8477c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
859b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
869b94acceSBarry Smith       }
879b94acceSBarry Smith     }
88c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
8982bf6240SBarry Smith     SNESGetType(snes,&method);
90c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
9119bcc07fSBarry Smith   }
929b94acceSBarry Smith   SNESGetSLES(snes,&sles);
939b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
943a40ed3dSBarry Smith   PetscFunctionReturn(0);
959b94acceSBarry Smith }
969b94acceSBarry Smith 
97639f9d9dSBarry Smith /*
98639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
99639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
100639f9d9dSBarry Smith */
101639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
102639f9d9dSBarry Smith static int numberofsetfromoptions;
103639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
104639f9d9dSBarry Smith 
1055615d1e5SSatish Balay #undef __FUNC__
106d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker"
107639f9d9dSBarry Smith /*@
108639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
109639f9d9dSBarry Smith 
110639f9d9dSBarry Smith     Input Parameter:
111639f9d9dSBarry Smith .   snescheck - function that checks for options
112639f9d9dSBarry Smith 
113639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
114639f9d9dSBarry Smith @*/
115639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
116639f9d9dSBarry Smith {
1173a40ed3dSBarry Smith   PetscFunctionBegin;
118639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
119a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed");
120639f9d9dSBarry Smith   }
121639f9d9dSBarry Smith 
122639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
1233a40ed3dSBarry Smith   PetscFunctionReturn(0);
124639f9d9dSBarry Smith }
125639f9d9dSBarry Smith 
1265615d1e5SSatish Balay #undef __FUNC__
1275615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1289b94acceSBarry Smith /*@
129682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1309b94acceSBarry Smith 
1319b94acceSBarry Smith    Input Parameter:
1329b94acceSBarry Smith .  snes - the SNES context
1339b94acceSBarry Smith 
13411ca99fdSLois Curfman McInnes    Notes:
13511ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
13611ca99fdSLois Curfman McInnes    the users manual.
13783e2fdc7SBarry Smith 
1389b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1399b94acceSBarry Smith 
140a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1419b94acceSBarry Smith @*/
1429b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1439b94acceSBarry Smith {
14482bf6240SBarry Smith   char     method[256];
1459b94acceSBarry Smith   double   tmp;
1469b94acceSBarry Smith   SLES     sles;
14781f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1483c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1499b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1509b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1519b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1529b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1539b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1549b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1559b94acceSBarry Smith 
1563a40ed3dSBarry Smith   PetscFunctionBegin;
15781f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
15881f57ec7SBarry Smith 
15977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16082bf6240SBarry Smith   if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp");
161ca161407SBarry Smith 
16282bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
16382bf6240SBarry Smith   ierr = OptionsGetString(snes->prefix,"-snes_type",method,256,&flg);
164052efed2SBarry Smith   if (flg) {
16582bf6240SBarry Smith     ierr = SNESSetType(snes,(SNESType) method); CHKERRQ(ierr);
1665334005bSBarry Smith   }
1673c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
168d64ed03dSBarry Smith   if (flg) {
169d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
170d64ed03dSBarry Smith   }
1713c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
172d64ed03dSBarry Smith   if (flg) {
173d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
174d64ed03dSBarry Smith   }
1753c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
176d64ed03dSBarry Smith   if (flg) {
177d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
178d64ed03dSBarry Smith   }
1793c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1803c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
181d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
182d64ed03dSBarry Smith   if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); }
183d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
184d64ed03dSBarry Smith   if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);  CHKERRQ(ierr);}
1853c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1863c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
187b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
188b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
189b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
190b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
191b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
192b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
193b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
19493c39befSBarry Smith 
19593c39befSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr);
19693c39befSBarry Smith   if (flg) {snes->converged = 0;}
19793c39befSBarry Smith 
1989b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
1999b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
2005bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
201*5cd90555SBarry Smith   if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
2023c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
203639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
2043c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
205d31a3109SLois Curfman McInnes   if (flg) {
206639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
207d31a3109SLois Curfman McInnes   }
20881f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
2093c7409f5SSatish Balay   if (flg) {
21017699dbbSLois Curfman McInnes     int    rank = 0;
211d7e8b826SBarry Smith     DrawLG lg;
21217699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
21317699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
21417699dbbSLois Curfman McInnes     if (!rank) {
21581f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2169b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
217f8c826e1SBarry Smith       snes->xmonitor = lg;
2189b94acceSBarry Smith       PLogObjectParent(snes,lg);
2199b94acceSBarry Smith     }
2209b94acceSBarry Smith   }
2213c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2223c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2239b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2249b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
225981c4779SBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
226981c4779SBarry Smith   } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
22731615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
22831615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
229d64ed03dSBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2309b94acceSBarry Smith   }
231639f9d9dSBarry Smith 
232639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
233639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
234639f9d9dSBarry Smith   }
235639f9d9dSBarry Smith 
2369b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2379b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
23882bf6240SBarry Smith   /*
23982bf6240SBarry Smith         Since the private setfromoptions requires the type to all ready have
24082bf6240SBarry Smith       been set we make sure a type is set by this time
24182bf6240SBarry Smith   */
24282bf6240SBarry Smith   if (!snes->type_name) {
24382bf6240SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
24482bf6240SBarry Smith       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
24582bf6240SBarry Smith     } else {
24682bf6240SBarry Smith       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
24782bf6240SBarry Smith     }
24882bf6240SBarry Smith   }
24982bf6240SBarry Smith 
25082bf6240SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
25182bf6240SBarry Smith   if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);}
25282bf6240SBarry Smith 
2533a40ed3dSBarry Smith   if (snes->setfromoptions) {
2543a40ed3dSBarry Smith     ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
2553a40ed3dSBarry Smith   }
2563a40ed3dSBarry Smith   PetscFunctionReturn(0);
2579b94acceSBarry Smith }
2589b94acceSBarry Smith 
259a847f771SSatish Balay 
2605615d1e5SSatish Balay #undef __FUNC__
261d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext"
2629b94acceSBarry Smith /*@
2639b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2649b94acceSBarry Smith    the nonlinear solvers.
2659b94acceSBarry Smith 
2669b94acceSBarry Smith    Input Parameters:
2679b94acceSBarry Smith .  snes - the SNES context
2689b94acceSBarry Smith .  usrP - optional user context
2699b94acceSBarry Smith 
2709b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2719b94acceSBarry Smith 
2729b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2739b94acceSBarry Smith @*/
2749b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
2759b94acceSBarry Smith {
2763a40ed3dSBarry Smith   PetscFunctionBegin;
27777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2789b94acceSBarry Smith   snes->user		= usrP;
2793a40ed3dSBarry Smith   PetscFunctionReturn(0);
2809b94acceSBarry Smith }
28174679c65SBarry Smith 
2825615d1e5SSatish Balay #undef __FUNC__
283d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext"
2849b94acceSBarry Smith /*@C
2859b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
2869b94acceSBarry Smith    nonlinear solvers.
2879b94acceSBarry Smith 
2889b94acceSBarry Smith    Input Parameter:
2899b94acceSBarry Smith .  snes - SNES context
2909b94acceSBarry Smith 
2919b94acceSBarry Smith    Output Parameter:
2929b94acceSBarry Smith .  usrP - user context
2939b94acceSBarry Smith 
2949b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
2959b94acceSBarry Smith 
2969b94acceSBarry Smith .seealso: SNESSetApplicationContext()
2979b94acceSBarry Smith @*/
2989b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
2999b94acceSBarry Smith {
3003a40ed3dSBarry Smith   PetscFunctionBegin;
30177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3029b94acceSBarry Smith   *usrP = snes->user;
3033a40ed3dSBarry Smith   PetscFunctionReturn(0);
3049b94acceSBarry Smith }
30574679c65SBarry Smith 
3065615d1e5SSatish Balay #undef __FUNC__
307d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber"
3089b94acceSBarry Smith /*@
3099b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3109b94acceSBarry Smith    nonlinear solver.
3119b94acceSBarry Smith 
3129b94acceSBarry Smith    Input Parameter:
3139b94acceSBarry Smith .  snes - SNES context
3149b94acceSBarry Smith 
3159b94acceSBarry Smith    Output Parameter:
3169b94acceSBarry Smith .  iter - iteration number
3179b94acceSBarry Smith 
3189b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3199b94acceSBarry Smith @*/
3209b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3219b94acceSBarry Smith {
3223a40ed3dSBarry Smith   PetscFunctionBegin;
32377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
32474679c65SBarry Smith   PetscValidIntPointer(iter);
3259b94acceSBarry Smith   *iter = snes->iter;
3263a40ed3dSBarry Smith   PetscFunctionReturn(0);
3279b94acceSBarry Smith }
32874679c65SBarry Smith 
3295615d1e5SSatish Balay #undef __FUNC__
3305615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
3319b94acceSBarry Smith /*@
3329b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3339b94acceSBarry Smith    with SNESSSetFunction().
3349b94acceSBarry Smith 
3359b94acceSBarry Smith    Input Parameter:
3369b94acceSBarry Smith .  snes - SNES context
3379b94acceSBarry Smith 
3389b94acceSBarry Smith    Output Parameter:
3399b94acceSBarry Smith .  fnorm - 2-norm of function
3409b94acceSBarry Smith 
3419b94acceSBarry Smith    Note:
3429b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
343a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
344a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3459b94acceSBarry Smith 
3469b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
347a86d99e1SLois Curfman McInnes 
348a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3499b94acceSBarry Smith @*/
3509b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3519b94acceSBarry Smith {
3523a40ed3dSBarry Smith   PetscFunctionBegin;
35377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
35474679c65SBarry Smith   PetscValidScalarPointer(fnorm);
35574679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
356d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only");
35774679c65SBarry Smith   }
3589b94acceSBarry Smith   *fnorm = snes->norm;
3593a40ed3dSBarry Smith   PetscFunctionReturn(0);
3609b94acceSBarry Smith }
36174679c65SBarry Smith 
3625615d1e5SSatish Balay #undef __FUNC__
3635615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
3649b94acceSBarry Smith /*@
3659b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
3669b94acceSBarry Smith    with SNESSSetGradient().
3679b94acceSBarry Smith 
3689b94acceSBarry Smith    Input Parameter:
3699b94acceSBarry Smith .  snes - SNES context
3709b94acceSBarry Smith 
3719b94acceSBarry Smith    Output Parameter:
3729b94acceSBarry Smith .  fnorm - 2-norm of gradient
3739b94acceSBarry Smith 
3749b94acceSBarry Smith    Note:
3759b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
376a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
377a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
3789b94acceSBarry Smith 
3799b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
380a86d99e1SLois Curfman McInnes 
381a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
3829b94acceSBarry Smith @*/
3839b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
3849b94acceSBarry Smith {
3853a40ed3dSBarry Smith   PetscFunctionBegin;
38677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
38774679c65SBarry Smith   PetscValidScalarPointer(gnorm);
38874679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
389d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
39074679c65SBarry Smith   }
3919b94acceSBarry Smith   *gnorm = snes->norm;
3923a40ed3dSBarry Smith   PetscFunctionReturn(0);
3939b94acceSBarry Smith }
39474679c65SBarry Smith 
3955615d1e5SSatish Balay #undef __FUNC__
396d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps"
3979b94acceSBarry Smith /*@
3989b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
3999b94acceSBarry Smith    attempted by the nonlinear solver.
4009b94acceSBarry Smith 
4019b94acceSBarry Smith    Input Parameter:
4029b94acceSBarry Smith .  snes - SNES context
4039b94acceSBarry Smith 
4049b94acceSBarry Smith    Output Parameter:
4059b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4069b94acceSBarry Smith 
407c96a6f78SLois Curfman McInnes    Notes:
408c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
409c96a6f78SLois Curfman McInnes 
4109b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4119b94acceSBarry Smith @*/
4129b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4139b94acceSBarry Smith {
4143a40ed3dSBarry Smith   PetscFunctionBegin;
41577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
41674679c65SBarry Smith   PetscValidIntPointer(nfails);
4179b94acceSBarry Smith   *nfails = snes->nfailures;
4183a40ed3dSBarry Smith   PetscFunctionReturn(0);
4199b94acceSBarry Smith }
420a847f771SSatish Balay 
4215615d1e5SSatish Balay #undef __FUNC__
422d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations"
423c96a6f78SLois Curfman McInnes /*@
424c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
425c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
426c96a6f78SLois Curfman McInnes 
427c96a6f78SLois Curfman McInnes    Input Parameter:
428c96a6f78SLois Curfman McInnes .  snes - SNES context
429c96a6f78SLois Curfman McInnes 
430c96a6f78SLois Curfman McInnes    Output Parameter:
431c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
432c96a6f78SLois Curfman McInnes 
433c96a6f78SLois Curfman McInnes    Notes:
434c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
435c96a6f78SLois Curfman McInnes 
436c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
437c96a6f78SLois Curfman McInnes @*/
43886bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
439c96a6f78SLois Curfman McInnes {
4403a40ed3dSBarry Smith   PetscFunctionBegin;
441c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
442c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
443c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
4443a40ed3dSBarry Smith   PetscFunctionReturn(0);
445c96a6f78SLois Curfman McInnes }
446c96a6f78SLois Curfman McInnes 
447c96a6f78SLois Curfman McInnes #undef __FUNC__
448d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES"
4499b94acceSBarry Smith /*@C
4509b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4519b94acceSBarry Smith 
4529b94acceSBarry Smith    Input Parameter:
4539b94acceSBarry Smith .  snes - the SNES context
4549b94acceSBarry Smith 
4559b94acceSBarry Smith    Output Parameter:
4569b94acceSBarry Smith .  sles - the SLES context
4579b94acceSBarry Smith 
4589b94acceSBarry Smith    Notes:
4599b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4609b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4619b94acceSBarry Smith    KSP and PC contexts as well.
4629b94acceSBarry Smith 
4639b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4649b94acceSBarry Smith 
4659b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4669b94acceSBarry Smith @*/
4679b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4689b94acceSBarry Smith {
4693a40ed3dSBarry Smith   PetscFunctionBegin;
47077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4719b94acceSBarry Smith   *sles = snes->sles;
4723a40ed3dSBarry Smith   PetscFunctionReturn(0);
4739b94acceSBarry Smith }
47482bf6240SBarry Smith 
4759b94acceSBarry Smith /* -----------------------------------------------------------*/
4765615d1e5SSatish Balay #undef __FUNC__
4775615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
4789b94acceSBarry Smith /*@C
4799b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
4809b94acceSBarry Smith 
4819b94acceSBarry Smith    Input Parameter:
4829b94acceSBarry Smith .  comm - MPI communicator
4839b94acceSBarry Smith .  type - type of method, one of
4849b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
4859b94acceSBarry Smith $      (for systems of nonlinear equations)
4869b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
4879b94acceSBarry Smith $      (for unconstrained minimization)
4889b94acceSBarry Smith 
4899b94acceSBarry Smith    Output Parameter:
4909b94acceSBarry Smith .  outsnes - the new SNES context
4919b94acceSBarry Smith 
492c1f60f51SBarry Smith    Options Database Key:
49319bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
49419bd6747SLois Curfman McInnes $              and no preconditioning matrix
49519bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
49619bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
49719bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
49819bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
499c1f60f51SBarry Smith 
5009b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5019b94acceSBarry Smith 
50263a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5039b94acceSBarry Smith @*/
5044b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5059b94acceSBarry Smith {
5069b94acceSBarry Smith   int                 ierr;
5079b94acceSBarry Smith   SNES                snes;
5089b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
50937fcc0dbSBarry Smith 
5103a40ed3dSBarry Smith   PetscFunctionBegin;
5119b94acceSBarry Smith   *outsnes = 0;
512d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
513d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
514d64ed03dSBarry Smith   }
515f830108cSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,comm,SNESDestroy,SNESView);
5169b94acceSBarry Smith   PLogObjectCreate(snes);
5179b94acceSBarry Smith   snes->max_its           = 50;
5189750a799SBarry Smith   snes->max_funcs	  = 10000;
5199b94acceSBarry Smith   snes->norm		  = 0.0;
5205a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
521b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
522b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5239b94acceSBarry Smith     snes->atol		  = 1.e-10;
524d64ed03dSBarry Smith   } else {
525b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
526b4874afaSBarry Smith     snes->ttol            = 0.0;
527b4874afaSBarry Smith     snes->atol		  = 1.e-50;
528b4874afaSBarry Smith   }
5299b94acceSBarry Smith   snes->xtol		  = 1.e-8;
530b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5319b94acceSBarry Smith   snes->nfuncs            = 0;
5329b94acceSBarry Smith   snes->nfailures         = 0;
5337a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
534639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5359b94acceSBarry Smith   snes->data              = 0;
5369b94acceSBarry Smith   snes->view              = 0;
5379b94acceSBarry Smith   snes->computeumfunction = 0;
5389b94acceSBarry Smith   snes->umfunP            = 0;
5399b94acceSBarry Smith   snes->fc                = 0;
5409b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5419b94acceSBarry Smith   snes->fmin              = -1.e30;
5429b94acceSBarry Smith   snes->method_class      = type;
5439b94acceSBarry Smith   snes->set_method_called = 0;
54482bf6240SBarry Smith   snes->setupcalled      = 0;
5459b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5466f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5476f24a144SLois Curfman McInnes   snes->vwork             = 0;
5486f24a144SLois Curfman McInnes   snes->nwork             = 0;
549c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
550c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
551c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
5529b94acceSBarry Smith 
5539b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5540452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
555eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
5569b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5579b94acceSBarry Smith   kctx->version     = 2;
5589b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5599b94acceSBarry Smith                              this was too large for some test cases */
5609b94acceSBarry Smith   kctx->rtol_last   = 0;
5619b94acceSBarry Smith   kctx->rtol_max    = .9;
5629b94acceSBarry Smith   kctx->gamma       = 1.0;
5639b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5649b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5659b94acceSBarry Smith   kctx->threshold   = .1;
5669b94acceSBarry Smith   kctx->lresid_last = 0;
5679b94acceSBarry Smith   kctx->norm_last   = 0;
5689b94acceSBarry Smith 
5699b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5709b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5715334005bSBarry Smith 
5729b94acceSBarry Smith   *outsnes = snes;
5733a40ed3dSBarry Smith   PetscFunctionReturn(0);
5749b94acceSBarry Smith }
5759b94acceSBarry Smith 
5769b94acceSBarry Smith /* --------------------------------------------------------------- */
5775615d1e5SSatish Balay #undef __FUNC__
578d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction"
5799b94acceSBarry Smith /*@C
5809b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
5819b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
5829b94acceSBarry Smith    equations.
5839b94acceSBarry Smith 
5849b94acceSBarry Smith    Input Parameters:
5859b94acceSBarry Smith .  snes - the SNES context
5869b94acceSBarry Smith .  func - function evaluation routine
5879b94acceSBarry Smith .  r - vector to store function value
5882cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
5892cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
5909b94acceSBarry Smith 
5919b94acceSBarry Smith    Calling sequence of func:
592313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
5939b94acceSBarry Smith 
5949b94acceSBarry Smith .  x - input vector
595313e4042SLois Curfman McInnes .  f - function vector
5962cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
5979b94acceSBarry Smith 
5989b94acceSBarry Smith    Notes:
5999b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6009b94acceSBarry Smith $      f'(x) x = -f(x),
6019b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
6029b94acceSBarry Smith 
6039b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6049b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6059b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6069b94acceSBarry Smith 
6079b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6089b94acceSBarry Smith 
609a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6109b94acceSBarry Smith @*/
6115334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6129b94acceSBarry Smith {
6133a40ed3dSBarry Smith   PetscFunctionBegin;
61477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
615a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
616a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
617a8c6a408SBarry Smith   }
6189b94acceSBarry Smith   snes->computefunction     = func;
6199b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6209b94acceSBarry Smith   snes->funP                = ctx;
6213a40ed3dSBarry Smith   PetscFunctionReturn(0);
6229b94acceSBarry Smith }
6239b94acceSBarry Smith 
6245615d1e5SSatish Balay #undef __FUNC__
6255615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
6269b94acceSBarry Smith /*@
6279b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6289b94acceSBarry Smith    SNESSetFunction().
6299b94acceSBarry Smith 
6309b94acceSBarry Smith    Input Parameters:
6319b94acceSBarry Smith .  snes - the SNES context
6329b94acceSBarry Smith .  x - input vector
6339b94acceSBarry Smith 
6349b94acceSBarry Smith    Output Parameter:
6353638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6369b94acceSBarry Smith 
6371bffabb2SLois Curfman McInnes    Notes:
6389b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6399b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6409b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6419b94acceSBarry Smith 
6429b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6439b94acceSBarry Smith 
644a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6459b94acceSBarry Smith @*/
6469b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6479b94acceSBarry Smith {
6489b94acceSBarry Smith   int    ierr;
6499b94acceSBarry Smith 
6503a40ed3dSBarry Smith   PetscFunctionBegin;
651d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
652a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
653d4bb536fSBarry Smith   }
6549b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
655d64ed03dSBarry Smith   PetscStackPush("SNES user function");
6569b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
657d64ed03dSBarry Smith   PetscStackPop;
658ae3c334cSLois Curfman McInnes   snes->nfuncs++;
6599b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6603a40ed3dSBarry Smith   PetscFunctionReturn(0);
6619b94acceSBarry Smith }
6629b94acceSBarry Smith 
6635615d1e5SSatish Balay #undef __FUNC__
664d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction"
6659b94acceSBarry Smith /*@C
6669b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6679b94acceSBarry Smith    unconstrained minimization.
6689b94acceSBarry Smith 
6699b94acceSBarry Smith    Input Parameters:
6709b94acceSBarry Smith .  snes - the SNES context
6719b94acceSBarry Smith .  func - function evaluation routine
672044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
673044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6749b94acceSBarry Smith 
6759b94acceSBarry Smith    Calling sequence of func:
6769b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
6779b94acceSBarry Smith 
6789b94acceSBarry Smith .  x - input vector
6799b94acceSBarry Smith .  f - function
680044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
6819b94acceSBarry Smith 
6829b94acceSBarry Smith    Notes:
6839b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
6849b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
6859b94acceSBarry Smith    SNESSetFunction().
6869b94acceSBarry Smith 
6879b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
6889b94acceSBarry Smith 
689a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
690a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
6919b94acceSBarry Smith @*/
6929b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
6939b94acceSBarry Smith                       void *ctx)
6949b94acceSBarry Smith {
6953a40ed3dSBarry Smith   PetscFunctionBegin;
69677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
697a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
698a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
699a8c6a408SBarry Smith   }
7009b94acceSBarry Smith   snes->computeumfunction   = func;
7019b94acceSBarry Smith   snes->umfunP              = ctx;
7023a40ed3dSBarry Smith   PetscFunctionReturn(0);
7039b94acceSBarry Smith }
7049b94acceSBarry Smith 
7055615d1e5SSatish Balay #undef __FUNC__
7065615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7079b94acceSBarry Smith /*@
7089b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7099b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7109b94acceSBarry Smith 
7119b94acceSBarry Smith    Input Parameters:
7129b94acceSBarry Smith .  snes - the SNES context
7139b94acceSBarry Smith .  x - input vector
7149b94acceSBarry Smith 
7159b94acceSBarry Smith    Output Parameter:
7169b94acceSBarry Smith .  y - function value
7179b94acceSBarry Smith 
7189b94acceSBarry Smith    Notes:
7199b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7209b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7219b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
722a86d99e1SLois Curfman McInnes 
723a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
724a86d99e1SLois Curfman McInnes 
725a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
726a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7279b94acceSBarry Smith @*/
7289b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7299b94acceSBarry Smith {
7309b94acceSBarry Smith   int    ierr;
7313a40ed3dSBarry Smith 
7323a40ed3dSBarry Smith   PetscFunctionBegin;
733a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
734a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
735a8c6a408SBarry Smith   }
7369b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
737d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
7389b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
739d64ed03dSBarry Smith   PetscStackPop;
740ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7419b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7423a40ed3dSBarry Smith   PetscFunctionReturn(0);
7439b94acceSBarry Smith }
7449b94acceSBarry Smith 
7455615d1e5SSatish Balay #undef __FUNC__
746d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient"
7479b94acceSBarry Smith /*@C
7489b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7499b94acceSBarry Smith    vector for use by the SNES routines.
7509b94acceSBarry Smith 
7519b94acceSBarry Smith    Input Parameters:
7529b94acceSBarry Smith .  snes - the SNES context
7539b94acceSBarry Smith .  func - function evaluation routine
754044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
755044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
7569b94acceSBarry Smith .  r - vector to store gradient value
7579b94acceSBarry Smith 
7589b94acceSBarry Smith    Calling sequence of func:
7599b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7609b94acceSBarry Smith 
7619b94acceSBarry Smith .  x - input vector
7629b94acceSBarry Smith .  g - gradient vector
763044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7649b94acceSBarry Smith 
7659b94acceSBarry Smith    Notes:
7669b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7679b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7689b94acceSBarry Smith    SNESSetFunction().
7699b94acceSBarry Smith 
7709b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7719b94acceSBarry Smith 
772a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
773a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
7749b94acceSBarry Smith @*/
77574679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
7769b94acceSBarry Smith {
7773a40ed3dSBarry Smith   PetscFunctionBegin;
77877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
779a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
780a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
781a8c6a408SBarry Smith   }
7829b94acceSBarry Smith   snes->computefunction     = func;
7839b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7849b94acceSBarry Smith   snes->funP                = ctx;
7853a40ed3dSBarry Smith   PetscFunctionReturn(0);
7869b94acceSBarry Smith }
7879b94acceSBarry Smith 
7885615d1e5SSatish Balay #undef __FUNC__
7895615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
7909b94acceSBarry Smith /*@
791a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
792a86d99e1SLois Curfman McInnes    SNESSetGradient().
7939b94acceSBarry Smith 
7949b94acceSBarry Smith    Input Parameters:
7959b94acceSBarry Smith .  snes - the SNES context
7969b94acceSBarry Smith .  x - input vector
7979b94acceSBarry Smith 
7989b94acceSBarry Smith    Output Parameter:
7999b94acceSBarry Smith .  y - gradient vector
8009b94acceSBarry Smith 
8019b94acceSBarry Smith    Notes:
8029b94acceSBarry Smith    SNESComputeGradient() is valid only for
8039b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8049b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
805a86d99e1SLois Curfman McInnes 
806a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
807a86d99e1SLois Curfman McInnes 
808a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
809a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8109b94acceSBarry Smith @*/
8119b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8129b94acceSBarry Smith {
8139b94acceSBarry Smith   int    ierr;
8143a40ed3dSBarry Smith 
8153a40ed3dSBarry Smith   PetscFunctionBegin;
8163a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
817a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8183a40ed3dSBarry Smith   }
8199b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
820d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
8219b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
822d64ed03dSBarry Smith   PetscStackPop;
8239b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8243a40ed3dSBarry Smith   PetscFunctionReturn(0);
8259b94acceSBarry Smith }
8269b94acceSBarry Smith 
8275615d1e5SSatish Balay #undef __FUNC__
8285615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
82962fef451SLois Curfman McInnes /*@
83062fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
83162fef451SLois Curfman McInnes    set with SNESSetJacobian().
83262fef451SLois Curfman McInnes 
83362fef451SLois Curfman McInnes    Input Parameters:
83462fef451SLois Curfman McInnes .  snes - the SNES context
83562fef451SLois Curfman McInnes .  x - input vector
83662fef451SLois Curfman McInnes 
83762fef451SLois Curfman McInnes    Output Parameters:
83862fef451SLois Curfman McInnes .  A - Jacobian matrix
83962fef451SLois Curfman McInnes .  B - optional preconditioning matrix
84062fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
84162fef451SLois Curfman McInnes 
84262fef451SLois Curfman McInnes    Notes:
84362fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
84462fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
84562fef451SLois Curfman McInnes 
846dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
847dc5a77f8SLois Curfman McInnes    flag parameter.
84862fef451SLois Curfman McInnes 
84962fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
85062fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
851005c665bSBarry Smith    methods is SNESComputeHessian().
85262fef451SLois Curfman McInnes 
85362fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
85462fef451SLois Curfman McInnes 
85562fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
85662fef451SLois Curfman McInnes @*/
8579b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8589b94acceSBarry Smith {
8599b94acceSBarry Smith   int    ierr;
8603a40ed3dSBarry Smith 
8613a40ed3dSBarry Smith   PetscFunctionBegin;
8623a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
863a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
8643a40ed3dSBarry Smith   }
8653a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
8669b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
867c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
868d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
8699b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
870d64ed03dSBarry Smith   PetscStackPop;
8719b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
8726d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
87377c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
87477c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
8753a40ed3dSBarry Smith   PetscFunctionReturn(0);
8769b94acceSBarry Smith }
8779b94acceSBarry Smith 
8785615d1e5SSatish Balay #undef __FUNC__
8795615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
88062fef451SLois Curfman McInnes /*@
88162fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
88262fef451SLois Curfman McInnes    set with SNESSetHessian().
88362fef451SLois Curfman McInnes 
88462fef451SLois Curfman McInnes    Input Parameters:
88562fef451SLois Curfman McInnes .  snes - the SNES context
88662fef451SLois Curfman McInnes .  x - input vector
88762fef451SLois Curfman McInnes 
88862fef451SLois Curfman McInnes    Output Parameters:
88962fef451SLois Curfman McInnes .  A - Hessian matrix
89062fef451SLois Curfman McInnes .  B - optional preconditioning matrix
89162fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
89262fef451SLois Curfman McInnes 
89362fef451SLois Curfman McInnes    Notes:
89462fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
89562fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
89662fef451SLois Curfman McInnes 
897dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
898dc5a77f8SLois Curfman McInnes    flag parameter.
89962fef451SLois Curfman McInnes 
90062fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
90162fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
90262fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
90362fef451SLois Curfman McInnes 
90462fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
90562fef451SLois Curfman McInnes 
906a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
907a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
90862fef451SLois Curfman McInnes @*/
90962fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
9109b94acceSBarry Smith {
9119b94acceSBarry Smith   int    ierr;
9123a40ed3dSBarry Smith 
9133a40ed3dSBarry Smith   PetscFunctionBegin;
9143a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
915a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9163a40ed3dSBarry Smith   }
9173a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
91862fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9190f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
920d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
92162fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
922d64ed03dSBarry Smith   PetscStackPop;
92362fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9240f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
92577c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
92677c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9273a40ed3dSBarry Smith   PetscFunctionReturn(0);
9289b94acceSBarry Smith }
9299b94acceSBarry Smith 
9305615d1e5SSatish Balay #undef __FUNC__
931d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian"
9329b94acceSBarry Smith /*@C
9339b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
934044dda88SLois Curfman McInnes    location to store the matrix.
9359b94acceSBarry Smith 
9369b94acceSBarry Smith    Input Parameters:
9379b94acceSBarry Smith .  snes - the SNES context
9389b94acceSBarry Smith .  A - Jacobian matrix
9399b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9409b94acceSBarry Smith .  func - Jacobian evaluation routine
9412cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9422cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9439b94acceSBarry Smith 
9449b94acceSBarry Smith    Calling sequence of func:
945313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9469b94acceSBarry Smith 
9479b94acceSBarry Smith .  x - input vector
9489b94acceSBarry Smith .  A - Jacobian matrix
9499b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
950ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
951ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9522cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
9539b94acceSBarry Smith 
9549b94acceSBarry Smith    Notes:
955dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9562cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
957ac21db08SLois Curfman McInnes 
958ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9599b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9609b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9619b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9629b94acceSBarry Smith    throughout the global iterations.
9639b94acceSBarry Smith 
9649b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9659b94acceSBarry Smith 
966ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
9679b94acceSBarry Smith @*/
9689b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
9699b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
9709b94acceSBarry Smith {
9713a40ed3dSBarry Smith   PetscFunctionBegin;
97277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
973a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
974a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
975a8c6a408SBarry Smith   }
9769b94acceSBarry Smith   snes->computejacobian = func;
9779b94acceSBarry Smith   snes->jacP            = ctx;
9789b94acceSBarry Smith   snes->jacobian        = A;
9799b94acceSBarry Smith   snes->jacobian_pre    = B;
9803a40ed3dSBarry Smith   PetscFunctionReturn(0);
9819b94acceSBarry Smith }
98262fef451SLois Curfman McInnes 
9835615d1e5SSatish Balay #undef __FUNC__
984d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian"
985b4fd4287SBarry Smith /*@
986b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
987b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
988b4fd4287SBarry Smith 
989b4fd4287SBarry Smith    Input Parameter:
990b4fd4287SBarry Smith .  snes - the nonlinear solver context
991b4fd4287SBarry Smith 
992b4fd4287SBarry Smith    Output Parameters:
993b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
994b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
995b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
996b4fd4287SBarry Smith 
997b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
998b4fd4287SBarry Smith @*/
999b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
1000b4fd4287SBarry Smith {
10013a40ed3dSBarry Smith   PetscFunctionBegin;
10023a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1003a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
10043a40ed3dSBarry Smith   }
1005b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
1006b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
1007b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
10083a40ed3dSBarry Smith   PetscFunctionReturn(0);
1009b4fd4287SBarry Smith }
1010b4fd4287SBarry Smith 
10115615d1e5SSatish Balay #undef __FUNC__
1012d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian"
10139b94acceSBarry Smith /*@C
10149b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1015044dda88SLois Curfman McInnes    location to store the matrix.
10169b94acceSBarry Smith 
10179b94acceSBarry Smith    Input Parameters:
10189b94acceSBarry Smith .  snes - the SNES context
10199b94acceSBarry Smith .  A - Hessian matrix
10209b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
10219b94acceSBarry Smith .  func - Jacobian evaluation routine
1022313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1023313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10249b94acceSBarry Smith 
10259b94acceSBarry Smith    Calling sequence of func:
1026313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10279b94acceSBarry Smith 
10289b94acceSBarry Smith .  x - input vector
10299b94acceSBarry Smith .  A - Hessian matrix
10309b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1031ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1032ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10332cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10349b94acceSBarry Smith 
10359b94acceSBarry Smith    Notes:
1036dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10372cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1038ac21db08SLois Curfman McInnes 
10399b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10409b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10419b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10429b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10439b94acceSBarry Smith    throughout the global iterations.
10449b94acceSBarry Smith 
10459b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10469b94acceSBarry Smith 
1047ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10489b94acceSBarry Smith @*/
10499b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10509b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10519b94acceSBarry Smith {
10523a40ed3dSBarry Smith   PetscFunctionBegin;
105377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1054d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1055a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1056d4bb536fSBarry Smith   }
10579b94acceSBarry Smith   snes->computejacobian = func;
10589b94acceSBarry Smith   snes->jacP            = ctx;
10599b94acceSBarry Smith   snes->jacobian        = A;
10609b94acceSBarry Smith   snes->jacobian_pre    = B;
10613a40ed3dSBarry Smith   PetscFunctionReturn(0);
10629b94acceSBarry Smith }
10639b94acceSBarry Smith 
10645615d1e5SSatish Balay #undef __FUNC__
1065d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian"
106662fef451SLois Curfman McInnes /*@
106762fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
106862fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
106962fef451SLois Curfman McInnes 
107062fef451SLois Curfman McInnes    Input Parameter:
107162fef451SLois Curfman McInnes .  snes - the nonlinear solver context
107262fef451SLois Curfman McInnes 
107362fef451SLois Curfman McInnes    Output Parameters:
107462fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
107562fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
107662fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
107762fef451SLois Curfman McInnes 
107862fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
107962fef451SLois Curfman McInnes @*/
108062fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
108162fef451SLois Curfman McInnes {
10823a40ed3dSBarry Smith   PetscFunctionBegin;
10833a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1084a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10853a40ed3dSBarry Smith   }
108662fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
108762fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
108862fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
10893a40ed3dSBarry Smith   PetscFunctionReturn(0);
109062fef451SLois Curfman McInnes }
109162fef451SLois Curfman McInnes 
10929b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
10939b94acceSBarry Smith 
10945615d1e5SSatish Balay #undef __FUNC__
10955615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
10969b94acceSBarry Smith /*@
10979b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1098272ac6f2SLois Curfman McInnes    of a nonlinear solver.
10999b94acceSBarry Smith 
11009b94acceSBarry Smith    Input Parameter:
11019b94acceSBarry Smith .  snes - the SNES context
11028ddd3da0SLois Curfman McInnes .  x - the solution vector
11039b94acceSBarry Smith 
1104272ac6f2SLois Curfman McInnes    Notes:
1105272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1106272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1107272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1108272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1109272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1110272ac6f2SLois Curfman McInnes 
11119b94acceSBarry Smith .keywords: SNES, nonlinear, setup
11129b94acceSBarry Smith 
11139b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
11149b94acceSBarry Smith @*/
11158ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
11169b94acceSBarry Smith {
1117272ac6f2SLois Curfman McInnes   int ierr, flg;
11183a40ed3dSBarry Smith 
11193a40ed3dSBarry Smith   PetscFunctionBegin;
112077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
112177c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
11228ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
11239b94acceSBarry Smith 
1124c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1125c1f60f51SBarry Smith   /*
1126c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1127dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1128c1f60f51SBarry Smith   */
1129c1f60f51SBarry Smith   if (flg) {
1130c1f60f51SBarry Smith     Mat J;
1131c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1132c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1133c1f60f51SBarry Smith     snes->mfshell = J;
1134c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1135c1f60f51SBarry Smith       snes->jacobian = J;
1136c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1137d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1138c1f60f51SBarry Smith       snes->jacobian = J;
1139c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1140d4bb536fSBarry Smith     } else {
1141a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option");
1142d4bb536fSBarry Smith     }
1143c1f60f51SBarry Smith   }
1144272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1145c1f60f51SBarry Smith   /*
1146dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1147c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1148c1f60f51SBarry Smith    */
114931615d04SLois Curfman McInnes   if (flg) {
1150272ac6f2SLois Curfman McInnes     Mat J;
1151272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1152272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1153272ac6f2SLois Curfman McInnes     snes->mfshell = J;
115483e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
115583e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
115694a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1157d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
115883e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
115994a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1160d4bb536fSBarry Smith     } else {
1161a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option");
1162d4bb536fSBarry Smith     }
1163272ac6f2SLois Curfman McInnes   }
11649b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1165a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1166a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1167a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first");
1168a8c6a408SBarry Smith     if (snes->vec_func == snes->vec_sol) {
1169a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector");
1170a8c6a408SBarry Smith     }
1171a703fe33SLois Curfman McInnes 
1172a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
117382bf6240SBarry Smith     if (snes->ksp_ewconv && PetscStrcmp(snes->type_name,SNES_EQ_TR)) {
1174a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1175a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1176a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1177a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1178a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1179a703fe33SLois Curfman McInnes     }
11809b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1181a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1182a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1183a8c6a408SBarry Smith     if (!snes->computeumfunction) {
1184a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first");
1185a8c6a408SBarry Smith     }
1186a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first");
1187d4bb536fSBarry Smith   } else {
1188a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
1189d4bb536fSBarry Smith   }
1190a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
119182bf6240SBarry Smith   snes->setupcalled = 1;
11923a40ed3dSBarry Smith   PetscFunctionReturn(0);
11939b94acceSBarry Smith }
11949b94acceSBarry Smith 
11955615d1e5SSatish Balay #undef __FUNC__
1196d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy"
11979b94acceSBarry Smith /*@C
11989b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11999b94acceSBarry Smith    with SNESCreate().
12009b94acceSBarry Smith 
12019b94acceSBarry Smith    Input Parameter:
12029b94acceSBarry Smith .  snes - the SNES context
12039b94acceSBarry Smith 
12049b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
12059b94acceSBarry Smith 
120663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
12079b94acceSBarry Smith @*/
12089b94acceSBarry Smith int SNESDestroy(SNES snes)
12099b94acceSBarry Smith {
12109b94acceSBarry Smith   int ierr;
12113a40ed3dSBarry Smith 
12123a40ed3dSBarry Smith   PetscFunctionBegin;
121377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12143a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1215d4bb536fSBarry Smith 
1216e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes); CHKERRQ(ierr);}
12170452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
1218522c5e43SBarry Smith   if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);}
12199b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
1220522c5e43SBarry Smith   if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);}
1221522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
12229b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
12230452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
12243a40ed3dSBarry Smith   PetscFunctionReturn(0);
12259b94acceSBarry Smith }
12269b94acceSBarry Smith 
12279b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12289b94acceSBarry Smith 
12295615d1e5SSatish Balay #undef __FUNC__
12305615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
12319b94acceSBarry Smith /*@
1232d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12339b94acceSBarry Smith 
12349b94acceSBarry Smith    Input Parameters:
12359b94acceSBarry Smith .  snes - the SNES context
123633174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
123733174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
123833174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
123933174efeSLois Curfman McInnes            of the change in the solution between steps
124033174efeSLois Curfman McInnes .  maxit - maximum number of iterations
124133174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
12429b94acceSBarry Smith 
124333174efeSLois Curfman McInnes    Options Database Keys:
124433174efeSLois Curfman McInnes $    -snes_atol <atol>
124533174efeSLois Curfman McInnes $    -snes_rtol <rtol>
124633174efeSLois Curfman McInnes $    -snes_stol <stol>
124733174efeSLois Curfman McInnes $    -snes_max_it <maxit>
124833174efeSLois Curfman McInnes $    -snes_max_funcs <maxf>
12499b94acceSBarry Smith 
1250d7a720efSLois Curfman McInnes    Notes:
12519b94acceSBarry Smith    The default maximum number of iterations is 50.
12529b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
12539b94acceSBarry Smith 
125433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
12559b94acceSBarry Smith 
1256d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
12579b94acceSBarry Smith @*/
1258d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
12599b94acceSBarry Smith {
12603a40ed3dSBarry Smith   PetscFunctionBegin;
126177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1262d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1263d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1264d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1265d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1266d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
12673a40ed3dSBarry Smith   PetscFunctionReturn(0);
12689b94acceSBarry Smith }
12699b94acceSBarry Smith 
12705615d1e5SSatish Balay #undef __FUNC__
12715615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
12729b94acceSBarry Smith /*@
127333174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
127433174efeSLois Curfman McInnes 
127533174efeSLois Curfman McInnes    Input Parameters:
127633174efeSLois Curfman McInnes .  snes - the SNES context
127733174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
127833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
127933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
128033174efeSLois Curfman McInnes            of the change in the solution between steps
128133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
128233174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
128333174efeSLois Curfman McInnes 
128433174efeSLois Curfman McInnes    Notes:
128533174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
128633174efeSLois Curfman McInnes 
128733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
128833174efeSLois Curfman McInnes 
128933174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
129033174efeSLois Curfman McInnes @*/
129133174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
129233174efeSLois Curfman McInnes {
12933a40ed3dSBarry Smith   PetscFunctionBegin;
129433174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
129533174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
129633174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
129733174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
129833174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
129933174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
13003a40ed3dSBarry Smith   PetscFunctionReturn(0);
130133174efeSLois Curfman McInnes }
130233174efeSLois Curfman McInnes 
13035615d1e5SSatish Balay #undef __FUNC__
13045615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
130533174efeSLois Curfman McInnes /*@
13069b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
13079b94acceSBarry Smith 
13089b94acceSBarry Smith    Input Parameters:
13099b94acceSBarry Smith .  snes - the SNES context
13109b94acceSBarry Smith .  tol - tolerance
13119b94acceSBarry Smith 
13129b94acceSBarry Smith    Options Database Key:
1313d7a720efSLois Curfman McInnes $    -snes_trtol <tol>
13149b94acceSBarry Smith 
13159b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
13169b94acceSBarry Smith 
1317d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
13189b94acceSBarry Smith @*/
13199b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
13209b94acceSBarry Smith {
13213a40ed3dSBarry Smith   PetscFunctionBegin;
132277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13239b94acceSBarry Smith   snes->deltatol = tol;
13243a40ed3dSBarry Smith   PetscFunctionReturn(0);
13259b94acceSBarry Smith }
13269b94acceSBarry Smith 
13275615d1e5SSatish Balay #undef __FUNC__
13285615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
13299b94acceSBarry Smith /*@
133077c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
13319b94acceSBarry Smith    for unconstrained minimization solvers.
13329b94acceSBarry Smith 
13339b94acceSBarry Smith    Input Parameters:
13349b94acceSBarry Smith .  snes - the SNES context
13359b94acceSBarry Smith .  ftol - minimum function tolerance
13369b94acceSBarry Smith 
13379b94acceSBarry Smith    Options Database Key:
1338d7a720efSLois Curfman McInnes $    -snes_fmin <ftol>
13399b94acceSBarry Smith 
13409b94acceSBarry Smith    Note:
134177c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
13429b94acceSBarry Smith    methods only.
13439b94acceSBarry Smith 
13449b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
13459b94acceSBarry Smith 
1346d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
13479b94acceSBarry Smith @*/
134877c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
13499b94acceSBarry Smith {
13503a40ed3dSBarry Smith   PetscFunctionBegin;
135177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13529b94acceSBarry Smith   snes->fmin = ftol;
13533a40ed3dSBarry Smith   PetscFunctionReturn(0);
13549b94acceSBarry Smith }
13559b94acceSBarry Smith 
13569b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
13579b94acceSBarry Smith 
13585615d1e5SSatish Balay #undef __FUNC__
1359d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor"
13609b94acceSBarry Smith /*@C
1361*5cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
13629b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
13639b94acceSBarry Smith    progress.
13649b94acceSBarry Smith 
13659b94acceSBarry Smith    Input Parameters:
13669b94acceSBarry Smith .  snes - the SNES context
13679b94acceSBarry Smith .  func - monitoring routine
1368044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1369044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
13709b94acceSBarry Smith 
13719b94acceSBarry Smith    Calling sequence of func:
1372682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
13739b94acceSBarry Smith 
13749b94acceSBarry Smith $    snes - the SNES context
13759b94acceSBarry Smith $    its - iteration number
1376044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
13779b94acceSBarry Smith $
13789b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
13799b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
13809b94acceSBarry Smith $
13819b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
13829b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
13839b94acceSBarry Smith 
13849665c990SLois Curfman McInnes    Options Database Keys:
13859665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
13869665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
13879665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
13889665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
13899665c990SLois Curfman McInnes $                           been hardwired into a code by
13909665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
13919665c990SLois Curfman McInnes $                           does not cancel those set via
13929665c990SLois Curfman McInnes $                           the options database.
13939665c990SLois Curfman McInnes 
13949665c990SLois Curfman McInnes 
1395639f9d9dSBarry Smith    Notes:
13966bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13976bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13986bc08f3fSLois Curfman McInnes    order in which they were set.
1399639f9d9dSBarry Smith 
14009b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
14019b94acceSBarry Smith 
1402*5cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
14039b94acceSBarry Smith @*/
140474679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
14059b94acceSBarry Smith {
14063a40ed3dSBarry Smith   PetscFunctionBegin;
1407639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1408a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
1409639f9d9dSBarry Smith   }
1410639f9d9dSBarry Smith 
1411639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1412639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
14133a40ed3dSBarry Smith   PetscFunctionReturn(0);
14149b94acceSBarry Smith }
14159b94acceSBarry Smith 
14165615d1e5SSatish Balay #undef __FUNC__
1417*5cd90555SBarry Smith #define __FUNC__ "SNESClearMonitor"
1418*5cd90555SBarry Smith /*@C
1419*5cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
1420*5cd90555SBarry Smith 
1421*5cd90555SBarry Smith    Input Parameters:
1422*5cd90555SBarry Smith .  snes - the SNES context
1423*5cd90555SBarry Smith 
1424*5cd90555SBarry Smith    Options Database:
1425*5cd90555SBarry Smith $    -snes_cancelmonitors : cancels all monitors that have
1426*5cd90555SBarry Smith $                           been hardwired into a code by
1427*5cd90555SBarry Smith $                           calls to SNESSetMonitor(), but
1428*5cd90555SBarry Smith $                           does not cancel those set via
1429*5cd90555SBarry Smith $                           the options database.
1430*5cd90555SBarry Smith 
1431*5cd90555SBarry Smith    Notes:
1432*5cd90555SBarry Smith      There is no way to clear one specific monitor from a SNES object.
1433*5cd90555SBarry Smith 
1434*5cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
1435*5cd90555SBarry Smith 
1436*5cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
1437*5cd90555SBarry Smith @*/
1438*5cd90555SBarry Smith int SNESClearMonitor( SNES snes )
1439*5cd90555SBarry Smith {
1440*5cd90555SBarry Smith   PetscFunctionBegin;
1441*5cd90555SBarry Smith   snes->numbermonitors = 0;
1442*5cd90555SBarry Smith   PetscFunctionReturn(0);
1443*5cd90555SBarry Smith }
1444*5cd90555SBarry Smith 
1445*5cd90555SBarry Smith #undef __FUNC__
1446d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest"
14479b94acceSBarry Smith /*@C
14489b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
14499b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
14509b94acceSBarry Smith 
14519b94acceSBarry Smith    Input Parameters:
14529b94acceSBarry Smith .  snes - the SNES context
14539b94acceSBarry Smith .  func - routine to test for convergence
1454044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1455044dda88SLois Curfman McInnes           (may be PETSC_NULL)
14569b94acceSBarry Smith 
14579b94acceSBarry Smith    Calling sequence of func:
14589b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
14599b94acceSBarry Smith              double f,void *cctx)
14609b94acceSBarry Smith 
14619b94acceSBarry Smith $    snes - the SNES context
1462044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
14639b94acceSBarry Smith $    xnorm - 2-norm of current iterate
14649b94acceSBarry Smith $
14659b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14669b94acceSBarry Smith $    gnorm - 2-norm of current step
14679b94acceSBarry Smith $    f - 2-norm of function
14689b94acceSBarry Smith $
14699b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14709b94acceSBarry Smith $    gnorm - 2-norm of current gradient
14719b94acceSBarry Smith $    f - function value
14729b94acceSBarry Smith 
14739b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14749b94acceSBarry Smith 
147540191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
147640191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
14779b94acceSBarry Smith @*/
147874679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
14799b94acceSBarry Smith {
14803a40ed3dSBarry Smith   PetscFunctionBegin;
14819b94acceSBarry Smith   (snes)->converged = func;
14829b94acceSBarry Smith   (snes)->cnvP      = cctx;
14833a40ed3dSBarry Smith   PetscFunctionReturn(0);
14849b94acceSBarry Smith }
14859b94acceSBarry Smith 
14865615d1e5SSatish Balay #undef __FUNC__
1487d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory"
1488c9005455SLois Curfman McInnes /*@
1489c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1490c9005455SLois Curfman McInnes 
1491c9005455SLois Curfman McInnes    Input Parameters:
1492c9005455SLois Curfman McInnes .  snes - iterative context obtained from SNESCreate()
1493c9005455SLois Curfman McInnes .  a   - array to hold history
1494c9005455SLois Curfman McInnes .  na  - size of a
1495c9005455SLois Curfman McInnes 
1496c9005455SLois Curfman McInnes    Notes:
1497c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1498c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1499c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1500c9005455SLois Curfman McInnes    at each step.
1501c9005455SLois Curfman McInnes 
1502c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1503c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1504c9005455SLois Curfman McInnes    during the section of code that is being timed.
1505c9005455SLois Curfman McInnes 
1506c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1507c9005455SLois Curfman McInnes @*/
1508c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1509c9005455SLois Curfman McInnes {
15103a40ed3dSBarry Smith   PetscFunctionBegin;
1511c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1512c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1513c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1514c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
15153a40ed3dSBarry Smith   PetscFunctionReturn(0);
1516c9005455SLois Curfman McInnes }
1517c9005455SLois Curfman McInnes 
1518c9005455SLois Curfman McInnes #undef __FUNC__
15195615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
15209b94acceSBarry Smith /*
15219b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
15229b94acceSBarry Smith    positive parameter delta.
15239b94acceSBarry Smith 
15249b94acceSBarry Smith     Input Parameters:
15259b94acceSBarry Smith .   snes - the SNES context
15269b94acceSBarry Smith .   y - approximate solution of linear system
15279b94acceSBarry Smith .   fnorm - 2-norm of current function
15289b94acceSBarry Smith .   delta - trust region size
15299b94acceSBarry Smith 
15309b94acceSBarry Smith     Output Parameters:
15319b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
15329b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
15339b94acceSBarry Smith     region, and exceeds zero otherwise.
15349b94acceSBarry Smith .   ynorm - 2-norm of the step
15359b94acceSBarry Smith 
15369b94acceSBarry Smith     Note:
153740191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
15389b94acceSBarry Smith     is set to be the maximum allowable step size.
15399b94acceSBarry Smith 
15409b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
15419b94acceSBarry Smith */
15429b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
15439b94acceSBarry Smith                           double *gpnorm,double *ynorm)
15449b94acceSBarry Smith {
15459b94acceSBarry Smith   double norm;
15469b94acceSBarry Smith   Scalar cnorm;
15473a40ed3dSBarry Smith   int    ierr;
15483a40ed3dSBarry Smith 
15493a40ed3dSBarry Smith   PetscFunctionBegin;
15503a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
15519b94acceSBarry Smith   if (norm > *delta) {
15529b94acceSBarry Smith      norm = *delta/norm;
15539b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
15549b94acceSBarry Smith      cnorm = norm;
15559b94acceSBarry Smith      VecScale( &cnorm, y );
15569b94acceSBarry Smith      *ynorm = *delta;
15579b94acceSBarry Smith   } else {
15589b94acceSBarry Smith      *gpnorm = 0.0;
15599b94acceSBarry Smith      *ynorm = norm;
15609b94acceSBarry Smith   }
15613a40ed3dSBarry Smith   PetscFunctionReturn(0);
15629b94acceSBarry Smith }
15639b94acceSBarry Smith 
15645615d1e5SSatish Balay #undef __FUNC__
15655615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
15669b94acceSBarry Smith /*@
15679b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
156863a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
15699b94acceSBarry Smith 
15709b94acceSBarry Smith    Input Parameter:
15719b94acceSBarry Smith .  snes - the SNES context
15728ddd3da0SLois Curfman McInnes .  x - the solution vector
15739b94acceSBarry Smith 
15749b94acceSBarry Smith    Output Parameter:
15759b94acceSBarry Smith    its - number of iterations until termination
15769b94acceSBarry Smith 
15778ddd3da0SLois Curfman McInnes    Note:
15788ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
15798ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
15808ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
15818ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
15828ddd3da0SLois Curfman McInnes 
15839b94acceSBarry Smith .keywords: SNES, nonlinear, solve
15849b94acceSBarry Smith 
158563a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
15869b94acceSBarry Smith @*/
15878ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
15889b94acceSBarry Smith {
15893c7409f5SSatish Balay   int ierr, flg;
1590052efed2SBarry Smith 
15913a40ed3dSBarry Smith   PetscFunctionBegin;
159277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
159374679c65SBarry Smith   PetscValidIntPointer(its);
159482bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1595c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
15969b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1597c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
15989b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
15999b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
16003c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
16016d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
16023a40ed3dSBarry Smith   PetscFunctionReturn(0);
16039b94acceSBarry Smith }
16049b94acceSBarry Smith 
16059b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
16069b94acceSBarry Smith 
16075615d1e5SSatish Balay #undef __FUNC__
16085615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
160982bf6240SBarry Smith /*@C
16104b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
16119b94acceSBarry Smith 
16129b94acceSBarry Smith    Input Parameters:
16139b94acceSBarry Smith .  snes - the SNES context
16149b94acceSBarry Smith .  method - a known method
16159b94acceSBarry Smith 
1616ae12b187SLois Curfman McInnes   Options Database Command:
1617ae12b187SLois Curfman McInnes $ -snes_type  <method>
1618ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1619ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1620ae12b187SLois Curfman McInnes 
16219b94acceSBarry Smith    Notes:
16229b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
16239b94acceSBarry Smith $  Systems of nonlinear equations:
162440191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
162540191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
16269b94acceSBarry Smith $  Unconstrained minimization:
162740191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
162840191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
16299b94acceSBarry Smith 
1630ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1631ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1632ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1633ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1634ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1635ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1636ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1637ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1638ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1639ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1640a703fe33SLois Curfman McInnes 
1641f536c699SSatish Balay .keywords: SNES, set, method
16429b94acceSBarry Smith @*/
16434b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
16449b94acceSBarry Smith {
164584cb2905SBarry Smith   int ierr;
16469b94acceSBarry Smith   int (*r)(SNES);
16473a40ed3dSBarry Smith 
16483a40ed3dSBarry Smith   PetscFunctionBegin;
164977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
165082bf6240SBarry Smith 
165182bf6240SBarry Smith   if (!PetscStrcmp(snes->type_name,method)) PetscFunctionReturn(0);
165282bf6240SBarry Smith 
165382bf6240SBarry Smith   if (snes->setupcalled) {
1654e1311b90SBarry Smith     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
165582bf6240SBarry Smith     snes->data = 0;
165682bf6240SBarry Smith   }
16577d1a2b2bSBarry Smith 
16589b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
165982bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);}
166082bf6240SBarry Smith 
1661ecf371e4SBarry Smith   ierr =  DLRegisterFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr);
166282bf6240SBarry Smith 
16630452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
166482bf6240SBarry Smith   snes->data = 0;
16653a40ed3dSBarry Smith   ierr = (*r)(snes); CHKERRQ(ierr);
166682bf6240SBarry Smith 
166782bf6240SBarry Smith   if (snes->type_name) PetscFree(snes->type_name);
166882bf6240SBarry Smith   snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name);
166982bf6240SBarry Smith   PetscStrcpy(snes->type_name,method);
167082bf6240SBarry Smith   snes->set_method_called = 1;
167182bf6240SBarry Smith 
16723a40ed3dSBarry Smith   PetscFunctionReturn(0);
16739b94acceSBarry Smith }
16749b94acceSBarry Smith 
1675a847f771SSatish Balay 
16769b94acceSBarry Smith /* --------------------------------------------------------------------- */
16775615d1e5SSatish Balay #undef __FUNC__
1678d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy"
16799b94acceSBarry Smith /*@C
16809b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
16819b94acceSBarry Smith    registered by SNESRegister().
16829b94acceSBarry Smith 
16839b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
16849b94acceSBarry Smith 
16859b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
16869b94acceSBarry Smith @*/
1687cf256101SBarry Smith int SNESRegisterDestroy(void)
16889b94acceSBarry Smith {
168982bf6240SBarry Smith   int ierr;
169082bf6240SBarry Smith 
16913a40ed3dSBarry Smith   PetscFunctionBegin;
169282bf6240SBarry Smith   if (SNESList) {
169382bf6240SBarry Smith     ierr = DLRegisterDestroy( SNESList );CHKERRQ(ierr);
169482bf6240SBarry Smith     SNESList = 0;
16959b94acceSBarry Smith   }
169684cb2905SBarry Smith   SNESRegisterAllCalled = 0;
16973a40ed3dSBarry Smith   PetscFunctionReturn(0);
16989b94acceSBarry Smith }
16999b94acceSBarry Smith 
17005615d1e5SSatish Balay #undef __FUNC__
1701d4bb536fSBarry Smith #define __FUNC__ "SNESGetType"
17029b94acceSBarry Smith /*@C
17039a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
17049b94acceSBarry Smith 
17059b94acceSBarry Smith    Input Parameter:
17064b0e389bSBarry Smith .  snes - nonlinear solver context
17079b94acceSBarry Smith 
17089b94acceSBarry Smith    Output Parameter:
170982bf6240SBarry Smith .  method - SNES method (a charactor string)
17109b94acceSBarry Smith 
17119b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
17129b94acceSBarry Smith @*/
171382bf6240SBarry Smith int SNESGetType(SNES snes, SNESType *method)
17149b94acceSBarry Smith {
17153a40ed3dSBarry Smith   PetscFunctionBegin;
171682bf6240SBarry Smith   *method = snes->type_name;
17173a40ed3dSBarry Smith   PetscFunctionReturn(0);
17189b94acceSBarry Smith }
17199b94acceSBarry Smith 
17205615d1e5SSatish Balay #undef __FUNC__
1721d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution"
17229b94acceSBarry Smith /*@C
17239b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
17249b94acceSBarry Smith    stored.
17259b94acceSBarry Smith 
17269b94acceSBarry Smith    Input Parameter:
17279b94acceSBarry Smith .  snes - the SNES context
17289b94acceSBarry Smith 
17299b94acceSBarry Smith    Output Parameter:
17309b94acceSBarry Smith .  x - the solution
17319b94acceSBarry Smith 
17329b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
17339b94acceSBarry Smith 
1734a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
17359b94acceSBarry Smith @*/
17369b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
17379b94acceSBarry Smith {
17383a40ed3dSBarry Smith   PetscFunctionBegin;
173977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17409b94acceSBarry Smith   *x = snes->vec_sol_always;
17413a40ed3dSBarry Smith   PetscFunctionReturn(0);
17429b94acceSBarry Smith }
17439b94acceSBarry Smith 
17445615d1e5SSatish Balay #undef __FUNC__
1745d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate"
17469b94acceSBarry Smith /*@C
17479b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
17489b94acceSBarry Smith    stored.
17499b94acceSBarry Smith 
17509b94acceSBarry Smith    Input Parameter:
17519b94acceSBarry Smith .  snes - the SNES context
17529b94acceSBarry Smith 
17539b94acceSBarry Smith    Output Parameter:
17549b94acceSBarry Smith .  x - the solution update
17559b94acceSBarry Smith 
17569b94acceSBarry Smith    Notes:
17579b94acceSBarry Smith    This vector is implementation dependent.
17589b94acceSBarry Smith 
17599b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
17609b94acceSBarry Smith 
17619b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
17629b94acceSBarry Smith @*/
17639b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
17649b94acceSBarry Smith {
17653a40ed3dSBarry Smith   PetscFunctionBegin;
176677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17679b94acceSBarry Smith   *x = snes->vec_sol_update_always;
17683a40ed3dSBarry Smith   PetscFunctionReturn(0);
17699b94acceSBarry Smith }
17709b94acceSBarry Smith 
17715615d1e5SSatish Balay #undef __FUNC__
1772d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction"
17739b94acceSBarry Smith /*@C
17743638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
17759b94acceSBarry Smith 
17769b94acceSBarry Smith    Input Parameter:
17779b94acceSBarry Smith .  snes - the SNES context
17789b94acceSBarry Smith 
17799b94acceSBarry Smith    Output Parameter:
17803638b69dSLois Curfman McInnes .  r - the function
17819b94acceSBarry Smith 
17829b94acceSBarry Smith    Notes:
17839b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
17849b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
17859b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
17869b94acceSBarry Smith 
1787a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
17889b94acceSBarry Smith 
178931615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
179031615d04SLois Curfman McInnes           SNESGetGradient()
17919b94acceSBarry Smith @*/
17929b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
17939b94acceSBarry Smith {
17943a40ed3dSBarry Smith   PetscFunctionBegin;
179577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1796a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1797a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1798a8c6a408SBarry Smith   }
17999b94acceSBarry Smith   *r = snes->vec_func_always;
18003a40ed3dSBarry Smith   PetscFunctionReturn(0);
18019b94acceSBarry Smith }
18029b94acceSBarry Smith 
18035615d1e5SSatish Balay #undef __FUNC__
1804d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient"
18059b94acceSBarry Smith /*@C
18063638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
18079b94acceSBarry Smith 
18089b94acceSBarry Smith    Input Parameter:
18099b94acceSBarry Smith .  snes - the SNES context
18109b94acceSBarry Smith 
18119b94acceSBarry Smith    Output Parameter:
18129b94acceSBarry Smith .  r - the gradient
18139b94acceSBarry Smith 
18149b94acceSBarry Smith    Notes:
18159b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
18169b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18179b94acceSBarry Smith    SNESGetFunction().
18189b94acceSBarry Smith 
18199b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
18209b94acceSBarry Smith 
182131615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
18229b94acceSBarry Smith @*/
18239b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
18249b94acceSBarry Smith {
18253a40ed3dSBarry Smith   PetscFunctionBegin;
182677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18273a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1828a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
18293a40ed3dSBarry Smith   }
18309b94acceSBarry Smith   *r = snes->vec_func_always;
18313a40ed3dSBarry Smith   PetscFunctionReturn(0);
18329b94acceSBarry Smith }
18339b94acceSBarry Smith 
18345615d1e5SSatish Balay #undef __FUNC__
1835d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction"
18369b94acceSBarry Smith /*@
18379b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
18389b94acceSBarry Smith    unconstrained minimization problems.
18399b94acceSBarry Smith 
18409b94acceSBarry Smith    Input Parameter:
18419b94acceSBarry Smith .  snes - the SNES context
18429b94acceSBarry Smith 
18439b94acceSBarry Smith    Output Parameter:
18449b94acceSBarry Smith .  r - the function
18459b94acceSBarry Smith 
18469b94acceSBarry Smith    Notes:
18479b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
18489b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18499b94acceSBarry Smith    SNESGetFunction().
18509b94acceSBarry Smith 
18519b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
18529b94acceSBarry Smith 
185331615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
18549b94acceSBarry Smith @*/
18559b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
18569b94acceSBarry Smith {
18573a40ed3dSBarry Smith   PetscFunctionBegin;
185877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
185974679c65SBarry Smith   PetscValidScalarPointer(r);
18603a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1861a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
18623a40ed3dSBarry Smith   }
18639b94acceSBarry Smith   *r = snes->fc;
18643a40ed3dSBarry Smith   PetscFunctionReturn(0);
18659b94acceSBarry Smith }
18669b94acceSBarry Smith 
18675615d1e5SSatish Balay #undef __FUNC__
1868d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix"
18693c7409f5SSatish Balay /*@C
18703c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1871d850072dSLois Curfman McInnes    SNES options in the database.
18723c7409f5SSatish Balay 
18733c7409f5SSatish Balay    Input Parameter:
18743c7409f5SSatish Balay .  snes - the SNES context
18753c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
18763c7409f5SSatish Balay 
1877d850072dSLois Curfman McInnes    Notes:
1878a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1879a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1880a83b1b31SSatish Balay    hyphen.
1881d850072dSLois Curfman McInnes 
18823c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1883a86d99e1SLois Curfman McInnes 
1884a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
18853c7409f5SSatish Balay @*/
18863c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
18873c7409f5SSatish Balay {
18883c7409f5SSatish Balay   int ierr;
18893c7409f5SSatish Balay 
18903a40ed3dSBarry Smith   PetscFunctionBegin;
189177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1892639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
18933c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
18943a40ed3dSBarry Smith   PetscFunctionReturn(0);
18953c7409f5SSatish Balay }
18963c7409f5SSatish Balay 
18975615d1e5SSatish Balay #undef __FUNC__
1898d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix"
18993c7409f5SSatish Balay /*@C
1900f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1901d850072dSLois Curfman McInnes    SNES options in the database.
19023c7409f5SSatish Balay 
19033c7409f5SSatish Balay    Input Parameter:
19043c7409f5SSatish Balay .  snes - the SNES context
19053c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
19063c7409f5SSatish Balay 
1907d850072dSLois Curfman McInnes    Notes:
1908a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1909a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1910a83b1b31SSatish Balay    hyphen.
1911d850072dSLois Curfman McInnes 
19123c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1913a86d99e1SLois Curfman McInnes 
1914a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
19153c7409f5SSatish Balay @*/
19163c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
19173c7409f5SSatish Balay {
19183c7409f5SSatish Balay   int ierr;
19193c7409f5SSatish Balay 
19203a40ed3dSBarry Smith   PetscFunctionBegin;
192177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1922639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19233c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19243a40ed3dSBarry Smith   PetscFunctionReturn(0);
19253c7409f5SSatish Balay }
19263c7409f5SSatish Balay 
19275615d1e5SSatish Balay #undef __FUNC__
1928d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix"
1929c83590e2SSatish Balay /*@
19303c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
19313c7409f5SSatish Balay    SNES options in the database.
19323c7409f5SSatish Balay 
19333c7409f5SSatish Balay    Input Parameter:
19343c7409f5SSatish Balay .  snes - the SNES context
19353c7409f5SSatish Balay 
19363c7409f5SSatish Balay    Output Parameter:
19373c7409f5SSatish Balay .  prefix - pointer to the prefix string used
19383c7409f5SSatish Balay 
19393c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
1940a86d99e1SLois Curfman McInnes 
1941a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
19423c7409f5SSatish Balay @*/
19433c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
19443c7409f5SSatish Balay {
19453c7409f5SSatish Balay   int ierr;
19463c7409f5SSatish Balay 
19473a40ed3dSBarry Smith   PetscFunctionBegin;
194877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1949639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19503a40ed3dSBarry Smith   PetscFunctionReturn(0);
19513c7409f5SSatish Balay }
19523c7409f5SSatish Balay 
1953ca161407SBarry Smith #undef __FUNC__
1954ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp"
1955ca161407SBarry Smith /*@
1956ca161407SBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
1957ca161407SBarry Smith 
1958ca161407SBarry Smith    Input Parameter:
1959ca161407SBarry Smith .  snes - the SNES context
1960ca161407SBarry Smith 
1961ca161407SBarry Smith    Options Database Keys:
1962ca161407SBarry Smith $  -help, -h
1963ca161407SBarry Smith 
1964ca161407SBarry Smith .keywords: SNES, nonlinear, help
1965ca161407SBarry Smith 
1966ca161407SBarry Smith .seealso: SNESSetFromOptions()
1967ca161407SBarry Smith @*/
1968ca161407SBarry Smith int SNESPrintHelp(SNES snes)
1969ca161407SBarry Smith {
1970ca161407SBarry Smith   char                p[64];
1971ca161407SBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
1972ca161407SBarry Smith   int                 ierr;
1973ca161407SBarry Smith 
1974ca161407SBarry Smith   PetscFunctionBegin;
1975ca161407SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1976ca161407SBarry Smith 
1977ca161407SBarry Smith   PetscStrcpy(p,"-");
1978ca161407SBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
1979ca161407SBarry Smith 
1980ca161407SBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
1981ca161407SBarry Smith 
198276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");
198382bf6240SBarry Smith   ierr = DLRegisterPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr);
198476be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
198576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
198676be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
198776be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
198876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
198976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
199076be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
199176be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");
199276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
199376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
1994ca161407SBarry Smith     residual norm at each iteration.\n",p);
199576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
1996ca161407SBarry Smith     residual norm for small residual norms. This is useful to conceal\n\
1997ca161407SBarry Smith     meaningless digits that may be different on different machines.\n",p);
199876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
1999ca161407SBarry Smith   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
200076be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2001ca161407SBarry Smith      " Options for solving systems of nonlinear equations only:\n");
200276be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
200376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
200476be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
200576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
200676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
200776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2008ca161407SBarry Smith      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
200976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2010ca161407SBarry Smith      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
201176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2012ca161407SBarry Smith      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
201376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2014ca161407SBarry Smith      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
201576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2016ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
201776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2018ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
201976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2020ca161407SBarry Smith      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
2021ca161407SBarry Smith   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
202276be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2023ca161407SBarry Smith      " Options for solving unconstrained minimization problems only:\n");
202476be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
202576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
202676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
2027ca161407SBarry Smith   }
202876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
202976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"a particular method\n");
2030ca161407SBarry Smith   if (snes->printhelp) {
2031ca161407SBarry Smith     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2032ca161407SBarry Smith   }
2033ca161407SBarry Smith   PetscFunctionReturn(0);
2034ca161407SBarry Smith }
20353c7409f5SSatish Balay 
20369b94acceSBarry Smith 
20379b94acceSBarry Smith 
20389b94acceSBarry Smith 
20393a40ed3dSBarry Smith 
2040