xref: /petsc/src/snes/interface/snes.c (revision 0ef38995692ee30aa448304a98abd5bc1b8605a7)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*0ef38995SBarry Smith static char vcid[] = "$Id: snes.c,v 1.164 1998/12/17 22:12:19 bsmith Exp bsmith $";
39b94acceSBarry Smith #endif
49b94acceSBarry Smith 
570f55243SBarry Smith #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
6f5eb4b81SSatish Balay #include "src/sys/nreg.h"
79b94acceSBarry Smith 
884cb2905SBarry Smith int SNESRegisterAllCalled = 0;
9488ecbafSBarry Smith FList SNESList = 0;
1082bf6240SBarry Smith 
115615d1e5SSatish Balay #undef __FUNC__
12d4bb536fSBarry Smith #define __FUNC__ "SNESView"
139b94acceSBarry Smith /*@
149b94acceSBarry Smith    SNESView - Prints the SNES data structure.
159b94acceSBarry Smith 
16fee21e36SBarry Smith    Collective on SNES, unless Viewer is VIEWER_STDOUT_SELF
17fee21e36SBarry Smith 
18c7afd0dbSLois Curfman McInnes    Input Parameters:
19c7afd0dbSLois Curfman McInnes +  SNES - the SNES context
20c7afd0dbSLois Curfman McInnes -  viewer - visualization context
21c7afd0dbSLois Curfman McInnes 
229b94acceSBarry Smith    Options Database Key:
23c8a8ba5cSLois Curfman McInnes .  -snes_view - Calls SNESView() at end of SNESSolve()
249b94acceSBarry Smith 
259b94acceSBarry Smith    Notes:
269b94acceSBarry Smith    The available visualization contexts include
27c7afd0dbSLois Curfman McInnes +     VIEWER_STDOUT_SELF - standard output (default)
28c7afd0dbSLois Curfman McInnes -     VIEWER_STDOUT_WORLD - synchronized standard
29c8a8ba5cSLois Curfman McInnes          output where only the first processor opens
30c8a8ba5cSLois Curfman McInnes          the file.  All other processors send their
31c8a8ba5cSLois Curfman McInnes          data to the first processor to print.
329b94acceSBarry Smith 
333e081fefSLois Curfman McInnes    The user can open an alternative visualization context with
3477ed5343SBarry Smith    ViewerASCIIOpen() - output to a specified file.
359b94acceSBarry Smith 
369b94acceSBarry Smith .keywords: SNES, view
379b94acceSBarry Smith 
3877ed5343SBarry Smith .seealso: ViewerASCIIOpen()
399b94acceSBarry Smith @*/
409b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
419b94acceSBarry Smith {
429b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
439b94acceSBarry Smith   int                 ierr;
449b94acceSBarry Smith   SLES                sles;
459b94acceSBarry Smith   char                *method;
4619bcc07fSBarry Smith   ViewerType          vtype;
479b94acceSBarry Smith 
483a40ed3dSBarry Smith   PetscFunctionBegin;
4974679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5074679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5174679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5274679c65SBarry Smith 
5319bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
543f1db9ecSBarry Smith   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
55*0ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"SNES Object:\n");
5682bf6240SBarry Smith     SNESGetType(snes,&method);
57*0ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  method: %s\n",method);
58*0ef38995SBarry Smith     if (snes->view) {
59*0ef38995SBarry Smith       ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
60*0ef38995SBarry Smith       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
61*0ef38995SBarry Smith       ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
62*0ef38995SBarry Smith     }
63*0ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);
64*0ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
659b94acceSBarry Smith                  snes->rtol, snes->atol, snes->trunctol, snes->xtol);
66*0ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",snes->linear_its);
67*0ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  total number of function evaluations=%d\n",snes->nfuncs);
68*0ef38995SBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
69*0ef38995SBarry Smith       ViewerASCIIPrintf(viewer,"  min function tolerance=%g\n",snes->fmin);
70*0ef38995SBarry Smith     }
719b94acceSBarry Smith     if (snes->ksp_ewconv) {
729b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
739b94acceSBarry Smith       if (kctx) {
74*0ef38995SBarry Smith         ViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);
75*0ef38995SBarry Smith         ViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);
76*0ef38995SBarry Smith         ViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);
779b94acceSBarry Smith       }
789b94acceSBarry Smith     }
793f1db9ecSBarry Smith   } else if (PetscTypeCompare(vtype,STRING_VIEWER)) {
80*0ef38995SBarry Smith     ierr = SNESGetType(snes,&method);CHKERRQ(ierr);
81*0ef38995SBarry Smith     ierr = ViewerStringSPrintf(viewer," %-3.3s",method);CHKERRQ(ierr);
8219bcc07fSBarry Smith   }
8377ed5343SBarry Smith   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
84*0ef38995SBarry Smith   ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
859b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
86*0ef38995SBarry Smith   ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
873a40ed3dSBarry Smith   PetscFunctionReturn(0);
889b94acceSBarry Smith }
899b94acceSBarry Smith 
90639f9d9dSBarry Smith /*
91639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
92639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
93639f9d9dSBarry Smith */
94639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
95639f9d9dSBarry Smith static int numberofsetfromoptions;
96639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
97639f9d9dSBarry Smith 
985615d1e5SSatish Balay #undef __FUNC__
99d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker"
100639f9d9dSBarry Smith /*@
101639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
102639f9d9dSBarry Smith 
103c7afd0dbSLois Curfman McInnes     Not Collective
104c7afd0dbSLois Curfman McInnes 
105639f9d9dSBarry Smith     Input Parameter:
106639f9d9dSBarry Smith .   snescheck - function that checks for options
107639f9d9dSBarry Smith 
108639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
109639f9d9dSBarry Smith @*/
110639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
111639f9d9dSBarry Smith {
1123a40ed3dSBarry Smith   PetscFunctionBegin;
113639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
114a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed");
115639f9d9dSBarry Smith   }
116639f9d9dSBarry Smith 
117639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
1183a40ed3dSBarry Smith   PetscFunctionReturn(0);
119639f9d9dSBarry Smith }
120639f9d9dSBarry Smith 
1215615d1e5SSatish Balay #undef __FUNC__
1225615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1239b94acceSBarry Smith /*@
124682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1259b94acceSBarry Smith 
126c7afd0dbSLois Curfman McInnes    Collective on SNES
127c7afd0dbSLois Curfman McInnes 
1289b94acceSBarry Smith    Input Parameter:
1299b94acceSBarry Smith .  snes - the SNES context
1309b94acceSBarry Smith 
13182738288SBarry Smith    Options Database:
13282738288SBarry Smith +  -snes_type type - SNES_EQ_LS, SNES_EQ_TR, SNES_UM_TR, SNES_UM_LS etc
13382738288SBarry Smith .  -snes_stol - convergence tolerance in terms of the norm
13482738288SBarry Smith                 of the change in the solution between steps
13582738288SBarry Smith .  -snes_atol atol - absolute tolerance of residual norm
13682738288SBarry Smith .  -snes_rtol rtol - relative decrease in tolerance norm from initial
13782738288SBarry Smith .  -snes_max_it max_it - maximum number of iterations
13882738288SBarry Smith .  -snes_max_funcs max_funcs - maximum number of function evaluations
13982738288SBarry Smith .  -snes_trtol trtol - trust region tolerance
14082738288SBarry Smith .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
14182738288SBarry Smith                                solver; hence iterations will continue until max_it
14282738288SBarry Smith                                or some other criteria is reached. Saves expense
14382738288SBarry Smith                                of convergence test
14482738288SBarry Smith .  -snes_monitor - prints residual norm at each iteration
14582738288SBarry Smith .  -snes_xmonitor - plots residual norm at each iteration
146e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
147e24b481bSBarry Smith -  -snes_mf_ksp_monitor - if using matrix-free multiply then prints h at each KSP iteration
14882738288SBarry Smith 
14982738288SBarry Smith     Options Database for Eisenstat-Walker method:
15082738288SBarry Smith +  -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence
15182738288SBarry Smith .  -snes_ksp_eq_version ver - version of  Eisenstat-Walker method
15282738288SBarry Smith .  -snes_ksp_ew_rtol0 rtol0 -
15382738288SBarry Smith .  -snes_ksp_ew_rtolmax rtolmax -
15482738288SBarry Smith .  -snes_ksp_ew_gamma gamma -
15582738288SBarry Smith .  -snes_ksp_ew_alpha alpha -
15682738288SBarry Smith .  -snes_ksp_ew_alpha2 alpha2 -
15782738288SBarry Smith -  -snes_ksp_ew_threshold threshold -
15882738288SBarry Smith 
15911ca99fdSLois Curfman McInnes    Notes:
16011ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
16111ca99fdSLois Curfman McInnes    the users manual.
16283e2fdc7SBarry Smith 
1639b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1649b94acceSBarry Smith 
165a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1669b94acceSBarry Smith @*/
1679b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1689b94acceSBarry Smith {
16982bf6240SBarry Smith   char     method[256];
1709b94acceSBarry Smith   double   tmp;
1719b94acceSBarry Smith   SLES     sles;
17281f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1733c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1749b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1759b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1769b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1779b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1789b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1799b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1809b94acceSBarry Smith 
1813a40ed3dSBarry Smith   PetscFunctionBegin;
18281f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
18381f57ec7SBarry Smith 
18477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18582bf6240SBarry Smith   if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp");
186ca161407SBarry Smith 
18782bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
18882bf6240SBarry Smith   ierr = OptionsGetString(snes->prefix,"-snes_type",method,256,&flg);
189052efed2SBarry Smith   if (flg) {
19082bf6240SBarry Smith     ierr = SNESSetType(snes,(SNESType) method); CHKERRQ(ierr);
1915334005bSBarry Smith   }
1923c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
193d64ed03dSBarry Smith   if (flg) {
194d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
195d64ed03dSBarry Smith   }
1963c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
197d64ed03dSBarry Smith   if (flg) {
198d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
199d64ed03dSBarry Smith   }
2003c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
201d64ed03dSBarry Smith   if (flg) {
202d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
203d64ed03dSBarry Smith   }
2043c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
2053c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
206d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
207d64ed03dSBarry Smith   if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); }
208d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
209d64ed03dSBarry Smith   if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);  CHKERRQ(ierr);}
2103c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
2113c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
212b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
213b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
214b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
215b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
216b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
217b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
218b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
21993c39befSBarry Smith 
22093c39befSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr);
22193c39befSBarry Smith   if (flg) {snes->converged = 0;}
22293c39befSBarry Smith 
2239b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
2249b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
2255bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
2265cd90555SBarry Smith   if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
2273c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
228639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
2293c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
2303f1db9ecSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);}
2313f1db9ecSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor",&flg); CHKERRQ(ierr);
2323f1db9ecSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0);  CHKERRQ(ierr);}
23381f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
2343c7409f5SSatish Balay   if (flg) {
23517699dbbSLois Curfman McInnes     int    rank = 0;
236d7e8b826SBarry Smith     DrawLG lg;
23717699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
23817699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
23917699dbbSLois Curfman McInnes     if (!rank) {
24081f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2419b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
242f8c826e1SBarry Smith       snes->xmonitor = lg;
2439b94acceSBarry Smith       PLogObjectParent(snes,lg);
2449b94acceSBarry Smith     }
2459b94acceSBarry Smith   }
246e24b481bSBarry Smith 
247e24b481bSBarry Smith 
2483c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2493c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2509b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2519b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
252981c4779SBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
253981c4779SBarry Smith   } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
25431615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
25531615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
256d64ed03dSBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2579b94acceSBarry Smith   }
258639f9d9dSBarry Smith 
259639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
260639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
261639f9d9dSBarry Smith   }
262639f9d9dSBarry Smith 
2639b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2649b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
26593993e2dSLois Curfman McInnes 
266e24b481bSBarry Smith   /* set the special KSP monitor for matrix-free application */
267e24b481bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_ksp_monitor",&flg); CHKERRQ(ierr);
268e24b481bSBarry Smith   if (flg) {
269e24b481bSBarry Smith     KSP ksp;
270e24b481bSBarry Smith     ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
271e24b481bSBarry Smith     ierr = KSPSetMonitor(ksp,MatSNESFDMFKSPMonitor,PETSC_NULL);CHKERRQ(ierr);
272e24b481bSBarry Smith   }
273e24b481bSBarry Smith 
27482bf6240SBarry Smith   /*
27593993e2dSLois Curfman McInnes       Since the private setfromoptions requires the type to have
27693993e2dSLois Curfman McInnes       been set already, we make sure a type is set by this time.
27782bf6240SBarry Smith   */
27882bf6240SBarry Smith   if (!snes->type_name) {
27982bf6240SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
28082bf6240SBarry Smith       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
28182bf6240SBarry Smith     } else {
28282bf6240SBarry Smith       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
28382bf6240SBarry Smith     }
28482bf6240SBarry Smith   }
28582bf6240SBarry Smith 
28682bf6240SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
28782bf6240SBarry Smith   if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);}
28882bf6240SBarry Smith 
2893a40ed3dSBarry Smith   if (snes->setfromoptions) {
2903a40ed3dSBarry Smith     ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
2913a40ed3dSBarry Smith   }
2923a40ed3dSBarry Smith   PetscFunctionReturn(0);
2939b94acceSBarry Smith }
2949b94acceSBarry Smith 
295a847f771SSatish Balay 
2965615d1e5SSatish Balay #undef __FUNC__
297d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext"
2989b94acceSBarry Smith /*@
2999b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3009b94acceSBarry Smith    the nonlinear solvers.
3019b94acceSBarry Smith 
302fee21e36SBarry Smith    Collective on SNES
303fee21e36SBarry Smith 
304c7afd0dbSLois Curfman McInnes    Input Parameters:
305c7afd0dbSLois Curfman McInnes +  snes - the SNES context
306c7afd0dbSLois Curfman McInnes -  usrP - optional user context
307c7afd0dbSLois Curfman McInnes 
3089b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3099b94acceSBarry Smith 
3109b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3119b94acceSBarry Smith @*/
3129b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3139b94acceSBarry Smith {
3143a40ed3dSBarry Smith   PetscFunctionBegin;
31577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3169b94acceSBarry Smith   snes->user		= usrP;
3173a40ed3dSBarry Smith   PetscFunctionReturn(0);
3189b94acceSBarry Smith }
31974679c65SBarry Smith 
3205615d1e5SSatish Balay #undef __FUNC__
321d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext"
3229b94acceSBarry Smith /*@C
3239b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3249b94acceSBarry Smith    nonlinear solvers.
3259b94acceSBarry Smith 
326c7afd0dbSLois Curfman McInnes    Not Collective
327c7afd0dbSLois Curfman McInnes 
3289b94acceSBarry Smith    Input Parameter:
3299b94acceSBarry Smith .  snes - SNES context
3309b94acceSBarry Smith 
3319b94acceSBarry Smith    Output Parameter:
3329b94acceSBarry Smith .  usrP - user context
3339b94acceSBarry Smith 
3349b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3359b94acceSBarry Smith 
3369b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3379b94acceSBarry Smith @*/
3389b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3399b94acceSBarry Smith {
3403a40ed3dSBarry Smith   PetscFunctionBegin;
34177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3429b94acceSBarry Smith   *usrP = snes->user;
3433a40ed3dSBarry Smith   PetscFunctionReturn(0);
3449b94acceSBarry Smith }
34574679c65SBarry Smith 
3465615d1e5SSatish Balay #undef __FUNC__
347d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber"
3489b94acceSBarry Smith /*@
3499b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3509b94acceSBarry Smith    nonlinear solver.
3519b94acceSBarry Smith 
352c7afd0dbSLois Curfman McInnes    Not Collective
353c7afd0dbSLois Curfman McInnes 
3549b94acceSBarry Smith    Input Parameter:
3559b94acceSBarry Smith .  snes - SNES context
3569b94acceSBarry Smith 
3579b94acceSBarry Smith    Output Parameter:
3589b94acceSBarry Smith .  iter - iteration number
3599b94acceSBarry Smith 
3609b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3619b94acceSBarry Smith @*/
3629b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3639b94acceSBarry Smith {
3643a40ed3dSBarry Smith   PetscFunctionBegin;
36577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
36674679c65SBarry Smith   PetscValidIntPointer(iter);
3679b94acceSBarry Smith   *iter = snes->iter;
3683a40ed3dSBarry Smith   PetscFunctionReturn(0);
3699b94acceSBarry Smith }
37074679c65SBarry Smith 
3715615d1e5SSatish Balay #undef __FUNC__
3725615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
3739b94acceSBarry Smith /*@
3749b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3759b94acceSBarry Smith    with SNESSSetFunction().
3769b94acceSBarry Smith 
377c7afd0dbSLois Curfman McInnes    Collective on SNES
378c7afd0dbSLois Curfman McInnes 
3799b94acceSBarry Smith    Input Parameter:
3809b94acceSBarry Smith .  snes - SNES context
3819b94acceSBarry Smith 
3829b94acceSBarry Smith    Output Parameter:
3839b94acceSBarry Smith .  fnorm - 2-norm of function
3849b94acceSBarry Smith 
3859b94acceSBarry Smith    Note:
3869b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
387a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
388a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3899b94acceSBarry Smith 
3909b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
391a86d99e1SLois Curfman McInnes 
392a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3939b94acceSBarry Smith @*/
3949b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3959b94acceSBarry Smith {
3963a40ed3dSBarry Smith   PetscFunctionBegin;
39777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
39874679c65SBarry Smith   PetscValidScalarPointer(fnorm);
39974679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
400d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only");
40174679c65SBarry Smith   }
4029b94acceSBarry Smith   *fnorm = snes->norm;
4033a40ed3dSBarry Smith   PetscFunctionReturn(0);
4049b94acceSBarry Smith }
40574679c65SBarry Smith 
4065615d1e5SSatish Balay #undef __FUNC__
4075615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
4089b94acceSBarry Smith /*@
4099b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
4109b94acceSBarry Smith    with SNESSSetGradient().
4119b94acceSBarry Smith 
412c7afd0dbSLois Curfman McInnes    Collective on SNES
413c7afd0dbSLois Curfman McInnes 
4149b94acceSBarry Smith    Input Parameter:
4159b94acceSBarry Smith .  snes - SNES context
4169b94acceSBarry Smith 
4179b94acceSBarry Smith    Output Parameter:
4189b94acceSBarry Smith .  fnorm - 2-norm of gradient
4199b94acceSBarry Smith 
4209b94acceSBarry Smith    Note:
4219b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
422a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
423a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4249b94acceSBarry Smith 
4259b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
426a86d99e1SLois Curfman McInnes 
427a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4289b94acceSBarry Smith @*/
4299b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
4309b94acceSBarry Smith {
4313a40ed3dSBarry Smith   PetscFunctionBegin;
43277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43374679c65SBarry Smith   PetscValidScalarPointer(gnorm);
43474679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
435d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
43674679c65SBarry Smith   }
4379b94acceSBarry Smith   *gnorm = snes->norm;
4383a40ed3dSBarry Smith   PetscFunctionReturn(0);
4399b94acceSBarry Smith }
44074679c65SBarry Smith 
4415615d1e5SSatish Balay #undef __FUNC__
442d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps"
4439b94acceSBarry Smith /*@
4449b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4459b94acceSBarry Smith    attempted by the nonlinear solver.
4469b94acceSBarry Smith 
447c7afd0dbSLois Curfman McInnes    Not Collective
448c7afd0dbSLois Curfman McInnes 
4499b94acceSBarry Smith    Input Parameter:
4509b94acceSBarry Smith .  snes - SNES context
4519b94acceSBarry Smith 
4529b94acceSBarry Smith    Output Parameter:
4539b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4549b94acceSBarry Smith 
455c96a6f78SLois Curfman McInnes    Notes:
456c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
457c96a6f78SLois Curfman McInnes 
4589b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4599b94acceSBarry Smith @*/
4609b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4619b94acceSBarry Smith {
4623a40ed3dSBarry Smith   PetscFunctionBegin;
46377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
46474679c65SBarry Smith   PetscValidIntPointer(nfails);
4659b94acceSBarry Smith   *nfails = snes->nfailures;
4663a40ed3dSBarry Smith   PetscFunctionReturn(0);
4679b94acceSBarry Smith }
468a847f771SSatish Balay 
4695615d1e5SSatish Balay #undef __FUNC__
470d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations"
471c96a6f78SLois Curfman McInnes /*@
472c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
473c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
474c96a6f78SLois Curfman McInnes 
475c7afd0dbSLois Curfman McInnes    Not Collective
476c7afd0dbSLois Curfman McInnes 
477c96a6f78SLois Curfman McInnes    Input Parameter:
478c96a6f78SLois Curfman McInnes .  snes - SNES context
479c96a6f78SLois Curfman McInnes 
480c96a6f78SLois Curfman McInnes    Output Parameter:
481c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
482c96a6f78SLois Curfman McInnes 
483c96a6f78SLois Curfman McInnes    Notes:
484c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
485c96a6f78SLois Curfman McInnes 
486c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
487c96a6f78SLois Curfman McInnes @*/
48886bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
489c96a6f78SLois Curfman McInnes {
4903a40ed3dSBarry Smith   PetscFunctionBegin;
491c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
492c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
493c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
4943a40ed3dSBarry Smith   PetscFunctionReturn(0);
495c96a6f78SLois Curfman McInnes }
496c96a6f78SLois Curfman McInnes 
497c96a6f78SLois Curfman McInnes #undef __FUNC__
498d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES"
4999b94acceSBarry Smith /*@C
5009b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
5019b94acceSBarry Smith 
502c7afd0dbSLois Curfman McInnes    Not Collective, but if SNES object is parallel, then SLES object is parallel
503c7afd0dbSLois Curfman McInnes 
5049b94acceSBarry Smith    Input Parameter:
5059b94acceSBarry Smith .  snes - the SNES context
5069b94acceSBarry Smith 
5079b94acceSBarry Smith    Output Parameter:
5089b94acceSBarry Smith .  sles - the SLES context
5099b94acceSBarry Smith 
5109b94acceSBarry Smith    Notes:
5119b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
5129b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5139b94acceSBarry Smith    KSP and PC contexts as well.
5149b94acceSBarry Smith 
5159b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
5169b94acceSBarry Smith 
5179b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
5189b94acceSBarry Smith @*/
5199b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
5209b94acceSBarry Smith {
5213a40ed3dSBarry Smith   PetscFunctionBegin;
52277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5239b94acceSBarry Smith   *sles = snes->sles;
5243a40ed3dSBarry Smith   PetscFunctionReturn(0);
5259b94acceSBarry Smith }
52682bf6240SBarry Smith 
527e24b481bSBarry Smith #undef __FUNC__
528e24b481bSBarry Smith #define __FUNC__ "SNESPublish_Petsc"
529e24b481bSBarry Smith static int SNESPublish_Petsc(PetscObject object)
530e24b481bSBarry Smith {
531e24b481bSBarry Smith #if defined(HAVE_AMS)
532e24b481bSBarry Smith   SNES          v = (SNES) object;
533e24b481bSBarry Smith   int          ierr;
534e24b481bSBarry Smith 
535e24b481bSBarry Smith   PetscFunctionBegin;
536e24b481bSBarry Smith 
537e24b481bSBarry Smith   /* if it is already published then return */
538e24b481bSBarry Smith   if (v->amem >=0 ) PetscFunctionReturn(0);
539e24b481bSBarry Smith 
5403f1db9ecSBarry Smith   ierr = PetscObjectPublishBaseBegin(object);CHKERRQ(ierr);
541e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
542e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
543e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
544e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
545e24b481bSBarry Smith   ierr = PetscObjectPublishBaseEnd(object);CHKERRQ(ierr);
546e24b481bSBarry Smith #else
547e24b481bSBarry Smith   PetscFunctionBegin;
548e24b481bSBarry Smith #endif
549e24b481bSBarry Smith   PetscFunctionReturn(0);
550e24b481bSBarry Smith }
551e24b481bSBarry Smith 
5529b94acceSBarry Smith /* -----------------------------------------------------------*/
5535615d1e5SSatish Balay #undef __FUNC__
5545615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
5559b94acceSBarry Smith /*@C
5569b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5579b94acceSBarry Smith 
558c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
559c7afd0dbSLois Curfman McInnes 
560c7afd0dbSLois Curfman McInnes    Input Parameters:
561c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
562c7afd0dbSLois Curfman McInnes -  type - type of method, either
563c7afd0dbSLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations)
564c7afd0dbSLois Curfman McInnes    or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization)
5659b94acceSBarry Smith 
5669b94acceSBarry Smith    Output Parameter:
5679b94acceSBarry Smith .  outsnes - the new SNES context
5689b94acceSBarry Smith 
569c7afd0dbSLois Curfman McInnes    Options Database Keys:
570c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
571c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
572c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
573c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
574c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
575c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
576c1f60f51SBarry Smith 
5779b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5789b94acceSBarry Smith 
57963a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5809b94acceSBarry Smith @*/
5814b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5829b94acceSBarry Smith {
5839b94acceSBarry Smith   int                 ierr;
5849b94acceSBarry Smith   SNES                snes;
5859b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
58637fcc0dbSBarry Smith 
5873a40ed3dSBarry Smith   PetscFunctionBegin;
5889b94acceSBarry Smith   *outsnes = 0;
589d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
590d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
591d64ed03dSBarry Smith   }
5923f1db9ecSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);
5939b94acceSBarry Smith   PLogObjectCreate(snes);
594e24b481bSBarry Smith   snes->bops->publish     = SNESPublish_Petsc;
5959b94acceSBarry Smith   snes->max_its           = 50;
5969750a799SBarry Smith   snes->max_funcs	  = 10000;
5979b94acceSBarry Smith   snes->norm		  = 0.0;
5985a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
599b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
600b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
6019b94acceSBarry Smith     snes->atol		  = 1.e-10;
602d64ed03dSBarry Smith   } else {
603b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
604b4874afaSBarry Smith     snes->ttol            = 0.0;
605b4874afaSBarry Smith     snes->atol		  = 1.e-50;
606b4874afaSBarry Smith   }
6079b94acceSBarry Smith   snes->xtol		  = 1.e-8;
608b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
6099b94acceSBarry Smith   snes->nfuncs            = 0;
6109b94acceSBarry Smith   snes->nfailures         = 0;
6117a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
612639f9d9dSBarry Smith   snes->numbermonitors    = 0;
6139b94acceSBarry Smith   snes->data              = 0;
6149b94acceSBarry Smith   snes->view              = 0;
6159b94acceSBarry Smith   snes->computeumfunction = 0;
6169b94acceSBarry Smith   snes->umfunP            = 0;
6179b94acceSBarry Smith   snes->fc                = 0;
6189b94acceSBarry Smith   snes->deltatol          = 1.e-12;
6199b94acceSBarry Smith   snes->fmin              = -1.e30;
6209b94acceSBarry Smith   snes->method_class      = type;
6219b94acceSBarry Smith   snes->set_method_called = 0;
62282bf6240SBarry Smith   snes->setupcalled      = 0;
6239b94acceSBarry Smith   snes->ksp_ewconv        = 0;
6246f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
6256f24a144SLois Curfman McInnes   snes->vwork             = 0;
6266f24a144SLois Curfman McInnes   snes->nwork             = 0;
627c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
628c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
629c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
6309b94acceSBarry Smith 
6319b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
6320452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
633eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
6349b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6359b94acceSBarry Smith   kctx->version     = 2;
6369b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6379b94acceSBarry Smith                              this was too large for some test cases */
6389b94acceSBarry Smith   kctx->rtol_last   = 0;
6399b94acceSBarry Smith   kctx->rtol_max    = .9;
6409b94acceSBarry Smith   kctx->gamma       = 1.0;
6419b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6429b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6439b94acceSBarry Smith   kctx->threshold   = .1;
6449b94acceSBarry Smith   kctx->lresid_last = 0;
6459b94acceSBarry Smith   kctx->norm_last   = 0;
6469b94acceSBarry Smith 
6479b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
6489b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
6495334005bSBarry Smith 
6509b94acceSBarry Smith   *outsnes = snes;
651e24b481bSBarry Smith   PetscPublishAll(snes);
6523a40ed3dSBarry Smith   PetscFunctionReturn(0);
6539b94acceSBarry Smith }
6549b94acceSBarry Smith 
6559b94acceSBarry Smith /* --------------------------------------------------------------- */
6565615d1e5SSatish Balay #undef __FUNC__
657d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction"
6589b94acceSBarry Smith /*@C
6599b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6609b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6619b94acceSBarry Smith    equations.
6629b94acceSBarry Smith 
663fee21e36SBarry Smith    Collective on SNES
664fee21e36SBarry Smith 
665c7afd0dbSLois Curfman McInnes    Input Parameters:
666c7afd0dbSLois Curfman McInnes +  snes - the SNES context
667c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
668c7afd0dbSLois Curfman McInnes .  r - vector to store function value
669c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
670c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6719b94acceSBarry Smith 
672c7afd0dbSLois Curfman McInnes    Calling sequence of func:
6738d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
674c7afd0dbSLois Curfman McInnes 
675313e4042SLois Curfman McInnes .  f - function vector
676c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
6779b94acceSBarry Smith 
6789b94acceSBarry Smith    Notes:
6799b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6809b94acceSBarry Smith $      f'(x) x = -f(x),
681c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
6829b94acceSBarry Smith 
6839b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6849b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6859b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6869b94acceSBarry Smith 
6879b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6889b94acceSBarry Smith 
689a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6909b94acceSBarry Smith @*/
6915334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6929b94acceSBarry Smith {
6933a40ed3dSBarry Smith   PetscFunctionBegin;
69477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
695a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
696a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
697a8c6a408SBarry Smith   }
6989b94acceSBarry Smith   snes->computefunction     = func;
6999b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7009b94acceSBarry Smith   snes->funP                = ctx;
7013a40ed3dSBarry Smith   PetscFunctionReturn(0);
7029b94acceSBarry Smith }
7039b94acceSBarry Smith 
7045615d1e5SSatish Balay #undef __FUNC__
7055615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
7069b94acceSBarry Smith /*@
7079b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
7089b94acceSBarry Smith    SNESSetFunction().
7099b94acceSBarry Smith 
710c7afd0dbSLois Curfman McInnes    Collective on SNES
711c7afd0dbSLois Curfman McInnes 
7129b94acceSBarry Smith    Input Parameters:
713c7afd0dbSLois Curfman McInnes +  snes - the SNES context
714c7afd0dbSLois Curfman McInnes -  x - input vector
7159b94acceSBarry Smith 
7169b94acceSBarry Smith    Output Parameter:
7173638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
7189b94acceSBarry Smith 
7191bffabb2SLois Curfman McInnes    Notes:
7209b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
7219b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
7229b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
7239b94acceSBarry Smith 
7249b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7259b94acceSBarry Smith 
726a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7279b94acceSBarry Smith @*/
7289b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
7299b94acceSBarry Smith {
7309b94acceSBarry Smith   int    ierr;
7319b94acceSBarry Smith 
7323a40ed3dSBarry Smith   PetscFunctionBegin;
733d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
734a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
735d4bb536fSBarry Smith   }
7369b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
737d64ed03dSBarry Smith   PetscStackPush("SNES user function");
7389b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
739d64ed03dSBarry Smith   PetscStackPop;
740ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7419b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
7423a40ed3dSBarry Smith   PetscFunctionReturn(0);
7439b94acceSBarry Smith }
7449b94acceSBarry Smith 
7455615d1e5SSatish Balay #undef __FUNC__
746d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction"
7479b94acceSBarry Smith /*@C
7489b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
7499b94acceSBarry Smith    unconstrained minimization.
7509b94acceSBarry Smith 
751fee21e36SBarry Smith    Collective on SNES
752fee21e36SBarry Smith 
753c7afd0dbSLois Curfman McInnes    Input Parameters:
754c7afd0dbSLois Curfman McInnes +  snes - the SNES context
755c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
756c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
757c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7589b94acceSBarry Smith 
759c7afd0dbSLois Curfman McInnes    Calling sequence of func:
7608d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,double *f,void *ctx);
761c7afd0dbSLois Curfman McInnes 
762c7afd0dbSLois Curfman McInnes +  x - input vector
7639b94acceSBarry Smith .  f - function
764c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined function context
7659b94acceSBarry Smith 
7669b94acceSBarry Smith    Notes:
7679b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7689b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7699b94acceSBarry Smith    SNESSetFunction().
7709b94acceSBarry Smith 
7719b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
7729b94acceSBarry Smith 
773a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
774a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
7759b94acceSBarry Smith @*/
7769b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
7779b94acceSBarry Smith                       void *ctx)
7789b94acceSBarry Smith {
7793a40ed3dSBarry Smith   PetscFunctionBegin;
78077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
781a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
782a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
783a8c6a408SBarry Smith   }
7849b94acceSBarry Smith   snes->computeumfunction   = func;
7859b94acceSBarry Smith   snes->umfunP              = ctx;
7863a40ed3dSBarry Smith   PetscFunctionReturn(0);
7879b94acceSBarry Smith }
7889b94acceSBarry Smith 
7895615d1e5SSatish Balay #undef __FUNC__
7905615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7919b94acceSBarry Smith /*@
7929b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7939b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7949b94acceSBarry Smith 
795c7afd0dbSLois Curfman McInnes    Collective on SNES
796c7afd0dbSLois Curfman McInnes 
7979b94acceSBarry Smith    Input Parameters:
798c7afd0dbSLois Curfman McInnes +  snes - the SNES context
799c7afd0dbSLois Curfman McInnes -  x - input vector
8009b94acceSBarry Smith 
8019b94acceSBarry Smith    Output Parameter:
8029b94acceSBarry Smith .  y - function value
8039b94acceSBarry Smith 
8049b94acceSBarry Smith    Notes:
8059b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
8069b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8079b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
808a86d99e1SLois Curfman McInnes 
809a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
810a86d99e1SLois Curfman McInnes 
811a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
812a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
8139b94acceSBarry Smith @*/
8149b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
8159b94acceSBarry Smith {
8169b94acceSBarry Smith   int    ierr;
8173a40ed3dSBarry Smith 
8183a40ed3dSBarry Smith   PetscFunctionBegin;
819a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
820a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
821a8c6a408SBarry Smith   }
8229b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
823d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
8249b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
825d64ed03dSBarry Smith   PetscStackPop;
826ae3c334cSLois Curfman McInnes   snes->nfuncs++;
8279b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
8283a40ed3dSBarry Smith   PetscFunctionReturn(0);
8299b94acceSBarry Smith }
8309b94acceSBarry Smith 
8315615d1e5SSatish Balay #undef __FUNC__
832d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient"
8339b94acceSBarry Smith /*@C
8349b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
8359b94acceSBarry Smith    vector for use by the SNES routines.
8369b94acceSBarry Smith 
837c7afd0dbSLois Curfman McInnes    Collective on SNES
838c7afd0dbSLois Curfman McInnes 
8399b94acceSBarry Smith    Input Parameters:
840c7afd0dbSLois Curfman McInnes +  snes - the SNES context
8419b94acceSBarry Smith .  func - function evaluation routine
842044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
843044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
844c7afd0dbSLois Curfman McInnes -  r - vector to store gradient value
845fee21e36SBarry Smith 
8469b94acceSBarry Smith    Calling sequence of func:
8478d76a1e5SLois Curfman McInnes $     func (SNES, Vec x, Vec g, void *ctx);
8489b94acceSBarry Smith 
849c7afd0dbSLois Curfman McInnes +  x - input vector
8509b94acceSBarry Smith .  g - gradient vector
851c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined gradient context
8529b94acceSBarry Smith 
8539b94acceSBarry Smith    Notes:
8549b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
8559b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
8569b94acceSBarry Smith    SNESSetFunction().
8579b94acceSBarry Smith 
8589b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
8599b94acceSBarry Smith 
860a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
861a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
8629b94acceSBarry Smith @*/
86374679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
8649b94acceSBarry Smith {
8653a40ed3dSBarry Smith   PetscFunctionBegin;
86677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
867a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
868a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
869a8c6a408SBarry Smith   }
8709b94acceSBarry Smith   snes->computefunction     = func;
8719b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
8729b94acceSBarry Smith   snes->funP                = ctx;
8733a40ed3dSBarry Smith   PetscFunctionReturn(0);
8749b94acceSBarry Smith }
8759b94acceSBarry Smith 
8765615d1e5SSatish Balay #undef __FUNC__
8775615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
8789b94acceSBarry Smith /*@
879a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
880a86d99e1SLois Curfman McInnes    SNESSetGradient().
8819b94acceSBarry Smith 
882c7afd0dbSLois Curfman McInnes    Collective on SNES
883c7afd0dbSLois Curfman McInnes 
8849b94acceSBarry Smith    Input Parameters:
885c7afd0dbSLois Curfman McInnes +  snes - the SNES context
886c7afd0dbSLois Curfman McInnes -  x - input vector
8879b94acceSBarry Smith 
8889b94acceSBarry Smith    Output Parameter:
8899b94acceSBarry Smith .  y - gradient vector
8909b94acceSBarry Smith 
8919b94acceSBarry Smith    Notes:
8929b94acceSBarry Smith    SNESComputeGradient() is valid only for
8939b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8949b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
895a86d99e1SLois Curfman McInnes 
896a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
897a86d99e1SLois Curfman McInnes 
898a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
899a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
9009b94acceSBarry Smith @*/
9019b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
9029b94acceSBarry Smith {
9039b94acceSBarry Smith   int    ierr;
9043a40ed3dSBarry Smith 
9053a40ed3dSBarry Smith   PetscFunctionBegin;
9063a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
907a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9083a40ed3dSBarry Smith   }
9099b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
910d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
9119b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
912d64ed03dSBarry Smith   PetscStackPop;
9139b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
9143a40ed3dSBarry Smith   PetscFunctionReturn(0);
9159b94acceSBarry Smith }
9169b94acceSBarry Smith 
9175615d1e5SSatish Balay #undef __FUNC__
9185615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
91962fef451SLois Curfman McInnes /*@
92062fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
92162fef451SLois Curfman McInnes    set with SNESSetJacobian().
92262fef451SLois Curfman McInnes 
923c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
924c7afd0dbSLois Curfman McInnes 
92562fef451SLois Curfman McInnes    Input Parameters:
926c7afd0dbSLois Curfman McInnes +  snes - the SNES context
927c7afd0dbSLois Curfman McInnes -  x - input vector
92862fef451SLois Curfman McInnes 
92962fef451SLois Curfman McInnes    Output Parameters:
930c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
93162fef451SLois Curfman McInnes .  B - optional preconditioning matrix
932c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
933fee21e36SBarry Smith 
93462fef451SLois Curfman McInnes    Notes:
93562fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
93662fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
93762fef451SLois Curfman McInnes 
938dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
939dc5a77f8SLois Curfman McInnes    flag parameter.
94062fef451SLois Curfman McInnes 
94162fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
94262fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
943005c665bSBarry Smith    methods is SNESComputeHessian().
94462fef451SLois Curfman McInnes 
94562fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
94662fef451SLois Curfman McInnes 
94762fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
94862fef451SLois Curfman McInnes @*/
9499b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
9509b94acceSBarry Smith {
9519b94acceSBarry Smith   int    ierr;
9523a40ed3dSBarry Smith 
9533a40ed3dSBarry Smith   PetscFunctionBegin;
9543a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
955a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
9563a40ed3dSBarry Smith   }
9573a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
9589b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
959c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
960d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
9619b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
962d64ed03dSBarry Smith   PetscStackPop;
9639b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
9646d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
96577c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
96677c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9673a40ed3dSBarry Smith   PetscFunctionReturn(0);
9689b94acceSBarry Smith }
9699b94acceSBarry Smith 
9705615d1e5SSatish Balay #undef __FUNC__
9715615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
97262fef451SLois Curfman McInnes /*@
97362fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
97462fef451SLois Curfman McInnes    set with SNESSetHessian().
97562fef451SLois Curfman McInnes 
976c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
977c7afd0dbSLois Curfman McInnes 
97862fef451SLois Curfman McInnes    Input Parameters:
979c7afd0dbSLois Curfman McInnes +  snes - the SNES context
980c7afd0dbSLois Curfman McInnes -  x - input vector
98162fef451SLois Curfman McInnes 
98262fef451SLois Curfman McInnes    Output Parameters:
983c7afd0dbSLois Curfman McInnes +  A - Hessian matrix
98462fef451SLois Curfman McInnes .  B - optional preconditioning matrix
985c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
986fee21e36SBarry Smith 
98762fef451SLois Curfman McInnes    Notes:
98862fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
98962fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
99062fef451SLois Curfman McInnes 
991dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
992dc5a77f8SLois Curfman McInnes    flag parameter.
99362fef451SLois Curfman McInnes 
99462fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
99562fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
99662fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
99762fef451SLois Curfman McInnes 
99862fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
99962fef451SLois Curfman McInnes 
1000a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
1001a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
100262fef451SLois Curfman McInnes @*/
100362fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
10049b94acceSBarry Smith {
10059b94acceSBarry Smith   int    ierr;
10063a40ed3dSBarry Smith 
10073a40ed3dSBarry Smith   PetscFunctionBegin;
10083a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1009a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10103a40ed3dSBarry Smith   }
10113a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
101262fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
10130f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
1014d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
101562fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
1016d64ed03dSBarry Smith   PetscStackPop;
101762fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
10180f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
101977c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
102077c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
10213a40ed3dSBarry Smith   PetscFunctionReturn(0);
10229b94acceSBarry Smith }
10239b94acceSBarry Smith 
10245615d1e5SSatish Balay #undef __FUNC__
1025d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian"
10269b94acceSBarry Smith /*@C
10279b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1028044dda88SLois Curfman McInnes    location to store the matrix.
10299b94acceSBarry Smith 
1030c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1031c7afd0dbSLois Curfman McInnes 
10329b94acceSBarry Smith    Input Parameters:
1033c7afd0dbSLois Curfman McInnes +  snes - the SNES context
10349b94acceSBarry Smith .  A - Jacobian matrix
10359b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
10369b94acceSBarry Smith .  func - Jacobian evaluation routine
1037c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
10382cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
10399b94acceSBarry Smith 
10409b94acceSBarry Smith    Calling sequence of func:
10418d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10429b94acceSBarry Smith 
1043c7afd0dbSLois Curfman McInnes +  x - input vector
10449b94acceSBarry Smith .  A - Jacobian matrix
10459b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1046ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1047ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
1048c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
10499b94acceSBarry Smith 
10509b94acceSBarry Smith    Notes:
1051dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10522cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1053ac21db08SLois Curfman McInnes 
1054ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
10559b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
10569b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10579b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10589b94acceSBarry Smith    throughout the global iterations.
10599b94acceSBarry Smith 
10609b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
10619b94acceSBarry Smith 
1062ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
10639b94acceSBarry Smith @*/
10649b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10659b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10669b94acceSBarry Smith {
10673a40ed3dSBarry Smith   PetscFunctionBegin;
106877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1069a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1070a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1071a8c6a408SBarry Smith   }
10729b94acceSBarry Smith   snes->computejacobian = func;
10739b94acceSBarry Smith   snes->jacP            = ctx;
10749b94acceSBarry Smith   snes->jacobian        = A;
10759b94acceSBarry Smith   snes->jacobian_pre    = B;
10763a40ed3dSBarry Smith   PetscFunctionReturn(0);
10779b94acceSBarry Smith }
107862fef451SLois Curfman McInnes 
10795615d1e5SSatish Balay #undef __FUNC__
1080d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian"
1081b4fd4287SBarry Smith /*@
1082b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1083b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1084b4fd4287SBarry Smith 
1085c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
1086c7afd0dbSLois Curfman McInnes 
1087b4fd4287SBarry Smith    Input Parameter:
1088b4fd4287SBarry Smith .  snes - the nonlinear solver context
1089b4fd4287SBarry Smith 
1090b4fd4287SBarry Smith    Output Parameters:
1091c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
1092b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
1093c7afd0dbSLois Curfman McInnes -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1094fee21e36SBarry Smith 
1095b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1096b4fd4287SBarry Smith @*/
1097b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
1098b4fd4287SBarry Smith {
10993a40ed3dSBarry Smith   PetscFunctionBegin;
11003a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1101a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
11023a40ed3dSBarry Smith   }
1103b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
1104b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
1105b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
11063a40ed3dSBarry Smith   PetscFunctionReturn(0);
1107b4fd4287SBarry Smith }
1108b4fd4287SBarry Smith 
11095615d1e5SSatish Balay #undef __FUNC__
1110d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian"
11119b94acceSBarry Smith /*@C
11129b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1113044dda88SLois Curfman McInnes    location to store the matrix.
11149b94acceSBarry Smith 
1115c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1116c7afd0dbSLois Curfman McInnes 
11179b94acceSBarry Smith    Input Parameters:
1118c7afd0dbSLois Curfman McInnes +  snes - the SNES context
11199b94acceSBarry Smith .  A - Hessian matrix
11209b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
11219b94acceSBarry Smith .  func - Jacobian evaluation routine
1122c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
1123313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
11249b94acceSBarry Smith 
11259b94acceSBarry Smith    Calling sequence of func:
11268d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
11279b94acceSBarry Smith 
1128c7afd0dbSLois Curfman McInnes +  x - input vector
11299b94acceSBarry Smith .  A - Hessian matrix
11309b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1131ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1132ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
1133c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Hessian context
11349b94acceSBarry Smith 
11359b94acceSBarry Smith    Notes:
1136dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
11372cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1138ac21db08SLois Curfman McInnes 
11399b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
11409b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
11419b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
11429b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
11439b94acceSBarry Smith    throughout the global iterations.
11449b94acceSBarry Smith 
11459b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
11469b94acceSBarry Smith 
1147ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
11489b94acceSBarry Smith @*/
11499b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
11509b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
11519b94acceSBarry Smith {
11523a40ed3dSBarry Smith   PetscFunctionBegin;
115377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1154d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1155a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1156d4bb536fSBarry Smith   }
11579b94acceSBarry Smith   snes->computejacobian = func;
11589b94acceSBarry Smith   snes->jacP            = ctx;
11599b94acceSBarry Smith   snes->jacobian        = A;
11609b94acceSBarry Smith   snes->jacobian_pre    = B;
11613a40ed3dSBarry Smith   PetscFunctionReturn(0);
11629b94acceSBarry Smith }
11639b94acceSBarry Smith 
11645615d1e5SSatish Balay #undef __FUNC__
1165d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian"
116662fef451SLois Curfman McInnes /*@
116762fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
116862fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
116962fef451SLois Curfman McInnes 
1170c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object is parallel if SNES object is parallel
1171c7afd0dbSLois Curfman McInnes 
117262fef451SLois Curfman McInnes    Input Parameter:
117362fef451SLois Curfman McInnes .  snes - the nonlinear solver context
117462fef451SLois Curfman McInnes 
117562fef451SLois Curfman McInnes    Output Parameters:
1176c7afd0dbSLois Curfman McInnes +  A - location to stash Hessian matrix (or PETSC_NULL)
117762fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
1178c7afd0dbSLois Curfman McInnes -  ctx - location to stash Hessian ctx (or PETSC_NULL)
1179fee21e36SBarry Smith 
118062fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
1181c7afd0dbSLois Curfman McInnes 
1182c7afd0dbSLois Curfman McInnes .keywords: SNES, get, Hessian
118362fef451SLois Curfman McInnes @*/
118462fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
118562fef451SLois Curfman McInnes {
11863a40ed3dSBarry Smith   PetscFunctionBegin;
11873a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1188a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
11893a40ed3dSBarry Smith   }
119062fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
119162fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
119262fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
11933a40ed3dSBarry Smith   PetscFunctionReturn(0);
119462fef451SLois Curfman McInnes }
119562fef451SLois Curfman McInnes 
11969b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
11979b94acceSBarry Smith 
11985615d1e5SSatish Balay #undef __FUNC__
11995615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
12009b94acceSBarry Smith /*@
12019b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1202272ac6f2SLois Curfman McInnes    of a nonlinear solver.
12039b94acceSBarry Smith 
1204fee21e36SBarry Smith    Collective on SNES
1205fee21e36SBarry Smith 
1206c7afd0dbSLois Curfman McInnes    Input Parameters:
1207c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1208c7afd0dbSLois Curfman McInnes -  x - the solution vector
1209c7afd0dbSLois Curfman McInnes 
1210272ac6f2SLois Curfman McInnes    Notes:
1211272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1212272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1213272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1214272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1215272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1216272ac6f2SLois Curfman McInnes 
12179b94acceSBarry Smith .keywords: SNES, nonlinear, setup
12189b94acceSBarry Smith 
12199b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
12209b94acceSBarry Smith @*/
12218ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
12229b94acceSBarry Smith {
1223272ac6f2SLois Curfman McInnes   int ierr, flg;
12243a40ed3dSBarry Smith 
12253a40ed3dSBarry Smith   PetscFunctionBegin;
122677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
122777c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
12288ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
12299b94acceSBarry Smith 
1230c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1231c1f60f51SBarry Smith   /*
1232c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1233dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1234c1f60f51SBarry Smith   */
1235c1f60f51SBarry Smith   if (flg) {
1236c1f60f51SBarry Smith     Mat J;
1237e24b481bSBarry Smith     ierr = MatCreateSNESFDMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1238c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1239c1f60f51SBarry Smith     snes->mfshell = J;
1240c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1241c1f60f51SBarry Smith       snes->jacobian = J;
1242c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1243d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1244c1f60f51SBarry Smith       snes->jacobian = J;
1245c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1246d4bb536fSBarry Smith     } else {
1247a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option");
1248d4bb536fSBarry Smith     }
1249e24b481bSBarry Smith     ierr = MatSNESFDMFSetFromOptions(J);CHKERRQ(ierr);
1250c1f60f51SBarry Smith   }
1251272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1252c1f60f51SBarry Smith   /*
1253dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1254c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1255c1f60f51SBarry Smith    */
125631615d04SLois Curfman McInnes   if (flg) {
1257272ac6f2SLois Curfman McInnes     Mat J;
1258e24b481bSBarry Smith     ierr = MatCreateSNESFDMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1259272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1260272ac6f2SLois Curfman McInnes     snes->mfshell = J;
126183e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
126283e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
126394a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1264d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
126583e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
126694a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1267d4bb536fSBarry Smith     } else {
1268a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option");
1269d4bb536fSBarry Smith     }
1270e24b481bSBarry Smith     ierr = MatSNESFDMFSetFromOptions(J);CHKERRQ(ierr);
1271272ac6f2SLois Curfman McInnes   }
12729b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1273a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1274a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1275a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first");
1276a8c6a408SBarry Smith     if (snes->vec_func == snes->vec_sol) {
1277a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector");
1278a8c6a408SBarry Smith     }
1279a703fe33SLois Curfman McInnes 
1280a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
128182bf6240SBarry Smith     if (snes->ksp_ewconv && PetscStrcmp(snes->type_name,SNES_EQ_TR)) {
1282a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1283a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1284a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1285a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1286a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1287a703fe33SLois Curfman McInnes     }
12889b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1289a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1290a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1291a8c6a408SBarry Smith     if (!snes->computeumfunction) {
1292a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first");
1293a8c6a408SBarry Smith     }
1294a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first");
1295d4bb536fSBarry Smith   } else {
1296a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
1297d4bb536fSBarry Smith   }
1298a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
129982bf6240SBarry Smith   snes->setupcalled = 1;
13003a40ed3dSBarry Smith   PetscFunctionReturn(0);
13019b94acceSBarry Smith }
13029b94acceSBarry Smith 
13035615d1e5SSatish Balay #undef __FUNC__
1304d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy"
13059b94acceSBarry Smith /*@C
13069b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
13079b94acceSBarry Smith    with SNESCreate().
13089b94acceSBarry Smith 
1309c7afd0dbSLois Curfman McInnes    Collective on SNES
1310c7afd0dbSLois Curfman McInnes 
13119b94acceSBarry Smith    Input Parameter:
13129b94acceSBarry Smith .  snes - the SNES context
13139b94acceSBarry Smith 
13149b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
13159b94acceSBarry Smith 
131663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
13179b94acceSBarry Smith @*/
13189b94acceSBarry Smith int SNESDestroy(SNES snes)
13199b94acceSBarry Smith {
13209b94acceSBarry Smith   int ierr;
13213a40ed3dSBarry Smith 
13223a40ed3dSBarry Smith   PetscFunctionBegin;
132377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13243a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1325d4bb536fSBarry Smith 
1326e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes); CHKERRQ(ierr);}
13270452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
1328522c5e43SBarry Smith   if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);}
13299b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
1330522c5e43SBarry Smith   if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);}
1331522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
13329b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
13330452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
13343a40ed3dSBarry Smith   PetscFunctionReturn(0);
13359b94acceSBarry Smith }
13369b94acceSBarry Smith 
13379b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
13389b94acceSBarry Smith 
13395615d1e5SSatish Balay #undef __FUNC__
13405615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
13419b94acceSBarry Smith /*@
1342d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
13439b94acceSBarry Smith 
1344c7afd0dbSLois Curfman McInnes    Collective on SNES
1345c7afd0dbSLois Curfman McInnes 
13469b94acceSBarry Smith    Input Parameters:
1347c7afd0dbSLois Curfman McInnes +  snes - the SNES context
134833174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
134933174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
135033174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
135133174efeSLois Curfman McInnes            of the change in the solution between steps
135233174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1353c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1354fee21e36SBarry Smith 
135533174efeSLois Curfman McInnes    Options Database Keys:
1356c7afd0dbSLois Curfman McInnes +    -snes_atol <atol> - Sets atol
1357c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1358c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1359c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1360c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
13619b94acceSBarry Smith 
1362d7a720efSLois Curfman McInnes    Notes:
13639b94acceSBarry Smith    The default maximum number of iterations is 50.
13649b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
13659b94acceSBarry Smith 
136633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
13679b94acceSBarry Smith 
1368d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
13699b94acceSBarry Smith @*/
1370d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
13719b94acceSBarry Smith {
13723a40ed3dSBarry Smith   PetscFunctionBegin;
137377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1374d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1375d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1376d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1377d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1378d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
13793a40ed3dSBarry Smith   PetscFunctionReturn(0);
13809b94acceSBarry Smith }
13819b94acceSBarry Smith 
13825615d1e5SSatish Balay #undef __FUNC__
13835615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
13849b94acceSBarry Smith /*@
138533174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
138633174efeSLois Curfman McInnes 
1387c7afd0dbSLois Curfman McInnes    Not Collective
1388c7afd0dbSLois Curfman McInnes 
138933174efeSLois Curfman McInnes    Input Parameters:
1390c7afd0dbSLois Curfman McInnes +  snes - the SNES context
139133174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
139233174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
139333174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
139433174efeSLois Curfman McInnes            of the change in the solution between steps
139533174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1396c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1397fee21e36SBarry Smith 
139833174efeSLois Curfman McInnes    Notes:
139933174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
140033174efeSLois Curfman McInnes 
140133174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
140233174efeSLois Curfman McInnes 
140333174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
140433174efeSLois Curfman McInnes @*/
140533174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
140633174efeSLois Curfman McInnes {
14073a40ed3dSBarry Smith   PetscFunctionBegin;
140833174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
140933174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
141033174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
141133174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
141233174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
141333174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
14143a40ed3dSBarry Smith   PetscFunctionReturn(0);
141533174efeSLois Curfman McInnes }
141633174efeSLois Curfman McInnes 
14175615d1e5SSatish Balay #undef __FUNC__
14185615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
141933174efeSLois Curfman McInnes /*@
14209b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
14219b94acceSBarry Smith 
1422fee21e36SBarry Smith    Collective on SNES
1423fee21e36SBarry Smith 
1424c7afd0dbSLois Curfman McInnes    Input Parameters:
1425c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1426c7afd0dbSLois Curfman McInnes -  tol - tolerance
1427c7afd0dbSLois Curfman McInnes 
14289b94acceSBarry Smith    Options Database Key:
1429c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
14309b94acceSBarry Smith 
14319b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
14329b94acceSBarry Smith 
1433d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
14349b94acceSBarry Smith @*/
14359b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
14369b94acceSBarry Smith {
14373a40ed3dSBarry Smith   PetscFunctionBegin;
143877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
14399b94acceSBarry Smith   snes->deltatol = tol;
14403a40ed3dSBarry Smith   PetscFunctionReturn(0);
14419b94acceSBarry Smith }
14429b94acceSBarry Smith 
14435615d1e5SSatish Balay #undef __FUNC__
14445615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
14459b94acceSBarry Smith /*@
144677c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
14479b94acceSBarry Smith    for unconstrained minimization solvers.
14489b94acceSBarry Smith 
1449fee21e36SBarry Smith    Collective on SNES
1450fee21e36SBarry Smith 
1451c7afd0dbSLois Curfman McInnes    Input Parameters:
1452c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1453c7afd0dbSLois Curfman McInnes -  ftol - minimum function tolerance
1454c7afd0dbSLois Curfman McInnes 
14559b94acceSBarry Smith    Options Database Key:
1456c7afd0dbSLois Curfman McInnes .  -snes_fmin <ftol> - Sets ftol
14579b94acceSBarry Smith 
14589b94acceSBarry Smith    Note:
145977c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
14609b94acceSBarry Smith    methods only.
14619b94acceSBarry Smith 
14629b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
14639b94acceSBarry Smith 
1464d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
14659b94acceSBarry Smith @*/
146677c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
14679b94acceSBarry Smith {
14683a40ed3dSBarry Smith   PetscFunctionBegin;
146977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
14709b94acceSBarry Smith   snes->fmin = ftol;
14713a40ed3dSBarry Smith   PetscFunctionReturn(0);
14729b94acceSBarry Smith }
14739b94acceSBarry Smith 
14749b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
14759b94acceSBarry Smith 
14765615d1e5SSatish Balay #undef __FUNC__
1477d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor"
14789b94acceSBarry Smith /*@C
14795cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
14809b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
14819b94acceSBarry Smith    progress.
14829b94acceSBarry Smith 
1483fee21e36SBarry Smith    Collective on SNES
1484fee21e36SBarry Smith 
1485c7afd0dbSLois Curfman McInnes    Input Parameters:
1486c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1487c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1488c7afd0dbSLois Curfman McInnes -  mctx - [optional] user-defined context for private data for the
1489c7afd0dbSLois Curfman McInnes           monitor routine (may be PETSC_NULL)
14909b94acceSBarry Smith 
1491c7afd0dbSLois Curfman McInnes    Calling sequence of func:
14928d76a1e5SLois Curfman McInnes $     int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
1493c7afd0dbSLois Curfman McInnes 
1494c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1495c7afd0dbSLois Curfman McInnes .    its - iteration number
1496c7afd0dbSLois Curfman McInnes .    mctx - [optional] monitoring context
1497c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
1498c7afd0dbSLois Curfman McInnes             for SNES_NONLINEAR_EQUATIONS methods
1499c7afd0dbSLois Curfman McInnes -    norm - 2-norm gradient value (may be estimated)
1500c7afd0dbSLois Curfman McInnes             for SNES_UNCONSTRAINED_MINIMIZATION methods
15019b94acceSBarry Smith 
15029665c990SLois Curfman McInnes    Options Database Keys:
1503c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1504c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1505c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1506c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1507c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1508c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1509c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1510c7afd0dbSLois Curfman McInnes                             the options database.
15119665c990SLois Curfman McInnes 
1512639f9d9dSBarry Smith    Notes:
15136bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
15146bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
15156bc08f3fSLois Curfman McInnes    order in which they were set.
1516639f9d9dSBarry Smith 
15179b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
15189b94acceSBarry Smith 
15195cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
15209b94acceSBarry Smith @*/
152174679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
15229b94acceSBarry Smith {
15233a40ed3dSBarry Smith   PetscFunctionBegin;
1524639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1525a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
1526639f9d9dSBarry Smith   }
1527639f9d9dSBarry Smith 
1528639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1529639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
15303a40ed3dSBarry Smith   PetscFunctionReturn(0);
15319b94acceSBarry Smith }
15329b94acceSBarry Smith 
15335615d1e5SSatish Balay #undef __FUNC__
15345cd90555SBarry Smith #define __FUNC__ "SNESClearMonitor"
15355cd90555SBarry Smith /*@C
15365cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
15375cd90555SBarry Smith 
1538c7afd0dbSLois Curfman McInnes    Collective on SNES
1539c7afd0dbSLois Curfman McInnes 
15405cd90555SBarry Smith    Input Parameters:
15415cd90555SBarry Smith .  snes - the SNES context
15425cd90555SBarry Smith 
15435cd90555SBarry Smith    Options Database:
1544c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1545c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1546c7afd0dbSLois Curfman McInnes     set via the options database
15475cd90555SBarry Smith 
15485cd90555SBarry Smith    Notes:
15495cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
15505cd90555SBarry Smith 
15515cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
15525cd90555SBarry Smith 
15535cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
15545cd90555SBarry Smith @*/
15555cd90555SBarry Smith int SNESClearMonitor( SNES snes )
15565cd90555SBarry Smith {
15575cd90555SBarry Smith   PetscFunctionBegin;
15585cd90555SBarry Smith   snes->numbermonitors = 0;
15595cd90555SBarry Smith   PetscFunctionReturn(0);
15605cd90555SBarry Smith }
15615cd90555SBarry Smith 
15625cd90555SBarry Smith #undef __FUNC__
1563d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest"
15649b94acceSBarry Smith /*@C
15659b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
15669b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
15679b94acceSBarry Smith 
1568fee21e36SBarry Smith    Collective on SNES
1569fee21e36SBarry Smith 
1570c7afd0dbSLois Curfman McInnes    Input Parameters:
1571c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1572c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1573c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1574c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
15759b94acceSBarry Smith 
1576c7afd0dbSLois Curfman McInnes    Calling sequence of func:
15778d76a1e5SLois Curfman McInnes $     int func (SNES snes,double xnorm,double gnorm,double f,void *cctx)
1578c7afd0dbSLois Curfman McInnes 
1579c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1580c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1581c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
1582c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods)
1583c7afd0dbSLois Curfman McInnes .    f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods)
1584c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods)
1585c7afd0dbSLois Curfman McInnes -    f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods)
15869b94acceSBarry Smith 
15879b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
15889b94acceSBarry Smith 
158940191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
159040191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
15919b94acceSBarry Smith @*/
159274679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
15939b94acceSBarry Smith {
15943a40ed3dSBarry Smith   PetscFunctionBegin;
15959b94acceSBarry Smith   (snes)->converged = func;
15969b94acceSBarry Smith   (snes)->cnvP      = cctx;
15973a40ed3dSBarry Smith   PetscFunctionReturn(0);
15989b94acceSBarry Smith }
15999b94acceSBarry Smith 
16005615d1e5SSatish Balay #undef __FUNC__
1601d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory"
1602c9005455SLois Curfman McInnes /*@
1603c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1604c9005455SLois Curfman McInnes 
1605fee21e36SBarry Smith    Collective on SNES
1606fee21e36SBarry Smith 
1607c7afd0dbSLois Curfman McInnes    Input Parameters:
1608c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1609c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1610c7afd0dbSLois Curfman McInnes -  na  - size of a
1611c7afd0dbSLois Curfman McInnes 
1612c9005455SLois Curfman McInnes    Notes:
1613c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1614c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1615c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1616c9005455SLois Curfman McInnes    at each step.
1617c9005455SLois Curfman McInnes 
1618c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1619c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1620c9005455SLois Curfman McInnes    during the section of code that is being timed.
1621c9005455SLois Curfman McInnes 
1622c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1623c9005455SLois Curfman McInnes @*/
1624c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1625c9005455SLois Curfman McInnes {
16263a40ed3dSBarry Smith   PetscFunctionBegin;
1627c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1628c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1629c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1630c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
16313a40ed3dSBarry Smith   PetscFunctionReturn(0);
1632c9005455SLois Curfman McInnes }
1633c9005455SLois Curfman McInnes 
1634c9005455SLois Curfman McInnes #undef __FUNC__
16355615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
16369b94acceSBarry Smith /*
16379b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
16389b94acceSBarry Smith    positive parameter delta.
16399b94acceSBarry Smith 
16409b94acceSBarry Smith     Input Parameters:
1641c7afd0dbSLois Curfman McInnes +   snes - the SNES context
16429b94acceSBarry Smith .   y - approximate solution of linear system
16439b94acceSBarry Smith .   fnorm - 2-norm of current function
1644c7afd0dbSLois Curfman McInnes -   delta - trust region size
16459b94acceSBarry Smith 
16469b94acceSBarry Smith     Output Parameters:
1647c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
16489b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
16499b94acceSBarry Smith     region, and exceeds zero otherwise.
1650c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
16519b94acceSBarry Smith 
16529b94acceSBarry Smith     Note:
165340191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
16549b94acceSBarry Smith     is set to be the maximum allowable step size.
16559b94acceSBarry Smith 
16569b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
16579b94acceSBarry Smith */
16589b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
16599b94acceSBarry Smith                           double *gpnorm,double *ynorm)
16609b94acceSBarry Smith {
16619b94acceSBarry Smith   double norm;
16629b94acceSBarry Smith   Scalar cnorm;
16633a40ed3dSBarry Smith   int    ierr;
16643a40ed3dSBarry Smith 
16653a40ed3dSBarry Smith   PetscFunctionBegin;
16663a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
16679b94acceSBarry Smith   if (norm > *delta) {
16689b94acceSBarry Smith      norm = *delta/norm;
16699b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
16709b94acceSBarry Smith      cnorm = norm;
16719b94acceSBarry Smith      VecScale( &cnorm, y );
16729b94acceSBarry Smith      *ynorm = *delta;
16739b94acceSBarry Smith   } else {
16749b94acceSBarry Smith      *gpnorm = 0.0;
16759b94acceSBarry Smith      *ynorm = norm;
16769b94acceSBarry Smith   }
16773a40ed3dSBarry Smith   PetscFunctionReturn(0);
16789b94acceSBarry Smith }
16799b94acceSBarry Smith 
16805615d1e5SSatish Balay #undef __FUNC__
16815615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
16829b94acceSBarry Smith /*@
16839b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
168463a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
16859b94acceSBarry Smith 
1686c7afd0dbSLois Curfman McInnes    Collective on SNES
1687c7afd0dbSLois Curfman McInnes 
1688b2002411SLois Curfman McInnes    Input Parameters:
1689c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1690c7afd0dbSLois Curfman McInnes -  x - the solution vector
16919b94acceSBarry Smith 
16929b94acceSBarry Smith    Output Parameter:
1693b2002411SLois Curfman McInnes .  its - number of iterations until termination
16949b94acceSBarry Smith 
1695b2002411SLois Curfman McInnes    Notes:
16968ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
16978ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
16988ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
16998ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
17008ddd3da0SLois Curfman McInnes 
17019b94acceSBarry Smith .keywords: SNES, nonlinear, solve
17029b94acceSBarry Smith 
170363a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
17049b94acceSBarry Smith @*/
17058ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
17069b94acceSBarry Smith {
17073c7409f5SSatish Balay   int ierr, flg;
1708052efed2SBarry Smith 
17093a40ed3dSBarry Smith   PetscFunctionBegin;
171077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
171174679c65SBarry Smith   PetscValidIntPointer(its);
171282bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1713c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
17149b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1715c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
17169b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
17179b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
17183c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
17196d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
17203a40ed3dSBarry Smith   PetscFunctionReturn(0);
17219b94acceSBarry Smith }
17229b94acceSBarry Smith 
17239b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
17249b94acceSBarry Smith 
17255615d1e5SSatish Balay #undef __FUNC__
17265615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
172782bf6240SBarry Smith /*@C
17284b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
17299b94acceSBarry Smith 
1730fee21e36SBarry Smith    Collective on SNES
1731fee21e36SBarry Smith 
1732c7afd0dbSLois Curfman McInnes    Input Parameters:
1733c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1734c7afd0dbSLois Curfman McInnes -  method - a known method
1735c7afd0dbSLois Curfman McInnes 
1736c7afd0dbSLois Curfman McInnes    Options Database Key:
1737c7afd0dbSLois Curfman McInnes .  -snes_type <method> - Sets the method; use -help for a list
1738c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1739ae12b187SLois Curfman McInnes 
17409b94acceSBarry Smith    Notes:
17419b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
1742c7afd0dbSLois Curfman McInnes .    SNES_EQ_LS - Newton's method with line search
1743c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
1744c7afd0dbSLois Curfman McInnes .    SNES_EQ_TR - Newton's method with trust region
1745c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
1746c7afd0dbSLois Curfman McInnes .    SNES_UM_TR - Newton's method with trust region
1747c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
1748c7afd0dbSLois Curfman McInnes .    SNES_UM_LS - Newton's method with line search
1749c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
17509b94acceSBarry Smith 
1751ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1752ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1753ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1754ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1755ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1756ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1757ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1758ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1759ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1760ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1761a703fe33SLois Curfman McInnes 
1762f536c699SSatish Balay .keywords: SNES, set, method
17639b94acceSBarry Smith @*/
17644b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
17659b94acceSBarry Smith {
176684cb2905SBarry Smith   int ierr;
17679b94acceSBarry Smith   int (*r)(SNES);
17683a40ed3dSBarry Smith 
17693a40ed3dSBarry Smith   PetscFunctionBegin;
177077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
177182bf6240SBarry Smith 
17723f1db9ecSBarry Smith   if (PetscTypeCompare(snes->type_name,method)) PetscFunctionReturn(0);
177382bf6240SBarry Smith 
177482bf6240SBarry Smith   if (snes->setupcalled) {
1775e1311b90SBarry Smith     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
177682bf6240SBarry Smith     snes->data = 0;
177782bf6240SBarry Smith   }
17787d1a2b2bSBarry Smith 
17799b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
178082bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);}
178182bf6240SBarry Smith 
1782488ecbafSBarry Smith   ierr =  FListFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr);
178382bf6240SBarry Smith 
17841548b14aSBarry Smith   if (!r) SETERRQ(1,1,"Unable to find requested SNES type");
17851548b14aSBarry Smith 
17860452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
178782bf6240SBarry Smith   snes->data = 0;
17883a40ed3dSBarry Smith   ierr = (*r)(snes); CHKERRQ(ierr);
178982bf6240SBarry Smith 
179082bf6240SBarry Smith   if (snes->type_name) PetscFree(snes->type_name);
179182bf6240SBarry Smith   snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name);
179282bf6240SBarry Smith   PetscStrcpy(snes->type_name,method);
179382bf6240SBarry Smith   snes->set_method_called = 1;
179482bf6240SBarry Smith 
17953a40ed3dSBarry Smith   PetscFunctionReturn(0);
17969b94acceSBarry Smith }
17979b94acceSBarry Smith 
1798a847f771SSatish Balay 
17999b94acceSBarry Smith /* --------------------------------------------------------------------- */
18005615d1e5SSatish Balay #undef __FUNC__
1801d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy"
18029b94acceSBarry Smith /*@C
18039b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
18049b94acceSBarry Smith    registered by SNESRegister().
18059b94acceSBarry Smith 
1806fee21e36SBarry Smith    Not Collective
1807fee21e36SBarry Smith 
18089b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
18099b94acceSBarry Smith 
18109b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
18119b94acceSBarry Smith @*/
1812cf256101SBarry Smith int SNESRegisterDestroy(void)
18139b94acceSBarry Smith {
181482bf6240SBarry Smith   int ierr;
181582bf6240SBarry Smith 
18163a40ed3dSBarry Smith   PetscFunctionBegin;
181782bf6240SBarry Smith   if (SNESList) {
1818488ecbafSBarry Smith     ierr = FListDestroy( SNESList );CHKERRQ(ierr);
181982bf6240SBarry Smith     SNESList = 0;
18209b94acceSBarry Smith   }
182184cb2905SBarry Smith   SNESRegisterAllCalled = 0;
18223a40ed3dSBarry Smith   PetscFunctionReturn(0);
18239b94acceSBarry Smith }
18249b94acceSBarry Smith 
18255615d1e5SSatish Balay #undef __FUNC__
1826d4bb536fSBarry Smith #define __FUNC__ "SNESGetType"
18279b94acceSBarry Smith /*@C
18289a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
18299b94acceSBarry Smith 
1830c7afd0dbSLois Curfman McInnes    Not Collective
1831c7afd0dbSLois Curfman McInnes 
18329b94acceSBarry Smith    Input Parameter:
18334b0e389bSBarry Smith .  snes - nonlinear solver context
18349b94acceSBarry Smith 
18359b94acceSBarry Smith    Output Parameter:
183682bf6240SBarry Smith .  method - SNES method (a charactor string)
18379b94acceSBarry Smith 
18389b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
18399b94acceSBarry Smith @*/
184082bf6240SBarry Smith int SNESGetType(SNES snes, SNESType *method)
18419b94acceSBarry Smith {
18423a40ed3dSBarry Smith   PetscFunctionBegin;
184382bf6240SBarry Smith   *method = snes->type_name;
18443a40ed3dSBarry Smith   PetscFunctionReturn(0);
18459b94acceSBarry Smith }
18469b94acceSBarry Smith 
18475615d1e5SSatish Balay #undef __FUNC__
1848d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution"
18499b94acceSBarry Smith /*@C
18509b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
18519b94acceSBarry Smith    stored.
18529b94acceSBarry Smith 
1853c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1854c7afd0dbSLois Curfman McInnes 
18559b94acceSBarry Smith    Input Parameter:
18569b94acceSBarry Smith .  snes - the SNES context
18579b94acceSBarry Smith 
18589b94acceSBarry Smith    Output Parameter:
18599b94acceSBarry Smith .  x - the solution
18609b94acceSBarry Smith 
18619b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
18629b94acceSBarry Smith 
1863a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
18649b94acceSBarry Smith @*/
18659b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
18669b94acceSBarry Smith {
18673a40ed3dSBarry Smith   PetscFunctionBegin;
186877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18699b94acceSBarry Smith   *x = snes->vec_sol_always;
18703a40ed3dSBarry Smith   PetscFunctionReturn(0);
18719b94acceSBarry Smith }
18729b94acceSBarry Smith 
18735615d1e5SSatish Balay #undef __FUNC__
1874d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate"
18759b94acceSBarry Smith /*@C
18769b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
18779b94acceSBarry Smith    stored.
18789b94acceSBarry Smith 
1879c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1880c7afd0dbSLois Curfman McInnes 
18819b94acceSBarry Smith    Input Parameter:
18829b94acceSBarry Smith .  snes - the SNES context
18839b94acceSBarry Smith 
18849b94acceSBarry Smith    Output Parameter:
18859b94acceSBarry Smith .  x - the solution update
18869b94acceSBarry Smith 
18879b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
18889b94acceSBarry Smith 
18899b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
18909b94acceSBarry Smith @*/
18919b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
18929b94acceSBarry Smith {
18933a40ed3dSBarry Smith   PetscFunctionBegin;
189477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18959b94acceSBarry Smith   *x = snes->vec_sol_update_always;
18963a40ed3dSBarry Smith   PetscFunctionReturn(0);
18979b94acceSBarry Smith }
18989b94acceSBarry Smith 
18995615d1e5SSatish Balay #undef __FUNC__
1900d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction"
19019b94acceSBarry Smith /*@C
19023638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
19039b94acceSBarry Smith 
1904c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1905c7afd0dbSLois Curfman McInnes 
19069b94acceSBarry Smith    Input Parameter:
19079b94acceSBarry Smith .  snes - the SNES context
19089b94acceSBarry Smith 
19099b94acceSBarry Smith    Output Parameter:
19103638b69dSLois Curfman McInnes .  r - the function
19119b94acceSBarry Smith 
19129b94acceSBarry Smith    Notes:
19139b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
19149b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
19159b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
19169b94acceSBarry Smith 
1917a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
19189b94acceSBarry Smith 
191931615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
192031615d04SLois Curfman McInnes           SNESGetGradient()
19219b94acceSBarry Smith @*/
19229b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
19239b94acceSBarry Smith {
19243a40ed3dSBarry Smith   PetscFunctionBegin;
192577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1926a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1927a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1928a8c6a408SBarry Smith   }
19299b94acceSBarry Smith   *r = snes->vec_func_always;
19303a40ed3dSBarry Smith   PetscFunctionReturn(0);
19319b94acceSBarry Smith }
19329b94acceSBarry Smith 
19335615d1e5SSatish Balay #undef __FUNC__
1934d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient"
19359b94acceSBarry Smith /*@C
19363638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
19379b94acceSBarry Smith 
1938c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1939c7afd0dbSLois Curfman McInnes 
19409b94acceSBarry Smith    Input Parameter:
19419b94acceSBarry Smith .  snes - the SNES context
19429b94acceSBarry Smith 
19439b94acceSBarry Smith    Output Parameter:
19449b94acceSBarry Smith .  r - the gradient
19459b94acceSBarry Smith 
19469b94acceSBarry Smith    Notes:
19479b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
19489b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
19499b94acceSBarry Smith    SNESGetFunction().
19509b94acceSBarry Smith 
19519b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
19529b94acceSBarry Smith 
195331615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
19549b94acceSBarry Smith @*/
19559b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
19569b94acceSBarry Smith {
19573a40ed3dSBarry Smith   PetscFunctionBegin;
195877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
19593a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1960a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
19613a40ed3dSBarry Smith   }
19629b94acceSBarry Smith   *r = snes->vec_func_always;
19633a40ed3dSBarry Smith   PetscFunctionReturn(0);
19649b94acceSBarry Smith }
19659b94acceSBarry Smith 
19665615d1e5SSatish Balay #undef __FUNC__
1967d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction"
19689b94acceSBarry Smith /*@
19699b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
19709b94acceSBarry Smith    unconstrained minimization problems.
19719b94acceSBarry Smith 
1972c7afd0dbSLois Curfman McInnes    Not Collective
1973c7afd0dbSLois Curfman McInnes 
19749b94acceSBarry Smith    Input Parameter:
19759b94acceSBarry Smith .  snes - the SNES context
19769b94acceSBarry Smith 
19779b94acceSBarry Smith    Output Parameter:
19789b94acceSBarry Smith .  r - the function
19799b94acceSBarry Smith 
19809b94acceSBarry Smith    Notes:
19819b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
19829b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
19839b94acceSBarry Smith    SNESGetFunction().
19849b94acceSBarry Smith 
19859b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
19869b94acceSBarry Smith 
198731615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
19889b94acceSBarry Smith @*/
19899b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
19909b94acceSBarry Smith {
19913a40ed3dSBarry Smith   PetscFunctionBegin;
199277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
199374679c65SBarry Smith   PetscValidScalarPointer(r);
19943a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1995a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
19963a40ed3dSBarry Smith   }
19979b94acceSBarry Smith   *r = snes->fc;
19983a40ed3dSBarry Smith   PetscFunctionReturn(0);
19999b94acceSBarry Smith }
20009b94acceSBarry Smith 
20015615d1e5SSatish Balay #undef __FUNC__
2002d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix"
20033c7409f5SSatish Balay /*@C
20043c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2005d850072dSLois Curfman McInnes    SNES options in the database.
20063c7409f5SSatish Balay 
2007fee21e36SBarry Smith    Collective on SNES
2008fee21e36SBarry Smith 
2009c7afd0dbSLois Curfman McInnes    Input Parameter:
2010c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2011c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2012c7afd0dbSLois Curfman McInnes 
2013d850072dSLois Curfman McInnes    Notes:
2014a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2015c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2016d850072dSLois Curfman McInnes 
20173c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2018a86d99e1SLois Curfman McInnes 
2019a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
20203c7409f5SSatish Balay @*/
20213c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
20223c7409f5SSatish Balay {
20233c7409f5SSatish Balay   int ierr;
20243c7409f5SSatish Balay 
20253a40ed3dSBarry Smith   PetscFunctionBegin;
202677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2027639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20283c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
20293a40ed3dSBarry Smith   PetscFunctionReturn(0);
20303c7409f5SSatish Balay }
20313c7409f5SSatish Balay 
20325615d1e5SSatish Balay #undef __FUNC__
2033d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix"
20343c7409f5SSatish Balay /*@C
2035f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2036d850072dSLois Curfman McInnes    SNES options in the database.
20373c7409f5SSatish Balay 
2038fee21e36SBarry Smith    Collective on SNES
2039fee21e36SBarry Smith 
2040c7afd0dbSLois Curfman McInnes    Input Parameters:
2041c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2042c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2043c7afd0dbSLois Curfman McInnes 
2044d850072dSLois Curfman McInnes    Notes:
2045a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2046c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2047d850072dSLois Curfman McInnes 
20483c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2049a86d99e1SLois Curfman McInnes 
2050a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20513c7409f5SSatish Balay @*/
20523c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
20533c7409f5SSatish Balay {
20543c7409f5SSatish Balay   int ierr;
20553c7409f5SSatish Balay 
20563a40ed3dSBarry Smith   PetscFunctionBegin;
205777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2058639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20593c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
20603a40ed3dSBarry Smith   PetscFunctionReturn(0);
20613c7409f5SSatish Balay }
20623c7409f5SSatish Balay 
20635615d1e5SSatish Balay #undef __FUNC__
2064d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix"
2065c83590e2SSatish Balay /*@
20663c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20673c7409f5SSatish Balay    SNES options in the database.
20683c7409f5SSatish Balay 
2069c7afd0dbSLois Curfman McInnes    Not Collective
2070c7afd0dbSLois Curfman McInnes 
20713c7409f5SSatish Balay    Input Parameter:
20723c7409f5SSatish Balay .  snes - the SNES context
20733c7409f5SSatish Balay 
20743c7409f5SSatish Balay    Output Parameter:
20753c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20763c7409f5SSatish Balay 
20773c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2078a86d99e1SLois Curfman McInnes 
2079a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20803c7409f5SSatish Balay @*/
20813c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
20823c7409f5SSatish Balay {
20833c7409f5SSatish Balay   int ierr;
20843c7409f5SSatish Balay 
20853a40ed3dSBarry Smith   PetscFunctionBegin;
208677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2087639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20883a40ed3dSBarry Smith   PetscFunctionReturn(0);
20893c7409f5SSatish Balay }
20903c7409f5SSatish Balay 
2091ca161407SBarry Smith #undef __FUNC__
2092ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp"
2093ca161407SBarry Smith /*@
2094ca161407SBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2095ca161407SBarry Smith 
2096c7afd0dbSLois Curfman McInnes    Collective on SNES
2097c7afd0dbSLois Curfman McInnes 
2098ca161407SBarry Smith    Input Parameter:
2099ca161407SBarry Smith .  snes - the SNES context
2100ca161407SBarry Smith 
2101ca161407SBarry Smith    Options Database Keys:
2102c7afd0dbSLois Curfman McInnes +  -help - Prints SNES options
2103c7afd0dbSLois Curfman McInnes -  -h - Prints SNES options
2104ca161407SBarry Smith 
2105ca161407SBarry Smith .keywords: SNES, nonlinear, help
2106ca161407SBarry Smith 
2107ca161407SBarry Smith .seealso: SNESSetFromOptions()
2108ca161407SBarry Smith @*/
2109ca161407SBarry Smith int SNESPrintHelp(SNES snes)
2110ca161407SBarry Smith {
2111ca161407SBarry Smith   char                p[64];
2112ca161407SBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2113ca161407SBarry Smith   int                 ierr;
2114ca161407SBarry Smith 
2115ca161407SBarry Smith   PetscFunctionBegin;
2116ca161407SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2117ca161407SBarry Smith 
2118ca161407SBarry Smith   PetscStrcpy(p,"-");
2119ca161407SBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2120ca161407SBarry Smith 
2121ca161407SBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2122ca161407SBarry Smith 
212376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");
2124488ecbafSBarry Smith   ierr = FListPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr);
212576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
212676be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
212776be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
212876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
212976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
213076be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
213176be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
213276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");
213376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
213476be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
2135ca161407SBarry Smith     residual norm at each iteration.\n",p);
213676be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
2137ca161407SBarry Smith     residual norm for small residual norms. This is useful to conceal\n\
2138ca161407SBarry Smith     meaningless digits that may be different on different machines.\n",p);
213976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
2140ca161407SBarry Smith   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
214176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2142ca161407SBarry Smith      " Options for solving systems of nonlinear equations only:\n");
214376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
214476be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
214576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
2146e24b481bSBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_ksp_monitor - if using matrix-free multiply then prints h at each KSP iteration\n",p);
214776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
214876be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
214976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2150ca161407SBarry Smith      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
215176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2152ca161407SBarry Smith      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
215376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2154ca161407SBarry Smith      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
215576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2156ca161407SBarry Smith      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
215776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2158ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
215976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2160ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
216176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2162ca161407SBarry Smith      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
2163ca161407SBarry Smith   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
216476be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2165ca161407SBarry Smith      " Options for solving unconstrained minimization problems only:\n");
216676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
216776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
216876be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
2169ca161407SBarry Smith   }
217076be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
217176be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"a particular method\n");
2172ca161407SBarry Smith   if (snes->printhelp) {
2173ca161407SBarry Smith     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2174ca161407SBarry Smith   }
2175ca161407SBarry Smith   PetscFunctionReturn(0);
2176ca161407SBarry Smith }
21773c7409f5SSatish Balay 
2178acb85ae6SSatish Balay /*MC
2179b2002411SLois Curfman McInnes    SNESRegister - Adds a method to the nonlinear solver package.
21809b94acceSBarry Smith 
2181b2002411SLois Curfman McInnes    Synopsis:
2182b2002411SLois Curfman McInnes    SNESRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
21839b94acceSBarry Smith 
21848d76a1e5SLois Curfman McInnes    Not collective
21858d76a1e5SLois Curfman McInnes 
2186b2002411SLois Curfman McInnes    Input Parameters:
2187c7afd0dbSLois Curfman McInnes +  name_solver - name of a new user-defined solver
2188b2002411SLois Curfman McInnes .  path - path (either absolute or relative) the library containing this solver
2189b2002411SLois Curfman McInnes .  name_create - name of routine to create method context
2190c7afd0dbSLois Curfman McInnes -  routine_create - routine to create method context
21919b94acceSBarry Smith 
2192b2002411SLois Curfman McInnes    Notes:
2193b2002411SLois Curfman McInnes    SNESRegister() may be called multiple times to add several user-defined solvers.
21943a40ed3dSBarry Smith 
2195b2002411SLois Curfman McInnes    If dynamic libraries are used, then the fourth input argument (routine_create)
2196b2002411SLois Curfman McInnes    is ignored.
2197b2002411SLois Curfman McInnes 
2198b2002411SLois Curfman McInnes    Sample usage:
219969225d78SLois Curfman McInnes .vb
2200b2002411SLois Curfman McInnes    SNESRegister("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2201b2002411SLois Curfman McInnes                 "MySolverCreate",MySolverCreate);
220269225d78SLois Curfman McInnes .ve
2203b2002411SLois Curfman McInnes 
2204b2002411SLois Curfman McInnes    Then, your solver can be chosen with the procedural interface via
2205b2002411SLois Curfman McInnes $     SNESSetType(snes,"my_solver")
2206b2002411SLois Curfman McInnes    or at runtime via the option
2207b2002411SLois Curfman McInnes $     -snes_type my_solver
2208b2002411SLois Curfman McInnes 
2209b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register
2210b2002411SLois Curfman McInnes 
2211b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2212b2002411SLois Curfman McInnes M*/
2213b2002411SLois Curfman McInnes 
2214b2002411SLois Curfman McInnes #undef __FUNC__
2215b2002411SLois Curfman McInnes #define __FUNC__ "SNESRegister_Private"
2216b2002411SLois Curfman McInnes int SNESRegister_Private(char *sname,char *path,char *name,int (*function)(SNES))
2217b2002411SLois Curfman McInnes {
2218b2002411SLois Curfman McInnes   char fullname[256];
2219b2002411SLois Curfman McInnes   int  ierr;
2220b2002411SLois Curfman McInnes 
2221b2002411SLois Curfman McInnes   PetscFunctionBegin;
2222b2002411SLois Curfman McInnes   PetscStrcpy(fullname,path); PetscStrcat(fullname,":");PetscStrcat(fullname,name);
2223488ecbafSBarry Smith   ierr = FListAdd_Private(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr);
2224b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2225b2002411SLois Curfman McInnes }
2226