xref: /petsc/src/snes/interface/snes.c (revision c7afd0db7e14276b7be0e9afcfcc158da35539de)
15cd90555SBarry Smith 
2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
3*c7afd0dbSLois Curfman McInnes static char vcid[] = "$Id: snes.c,v 1.148 1998/04/22 14:14:26 curfman Exp curfman $";
49b94acceSBarry Smith #endif
59b94acceSBarry Smith 
670f55243SBarry Smith #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
7f5eb4b81SSatish Balay #include "src/sys/nreg.h"
89b94acceSBarry Smith #include "pinclude/pviewer.h"
99b94acceSBarry Smith #include <math.h>
109b94acceSBarry Smith 
1184cb2905SBarry Smith int SNESRegisterAllCalled = 0;
1282bf6240SBarry Smith DLList SNESList = 0;
1382bf6240SBarry Smith 
145615d1e5SSatish Balay #undef __FUNC__
15d4bb536fSBarry Smith #define __FUNC__ "SNESView"
169b94acceSBarry Smith /*@
179b94acceSBarry Smith    SNESView - Prints the SNES data structure.
189b94acceSBarry Smith 
19fee21e36SBarry Smith    Collective on SNES, unless Viewer is VIEWER_STDOUT_SELF
20fee21e36SBarry Smith 
21*c7afd0dbSLois Curfman McInnes    Input Parameters:
22*c7afd0dbSLois Curfman McInnes +  SNES - the SNES context
23*c7afd0dbSLois Curfman McInnes -  viewer - visualization context
24*c7afd0dbSLois Curfman McInnes 
259b94acceSBarry Smith    Options Database Key:
26c8a8ba5cSLois Curfman McInnes .  -snes_view - Calls SNESView() at end of SNESSolve()
279b94acceSBarry Smith 
289b94acceSBarry Smith    Notes:
299b94acceSBarry Smith    The available visualization contexts include
30*c7afd0dbSLois Curfman McInnes +     VIEWER_STDOUT_SELF - standard output (default)
31*c7afd0dbSLois Curfman McInnes -     VIEWER_STDOUT_WORLD - synchronized standard
32c8a8ba5cSLois Curfman McInnes          output where only the first processor opens
33c8a8ba5cSLois Curfman McInnes          the file.  All other processors send their
34c8a8ba5cSLois Curfman McInnes          data to the first processor to print.
359b94acceSBarry Smith 
36*c7afd0dbSLois Curfman McInnes    The user can open an alternative vistualization context with
37*c7afd0dbSLois Curfman McInnes    ViewerFileOpenASCII() - output to a specified file
389b94acceSBarry Smith 
399b94acceSBarry Smith .keywords: SNES, view
409b94acceSBarry Smith 
41dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
429b94acceSBarry Smith @*/
439b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
449b94acceSBarry Smith {
459b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
469b94acceSBarry Smith   FILE                *fd;
479b94acceSBarry Smith   int                 ierr;
489b94acceSBarry Smith   SLES                sles;
499b94acceSBarry Smith   char                *method;
5019bcc07fSBarry Smith   ViewerType          vtype;
519b94acceSBarry Smith 
523a40ed3dSBarry Smith   PetscFunctionBegin;
5374679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5474679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5574679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5674679c65SBarry Smith 
5719bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5819bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
5990ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
6077c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
6182bf6240SBarry Smith     SNESGetType(snes,&method);
6277c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
63e1311b90SBarry Smith     if (snes->view) (*snes->view)(snes,viewer);
6477c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
659b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
669b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6777c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
689b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
699b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
707a00f4a9SLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
717a00f4a9SLois Curfman McInnes     "  total number of linear solver iterations=%d\n",snes->linear_its);
72ae3c334cSLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
7368632166SLois Curfman McInnes      "  total number of function evaluations=%d\n",snes->nfuncs);
749b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
7577c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
769b94acceSBarry Smith     if (snes->ksp_ewconv) {
779b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
789b94acceSBarry Smith       if (kctx) {
7977c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
809b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
819b94acceSBarry Smith         kctx->version);
8277c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
839b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
849b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
8577c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
869b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
879b94acceSBarry Smith       }
889b94acceSBarry Smith     }
89c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
9082bf6240SBarry Smith     SNESGetType(snes,&method);
91c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
9219bcc07fSBarry Smith   }
939b94acceSBarry Smith   SNESGetSLES(snes,&sles);
949b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
953a40ed3dSBarry Smith   PetscFunctionReturn(0);
969b94acceSBarry Smith }
979b94acceSBarry Smith 
98639f9d9dSBarry Smith /*
99639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
100639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
101639f9d9dSBarry Smith */
102639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
103639f9d9dSBarry Smith static int numberofsetfromoptions;
104639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
105639f9d9dSBarry Smith 
1065615d1e5SSatish Balay #undef __FUNC__
107d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker"
108639f9d9dSBarry Smith /*@
109639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
110639f9d9dSBarry Smith 
111*c7afd0dbSLois Curfman McInnes     Not Collective
112*c7afd0dbSLois Curfman McInnes 
113639f9d9dSBarry Smith     Input Parameter:
114639f9d9dSBarry Smith .   snescheck - function that checks for options
115639f9d9dSBarry Smith 
116639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
117639f9d9dSBarry Smith @*/
118639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
119639f9d9dSBarry Smith {
1203a40ed3dSBarry Smith   PetscFunctionBegin;
121639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
122a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed");
123639f9d9dSBarry Smith   }
124639f9d9dSBarry Smith 
125639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
1263a40ed3dSBarry Smith   PetscFunctionReturn(0);
127639f9d9dSBarry Smith }
128639f9d9dSBarry Smith 
1295615d1e5SSatish Balay #undef __FUNC__
1305615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1319b94acceSBarry Smith /*@
132682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1339b94acceSBarry Smith 
134*c7afd0dbSLois Curfman McInnes    Collective on SNES
135*c7afd0dbSLois Curfman McInnes 
1369b94acceSBarry Smith    Input Parameter:
1379b94acceSBarry Smith .  snes - the SNES context
1389b94acceSBarry Smith 
13911ca99fdSLois Curfman McInnes    Notes:
14011ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
14111ca99fdSLois Curfman McInnes    the users manual.
14283e2fdc7SBarry Smith 
1439b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1449b94acceSBarry Smith 
145a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix()
1469b94acceSBarry Smith @*/
1479b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1489b94acceSBarry Smith {
14982bf6240SBarry Smith   char     method[256];
1509b94acceSBarry Smith   double   tmp;
1519b94acceSBarry Smith   SLES     sles;
15281f57ec7SBarry Smith   int      ierr, flg,i,loc[4],nmax = 4;
1533c7409f5SSatish Balay   int      version   = PETSC_DEFAULT;
1549b94acceSBarry Smith   double   rtol_0    = PETSC_DEFAULT;
1559b94acceSBarry Smith   double   rtol_max  = PETSC_DEFAULT;
1569b94acceSBarry Smith   double   gamma2    = PETSC_DEFAULT;
1579b94acceSBarry Smith   double   alpha     = PETSC_DEFAULT;
1589b94acceSBarry Smith   double   alpha2    = PETSC_DEFAULT;
1599b94acceSBarry Smith   double   threshold = PETSC_DEFAULT;
1609b94acceSBarry Smith 
1613a40ed3dSBarry Smith   PetscFunctionBegin;
16281f57ec7SBarry Smith   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
16381f57ec7SBarry Smith 
16477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16582bf6240SBarry Smith   if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp");
166ca161407SBarry Smith 
16782bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
16882bf6240SBarry Smith   ierr = OptionsGetString(snes->prefix,"-snes_type",method,256,&flg);
169052efed2SBarry Smith   if (flg) {
17082bf6240SBarry Smith     ierr = SNESSetType(snes,(SNESType) method); CHKERRQ(ierr);
1715334005bSBarry Smith   }
1723c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr);
173d64ed03dSBarry Smith   if (flg) {
174d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
175d64ed03dSBarry Smith   }
1763c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr);
177d64ed03dSBarry Smith   if (flg) {
178d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
179d64ed03dSBarry Smith   }
1803c7409f5SSatish Balay   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);  CHKERRQ(ierr);
181d64ed03dSBarry Smith   if (flg) {
182d64ed03dSBarry Smith     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
183d64ed03dSBarry Smith   }
1843c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr);
1853c7409f5SSatish Balay   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr);
186d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr);
187d64ed03dSBarry Smith   if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); }
188d7a720efSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr);
189d64ed03dSBarry Smith   if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);  CHKERRQ(ierr);}
1903c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr);
1913c7409f5SSatish Balay   if (flg) { snes->ksp_ewconv = 1; }
192b18e04deSLois Curfman McInnes   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr);
193b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr);
194b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr);
195b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr);
196b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr);
197b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr);
198b18e04deSLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr);
19993c39befSBarry Smith 
20093c39befSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr);
20193c39befSBarry Smith   if (flg) {snes->converged = 0;}
20293c39befSBarry Smith 
2039b94acceSBarry Smith   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
2049b94acceSBarry Smith                                   alpha2,threshold); CHKERRQ(ierr);
2055bbad29bSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr);
2065cd90555SBarry Smith   if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
2073c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr);
208639f9d9dSBarry Smith   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);}
2093c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr);
210d31a3109SLois Curfman McInnes   if (flg) {
211639f9d9dSBarry Smith    ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0);  CHKERRQ(ierr);
212d31a3109SLois Curfman McInnes   }
21381f57ec7SBarry Smith   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
2143c7409f5SSatish Balay   if (flg) {
21517699dbbSLois Curfman McInnes     int    rank = 0;
216d7e8b826SBarry Smith     DrawLG lg;
21717699dbbSLois Curfman McInnes     MPI_Initialized(&rank);
21817699dbbSLois Curfman McInnes     if (rank) MPI_Comm_rank(snes->comm,&rank);
21917699dbbSLois Curfman McInnes     if (!rank) {
22081f57ec7SBarry Smith       ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr);
2219b94acceSBarry Smith       ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr);
222f8c826e1SBarry Smith       snes->xmonitor = lg;
2239b94acceSBarry Smith       PLogObjectParent(snes,lg);
2249b94acceSBarry Smith     }
2259b94acceSBarry Smith   }
2263c7409f5SSatish Balay   ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);  CHKERRQ(ierr);
2273c7409f5SSatish Balay   if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
2289b94acceSBarry Smith     ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,
2299b94acceSBarry Smith                  SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr);
230981c4779SBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
231981c4779SBarry Smith   } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
23231615d04SLois Curfman McInnes     ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,
23331615d04SLois Curfman McInnes                  SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr);
234d64ed03dSBarry Smith     PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2359b94acceSBarry Smith   }
236639f9d9dSBarry Smith 
237639f9d9dSBarry Smith   for ( i=0; i<numberofsetfromoptions; i++ ) {
238639f9d9dSBarry Smith     ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr);
239639f9d9dSBarry Smith   }
240639f9d9dSBarry Smith 
2419b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
2429b94acceSBarry Smith   ierr = SLESSetFromOptions(sles); CHKERRQ(ierr);
24382bf6240SBarry Smith   /*
24482bf6240SBarry Smith         Since the private setfromoptions requires the type to all ready have
24582bf6240SBarry Smith       been set we make sure a type is set by this time
24682bf6240SBarry Smith   */
24782bf6240SBarry Smith   if (!snes->type_name) {
24882bf6240SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
24982bf6240SBarry Smith       ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr);
25082bf6240SBarry Smith     } else {
25182bf6240SBarry Smith       ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr);
25282bf6240SBarry Smith     }
25382bf6240SBarry Smith   }
25482bf6240SBarry Smith 
25582bf6240SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
25682bf6240SBarry Smith   if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);}
25782bf6240SBarry Smith 
2583a40ed3dSBarry Smith   if (snes->setfromoptions) {
2593a40ed3dSBarry Smith     ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
2603a40ed3dSBarry Smith   }
2613a40ed3dSBarry Smith   PetscFunctionReturn(0);
2629b94acceSBarry Smith }
2639b94acceSBarry Smith 
264a847f771SSatish Balay 
2655615d1e5SSatish Balay #undef __FUNC__
266d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext"
2679b94acceSBarry Smith /*@
2689b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2699b94acceSBarry Smith    the nonlinear solvers.
2709b94acceSBarry Smith 
271fee21e36SBarry Smith    Collective on SNES
272fee21e36SBarry Smith 
273*c7afd0dbSLois Curfman McInnes    Input Parameters:
274*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
275*c7afd0dbSLois Curfman McInnes -  usrP - optional user context
276*c7afd0dbSLois Curfman McInnes 
2779b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2789b94acceSBarry Smith 
2799b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2809b94acceSBarry Smith @*/
2819b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
2829b94acceSBarry Smith {
2833a40ed3dSBarry Smith   PetscFunctionBegin;
28477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2859b94acceSBarry Smith   snes->user		= usrP;
2863a40ed3dSBarry Smith   PetscFunctionReturn(0);
2879b94acceSBarry Smith }
28874679c65SBarry Smith 
2895615d1e5SSatish Balay #undef __FUNC__
290d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext"
2919b94acceSBarry Smith /*@C
2929b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
2939b94acceSBarry Smith    nonlinear solvers.
2949b94acceSBarry Smith 
295*c7afd0dbSLois Curfman McInnes    Not Collective
296*c7afd0dbSLois Curfman McInnes 
2979b94acceSBarry Smith    Input Parameter:
2989b94acceSBarry Smith .  snes - SNES context
2999b94acceSBarry Smith 
3009b94acceSBarry Smith    Output Parameter:
3019b94acceSBarry Smith .  usrP - user context
3029b94acceSBarry Smith 
3039b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3049b94acceSBarry Smith 
3059b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3069b94acceSBarry Smith @*/
3079b94acceSBarry Smith int SNESGetApplicationContext( SNES snes,  void **usrP )
3089b94acceSBarry Smith {
3093a40ed3dSBarry Smith   PetscFunctionBegin;
31077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3119b94acceSBarry Smith   *usrP = snes->user;
3123a40ed3dSBarry Smith   PetscFunctionReturn(0);
3139b94acceSBarry Smith }
31474679c65SBarry Smith 
3155615d1e5SSatish Balay #undef __FUNC__
316d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber"
3179b94acceSBarry Smith /*@
3189b94acceSBarry Smith    SNESGetIterationNumber - Gets the current iteration number of the
3199b94acceSBarry Smith    nonlinear solver.
3209b94acceSBarry Smith 
321*c7afd0dbSLois Curfman McInnes    Not Collective
322*c7afd0dbSLois Curfman McInnes 
3239b94acceSBarry Smith    Input Parameter:
3249b94acceSBarry Smith .  snes - SNES context
3259b94acceSBarry Smith 
3269b94acceSBarry Smith    Output Parameter:
3279b94acceSBarry Smith .  iter - iteration number
3289b94acceSBarry Smith 
3299b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3309b94acceSBarry Smith @*/
3319b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3329b94acceSBarry Smith {
3333a40ed3dSBarry Smith   PetscFunctionBegin;
33477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
33574679c65SBarry Smith   PetscValidIntPointer(iter);
3369b94acceSBarry Smith   *iter = snes->iter;
3373a40ed3dSBarry Smith   PetscFunctionReturn(0);
3389b94acceSBarry Smith }
33974679c65SBarry Smith 
3405615d1e5SSatish Balay #undef __FUNC__
3415615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm"
3429b94acceSBarry Smith /*@
3439b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3449b94acceSBarry Smith    with SNESSSetFunction().
3459b94acceSBarry Smith 
346*c7afd0dbSLois Curfman McInnes    Collective on SNES
347*c7afd0dbSLois Curfman McInnes 
3489b94acceSBarry Smith    Input Parameter:
3499b94acceSBarry Smith .  snes - SNES context
3509b94acceSBarry Smith 
3519b94acceSBarry Smith    Output Parameter:
3529b94acceSBarry Smith .  fnorm - 2-norm of function
3539b94acceSBarry Smith 
3549b94acceSBarry Smith    Note:
3559b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
356a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
357a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3589b94acceSBarry Smith 
3599b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
360a86d99e1SLois Curfman McInnes 
361a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction()
3629b94acceSBarry Smith @*/
3639b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3649b94acceSBarry Smith {
3653a40ed3dSBarry Smith   PetscFunctionBegin;
36677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
36774679c65SBarry Smith   PetscValidScalarPointer(fnorm);
36874679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
369d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only");
37074679c65SBarry Smith   }
3719b94acceSBarry Smith   *fnorm = snes->norm;
3723a40ed3dSBarry Smith   PetscFunctionReturn(0);
3739b94acceSBarry Smith }
37474679c65SBarry Smith 
3755615d1e5SSatish Balay #undef __FUNC__
3765615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm"
3779b94acceSBarry Smith /*@
3789b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
3799b94acceSBarry Smith    with SNESSSetGradient().
3809b94acceSBarry Smith 
381*c7afd0dbSLois Curfman McInnes    Collective on SNES
382*c7afd0dbSLois Curfman McInnes 
3839b94acceSBarry Smith    Input Parameter:
3849b94acceSBarry Smith .  snes - SNES context
3859b94acceSBarry Smith 
3869b94acceSBarry Smith    Output Parameter:
3879b94acceSBarry Smith .  fnorm - 2-norm of gradient
3889b94acceSBarry Smith 
3899b94acceSBarry Smith    Note:
3909b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
391a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
392a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
3939b94acceSBarry Smith 
3949b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
395a86d99e1SLois Curfman McInnes 
396a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
3979b94acceSBarry Smith @*/
3989b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
3999b94acceSBarry Smith {
4003a40ed3dSBarry Smith   PetscFunctionBegin;
40177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
40274679c65SBarry Smith   PetscValidScalarPointer(gnorm);
40374679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
404d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
40574679c65SBarry Smith   }
4069b94acceSBarry Smith   *gnorm = snes->norm;
4073a40ed3dSBarry Smith   PetscFunctionReturn(0);
4089b94acceSBarry Smith }
40974679c65SBarry Smith 
4105615d1e5SSatish Balay #undef __FUNC__
411d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps"
4129b94acceSBarry Smith /*@
4139b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4149b94acceSBarry Smith    attempted by the nonlinear solver.
4159b94acceSBarry Smith 
416*c7afd0dbSLois Curfman McInnes    Not Collective
417*c7afd0dbSLois Curfman McInnes 
4189b94acceSBarry Smith    Input Parameter:
4199b94acceSBarry Smith .  snes - SNES context
4209b94acceSBarry Smith 
4219b94acceSBarry Smith    Output Parameter:
4229b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4239b94acceSBarry Smith 
424c96a6f78SLois Curfman McInnes    Notes:
425c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
426c96a6f78SLois Curfman McInnes 
4279b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4289b94acceSBarry Smith @*/
4299b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4309b94acceSBarry Smith {
4313a40ed3dSBarry Smith   PetscFunctionBegin;
43277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43374679c65SBarry Smith   PetscValidIntPointer(nfails);
4349b94acceSBarry Smith   *nfails = snes->nfailures;
4353a40ed3dSBarry Smith   PetscFunctionReturn(0);
4369b94acceSBarry Smith }
437a847f771SSatish Balay 
4385615d1e5SSatish Balay #undef __FUNC__
439d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations"
440c96a6f78SLois Curfman McInnes /*@
441c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
442c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
443c96a6f78SLois Curfman McInnes 
444*c7afd0dbSLois Curfman McInnes    Not Collective
445*c7afd0dbSLois Curfman McInnes 
446c96a6f78SLois Curfman McInnes    Input Parameter:
447c96a6f78SLois Curfman McInnes .  snes - SNES context
448c96a6f78SLois Curfman McInnes 
449c96a6f78SLois Curfman McInnes    Output Parameter:
450c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
451c96a6f78SLois Curfman McInnes 
452c96a6f78SLois Curfman McInnes    Notes:
453c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
454c96a6f78SLois Curfman McInnes 
455c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
456c96a6f78SLois Curfman McInnes @*/
45786bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
458c96a6f78SLois Curfman McInnes {
4593a40ed3dSBarry Smith   PetscFunctionBegin;
460c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
461c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
462c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
4633a40ed3dSBarry Smith   PetscFunctionReturn(0);
464c96a6f78SLois Curfman McInnes }
465c96a6f78SLois Curfman McInnes 
466c96a6f78SLois Curfman McInnes #undef __FUNC__
467d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES"
4689b94acceSBarry Smith /*@C
4699b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4709b94acceSBarry Smith 
471*c7afd0dbSLois Curfman McInnes    Not Collective, but if SNES object is parallel, then SLES object is parallel
472*c7afd0dbSLois Curfman McInnes 
4739b94acceSBarry Smith    Input Parameter:
4749b94acceSBarry Smith .  snes - the SNES context
4759b94acceSBarry Smith 
4769b94acceSBarry Smith    Output Parameter:
4779b94acceSBarry Smith .  sles - the SLES context
4789b94acceSBarry Smith 
4799b94acceSBarry Smith    Notes:
4809b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4819b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4829b94acceSBarry Smith    KSP and PC contexts as well.
4839b94acceSBarry Smith 
4849b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4859b94acceSBarry Smith 
4869b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4879b94acceSBarry Smith @*/
4889b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4899b94acceSBarry Smith {
4903a40ed3dSBarry Smith   PetscFunctionBegin;
49177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4929b94acceSBarry Smith   *sles = snes->sles;
4933a40ed3dSBarry Smith   PetscFunctionReturn(0);
4949b94acceSBarry Smith }
49582bf6240SBarry Smith 
4969b94acceSBarry Smith /* -----------------------------------------------------------*/
4975615d1e5SSatish Balay #undef __FUNC__
4985615d1e5SSatish Balay #define __FUNC__ "SNESCreate"
4999b94acceSBarry Smith /*@C
5009b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5019b94acceSBarry Smith 
502*c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
503*c7afd0dbSLois Curfman McInnes 
504*c7afd0dbSLois Curfman McInnes    Input Parameters:
505*c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
506*c7afd0dbSLois Curfman McInnes -  type - type of method, either
507*c7afd0dbSLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations)
508*c7afd0dbSLois Curfman McInnes    or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization)
5099b94acceSBarry Smith 
5109b94acceSBarry Smith    Output Parameter:
5119b94acceSBarry Smith .  outsnes - the new SNES context
5129b94acceSBarry Smith 
513*c7afd0dbSLois Curfman McInnes    Options Database Keys:
514*c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
515*c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
516*c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
517*c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
518*c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
519*c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
520c1f60f51SBarry Smith 
5219b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5229b94acceSBarry Smith 
52363a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5249b94acceSBarry Smith @*/
5254b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5269b94acceSBarry Smith {
5279b94acceSBarry Smith   int                 ierr;
5289b94acceSBarry Smith   SNES                snes;
5299b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
53037fcc0dbSBarry Smith 
5313a40ed3dSBarry Smith   PetscFunctionBegin;
5329b94acceSBarry Smith   *outsnes = 0;
533d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
534d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
535d64ed03dSBarry Smith   }
536f830108cSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,comm,SNESDestroy,SNESView);
5379b94acceSBarry Smith   PLogObjectCreate(snes);
5389b94acceSBarry Smith   snes->max_its           = 50;
5399750a799SBarry Smith   snes->max_funcs	  = 10000;
5409b94acceSBarry Smith   snes->norm		  = 0.0;
5415a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
542b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
543b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5449b94acceSBarry Smith     snes->atol		  = 1.e-10;
545d64ed03dSBarry Smith   } else {
546b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
547b4874afaSBarry Smith     snes->ttol            = 0.0;
548b4874afaSBarry Smith     snes->atol		  = 1.e-50;
549b4874afaSBarry Smith   }
5509b94acceSBarry Smith   snes->xtol		  = 1.e-8;
551b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5529b94acceSBarry Smith   snes->nfuncs            = 0;
5539b94acceSBarry Smith   snes->nfailures         = 0;
5547a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
555639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5569b94acceSBarry Smith   snes->data              = 0;
5579b94acceSBarry Smith   snes->view              = 0;
5589b94acceSBarry Smith   snes->computeumfunction = 0;
5599b94acceSBarry Smith   snes->umfunP            = 0;
5609b94acceSBarry Smith   snes->fc                = 0;
5619b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5629b94acceSBarry Smith   snes->fmin              = -1.e30;
5639b94acceSBarry Smith   snes->method_class      = type;
5649b94acceSBarry Smith   snes->set_method_called = 0;
56582bf6240SBarry Smith   snes->setupcalled      = 0;
5669b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5676f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5686f24a144SLois Curfman McInnes   snes->vwork             = 0;
5696f24a144SLois Curfman McInnes   snes->nwork             = 0;
570c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
571c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
572c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
5739b94acceSBarry Smith 
5749b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5750452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
576eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
5779b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5789b94acceSBarry Smith   kctx->version     = 2;
5799b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5809b94acceSBarry Smith                              this was too large for some test cases */
5819b94acceSBarry Smith   kctx->rtol_last   = 0;
5829b94acceSBarry Smith   kctx->rtol_max    = .9;
5839b94acceSBarry Smith   kctx->gamma       = 1.0;
5849b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5859b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5869b94acceSBarry Smith   kctx->threshold   = .1;
5879b94acceSBarry Smith   kctx->lresid_last = 0;
5889b94acceSBarry Smith   kctx->norm_last   = 0;
5899b94acceSBarry Smith 
5909b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5919b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5925334005bSBarry Smith 
5939b94acceSBarry Smith   *outsnes = snes;
5943a40ed3dSBarry Smith   PetscFunctionReturn(0);
5959b94acceSBarry Smith }
5969b94acceSBarry Smith 
5979b94acceSBarry Smith /* --------------------------------------------------------------- */
5985615d1e5SSatish Balay #undef __FUNC__
599d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction"
6009b94acceSBarry Smith /*@C
6019b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6029b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6039b94acceSBarry Smith    equations.
6049b94acceSBarry Smith 
605fee21e36SBarry Smith    Collective on SNES
606fee21e36SBarry Smith 
607*c7afd0dbSLois Curfman McInnes    Input Parameters:
608*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
609*c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
610*c7afd0dbSLois Curfman McInnes .  r - vector to store function value
611*c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
612*c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6139b94acceSBarry Smith 
614*c7afd0dbSLois Curfman McInnes    Calling sequence of func:
615*c7afd0dbSLois Curfman McInnes $  func (SNES snes,Vec x,Vec f,void *ctx);
616*c7afd0dbSLois Curfman McInnes 
617*c7afd0dbSLois Curfman McInnes +  x - input vector
618313e4042SLois Curfman McInnes .  f - function vector
619*c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
6209b94acceSBarry Smith 
6219b94acceSBarry Smith    Notes:
6229b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6239b94acceSBarry Smith $      f'(x) x = -f(x),
624*c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
6259b94acceSBarry Smith 
6269b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6279b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6289b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6299b94acceSBarry Smith 
6309b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6319b94acceSBarry Smith 
632a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6339b94acceSBarry Smith @*/
6345334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6359b94acceSBarry Smith {
6363a40ed3dSBarry Smith   PetscFunctionBegin;
63777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
638a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
639a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
640a8c6a408SBarry Smith   }
6419b94acceSBarry Smith   snes->computefunction     = func;
6429b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6439b94acceSBarry Smith   snes->funP                = ctx;
6443a40ed3dSBarry Smith   PetscFunctionReturn(0);
6459b94acceSBarry Smith }
6469b94acceSBarry Smith 
6475615d1e5SSatish Balay #undef __FUNC__
6485615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
6499b94acceSBarry Smith /*@
6509b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6519b94acceSBarry Smith    SNESSetFunction().
6529b94acceSBarry Smith 
653*c7afd0dbSLois Curfman McInnes    Collective on SNES
654*c7afd0dbSLois Curfman McInnes 
6559b94acceSBarry Smith    Input Parameters:
656*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
657*c7afd0dbSLois Curfman McInnes -  x - input vector
6589b94acceSBarry Smith 
6599b94acceSBarry Smith    Output Parameter:
6603638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6619b94acceSBarry Smith 
6621bffabb2SLois Curfman McInnes    Notes:
6639b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6649b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6659b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6669b94acceSBarry Smith 
6679b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6689b94acceSBarry Smith 
669a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6709b94acceSBarry Smith @*/
6719b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6729b94acceSBarry Smith {
6739b94acceSBarry Smith   int    ierr;
6749b94acceSBarry Smith 
6753a40ed3dSBarry Smith   PetscFunctionBegin;
676d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
677a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
678d4bb536fSBarry Smith   }
6799b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
680d64ed03dSBarry Smith   PetscStackPush("SNES user function");
6819b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
682d64ed03dSBarry Smith   PetscStackPop;
683ae3c334cSLois Curfman McInnes   snes->nfuncs++;
6849b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6853a40ed3dSBarry Smith   PetscFunctionReturn(0);
6869b94acceSBarry Smith }
6879b94acceSBarry Smith 
6885615d1e5SSatish Balay #undef __FUNC__
689d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction"
6909b94acceSBarry Smith /*@C
6919b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6929b94acceSBarry Smith    unconstrained minimization.
6939b94acceSBarry Smith 
694fee21e36SBarry Smith    Collective on SNES
695fee21e36SBarry Smith 
696*c7afd0dbSLois Curfman McInnes    Input Parameters:
697*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
698*c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
699*c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
700*c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7019b94acceSBarry Smith 
702*c7afd0dbSLois Curfman McInnes    Calling sequence of func:
703*c7afd0dbSLois Curfman McInnes $  func (SNES snes,Vec x,double *f,void *ctx);
704*c7afd0dbSLois Curfman McInnes 
705*c7afd0dbSLois Curfman McInnes +  x - input vector
7069b94acceSBarry Smith .  f - function
707*c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined function context
7089b94acceSBarry Smith 
7099b94acceSBarry Smith    Notes:
7109b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7119b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7129b94acceSBarry Smith    SNESSetFunction().
7139b94acceSBarry Smith 
7149b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
7159b94acceSBarry Smith 
716a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
717a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
7189b94acceSBarry Smith @*/
7199b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
7209b94acceSBarry Smith                       void *ctx)
7219b94acceSBarry Smith {
7223a40ed3dSBarry Smith   PetscFunctionBegin;
72377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
724a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
725a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
726a8c6a408SBarry Smith   }
7279b94acceSBarry Smith   snes->computeumfunction   = func;
7289b94acceSBarry Smith   snes->umfunP              = ctx;
7293a40ed3dSBarry Smith   PetscFunctionReturn(0);
7309b94acceSBarry Smith }
7319b94acceSBarry Smith 
7325615d1e5SSatish Balay #undef __FUNC__
7335615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7349b94acceSBarry Smith /*@
7359b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7369b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7379b94acceSBarry Smith 
738*c7afd0dbSLois Curfman McInnes    Collective on SNES
739*c7afd0dbSLois Curfman McInnes 
7409b94acceSBarry Smith    Input Parameters:
741*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
742*c7afd0dbSLois Curfman McInnes -  x - input vector
7439b94acceSBarry Smith 
7449b94acceSBarry Smith    Output Parameter:
7459b94acceSBarry Smith .  y - function value
7469b94acceSBarry Smith 
7479b94acceSBarry Smith    Notes:
7489b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7499b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7509b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
751a86d99e1SLois Curfman McInnes 
752a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
753a86d99e1SLois Curfman McInnes 
754a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
755a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7569b94acceSBarry Smith @*/
7579b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7589b94acceSBarry Smith {
7599b94acceSBarry Smith   int    ierr;
7603a40ed3dSBarry Smith 
7613a40ed3dSBarry Smith   PetscFunctionBegin;
762a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
763a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
764a8c6a408SBarry Smith   }
7659b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
766d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
7679b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
768d64ed03dSBarry Smith   PetscStackPop;
769ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7709b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7713a40ed3dSBarry Smith   PetscFunctionReturn(0);
7729b94acceSBarry Smith }
7739b94acceSBarry Smith 
7745615d1e5SSatish Balay #undef __FUNC__
775d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient"
7769b94acceSBarry Smith /*@C
7779b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7789b94acceSBarry Smith    vector for use by the SNES routines.
7799b94acceSBarry Smith 
780*c7afd0dbSLois Curfman McInnes    Collective on SNES
781*c7afd0dbSLois Curfman McInnes 
7829b94acceSBarry Smith    Input Parameters:
783*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
7849b94acceSBarry Smith .  func - function evaluation routine
785044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
786044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
787*c7afd0dbSLois Curfman McInnes -  r - vector to store gradient value
788fee21e36SBarry Smith 
7899b94acceSBarry Smith    Calling sequence of func:
790*c7afd0dbSLois Curfman McInnes $  func (SNES, Vec x, Vec g, void *ctx);
7919b94acceSBarry Smith 
792*c7afd0dbSLois Curfman McInnes +  x - input vector
7939b94acceSBarry Smith .  g - gradient vector
794*c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined gradient context
7959b94acceSBarry Smith 
7969b94acceSBarry Smith    Notes:
7979b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7989b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7999b94acceSBarry Smith    SNESSetFunction().
8009b94acceSBarry Smith 
8019b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
8029b94acceSBarry Smith 
803a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
804a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
8059b94acceSBarry Smith @*/
80674679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
8079b94acceSBarry Smith {
8083a40ed3dSBarry Smith   PetscFunctionBegin;
80977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
810a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
811a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
812a8c6a408SBarry Smith   }
8139b94acceSBarry Smith   snes->computefunction     = func;
8149b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
8159b94acceSBarry Smith   snes->funP                = ctx;
8163a40ed3dSBarry Smith   PetscFunctionReturn(0);
8179b94acceSBarry Smith }
8189b94acceSBarry Smith 
8195615d1e5SSatish Balay #undef __FUNC__
8205615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
8219b94acceSBarry Smith /*@
822a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
823a86d99e1SLois Curfman McInnes    SNESSetGradient().
8249b94acceSBarry Smith 
825*c7afd0dbSLois Curfman McInnes    Collective on SNES
826*c7afd0dbSLois Curfman McInnes 
8279b94acceSBarry Smith    Input Parameters:
828*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
829*c7afd0dbSLois Curfman McInnes -  x - input vector
8309b94acceSBarry Smith 
8319b94acceSBarry Smith    Output Parameter:
8329b94acceSBarry Smith .  y - gradient vector
8339b94acceSBarry Smith 
8349b94acceSBarry Smith    Notes:
8359b94acceSBarry Smith    SNESComputeGradient() is valid only for
8369b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8379b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
838a86d99e1SLois Curfman McInnes 
839a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
840a86d99e1SLois Curfman McInnes 
841a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
842a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8439b94acceSBarry Smith @*/
8449b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8459b94acceSBarry Smith {
8469b94acceSBarry Smith   int    ierr;
8473a40ed3dSBarry Smith 
8483a40ed3dSBarry Smith   PetscFunctionBegin;
8493a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
850a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8513a40ed3dSBarry Smith   }
8529b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
853d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
8549b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
855d64ed03dSBarry Smith   PetscStackPop;
8569b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8573a40ed3dSBarry Smith   PetscFunctionReturn(0);
8589b94acceSBarry Smith }
8599b94acceSBarry Smith 
8605615d1e5SSatish Balay #undef __FUNC__
8615615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
86262fef451SLois Curfman McInnes /*@
86362fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
86462fef451SLois Curfman McInnes    set with SNESSetJacobian().
86562fef451SLois Curfman McInnes 
866*c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
867*c7afd0dbSLois Curfman McInnes 
86862fef451SLois Curfman McInnes    Input Parameters:
869*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
870*c7afd0dbSLois Curfman McInnes -  x - input vector
87162fef451SLois Curfman McInnes 
87262fef451SLois Curfman McInnes    Output Parameters:
873*c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
87462fef451SLois Curfman McInnes .  B - optional preconditioning matrix
875*c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
876fee21e36SBarry Smith 
87762fef451SLois Curfman McInnes    Notes:
87862fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
87962fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
88062fef451SLois Curfman McInnes 
881dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
882dc5a77f8SLois Curfman McInnes    flag parameter.
88362fef451SLois Curfman McInnes 
88462fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
88562fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
886005c665bSBarry Smith    methods is SNESComputeHessian().
88762fef451SLois Curfman McInnes 
88862fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
88962fef451SLois Curfman McInnes 
89062fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
89162fef451SLois Curfman McInnes @*/
8929b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8939b94acceSBarry Smith {
8949b94acceSBarry Smith   int    ierr;
8953a40ed3dSBarry Smith 
8963a40ed3dSBarry Smith   PetscFunctionBegin;
8973a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
898a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
8993a40ed3dSBarry Smith   }
9003a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
9019b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
902c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
903d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
9049b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
905d64ed03dSBarry Smith   PetscStackPop;
9069b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
9076d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
90877c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
90977c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9103a40ed3dSBarry Smith   PetscFunctionReturn(0);
9119b94acceSBarry Smith }
9129b94acceSBarry Smith 
9135615d1e5SSatish Balay #undef __FUNC__
9145615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
91562fef451SLois Curfman McInnes /*@
91662fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
91762fef451SLois Curfman McInnes    set with SNESSetHessian().
91862fef451SLois Curfman McInnes 
919*c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
920*c7afd0dbSLois Curfman McInnes 
92162fef451SLois Curfman McInnes    Input Parameters:
922*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
923*c7afd0dbSLois Curfman McInnes -  x - input vector
92462fef451SLois Curfman McInnes 
92562fef451SLois Curfman McInnes    Output Parameters:
926*c7afd0dbSLois Curfman McInnes +  A - Hessian matrix
92762fef451SLois Curfman McInnes .  B - optional preconditioning matrix
928*c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
929fee21e36SBarry Smith 
93062fef451SLois Curfman McInnes    Notes:
93162fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
93262fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
93362fef451SLois Curfman McInnes 
934dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
935dc5a77f8SLois Curfman McInnes    flag parameter.
93662fef451SLois Curfman McInnes 
93762fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
93862fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
93962fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
94062fef451SLois Curfman McInnes 
94162fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
94262fef451SLois Curfman McInnes 
943a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
944a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
94562fef451SLois Curfman McInnes @*/
94662fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
9479b94acceSBarry Smith {
9489b94acceSBarry Smith   int    ierr;
9493a40ed3dSBarry Smith 
9503a40ed3dSBarry Smith   PetscFunctionBegin;
9513a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
952a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9533a40ed3dSBarry Smith   }
9543a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
95562fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9560f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
957d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
95862fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
959d64ed03dSBarry Smith   PetscStackPop;
96062fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9610f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
96277c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
96377c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9643a40ed3dSBarry Smith   PetscFunctionReturn(0);
9659b94acceSBarry Smith }
9669b94acceSBarry Smith 
9675615d1e5SSatish Balay #undef __FUNC__
968d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian"
9699b94acceSBarry Smith /*@C
9709b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
971044dda88SLois Curfman McInnes    location to store the matrix.
9729b94acceSBarry Smith 
973*c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
974*c7afd0dbSLois Curfman McInnes 
9759b94acceSBarry Smith    Input Parameters:
976*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
9779b94acceSBarry Smith .  A - Jacobian matrix
9789b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
9799b94acceSBarry Smith .  func - Jacobian evaluation routine
980*c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
9812cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9829b94acceSBarry Smith 
9839b94acceSBarry Smith    Calling sequence of func:
984*c7afd0dbSLois Curfman McInnes $  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9859b94acceSBarry Smith 
986*c7afd0dbSLois Curfman McInnes +  x - input vector
9879b94acceSBarry Smith .  A - Jacobian matrix
9889b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
989ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
990ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
991*c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
9929b94acceSBarry Smith 
9939b94acceSBarry Smith    Notes:
994dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9952cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
996ac21db08SLois Curfman McInnes 
997ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9989b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9999b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10009b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10019b94acceSBarry Smith    throughout the global iterations.
10029b94acceSBarry Smith 
10039b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
10049b94acceSBarry Smith 
1005ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
10069b94acceSBarry Smith @*/
10079b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10089b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10099b94acceSBarry Smith {
10103a40ed3dSBarry Smith   PetscFunctionBegin;
101177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1012a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1013a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1014a8c6a408SBarry Smith   }
10159b94acceSBarry Smith   snes->computejacobian = func;
10169b94acceSBarry Smith   snes->jacP            = ctx;
10179b94acceSBarry Smith   snes->jacobian        = A;
10189b94acceSBarry Smith   snes->jacobian_pre    = B;
10193a40ed3dSBarry Smith   PetscFunctionReturn(0);
10209b94acceSBarry Smith }
102162fef451SLois Curfman McInnes 
10225615d1e5SSatish Balay #undef __FUNC__
1023d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian"
1024b4fd4287SBarry Smith /*@
1025b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1026b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1027b4fd4287SBarry Smith 
1028*c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
1029*c7afd0dbSLois Curfman McInnes 
1030b4fd4287SBarry Smith    Input Parameter:
1031b4fd4287SBarry Smith .  snes - the nonlinear solver context
1032b4fd4287SBarry Smith 
1033b4fd4287SBarry Smith    Output Parameters:
1034*c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
1035b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
1036*c7afd0dbSLois Curfman McInnes -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1037fee21e36SBarry Smith 
1038b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1039b4fd4287SBarry Smith @*/
1040b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
1041b4fd4287SBarry Smith {
10423a40ed3dSBarry Smith   PetscFunctionBegin;
10433a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1044a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
10453a40ed3dSBarry Smith   }
1046b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
1047b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
1048b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
10493a40ed3dSBarry Smith   PetscFunctionReturn(0);
1050b4fd4287SBarry Smith }
1051b4fd4287SBarry Smith 
10525615d1e5SSatish Balay #undef __FUNC__
1053d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian"
10549b94acceSBarry Smith /*@C
10559b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1056044dda88SLois Curfman McInnes    location to store the matrix.
10579b94acceSBarry Smith 
1058*c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1059*c7afd0dbSLois Curfman McInnes 
10609b94acceSBarry Smith    Input Parameters:
1061*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
10629b94acceSBarry Smith .  A - Hessian matrix
10639b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
10649b94acceSBarry Smith .  func - Jacobian evaluation routine
1065*c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
1066313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10679b94acceSBarry Smith 
10689b94acceSBarry Smith    Calling sequence of func:
1069*c7afd0dbSLois Curfman McInnes $  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10709b94acceSBarry Smith 
1071*c7afd0dbSLois Curfman McInnes +  x - input vector
10729b94acceSBarry Smith .  A - Hessian matrix
10739b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1074ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1075ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
1076*c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Hessian context
10779b94acceSBarry Smith 
10789b94acceSBarry Smith    Notes:
1079dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10802cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1081ac21db08SLois Curfman McInnes 
10829b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10839b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10849b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10859b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10869b94acceSBarry Smith    throughout the global iterations.
10879b94acceSBarry Smith 
10889b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10899b94acceSBarry Smith 
1090ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10919b94acceSBarry Smith @*/
10929b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10939b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10949b94acceSBarry Smith {
10953a40ed3dSBarry Smith   PetscFunctionBegin;
109677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1097d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1098a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1099d4bb536fSBarry Smith   }
11009b94acceSBarry Smith   snes->computejacobian = func;
11019b94acceSBarry Smith   snes->jacP            = ctx;
11029b94acceSBarry Smith   snes->jacobian        = A;
11039b94acceSBarry Smith   snes->jacobian_pre    = B;
11043a40ed3dSBarry Smith   PetscFunctionReturn(0);
11059b94acceSBarry Smith }
11069b94acceSBarry Smith 
11075615d1e5SSatish Balay #undef __FUNC__
1108d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian"
110962fef451SLois Curfman McInnes /*@
111062fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
111162fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
111262fef451SLois Curfman McInnes 
1113*c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object is parallel if SNES object is parallel
1114*c7afd0dbSLois Curfman McInnes 
111562fef451SLois Curfman McInnes    Input Parameter:
111662fef451SLois Curfman McInnes .  snes - the nonlinear solver context
111762fef451SLois Curfman McInnes 
111862fef451SLois Curfman McInnes    Output Parameters:
1119*c7afd0dbSLois Curfman McInnes +  A - location to stash Hessian matrix (or PETSC_NULL)
112062fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
1121*c7afd0dbSLois Curfman McInnes -  ctx - location to stash Hessian ctx (or PETSC_NULL)
1122fee21e36SBarry Smith 
112362fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
1124*c7afd0dbSLois Curfman McInnes 
1125*c7afd0dbSLois Curfman McInnes .keywords: SNES, get, Hessian
112662fef451SLois Curfman McInnes @*/
112762fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx)
112862fef451SLois Curfman McInnes {
11293a40ed3dSBarry Smith   PetscFunctionBegin;
11303a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1131a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
11323a40ed3dSBarry Smith   }
113362fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
113462fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
113562fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
11363a40ed3dSBarry Smith   PetscFunctionReturn(0);
113762fef451SLois Curfman McInnes }
113862fef451SLois Curfman McInnes 
11399b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
11409b94acceSBarry Smith 
11415615d1e5SSatish Balay #undef __FUNC__
11425615d1e5SSatish Balay #define __FUNC__ "SNESSetUp"
11439b94acceSBarry Smith /*@
11449b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1145272ac6f2SLois Curfman McInnes    of a nonlinear solver.
11469b94acceSBarry Smith 
1147fee21e36SBarry Smith    Collective on SNES
1148fee21e36SBarry Smith 
1149*c7afd0dbSLois Curfman McInnes    Input Parameters:
1150*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1151*c7afd0dbSLois Curfman McInnes -  x - the solution vector
1152*c7afd0dbSLois Curfman McInnes 
1153272ac6f2SLois Curfman McInnes    Notes:
1154272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1155272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1156272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1157272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1158272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1159272ac6f2SLois Curfman McInnes 
11609b94acceSBarry Smith .keywords: SNES, nonlinear, setup
11619b94acceSBarry Smith 
11629b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
11639b94acceSBarry Smith @*/
11648ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
11659b94acceSBarry Smith {
1166272ac6f2SLois Curfman McInnes   int ierr, flg;
11673a40ed3dSBarry Smith 
11683a40ed3dSBarry Smith   PetscFunctionBegin;
116977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
117077c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
11718ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
11729b94acceSBarry Smith 
1173c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);  CHKERRQ(ierr);
1174c1f60f51SBarry Smith   /*
1175c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1176dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1177c1f60f51SBarry Smith   */
1178c1f60f51SBarry Smith   if (flg) {
1179c1f60f51SBarry Smith     Mat J;
1180c1f60f51SBarry Smith     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1181c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1182c1f60f51SBarry Smith     snes->mfshell = J;
1183c1f60f51SBarry Smith     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1184c1f60f51SBarry Smith       snes->jacobian = J;
1185c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1186d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1187c1f60f51SBarry Smith       snes->jacobian = J;
1188c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1189d4bb536fSBarry Smith     } else {
1190a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option");
1191d4bb536fSBarry Smith     }
1192c1f60f51SBarry Smith   }
1193272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);  CHKERRQ(ierr);
1194c1f60f51SBarry Smith   /*
1195dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1196c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1197c1f60f51SBarry Smith    */
119831615d04SLois Curfman McInnes   if (flg) {
1199272ac6f2SLois Curfman McInnes     Mat J;
1200272ac6f2SLois Curfman McInnes     ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1201272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1202272ac6f2SLois Curfman McInnes     snes->mfshell = J;
120383e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
120483e56afdSLois Curfman McInnes       ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
120594a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1206d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
120783e56afdSLois Curfman McInnes       ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr);
120894a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1209d4bb536fSBarry Smith     } else {
1210a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option");
1211d4bb536fSBarry Smith     }
1212272ac6f2SLois Curfman McInnes   }
12139b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1214a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1215a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1216a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first");
1217a8c6a408SBarry Smith     if (snes->vec_func == snes->vec_sol) {
1218a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector");
1219a8c6a408SBarry Smith     }
1220a703fe33SLois Curfman McInnes 
1221a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
122282bf6240SBarry Smith     if (snes->ksp_ewconv && PetscStrcmp(snes->type_name,SNES_EQ_TR)) {
1223a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1224a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr);
1225a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr);
1226a703fe33SLois Curfman McInnes       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,
1227a703fe33SLois Curfman McInnes              (void *)snes); CHKERRQ(ierr);
1228a703fe33SLois Curfman McInnes     }
12299b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1230a8c6a408SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1231a8c6a408SBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1232a8c6a408SBarry Smith     if (!snes->computeumfunction) {
1233a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first");
1234a8c6a408SBarry Smith     }
1235a8c6a408SBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first");
1236d4bb536fSBarry Smith   } else {
1237a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
1238d4bb536fSBarry Smith   }
1239a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);}
124082bf6240SBarry Smith   snes->setupcalled = 1;
12413a40ed3dSBarry Smith   PetscFunctionReturn(0);
12429b94acceSBarry Smith }
12439b94acceSBarry Smith 
12445615d1e5SSatish Balay #undef __FUNC__
1245d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy"
12469b94acceSBarry Smith /*@C
12479b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
12489b94acceSBarry Smith    with SNESCreate().
12499b94acceSBarry Smith 
1250*c7afd0dbSLois Curfman McInnes    Collective on SNES
1251*c7afd0dbSLois Curfman McInnes 
12529b94acceSBarry Smith    Input Parameter:
12539b94acceSBarry Smith .  snes - the SNES context
12549b94acceSBarry Smith 
12559b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
12569b94acceSBarry Smith 
125763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
12589b94acceSBarry Smith @*/
12599b94acceSBarry Smith int SNESDestroy(SNES snes)
12609b94acceSBarry Smith {
12619b94acceSBarry Smith   int ierr;
12623a40ed3dSBarry Smith 
12633a40ed3dSBarry Smith   PetscFunctionBegin;
126477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12653a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1266d4bb536fSBarry Smith 
1267e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes); CHKERRQ(ierr);}
12680452661fSBarry Smith   if (snes->kspconvctx) PetscFree(snes->kspconvctx);
1269522c5e43SBarry Smith   if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);}
12709b94acceSBarry Smith   ierr = SLESDestroy(snes->sles); CHKERRQ(ierr);
1271522c5e43SBarry Smith   if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);}
1272522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
12739b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
12740452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
12753a40ed3dSBarry Smith   PetscFunctionReturn(0);
12769b94acceSBarry Smith }
12779b94acceSBarry Smith 
12789b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12799b94acceSBarry Smith 
12805615d1e5SSatish Balay #undef __FUNC__
12815615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances"
12829b94acceSBarry Smith /*@
1283d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12849b94acceSBarry Smith 
1285*c7afd0dbSLois Curfman McInnes    Collective on SNES
1286*c7afd0dbSLois Curfman McInnes 
12879b94acceSBarry Smith    Input Parameters:
1288*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
128933174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
129033174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
129133174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
129233174efeSLois Curfman McInnes            of the change in the solution between steps
129333174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1294*c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1295fee21e36SBarry Smith 
129633174efeSLois Curfman McInnes    Options Database Keys:
1297*c7afd0dbSLois Curfman McInnes +    -snes_atol <atol> - Sets atol
1298*c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1299*c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1300*c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1301*c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
13029b94acceSBarry Smith 
1303d7a720efSLois Curfman McInnes    Notes:
13049b94acceSBarry Smith    The default maximum number of iterations is 50.
13059b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
13069b94acceSBarry Smith 
130733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
13089b94acceSBarry Smith 
1309d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
13109b94acceSBarry Smith @*/
1311d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
13129b94acceSBarry Smith {
13133a40ed3dSBarry Smith   PetscFunctionBegin;
131477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1315d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1316d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1317d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1318d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1319d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
13203a40ed3dSBarry Smith   PetscFunctionReturn(0);
13219b94acceSBarry Smith }
13229b94acceSBarry Smith 
13235615d1e5SSatish Balay #undef __FUNC__
13245615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances"
13259b94acceSBarry Smith /*@
132633174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
132733174efeSLois Curfman McInnes 
1328*c7afd0dbSLois Curfman McInnes    Not Collective
1329*c7afd0dbSLois Curfman McInnes 
133033174efeSLois Curfman McInnes    Input Parameters:
1331*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
133233174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
133333174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
133433174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
133533174efeSLois Curfman McInnes            of the change in the solution between steps
133633174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1337*c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1338fee21e36SBarry Smith 
133933174efeSLois Curfman McInnes    Notes:
134033174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
134133174efeSLois Curfman McInnes 
134233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
134333174efeSLois Curfman McInnes 
134433174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
134533174efeSLois Curfman McInnes @*/
134633174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
134733174efeSLois Curfman McInnes {
13483a40ed3dSBarry Smith   PetscFunctionBegin;
134933174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
135033174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
135133174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
135233174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
135333174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
135433174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
13553a40ed3dSBarry Smith   PetscFunctionReturn(0);
135633174efeSLois Curfman McInnes }
135733174efeSLois Curfman McInnes 
13585615d1e5SSatish Balay #undef __FUNC__
13595615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance"
136033174efeSLois Curfman McInnes /*@
13619b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
13629b94acceSBarry Smith 
1363fee21e36SBarry Smith    Collective on SNES
1364fee21e36SBarry Smith 
1365*c7afd0dbSLois Curfman McInnes    Input Parameters:
1366*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1367*c7afd0dbSLois Curfman McInnes -  tol - tolerance
1368*c7afd0dbSLois Curfman McInnes 
13699b94acceSBarry Smith    Options Database Key:
1370*c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
13719b94acceSBarry Smith 
13729b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
13739b94acceSBarry Smith 
1374d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
13759b94acceSBarry Smith @*/
13769b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol)
13779b94acceSBarry Smith {
13783a40ed3dSBarry Smith   PetscFunctionBegin;
137977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13809b94acceSBarry Smith   snes->deltatol = tol;
13813a40ed3dSBarry Smith   PetscFunctionReturn(0);
13829b94acceSBarry Smith }
13839b94acceSBarry Smith 
13845615d1e5SSatish Balay #undef __FUNC__
13855615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
13869b94acceSBarry Smith /*@
138777c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
13889b94acceSBarry Smith    for unconstrained minimization solvers.
13899b94acceSBarry Smith 
1390fee21e36SBarry Smith    Collective on SNES
1391fee21e36SBarry Smith 
1392*c7afd0dbSLois Curfman McInnes    Input Parameters:
1393*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1394*c7afd0dbSLois Curfman McInnes -  ftol - minimum function tolerance
1395*c7afd0dbSLois Curfman McInnes 
13969b94acceSBarry Smith    Options Database Key:
1397*c7afd0dbSLois Curfman McInnes .  -snes_fmin <ftol> - Sets ftol
13989b94acceSBarry Smith 
13999b94acceSBarry Smith    Note:
140077c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
14019b94acceSBarry Smith    methods only.
14029b94acceSBarry Smith 
14039b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
14049b94acceSBarry Smith 
1405d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
14069b94acceSBarry Smith @*/
140777c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
14089b94acceSBarry Smith {
14093a40ed3dSBarry Smith   PetscFunctionBegin;
141077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
14119b94acceSBarry Smith   snes->fmin = ftol;
14123a40ed3dSBarry Smith   PetscFunctionReturn(0);
14139b94acceSBarry Smith }
14149b94acceSBarry Smith 
14159b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
14169b94acceSBarry Smith 
14175615d1e5SSatish Balay #undef __FUNC__
1418d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor"
14199b94acceSBarry Smith /*@C
14205cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
14219b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
14229b94acceSBarry Smith    progress.
14239b94acceSBarry Smith 
1424fee21e36SBarry Smith    Collective on SNES
1425fee21e36SBarry Smith 
1426*c7afd0dbSLois Curfman McInnes    Input Parameters:
1427*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1428*c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1429*c7afd0dbSLois Curfman McInnes -  mctx - [optional] user-defined context for private data for the
1430*c7afd0dbSLois Curfman McInnes           monitor routine (may be PETSC_NULL)
14319b94acceSBarry Smith 
1432*c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1433*c7afd0dbSLois Curfman McInnes $  int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
1434*c7afd0dbSLois Curfman McInnes 
1435*c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1436*c7afd0dbSLois Curfman McInnes .    its - iteration number
1437*c7afd0dbSLois Curfman McInnes .    mctx - [optional] monitoring context
1438*c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
1439*c7afd0dbSLois Curfman McInnes             for SNES_NONLINEAR_EQUATIONS methods
1440*c7afd0dbSLois Curfman McInnes -    norm - 2-norm gradient value (may be estimated)
1441*c7afd0dbSLois Curfman McInnes             for SNES_UNCONSTRAINED_MINIMIZATION methods
14429b94acceSBarry Smith 
14439665c990SLois Curfman McInnes    Options Database Keys:
1444*c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1445*c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1446*c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1447*c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1448*c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1449*c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1450*c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1451*c7afd0dbSLois Curfman McInnes                             the options database.
14529665c990SLois Curfman McInnes 
1453639f9d9dSBarry Smith    Notes:
14546bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
14556bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
14566bc08f3fSLois Curfman McInnes    order in which they were set.
1457639f9d9dSBarry Smith 
14589b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
14599b94acceSBarry Smith 
14605cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
14619b94acceSBarry Smith @*/
146274679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
14639b94acceSBarry Smith {
14643a40ed3dSBarry Smith   PetscFunctionBegin;
1465639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1466a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
1467639f9d9dSBarry Smith   }
1468639f9d9dSBarry Smith 
1469639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1470639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
14713a40ed3dSBarry Smith   PetscFunctionReturn(0);
14729b94acceSBarry Smith }
14739b94acceSBarry Smith 
14745615d1e5SSatish Balay #undef __FUNC__
14755cd90555SBarry Smith #define __FUNC__ "SNESClearMonitor"
14765cd90555SBarry Smith /*@C
14775cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
14785cd90555SBarry Smith 
1479*c7afd0dbSLois Curfman McInnes    Collective on SNES
1480*c7afd0dbSLois Curfman McInnes 
14815cd90555SBarry Smith    Input Parameters:
14825cd90555SBarry Smith .  snes - the SNES context
14835cd90555SBarry Smith 
14845cd90555SBarry Smith    Options Database:
1485*c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1486*c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1487*c7afd0dbSLois Curfman McInnes     set via the options database
14885cd90555SBarry Smith 
14895cd90555SBarry Smith    Notes:
14905cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
14915cd90555SBarry Smith 
14925cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
14935cd90555SBarry Smith 
14945cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
14955cd90555SBarry Smith @*/
14965cd90555SBarry Smith int SNESClearMonitor( SNES snes )
14975cd90555SBarry Smith {
14985cd90555SBarry Smith   PetscFunctionBegin;
14995cd90555SBarry Smith   snes->numbermonitors = 0;
15005cd90555SBarry Smith   PetscFunctionReturn(0);
15015cd90555SBarry Smith }
15025cd90555SBarry Smith 
15035cd90555SBarry Smith #undef __FUNC__
1504d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest"
15059b94acceSBarry Smith /*@C
15069b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
15079b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
15089b94acceSBarry Smith 
1509fee21e36SBarry Smith    Collective on SNES
1510fee21e36SBarry Smith 
1511*c7afd0dbSLois Curfman McInnes    Input Parameters:
1512*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1513*c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1514*c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1515*c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
15169b94acceSBarry Smith 
1517*c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1518*c7afd0dbSLois Curfman McInnes $  int func (SNES snes,double xnorm,double gnorm,double f,void *cctx)
1519*c7afd0dbSLois Curfman McInnes 
1520*c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1521*c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1522*c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
1523*c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods)
1524*c7afd0dbSLois Curfman McInnes .    f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods)
1525*c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods)
1526*c7afd0dbSLois Curfman McInnes -    f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods)
15279b94acceSBarry Smith 
15289b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
15299b94acceSBarry Smith 
153040191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
153140191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
15329b94acceSBarry Smith @*/
153374679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
15349b94acceSBarry Smith {
15353a40ed3dSBarry Smith   PetscFunctionBegin;
15369b94acceSBarry Smith   (snes)->converged = func;
15379b94acceSBarry Smith   (snes)->cnvP      = cctx;
15383a40ed3dSBarry Smith   PetscFunctionReturn(0);
15399b94acceSBarry Smith }
15409b94acceSBarry Smith 
15415615d1e5SSatish Balay #undef __FUNC__
1542d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory"
1543c9005455SLois Curfman McInnes /*@
1544c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1545c9005455SLois Curfman McInnes 
1546fee21e36SBarry Smith    Collective on SNES
1547fee21e36SBarry Smith 
1548*c7afd0dbSLois Curfman McInnes    Input Parameters:
1549*c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1550*c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1551*c7afd0dbSLois Curfman McInnes -  na  - size of a
1552*c7afd0dbSLois Curfman McInnes 
1553c9005455SLois Curfman McInnes    Notes:
1554c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1555c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1556c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1557c9005455SLois Curfman McInnes    at each step.
1558c9005455SLois Curfman McInnes 
1559c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1560c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1561c9005455SLois Curfman McInnes    during the section of code that is being timed.
1562c9005455SLois Curfman McInnes 
1563c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1564c9005455SLois Curfman McInnes @*/
1565c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1566c9005455SLois Curfman McInnes {
15673a40ed3dSBarry Smith   PetscFunctionBegin;
1568c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1569c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1570c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1571c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
15723a40ed3dSBarry Smith   PetscFunctionReturn(0);
1573c9005455SLois Curfman McInnes }
1574c9005455SLois Curfman McInnes 
1575c9005455SLois Curfman McInnes #undef __FUNC__
15765615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
15779b94acceSBarry Smith /*
15789b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
15799b94acceSBarry Smith    positive parameter delta.
15809b94acceSBarry Smith 
15819b94acceSBarry Smith     Input Parameters:
1582*c7afd0dbSLois Curfman McInnes +   snes - the SNES context
15839b94acceSBarry Smith .   y - approximate solution of linear system
15849b94acceSBarry Smith .   fnorm - 2-norm of current function
1585*c7afd0dbSLois Curfman McInnes -   delta - trust region size
15869b94acceSBarry Smith 
15879b94acceSBarry Smith     Output Parameters:
1588*c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
15899b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
15909b94acceSBarry Smith     region, and exceeds zero otherwise.
1591*c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
15929b94acceSBarry Smith 
15939b94acceSBarry Smith     Note:
159440191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
15959b94acceSBarry Smith     is set to be the maximum allowable step size.
15969b94acceSBarry Smith 
15979b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
15989b94acceSBarry Smith */
15999b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
16009b94acceSBarry Smith                           double *gpnorm,double *ynorm)
16019b94acceSBarry Smith {
16029b94acceSBarry Smith   double norm;
16039b94acceSBarry Smith   Scalar cnorm;
16043a40ed3dSBarry Smith   int    ierr;
16053a40ed3dSBarry Smith 
16063a40ed3dSBarry Smith   PetscFunctionBegin;
16073a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
16089b94acceSBarry Smith   if (norm > *delta) {
16099b94acceSBarry Smith      norm = *delta/norm;
16109b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
16119b94acceSBarry Smith      cnorm = norm;
16129b94acceSBarry Smith      VecScale( &cnorm, y );
16139b94acceSBarry Smith      *ynorm = *delta;
16149b94acceSBarry Smith   } else {
16159b94acceSBarry Smith      *gpnorm = 0.0;
16169b94acceSBarry Smith      *ynorm = norm;
16179b94acceSBarry Smith   }
16183a40ed3dSBarry Smith   PetscFunctionReturn(0);
16199b94acceSBarry Smith }
16209b94acceSBarry Smith 
16215615d1e5SSatish Balay #undef __FUNC__
16225615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
16239b94acceSBarry Smith /*@
16249b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
162563a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
16269b94acceSBarry Smith 
1627*c7afd0dbSLois Curfman McInnes    Collective on SNES
1628*c7afd0dbSLois Curfman McInnes 
1629b2002411SLois Curfman McInnes    Input Parameters:
1630*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1631*c7afd0dbSLois Curfman McInnes -  x - the solution vector
16329b94acceSBarry Smith 
16339b94acceSBarry Smith    Output Parameter:
1634b2002411SLois Curfman McInnes .  its - number of iterations until termination
16359b94acceSBarry Smith 
1636b2002411SLois Curfman McInnes    Notes:
16378ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
16388ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
16398ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
16408ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
16418ddd3da0SLois Curfman McInnes 
16429b94acceSBarry Smith .keywords: SNES, nonlinear, solve
16439b94acceSBarry Smith 
164463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
16459b94acceSBarry Smith @*/
16468ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
16479b94acceSBarry Smith {
16483c7409f5SSatish Balay   int ierr, flg;
1649052efed2SBarry Smith 
16503a40ed3dSBarry Smith   PetscFunctionBegin;
165177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
165274679c65SBarry Smith   PetscValidIntPointer(its);
165382bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1654c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
16559b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1656c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
16579b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
16589b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
16593c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
16606d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
16613a40ed3dSBarry Smith   PetscFunctionReturn(0);
16629b94acceSBarry Smith }
16639b94acceSBarry Smith 
16649b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
16659b94acceSBarry Smith 
16665615d1e5SSatish Balay #undef __FUNC__
16675615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
166882bf6240SBarry Smith /*@C
16694b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
16709b94acceSBarry Smith 
1671fee21e36SBarry Smith    Collective on SNES
1672fee21e36SBarry Smith 
1673*c7afd0dbSLois Curfman McInnes    Input Parameters:
1674*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1675*c7afd0dbSLois Curfman McInnes -  method - a known method
1676*c7afd0dbSLois Curfman McInnes 
1677*c7afd0dbSLois Curfman McInnes    Options Database Key:
1678*c7afd0dbSLois Curfman McInnes .  -snes_type <method> - Sets the method; use -help for a list
1679*c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1680ae12b187SLois Curfman McInnes 
16819b94acceSBarry Smith    Notes:
16829b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
1683*c7afd0dbSLois Curfman McInnes .    SNES_EQ_LS - Newton's method with line search
1684*c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
1685*c7afd0dbSLois Curfman McInnes .    SNES_EQ_TR - Newton's method with trust region
1686*c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
1687*c7afd0dbSLois Curfman McInnes .    SNES_UM_TR - Newton's method with trust region
1688*c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
1689*c7afd0dbSLois Curfman McInnes .    SNES_UM_LS - Newton's method with line search
1690*c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
16919b94acceSBarry Smith 
1692ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1693ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1694ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1695ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1696ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1697ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1698ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1699ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1700ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1701ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1702a703fe33SLois Curfman McInnes 
1703f536c699SSatish Balay .keywords: SNES, set, method
17049b94acceSBarry Smith @*/
17054b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
17069b94acceSBarry Smith {
170784cb2905SBarry Smith   int ierr;
17089b94acceSBarry Smith   int (*r)(SNES);
17093a40ed3dSBarry Smith 
17103a40ed3dSBarry Smith   PetscFunctionBegin;
171177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
171282bf6240SBarry Smith 
171382bf6240SBarry Smith   if (!PetscStrcmp(snes->type_name,method)) PetscFunctionReturn(0);
171482bf6240SBarry Smith 
171582bf6240SBarry Smith   if (snes->setupcalled) {
1716e1311b90SBarry Smith     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
171782bf6240SBarry Smith     snes->data = 0;
171882bf6240SBarry Smith   }
17197d1a2b2bSBarry Smith 
17209b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
172182bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);}
172282bf6240SBarry Smith 
1723ecf371e4SBarry Smith   ierr =  DLRegisterFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr);
172482bf6240SBarry Smith 
17250452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
172682bf6240SBarry Smith   snes->data = 0;
17273a40ed3dSBarry Smith   ierr = (*r)(snes); CHKERRQ(ierr);
172882bf6240SBarry Smith 
172982bf6240SBarry Smith   if (snes->type_name) PetscFree(snes->type_name);
173082bf6240SBarry Smith   snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name);
173182bf6240SBarry Smith   PetscStrcpy(snes->type_name,method);
173282bf6240SBarry Smith   snes->set_method_called = 1;
173382bf6240SBarry Smith 
17343a40ed3dSBarry Smith   PetscFunctionReturn(0);
17359b94acceSBarry Smith }
17369b94acceSBarry Smith 
1737a847f771SSatish Balay 
17389b94acceSBarry Smith /* --------------------------------------------------------------------- */
17395615d1e5SSatish Balay #undef __FUNC__
1740d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy"
17419b94acceSBarry Smith /*@C
17429b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
17439b94acceSBarry Smith    registered by SNESRegister().
17449b94acceSBarry Smith 
1745fee21e36SBarry Smith    Not Collective
1746fee21e36SBarry Smith 
17479b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
17489b94acceSBarry Smith 
17499b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
17509b94acceSBarry Smith @*/
1751cf256101SBarry Smith int SNESRegisterDestroy(void)
17529b94acceSBarry Smith {
175382bf6240SBarry Smith   int ierr;
175482bf6240SBarry Smith 
17553a40ed3dSBarry Smith   PetscFunctionBegin;
175682bf6240SBarry Smith   if (SNESList) {
175782bf6240SBarry Smith     ierr = DLRegisterDestroy( SNESList );CHKERRQ(ierr);
175882bf6240SBarry Smith     SNESList = 0;
17599b94acceSBarry Smith   }
176084cb2905SBarry Smith   SNESRegisterAllCalled = 0;
17613a40ed3dSBarry Smith   PetscFunctionReturn(0);
17629b94acceSBarry Smith }
17639b94acceSBarry Smith 
17645615d1e5SSatish Balay #undef __FUNC__
1765d4bb536fSBarry Smith #define __FUNC__ "SNESGetType"
17669b94acceSBarry Smith /*@C
17679a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
17689b94acceSBarry Smith 
1769*c7afd0dbSLois Curfman McInnes    Not Collective
1770*c7afd0dbSLois Curfman McInnes 
17719b94acceSBarry Smith    Input Parameter:
17724b0e389bSBarry Smith .  snes - nonlinear solver context
17739b94acceSBarry Smith 
17749b94acceSBarry Smith    Output Parameter:
177582bf6240SBarry Smith .  method - SNES method (a charactor string)
17769b94acceSBarry Smith 
17779b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
17789b94acceSBarry Smith @*/
177982bf6240SBarry Smith int SNESGetType(SNES snes, SNESType *method)
17809b94acceSBarry Smith {
17813a40ed3dSBarry Smith   PetscFunctionBegin;
178282bf6240SBarry Smith   *method = snes->type_name;
17833a40ed3dSBarry Smith   PetscFunctionReturn(0);
17849b94acceSBarry Smith }
17859b94acceSBarry Smith 
17865615d1e5SSatish Balay #undef __FUNC__
1787d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution"
17889b94acceSBarry Smith /*@C
17899b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
17909b94acceSBarry Smith    stored.
17919b94acceSBarry Smith 
1792*c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1793*c7afd0dbSLois Curfman McInnes 
17949b94acceSBarry Smith    Input Parameter:
17959b94acceSBarry Smith .  snes - the SNES context
17969b94acceSBarry Smith 
17979b94acceSBarry Smith    Output Parameter:
17989b94acceSBarry Smith .  x - the solution
17999b94acceSBarry Smith 
18009b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
18019b94acceSBarry Smith 
1802a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
18039b94acceSBarry Smith @*/
18049b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
18059b94acceSBarry Smith {
18063a40ed3dSBarry Smith   PetscFunctionBegin;
180777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18089b94acceSBarry Smith   *x = snes->vec_sol_always;
18093a40ed3dSBarry Smith   PetscFunctionReturn(0);
18109b94acceSBarry Smith }
18119b94acceSBarry Smith 
18125615d1e5SSatish Balay #undef __FUNC__
1813d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate"
18149b94acceSBarry Smith /*@C
18159b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
18169b94acceSBarry Smith    stored.
18179b94acceSBarry Smith 
1818*c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1819*c7afd0dbSLois Curfman McInnes 
18209b94acceSBarry Smith    Input Parameter:
18219b94acceSBarry Smith .  snes - the SNES context
18229b94acceSBarry Smith 
18239b94acceSBarry Smith    Output Parameter:
18249b94acceSBarry Smith .  x - the solution update
18259b94acceSBarry Smith 
18269b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
18279b94acceSBarry Smith 
18289b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
18299b94acceSBarry Smith @*/
18309b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
18319b94acceSBarry Smith {
18323a40ed3dSBarry Smith   PetscFunctionBegin;
183377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18349b94acceSBarry Smith   *x = snes->vec_sol_update_always;
18353a40ed3dSBarry Smith   PetscFunctionReturn(0);
18369b94acceSBarry Smith }
18379b94acceSBarry Smith 
18385615d1e5SSatish Balay #undef __FUNC__
1839d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction"
18409b94acceSBarry Smith /*@C
18413638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
18429b94acceSBarry Smith 
1843*c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1844*c7afd0dbSLois Curfman McInnes 
18459b94acceSBarry Smith    Input Parameter:
18469b94acceSBarry Smith .  snes - the SNES context
18479b94acceSBarry Smith 
18489b94acceSBarry Smith    Output Parameter:
18493638b69dSLois Curfman McInnes .  r - the function
18509b94acceSBarry Smith 
18519b94acceSBarry Smith    Notes:
18529b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
18539b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
18549b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
18559b94acceSBarry Smith 
1856a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
18579b94acceSBarry Smith 
185831615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
185931615d04SLois Curfman McInnes           SNESGetGradient()
18609b94acceSBarry Smith @*/
18619b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
18629b94acceSBarry Smith {
18633a40ed3dSBarry Smith   PetscFunctionBegin;
186477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1865a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1866a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1867a8c6a408SBarry Smith   }
18689b94acceSBarry Smith   *r = snes->vec_func_always;
18693a40ed3dSBarry Smith   PetscFunctionReturn(0);
18709b94acceSBarry Smith }
18719b94acceSBarry Smith 
18725615d1e5SSatish Balay #undef __FUNC__
1873d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient"
18749b94acceSBarry Smith /*@C
18753638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
18769b94acceSBarry Smith 
1877*c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1878*c7afd0dbSLois Curfman McInnes 
18799b94acceSBarry Smith    Input Parameter:
18809b94acceSBarry Smith .  snes - the SNES context
18819b94acceSBarry Smith 
18829b94acceSBarry Smith    Output Parameter:
18839b94acceSBarry Smith .  r - the gradient
18849b94acceSBarry Smith 
18859b94acceSBarry Smith    Notes:
18869b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
18879b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18889b94acceSBarry Smith    SNESGetFunction().
18899b94acceSBarry Smith 
18909b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
18919b94acceSBarry Smith 
189231615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
18939b94acceSBarry Smith @*/
18949b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
18959b94acceSBarry Smith {
18963a40ed3dSBarry Smith   PetscFunctionBegin;
189777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18983a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1899a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
19003a40ed3dSBarry Smith   }
19019b94acceSBarry Smith   *r = snes->vec_func_always;
19023a40ed3dSBarry Smith   PetscFunctionReturn(0);
19039b94acceSBarry Smith }
19049b94acceSBarry Smith 
19055615d1e5SSatish Balay #undef __FUNC__
1906d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction"
19079b94acceSBarry Smith /*@
19089b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
19099b94acceSBarry Smith    unconstrained minimization problems.
19109b94acceSBarry Smith 
1911*c7afd0dbSLois Curfman McInnes    Not Collective
1912*c7afd0dbSLois Curfman McInnes 
19139b94acceSBarry Smith    Input Parameter:
19149b94acceSBarry Smith .  snes - the SNES context
19159b94acceSBarry Smith 
19169b94acceSBarry Smith    Output Parameter:
19179b94acceSBarry Smith .  r - the function
19189b94acceSBarry Smith 
19199b94acceSBarry Smith    Notes:
19209b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
19219b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
19229b94acceSBarry Smith    SNESGetFunction().
19239b94acceSBarry Smith 
19249b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
19259b94acceSBarry Smith 
192631615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
19279b94acceSBarry Smith @*/
19289b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
19299b94acceSBarry Smith {
19303a40ed3dSBarry Smith   PetscFunctionBegin;
193177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
193274679c65SBarry Smith   PetscValidScalarPointer(r);
19333a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1934a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
19353a40ed3dSBarry Smith   }
19369b94acceSBarry Smith   *r = snes->fc;
19373a40ed3dSBarry Smith   PetscFunctionReturn(0);
19389b94acceSBarry Smith }
19399b94acceSBarry Smith 
19405615d1e5SSatish Balay #undef __FUNC__
1941d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix"
19423c7409f5SSatish Balay /*@C
19433c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1944d850072dSLois Curfman McInnes    SNES options in the database.
19453c7409f5SSatish Balay 
1946fee21e36SBarry Smith    Collective on SNES
1947fee21e36SBarry Smith 
1948*c7afd0dbSLois Curfman McInnes    Input Parameter:
1949*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1950*c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
1951*c7afd0dbSLois Curfman McInnes 
1952d850072dSLois Curfman McInnes    Notes:
1953a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1954*c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
1955d850072dSLois Curfman McInnes 
19563c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1957a86d99e1SLois Curfman McInnes 
1958a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
19593c7409f5SSatish Balay @*/
19603c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
19613c7409f5SSatish Balay {
19623c7409f5SSatish Balay   int ierr;
19633c7409f5SSatish Balay 
19643a40ed3dSBarry Smith   PetscFunctionBegin;
196577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1966639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19673c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19683a40ed3dSBarry Smith   PetscFunctionReturn(0);
19693c7409f5SSatish Balay }
19703c7409f5SSatish Balay 
19715615d1e5SSatish Balay #undef __FUNC__
1972d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix"
19733c7409f5SSatish Balay /*@C
1974f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1975d850072dSLois Curfman McInnes    SNES options in the database.
19763c7409f5SSatish Balay 
1977fee21e36SBarry Smith    Collective on SNES
1978fee21e36SBarry Smith 
1979*c7afd0dbSLois Curfman McInnes    Input Parameters:
1980*c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1981*c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
1982*c7afd0dbSLois Curfman McInnes 
1983d850072dSLois Curfman McInnes    Notes:
1984a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1985*c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
1986d850072dSLois Curfman McInnes 
19873c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1988a86d99e1SLois Curfman McInnes 
1989a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
19903c7409f5SSatish Balay @*/
19913c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
19923c7409f5SSatish Balay {
19933c7409f5SSatish Balay   int ierr;
19943c7409f5SSatish Balay 
19953a40ed3dSBarry Smith   PetscFunctionBegin;
199677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1997639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19983c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19993a40ed3dSBarry Smith   PetscFunctionReturn(0);
20003c7409f5SSatish Balay }
20013c7409f5SSatish Balay 
20025615d1e5SSatish Balay #undef __FUNC__
2003d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix"
2004c83590e2SSatish Balay /*@
20053c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20063c7409f5SSatish Balay    SNES options in the database.
20073c7409f5SSatish Balay 
2008*c7afd0dbSLois Curfman McInnes    Not Collective
2009*c7afd0dbSLois Curfman McInnes 
20103c7409f5SSatish Balay    Input Parameter:
20113c7409f5SSatish Balay .  snes - the SNES context
20123c7409f5SSatish Balay 
20133c7409f5SSatish Balay    Output Parameter:
20143c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20153c7409f5SSatish Balay 
20163c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2017a86d99e1SLois Curfman McInnes 
2018a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20193c7409f5SSatish Balay @*/
20203c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
20213c7409f5SSatish Balay {
20223c7409f5SSatish Balay   int ierr;
20233c7409f5SSatish Balay 
20243a40ed3dSBarry Smith   PetscFunctionBegin;
202577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2026639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20273a40ed3dSBarry Smith   PetscFunctionReturn(0);
20283c7409f5SSatish Balay }
20293c7409f5SSatish Balay 
2030ca161407SBarry Smith #undef __FUNC__
2031ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp"
2032ca161407SBarry Smith /*@
2033ca161407SBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2034ca161407SBarry Smith 
2035*c7afd0dbSLois Curfman McInnes    Collective on SNES
2036*c7afd0dbSLois Curfman McInnes 
2037ca161407SBarry Smith    Input Parameter:
2038ca161407SBarry Smith .  snes - the SNES context
2039ca161407SBarry Smith 
2040ca161407SBarry Smith    Options Database Keys:
2041*c7afd0dbSLois Curfman McInnes +  -help - Prints SNES options
2042*c7afd0dbSLois Curfman McInnes -  -h - Prints SNES options
2043ca161407SBarry Smith 
2044ca161407SBarry Smith .keywords: SNES, nonlinear, help
2045ca161407SBarry Smith 
2046ca161407SBarry Smith .seealso: SNESSetFromOptions()
2047ca161407SBarry Smith @*/
2048ca161407SBarry Smith int SNESPrintHelp(SNES snes)
2049ca161407SBarry Smith {
2050ca161407SBarry Smith   char                p[64];
2051ca161407SBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2052ca161407SBarry Smith   int                 ierr;
2053ca161407SBarry Smith 
2054ca161407SBarry Smith   PetscFunctionBegin;
2055ca161407SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2056ca161407SBarry Smith 
2057ca161407SBarry Smith   PetscStrcpy(p,"-");
2058ca161407SBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2059ca161407SBarry Smith 
2060ca161407SBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2061ca161407SBarry Smith 
206276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");
206382bf6240SBarry Smith   ierr = DLRegisterPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr);
206476be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
206576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
206676be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
206776be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
206876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
206976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
207076be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
207176be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");
207276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
207376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
2074ca161407SBarry Smith     residual norm at each iteration.\n",p);
207576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
2076ca161407SBarry Smith     residual norm for small residual norms. This is useful to conceal\n\
2077ca161407SBarry Smith     meaningless digits that may be different on different machines.\n",p);
207876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
2079ca161407SBarry Smith   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
208076be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2081ca161407SBarry Smith      " Options for solving systems of nonlinear equations only:\n");
208276be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
208376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
208476be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
208576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
208676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
208776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2088ca161407SBarry Smith      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
208976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2090ca161407SBarry Smith      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
209176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2092ca161407SBarry Smith      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
209376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2094ca161407SBarry Smith      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
209576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2096ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
209776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2098ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
209976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2100ca161407SBarry Smith      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
2101ca161407SBarry Smith   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
210276be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2103ca161407SBarry Smith      " Options for solving unconstrained minimization problems only:\n");
210476be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
210576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
210676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
2107ca161407SBarry Smith   }
210876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
210976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"a particular method\n");
2110ca161407SBarry Smith   if (snes->printhelp) {
2111ca161407SBarry Smith     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2112ca161407SBarry Smith   }
2113ca161407SBarry Smith   PetscFunctionReturn(0);
2114ca161407SBarry Smith }
21153c7409f5SSatish Balay 
2116acb85ae6SSatish Balay /*MC
2117b2002411SLois Curfman McInnes    SNESRegister - Adds a method to the nonlinear solver package.
21189b94acceSBarry Smith 
2119b2002411SLois Curfman McInnes    Synopsis:
2120b2002411SLois Curfman McInnes    SNESRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
21219b94acceSBarry Smith 
2122*c7afd0dbSLois Curfman McInnes    Not collective
2123*c7afd0dbSLois Curfman McInnes 
2124b2002411SLois Curfman McInnes    Input Parameters:
2125*c7afd0dbSLois Curfman McInnes +  name_solver - name of a new user-defined solver
2126b2002411SLois Curfman McInnes .  path - path (either absolute or relative) the library containing this solver
2127b2002411SLois Curfman McInnes .  name_create - name of routine to create method context
2128*c7afd0dbSLois Curfman McInnes -  routine_create - routine to create method context
21299b94acceSBarry Smith 
2130b2002411SLois Curfman McInnes    Notes:
2131b2002411SLois Curfman McInnes    SNESRegister() may be called multiple times to add several user-defined solvers.
21323a40ed3dSBarry Smith 
2133b2002411SLois Curfman McInnes    If dynamic libraries are used, then the fourth input argument (routine_create)
2134b2002411SLois Curfman McInnes    is ignored.
2135b2002411SLois Curfman McInnes 
2136b2002411SLois Curfman McInnes    Sample usage:
2137b2002411SLois Curfman McInnes    SNESRegister("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2138b2002411SLois Curfman McInnes                 "MySolverCreate",MySolverCreate);
2139b2002411SLois Curfman McInnes 
2140b2002411SLois Curfman McInnes    Then, your solver can be chosen with the procedural interface via
2141b2002411SLois Curfman McInnes $     SNESSetType(snes,"my_solver")
2142b2002411SLois Curfman McInnes    or at runtime via the option
2143b2002411SLois Curfman McInnes $     -snes_type my_solver
2144b2002411SLois Curfman McInnes 
2145b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register
2146b2002411SLois Curfman McInnes 
2147b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2148b2002411SLois Curfman McInnes M*/
2149b2002411SLois Curfman McInnes 
2150b2002411SLois Curfman McInnes #undef __FUNC__
2151b2002411SLois Curfman McInnes #define __FUNC__ "SNESRegister_Private"
2152b2002411SLois Curfman McInnes int SNESRegister_Private(char *sname,char *path,char *name,int (*function)(SNES))
2153b2002411SLois Curfman McInnes {
2154b2002411SLois Curfman McInnes   char fullname[256];
2155b2002411SLois Curfman McInnes   int  ierr;
2156b2002411SLois Curfman McInnes 
2157b2002411SLois Curfman McInnes   PetscFunctionBegin;
2158b2002411SLois Curfman McInnes   PetscStrcpy(fullname,path); PetscStrcat(fullname,":");PetscStrcat(fullname,name);
2159b2002411SLois Curfman McInnes   ierr = DLRegister_Private(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr);
2160b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2161b2002411SLois Curfman McInnes }
2162