xref: /petsc/src/snes/interface/snes.c (revision b2002411cd4bf8262b8aeca5b64ca630feca1e56)
15cd90555SBarry Smith 
2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
3*b2002411SLois Curfman McInnes static char vcid[] = "$Id: snes.c,v 1.145 1998/04/21 19:38:21 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 
1484cb2905SBarry Smith 
155615d1e5SSatish Balay #undef __FUNC__
16d4bb536fSBarry Smith #define __FUNC__ "SNESView"
179b94acceSBarry Smith /*@
189b94acceSBarry Smith    SNESView - Prints the SNES data structure.
199b94acceSBarry Smith 
209b94acceSBarry Smith    Input Parameters:
219b94acceSBarry Smith .  SNES - the SNES context
229b94acceSBarry Smith .  viewer - visualization context
239b94acceSBarry Smith 
24fee21e36SBarry Smith    Collective on SNES, unless Viewer is VIEWER_STDOUT_SELF
25fee21e36SBarry Smith 
269b94acceSBarry Smith    Options Database Key:
27c8a8ba5cSLois Curfman McInnes .  -snes_view - Calls SNESView() at end of SNESSolve()
289b94acceSBarry Smith 
299b94acceSBarry Smith    Notes:
309b94acceSBarry Smith    The available visualization contexts include
31c8a8ba5cSLois Curfman McInnes .     VIEWER_STDOUT_SELF - standard output (default)
32c8a8ba5cSLois Curfman McInnes .     VIEWER_STDOUT_WORLD - synchronized standard
33c8a8ba5cSLois Curfman McInnes          output where only the first processor opens
34c8a8ba5cSLois Curfman McInnes          the file.  All other processors send their
35c8a8ba5cSLois Curfman McInnes          data to the first processor to print.
369b94acceSBarry Smith 
379b94acceSBarry Smith    The user can open alternative vistualization contexts with
38dbb450caSBarry Smith $    ViewerFileOpenASCII() - output to a specified file
399b94acceSBarry Smith 
409b94acceSBarry Smith .keywords: SNES, view
419b94acceSBarry Smith 
42dbb450caSBarry Smith .seealso: ViewerFileOpenASCII()
439b94acceSBarry Smith @*/
449b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
459b94acceSBarry Smith {
469b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
479b94acceSBarry Smith   FILE                *fd;
489b94acceSBarry Smith   int                 ierr;
499b94acceSBarry Smith   SLES                sles;
509b94acceSBarry Smith   char                *method;
5119bcc07fSBarry Smith   ViewerType          vtype;
529b94acceSBarry Smith 
533a40ed3dSBarry Smith   PetscFunctionBegin;
5474679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5574679c65SBarry Smith   if (viewer) {PetscValidHeader(viewer);}
5674679c65SBarry Smith   else { viewer = VIEWER_STDOUT_SELF; }
5774679c65SBarry Smith 
5819bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
5919bcc07fSBarry Smith   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
6090ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
6177c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"SNES Object:\n");
6282bf6240SBarry Smith     SNESGetType(snes,&method);
6377c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"  method: %s\n",method);
64e1311b90SBarry Smith     if (snes->view) (*snes->view)(snes,viewer);
6577c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
669b94acceSBarry Smith       "  maximum iterations=%d, maximum function evaluations=%d\n",
679b94acceSBarry Smith       snes->max_its,snes->max_funcs);
6877c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,
699b94acceSBarry Smith     "  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
709b94acceSBarry Smith       snes->rtol, snes->atol, snes->trunctol, snes->xtol);
717a00f4a9SLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
727a00f4a9SLois Curfman McInnes     "  total number of linear solver iterations=%d\n",snes->linear_its);
73ae3c334cSLois Curfman McInnes     PetscFPrintf(snes->comm,fd,
7468632166SLois Curfman McInnes      "  total number of function evaluations=%d\n",snes->nfuncs);
759b94acceSBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)
7677c4ece6SBarry Smith       PetscFPrintf(snes->comm,fd,"  min function tolerance=%g\n",snes->fmin);
779b94acceSBarry Smith     if (snes->ksp_ewconv) {
789b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
799b94acceSBarry Smith       if (kctx) {
8077c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
819b94acceSBarry Smith      "  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",
829b94acceSBarry Smith         kctx->version);
8377c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,
849b94acceSBarry Smith           "    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,
859b94acceSBarry Smith           kctx->rtol_max,kctx->threshold);
8677c4ece6SBarry Smith         PetscFPrintf(snes->comm,fd,"    gamma=%g, alpha=%g, alpha2=%g\n",
879b94acceSBarry Smith           kctx->gamma,kctx->alpha,kctx->alpha2);
889b94acceSBarry Smith       }
899b94acceSBarry Smith     }
90c30f285eSLois Curfman McInnes   } else if (vtype == STRING_VIEWER) {
9182bf6240SBarry Smith     SNESGetType(snes,&method);
92c30f285eSLois Curfman McInnes     ViewerStringSPrintf(viewer," %-3.3s",method);
9319bcc07fSBarry Smith   }
949b94acceSBarry Smith   SNESGetSLES(snes,&sles);
959b94acceSBarry Smith   ierr = SLESView(sles,viewer); CHKERRQ(ierr);
963a40ed3dSBarry Smith   PetscFunctionReturn(0);
979b94acceSBarry Smith }
989b94acceSBarry Smith 
99639f9d9dSBarry Smith /*
100639f9d9dSBarry Smith        We retain a list of functions that also take SNES command
101639f9d9dSBarry Smith     line options. These are called at the end SNESSetFromOptions()
102639f9d9dSBarry Smith */
103639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5
104639f9d9dSBarry Smith static int numberofsetfromoptions;
105639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
106639f9d9dSBarry Smith 
1075615d1e5SSatish Balay #undef __FUNC__
108d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker"
109639f9d9dSBarry Smith /*@
110639f9d9dSBarry Smith     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
111639f9d9dSBarry Smith 
112639f9d9dSBarry Smith     Input Parameter:
113639f9d9dSBarry Smith .   snescheck - function that checks for options
114639f9d9dSBarry Smith 
115fee21e36SBarry Smith     Not Collective
116fee21e36SBarry Smith 
117639f9d9dSBarry Smith .seealso: SNESSetFromOptions()
118639f9d9dSBarry Smith @*/
119639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) )
120639f9d9dSBarry Smith {
1213a40ed3dSBarry Smith   PetscFunctionBegin;
122639f9d9dSBarry Smith   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
123a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed");
124639f9d9dSBarry Smith   }
125639f9d9dSBarry Smith 
126639f9d9dSBarry Smith   othersetfromoptions[numberofsetfromoptions++] = snescheck;
1273a40ed3dSBarry Smith   PetscFunctionReturn(0);
128639f9d9dSBarry Smith }
129639f9d9dSBarry Smith 
1305615d1e5SSatish Balay #undef __FUNC__
1315615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions"
1329b94acceSBarry Smith /*@
133682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1349b94acceSBarry Smith 
1359b94acceSBarry Smith    Input Parameter:
1369b94acceSBarry Smith .  snes - the SNES context
1379b94acceSBarry Smith 
138fee21e36SBarry Smith    Collective on SNES
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 
2719b94acceSBarry Smith    Input Parameters:
2729b94acceSBarry Smith .  snes - the SNES context
2739b94acceSBarry Smith .  usrP - optional user context
2749b94acceSBarry Smith 
275fee21e36SBarry Smith    Collective on SNES
276fee21e36SBarry Smith 
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 
2959b94acceSBarry Smith    Input Parameter:
2969b94acceSBarry Smith .  snes - SNES context
2979b94acceSBarry Smith 
2989b94acceSBarry Smith    Output Parameter:
2999b94acceSBarry Smith .  usrP - user context
3009b94acceSBarry Smith 
301fee21e36SBarry Smith    Not Collective
302fee21e36SBarry 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 
3219b94acceSBarry Smith    Input Parameter:
3229b94acceSBarry Smith .  snes - SNES context
3239b94acceSBarry Smith 
3249b94acceSBarry Smith    Output Parameter:
3259b94acceSBarry Smith .  iter - iteration number
3269b94acceSBarry Smith 
327fee21e36SBarry Smith    Not Collective
328fee21e36SBarry 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 
3469b94acceSBarry Smith    Input Parameter:
3479b94acceSBarry Smith .  snes - SNES context
3489b94acceSBarry Smith 
3499b94acceSBarry Smith    Output Parameter:
3509b94acceSBarry Smith .  fnorm - 2-norm of function
3519b94acceSBarry Smith 
352fee21e36SBarry Smith    Collective on SNES
353fee21e36SBarry 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 
3819b94acceSBarry Smith    Input Parameter:
3829b94acceSBarry Smith .  snes - SNES context
3839b94acceSBarry Smith 
3849b94acceSBarry Smith    Output Parameter:
3859b94acceSBarry Smith .  fnorm - 2-norm of gradient
3869b94acceSBarry Smith 
387fee21e36SBarry Smith    Collective on SNES
388fee21e36SBarry 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 
4169b94acceSBarry Smith    Input Parameter:
4179b94acceSBarry Smith .  snes - SNES context
4189b94acceSBarry Smith 
4199b94acceSBarry Smith    Output Parameter:
4209b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4219b94acceSBarry Smith 
422fee21e36SBarry Smith    Not Collective
423fee21e36SBarry 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 
444c96a6f78SLois Curfman McInnes    Input Parameter:
445c96a6f78SLois Curfman McInnes .  snes - SNES context
446c96a6f78SLois Curfman McInnes 
447c96a6f78SLois Curfman McInnes    Output Parameter:
448c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
449c96a6f78SLois Curfman McInnes 
450fee21e36SBarry Smith    Not Collective
451fee21e36SBarry Smith 
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 
4719b94acceSBarry Smith    Input Parameter:
4729b94acceSBarry Smith .  snes - the SNES context
4739b94acceSBarry Smith 
4749b94acceSBarry Smith    Output Parameter:
4759b94acceSBarry Smith .  sles - the SLES context
4769b94acceSBarry Smith 
477fee21e36SBarry Smith    Not Collective, but if SNES object is parallel, then SLES object is parallel
478fee21e36SBarry 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 
5029b94acceSBarry Smith    Input Parameter:
5039b94acceSBarry Smith .  comm - MPI communicator
5049b94acceSBarry Smith .  type - type of method, one of
5059b94acceSBarry Smith $    SNES_NONLINEAR_EQUATIONS
5069b94acceSBarry Smith $      (for systems of nonlinear equations)
5079b94acceSBarry Smith $    SNES_UNCONSTRAINED_MINIMIZATION
5089b94acceSBarry Smith $      (for unconstrained minimization)
5099b94acceSBarry Smith 
5109b94acceSBarry Smith    Output Parameter:
5119b94acceSBarry Smith .  outsnes - the new SNES context
5129b94acceSBarry Smith 
513fee21e36SBarry Smith    Collective on MPI_Comm
514fee21e36SBarry Smith 
515c1f60f51SBarry Smith    Options Database Key:
51619bd6747SLois Curfman McInnes $   -snes_mf - use default matrix-free Jacobian-vector products,
51719bd6747SLois Curfman McInnes $              and no preconditioning matrix
51819bd6747SLois Curfman McInnes $   -snes_mf_operator - use default matrix-free Jacobian-vector
51919bd6747SLois Curfman McInnes $             products, and a user-provided preconditioning matrix
52019bd6747SLois Curfman McInnes $             as set by SNESSetJacobian()
52119bd6747SLois Curfman McInnes $   -snes_fd - use (slow!) finite differences to compute Jacobian
522c1f60f51SBarry Smith 
5239b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5249b94acceSBarry Smith 
52563a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5269b94acceSBarry Smith @*/
5274b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5289b94acceSBarry Smith {
5299b94acceSBarry Smith   int                 ierr;
5309b94acceSBarry Smith   SNES                snes;
5319b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
53237fcc0dbSBarry Smith 
5333a40ed3dSBarry Smith   PetscFunctionBegin;
5349b94acceSBarry Smith   *outsnes = 0;
535d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
536d252947aSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type");
537d64ed03dSBarry Smith   }
538f830108cSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,comm,SNESDestroy,SNESView);
5399b94acceSBarry Smith   PLogObjectCreate(snes);
5409b94acceSBarry Smith   snes->max_its           = 50;
5419750a799SBarry Smith   snes->max_funcs	  = 10000;
5429b94acceSBarry Smith   snes->norm		  = 0.0;
5435a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
544b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
545b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5469b94acceSBarry Smith     snes->atol		  = 1.e-10;
547d64ed03dSBarry Smith   } else {
548b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
549b4874afaSBarry Smith     snes->ttol            = 0.0;
550b4874afaSBarry Smith     snes->atol		  = 1.e-50;
551b4874afaSBarry Smith   }
5529b94acceSBarry Smith   snes->xtol		  = 1.e-8;
553b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5549b94acceSBarry Smith   snes->nfuncs            = 0;
5559b94acceSBarry Smith   snes->nfailures         = 0;
5567a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
557639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5589b94acceSBarry Smith   snes->data              = 0;
5599b94acceSBarry Smith   snes->view              = 0;
5609b94acceSBarry Smith   snes->computeumfunction = 0;
5619b94acceSBarry Smith   snes->umfunP            = 0;
5629b94acceSBarry Smith   snes->fc                = 0;
5639b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5649b94acceSBarry Smith   snes->fmin              = -1.e30;
5659b94acceSBarry Smith   snes->method_class      = type;
5669b94acceSBarry Smith   snes->set_method_called = 0;
56782bf6240SBarry Smith   snes->setupcalled      = 0;
5689b94acceSBarry Smith   snes->ksp_ewconv        = 0;
5696f24a144SLois Curfman McInnes   snes->xmonitor          = 0;
5706f24a144SLois Curfman McInnes   snes->vwork             = 0;
5716f24a144SLois Curfman McInnes   snes->nwork             = 0;
572c9005455SLois Curfman McInnes   snes->conv_hist_size    = 0;
573c9005455SLois Curfman McInnes   snes->conv_act_size     = 0;
574c9005455SLois Curfman McInnes   snes->conv_hist         = 0;
5759b94acceSBarry Smith 
5769b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
5770452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx);
578eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
5799b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
5809b94acceSBarry Smith   kctx->version     = 2;
5819b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
5829b94acceSBarry Smith                              this was too large for some test cases */
5839b94acceSBarry Smith   kctx->rtol_last   = 0;
5849b94acceSBarry Smith   kctx->rtol_max    = .9;
5859b94acceSBarry Smith   kctx->gamma       = 1.0;
5869b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
5879b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
5889b94acceSBarry Smith   kctx->threshold   = .1;
5899b94acceSBarry Smith   kctx->lresid_last = 0;
5909b94acceSBarry Smith   kctx->norm_last   = 0;
5919b94acceSBarry Smith 
5929b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr);
5939b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
5945334005bSBarry Smith 
5959b94acceSBarry Smith   *outsnes = snes;
5963a40ed3dSBarry Smith   PetscFunctionReturn(0);
5979b94acceSBarry Smith }
5989b94acceSBarry Smith 
5999b94acceSBarry Smith /* --------------------------------------------------------------- */
6005615d1e5SSatish Balay #undef __FUNC__
601d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction"
6029b94acceSBarry Smith /*@C
6039b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6049b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6059b94acceSBarry Smith    equations.
6069b94acceSBarry Smith 
6079b94acceSBarry Smith    Input Parameters:
6089b94acceSBarry Smith .  snes - the SNES context
6099b94acceSBarry Smith .  func - function evaluation routine
6109b94acceSBarry Smith .  r - vector to store function value
6112cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
6122cd2dfdcSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6139b94acceSBarry Smith 
614fee21e36SBarry Smith    Collective on SNES
615fee21e36SBarry Smith 
6169b94acceSBarry Smith    Calling sequence of func:
617313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Vec f,void *ctx);
6189b94acceSBarry Smith 
6199b94acceSBarry Smith .  x - input vector
620313e4042SLois Curfman McInnes .  f - function vector
6212cd2dfdcSLois Curfman McInnes .  ctx - optional user-defined function context
6229b94acceSBarry Smith 
6239b94acceSBarry Smith    Notes:
6249b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6259b94acceSBarry Smith $      f'(x) x = -f(x),
6269b94acceSBarry Smith $  where f'(x) denotes the Jacobian matrix and f(x) is the function.
6279b94acceSBarry Smith 
6289b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6299b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6309b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6319b94acceSBarry Smith 
6329b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6339b94acceSBarry Smith 
634a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6359b94acceSBarry Smith @*/
6365334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx)
6379b94acceSBarry Smith {
6383a40ed3dSBarry Smith   PetscFunctionBegin;
63977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
640a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
641a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
642a8c6a408SBarry Smith   }
6439b94acceSBarry Smith   snes->computefunction     = func;
6449b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6459b94acceSBarry Smith   snes->funP                = ctx;
6463a40ed3dSBarry Smith   PetscFunctionReturn(0);
6479b94acceSBarry Smith }
6489b94acceSBarry Smith 
6495615d1e5SSatish Balay #undef __FUNC__
6505615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction"
6519b94acceSBarry Smith /*@
6529b94acceSBarry Smith    SNESComputeFunction - Computes the function that has been set with
6539b94acceSBarry Smith    SNESSetFunction().
6549b94acceSBarry Smith 
6559b94acceSBarry Smith    Input Parameters:
6569b94acceSBarry Smith .  snes - the SNES context
6579b94acceSBarry Smith .  x - input vector
6589b94acceSBarry Smith 
6599b94acceSBarry Smith    Output Parameter:
6603638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6619b94acceSBarry Smith 
662fee21e36SBarry Smith    Collective on SNES
663fee21e36SBarry Smith 
6641bffabb2SLois Curfman McInnes    Notes:
6659b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6669b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6679b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6689b94acceSBarry Smith 
6699b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
6709b94acceSBarry Smith 
671a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
6729b94acceSBarry Smith @*/
6739b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y)
6749b94acceSBarry Smith {
6759b94acceSBarry Smith   int    ierr;
6769b94acceSBarry Smith 
6773a40ed3dSBarry Smith   PetscFunctionBegin;
678d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
679a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
680d4bb536fSBarry Smith   }
6819b94acceSBarry Smith   PLogEventBegin(SNES_FunctionEval,snes,x,y,0);
682d64ed03dSBarry Smith   PetscStackPush("SNES user function");
6839b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
684d64ed03dSBarry Smith   PetscStackPop;
685ae3c334cSLois Curfman McInnes   snes->nfuncs++;
6869b94acceSBarry Smith   PLogEventEnd(SNES_FunctionEval,snes,x,y,0);
6873a40ed3dSBarry Smith   PetscFunctionReturn(0);
6889b94acceSBarry Smith }
6899b94acceSBarry Smith 
6905615d1e5SSatish Balay #undef __FUNC__
691d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction"
6929b94acceSBarry Smith /*@C
6939b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
6949b94acceSBarry Smith    unconstrained minimization.
6959b94acceSBarry Smith 
6969b94acceSBarry Smith    Input Parameters:
6979b94acceSBarry Smith .  snes - the SNES context
6989b94acceSBarry Smith .  func - function evaluation routine
699044dda88SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
700044dda88SLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7019b94acceSBarry Smith 
702fee21e36SBarry Smith    Collective on SNES
703fee21e36SBarry Smith 
7049b94acceSBarry Smith    Calling sequence of func:
7059b94acceSBarry Smith .  func (SNES snes,Vec x,double *f,void *ctx);
7069b94acceSBarry Smith 
7079b94acceSBarry Smith .  x - input vector
7089b94acceSBarry Smith .  f - function
709044dda88SLois Curfman McInnes .  ctx - [optional] user-defined function context
7109b94acceSBarry Smith 
7119b94acceSBarry Smith    Notes:
7129b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7139b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7149b94acceSBarry Smith    SNESSetFunction().
7159b94acceSBarry Smith 
7169b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
7179b94acceSBarry Smith 
718a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
719a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
7209b94acceSBarry Smith @*/
7219b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),
7229b94acceSBarry Smith                       void *ctx)
7239b94acceSBarry Smith {
7243a40ed3dSBarry Smith   PetscFunctionBegin;
72577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
726a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
727a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
728a8c6a408SBarry Smith   }
7299b94acceSBarry Smith   snes->computeumfunction   = func;
7309b94acceSBarry Smith   snes->umfunP              = ctx;
7313a40ed3dSBarry Smith   PetscFunctionReturn(0);
7329b94acceSBarry Smith }
7339b94acceSBarry Smith 
7345615d1e5SSatish Balay #undef __FUNC__
7355615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction"
7369b94acceSBarry Smith /*@
7379b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7389b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7399b94acceSBarry Smith 
7409b94acceSBarry Smith    Input Parameters:
7419b94acceSBarry Smith .  snes - the SNES context
7429b94acceSBarry Smith .  x - input vector
7439b94acceSBarry Smith 
7449b94acceSBarry Smith    Output Parameter:
7459b94acceSBarry Smith .  y - function value
7469b94acceSBarry Smith 
747fee21e36SBarry Smith    Collective on SNES
748fee21e36SBarry Smith 
7499b94acceSBarry Smith    Notes:
7509b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7519b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7529b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
753a86d99e1SLois Curfman McInnes 
754a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
755a86d99e1SLois Curfman McInnes 
756a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
757a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
7589b94acceSBarry Smith @*/
7599b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y)
7609b94acceSBarry Smith {
7619b94acceSBarry Smith   int    ierr;
7623a40ed3dSBarry Smith 
7633a40ed3dSBarry Smith   PetscFunctionBegin;
764a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
765a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
766a8c6a408SBarry Smith   }
7679b94acceSBarry Smith   PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);
768d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
7699b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr);
770d64ed03dSBarry Smith   PetscStackPop;
771ae3c334cSLois Curfman McInnes   snes->nfuncs++;
7729b94acceSBarry Smith   PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);
7733a40ed3dSBarry Smith   PetscFunctionReturn(0);
7749b94acceSBarry Smith }
7759b94acceSBarry Smith 
7765615d1e5SSatish Balay #undef __FUNC__
777d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient"
7789b94acceSBarry Smith /*@C
7799b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
7809b94acceSBarry Smith    vector for use by the SNES routines.
7819b94acceSBarry Smith 
7829b94acceSBarry Smith    Input Parameters:
7839b94acceSBarry Smith .  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)
7879b94acceSBarry Smith .  r - vector to store gradient value
7889b94acceSBarry Smith 
789fee21e36SBarry Smith    Collective on SNES
790fee21e36SBarry Smith 
7919b94acceSBarry Smith    Calling sequence of func:
7929b94acceSBarry Smith .  func (SNES, Vec x, Vec g, void *ctx);
7939b94acceSBarry Smith 
7949b94acceSBarry Smith .  x - input vector
7959b94acceSBarry Smith .  g - gradient vector
796044dda88SLois Curfman McInnes .  ctx - optional user-defined gradient context
7979b94acceSBarry Smith 
7989b94acceSBarry Smith    Notes:
7999b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
8009b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
8019b94acceSBarry Smith    SNESSetFunction().
8029b94acceSBarry Smith 
8039b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
8049b94acceSBarry Smith 
805a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
806a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
8079b94acceSBarry Smith @*/
80874679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
8099b94acceSBarry Smith {
8103a40ed3dSBarry Smith   PetscFunctionBegin;
81177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
812a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
813a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
814a8c6a408SBarry Smith   }
8159b94acceSBarry Smith   snes->computefunction     = func;
8169b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
8179b94acceSBarry Smith   snes->funP                = ctx;
8183a40ed3dSBarry Smith   PetscFunctionReturn(0);
8199b94acceSBarry Smith }
8209b94acceSBarry Smith 
8215615d1e5SSatish Balay #undef __FUNC__
8225615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient"
8239b94acceSBarry Smith /*@
824a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
825a86d99e1SLois Curfman McInnes    SNESSetGradient().
8269b94acceSBarry Smith 
8279b94acceSBarry Smith    Input Parameters:
8289b94acceSBarry Smith .  snes - the SNES context
8299b94acceSBarry Smith .  x - input vector
8309b94acceSBarry Smith 
8319b94acceSBarry Smith    Output Parameter:
8329b94acceSBarry Smith .  y - gradient vector
8339b94acceSBarry Smith 
834fee21e36SBarry Smith    Collective on SNES
835fee21e36SBarry Smith 
8369b94acceSBarry Smith    Notes:
8379b94acceSBarry Smith    SNESComputeGradient() is valid only for
8389b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8399b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
840a86d99e1SLois Curfman McInnes 
841a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
842a86d99e1SLois Curfman McInnes 
843a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
844a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
8459b94acceSBarry Smith @*/
8469b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y)
8479b94acceSBarry Smith {
8489b94acceSBarry Smith   int    ierr;
8493a40ed3dSBarry Smith 
8503a40ed3dSBarry Smith   PetscFunctionBegin;
8513a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
852a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
8533a40ed3dSBarry Smith   }
8549b94acceSBarry Smith   PLogEventBegin(SNES_GradientEval,snes,x,y,0);
855d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
8569b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr);
857d64ed03dSBarry Smith   PetscStackPop;
8589b94acceSBarry Smith   PLogEventEnd(SNES_GradientEval,snes,x,y,0);
8593a40ed3dSBarry Smith   PetscFunctionReturn(0);
8609b94acceSBarry Smith }
8619b94acceSBarry Smith 
8625615d1e5SSatish Balay #undef __FUNC__
8635615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian"
86462fef451SLois Curfman McInnes /*@
86562fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
86662fef451SLois Curfman McInnes    set with SNESSetJacobian().
86762fef451SLois Curfman McInnes 
86862fef451SLois Curfman McInnes    Input Parameters:
86962fef451SLois Curfman McInnes .  snes - the SNES context
87062fef451SLois Curfman McInnes .  x - input vector
87162fef451SLois Curfman McInnes 
87262fef451SLois Curfman McInnes    Output Parameters:
87362fef451SLois Curfman McInnes .  A - Jacobian matrix
87462fef451SLois Curfman McInnes .  B - optional preconditioning matrix
87562fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
87662fef451SLois Curfman McInnes 
877fee21e36SBarry Smith    Collective on SNES and Mat
878fee21e36SBarry Smith 
87962fef451SLois Curfman McInnes    Notes:
88062fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
88162fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
88262fef451SLois Curfman McInnes 
883dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
884dc5a77f8SLois Curfman McInnes    flag parameter.
88562fef451SLois Curfman McInnes 
88662fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
88762fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
888005c665bSBarry Smith    methods is SNESComputeHessian().
88962fef451SLois Curfman McInnes 
89062fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
89162fef451SLois Curfman McInnes 
89262fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
89362fef451SLois Curfman McInnes @*/
8949b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8959b94acceSBarry Smith {
8969b94acceSBarry Smith   int    ierr;
8973a40ed3dSBarry Smith 
8983a40ed3dSBarry Smith   PetscFunctionBegin;
8993a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
900a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
9013a40ed3dSBarry Smith   }
9023a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
9039b94acceSBarry Smith   PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);
904c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
905d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
9069b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr);
907d64ed03dSBarry Smith   PetscStackPop;
9089b94acceSBarry Smith   PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);
9096d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
91077c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
91177c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9123a40ed3dSBarry Smith   PetscFunctionReturn(0);
9139b94acceSBarry Smith }
9149b94acceSBarry Smith 
9155615d1e5SSatish Balay #undef __FUNC__
9165615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian"
91762fef451SLois Curfman McInnes /*@
91862fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
91962fef451SLois Curfman McInnes    set with SNESSetHessian().
92062fef451SLois Curfman McInnes 
92162fef451SLois Curfman McInnes    Input Parameters:
92262fef451SLois Curfman McInnes .  snes - the SNES context
92362fef451SLois Curfman McInnes .  x - input vector
92462fef451SLois Curfman McInnes 
92562fef451SLois Curfman McInnes    Output Parameters:
92662fef451SLois Curfman McInnes .  A - Hessian matrix
92762fef451SLois Curfman McInnes .  B - optional preconditioning matrix
92862fef451SLois Curfman McInnes .  flag - flag indicating matrix structure
92962fef451SLois Curfman McInnes 
930fee21e36SBarry Smith    Collective on SNES and Mat
931fee21e36SBarry Smith 
93262fef451SLois Curfman McInnes    Notes:
93362fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
93462fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
93562fef451SLois Curfman McInnes 
936dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
937dc5a77f8SLois Curfman McInnes    flag parameter.
93862fef451SLois Curfman McInnes 
93962fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
94062fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
94162fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
94262fef451SLois Curfman McInnes 
94362fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
94462fef451SLois Curfman McInnes 
945a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
946a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
94762fef451SLois Curfman McInnes @*/
94862fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
9499b94acceSBarry Smith {
9509b94acceSBarry Smith   int    ierr;
9513a40ed3dSBarry Smith 
9523a40ed3dSBarry Smith   PetscFunctionBegin;
9533a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
954a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9553a40ed3dSBarry Smith   }
9563a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
95762fef451SLois Curfman McInnes   PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);
9580f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
959d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
96062fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr);
961d64ed03dSBarry Smith   PetscStackPop;
96262fef451SLois Curfman McInnes   PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);
9630f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
96477c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
96577c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9663a40ed3dSBarry Smith   PetscFunctionReturn(0);
9679b94acceSBarry Smith }
9689b94acceSBarry Smith 
9695615d1e5SSatish Balay #undef __FUNC__
970d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian"
9719b94acceSBarry Smith /*@C
9729b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
973044dda88SLois Curfman McInnes    location to store the matrix.
9749b94acceSBarry Smith 
9759b94acceSBarry Smith    Input Parameters:
9769b94acceSBarry Smith .  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
9802cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
9812cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9829b94acceSBarry Smith 
983fee21e36SBarry Smith    Collective on SNES and Mat
984fee21e36SBarry Smith 
9859b94acceSBarry Smith    Calling sequence of func:
986313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9879b94acceSBarry Smith 
9889b94acceSBarry Smith .  x - input vector
9899b94acceSBarry Smith .  A - Jacobian matrix
9909b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
991ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
992ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
9932cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Jacobian context
9949b94acceSBarry Smith 
9959b94acceSBarry Smith    Notes:
996dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
9972cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
998ac21db08SLois Curfman McInnes 
999ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
10009b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
10019b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10029b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10039b94acceSBarry Smith    throughout the global iterations.
10049b94acceSBarry Smith 
10059b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
10069b94acceSBarry Smith 
1007ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
10089b94acceSBarry Smith @*/
10099b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10109b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10119b94acceSBarry Smith {
10123a40ed3dSBarry Smith   PetscFunctionBegin;
101377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1014a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1015a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1016a8c6a408SBarry Smith   }
10179b94acceSBarry Smith   snes->computejacobian = func;
10189b94acceSBarry Smith   snes->jacP            = ctx;
10199b94acceSBarry Smith   snes->jacobian        = A;
10209b94acceSBarry Smith   snes->jacobian_pre    = B;
10213a40ed3dSBarry Smith   PetscFunctionReturn(0);
10229b94acceSBarry Smith }
102362fef451SLois Curfman McInnes 
10245615d1e5SSatish Balay #undef __FUNC__
1025d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian"
1026b4fd4287SBarry Smith /*@
1027b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1028b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1029b4fd4287SBarry Smith 
1030b4fd4287SBarry Smith    Input Parameter:
1031b4fd4287SBarry Smith .  snes - the nonlinear solver context
1032b4fd4287SBarry Smith 
1033b4fd4287SBarry Smith    Output Parameters:
1034b4fd4287SBarry Smith .  A - location to stash Jacobian matrix (or PETSC_NULL)
1035b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
1036b4fd4287SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1037b4fd4287SBarry Smith 
1038fee21e36SBarry Smith    Not Collective, but Mat object will be parallel if SNES object is
1039fee21e36SBarry Smith 
1040b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1041b4fd4287SBarry Smith @*/
1042b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx)
1043b4fd4287SBarry Smith {
10443a40ed3dSBarry Smith   PetscFunctionBegin;
10453a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1046a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
10473a40ed3dSBarry Smith   }
1048b4fd4287SBarry Smith   if (A)   *A = snes->jacobian;
1049b4fd4287SBarry Smith   if (B)   *B = snes->jacobian_pre;
1050b4fd4287SBarry Smith   if (ctx) *ctx = snes->jacP;
10513a40ed3dSBarry Smith   PetscFunctionReturn(0);
1052b4fd4287SBarry Smith }
1053b4fd4287SBarry Smith 
10545615d1e5SSatish Balay #undef __FUNC__
1055d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian"
10569b94acceSBarry Smith /*@C
10579b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1058044dda88SLois Curfman McInnes    location to store the matrix.
10599b94acceSBarry Smith 
10609b94acceSBarry Smith    Input Parameters:
10619b94acceSBarry Smith .  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
1065313e4042SLois Curfman McInnes .  ctx - [optional] user-defined context for private data for the
1066313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
10679b94acceSBarry Smith 
1068fee21e36SBarry Smith    Collective on SNES and Mat
1069fee21e36SBarry Smith 
10709b94acceSBarry Smith    Calling sequence of func:
1071313e4042SLois Curfman McInnes .  func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10729b94acceSBarry Smith 
10739b94acceSBarry Smith .  x - input vector
10749b94acceSBarry Smith .  A - Hessian matrix
10759b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1076ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1077ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
10782cd2dfdcSLois Curfman McInnes .  ctx - [optional] user-defined Hessian context
10799b94acceSBarry Smith 
10809b94acceSBarry Smith    Notes:
1081dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10822cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1083ac21db08SLois Curfman McInnes 
10849b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
10859b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
10869b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10879b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10889b94acceSBarry Smith    throughout the global iterations.
10899b94acceSBarry Smith 
10909b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
10919b94acceSBarry Smith 
1092ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
10939b94acceSBarry Smith @*/
10949b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,
10959b94acceSBarry Smith                     MatStructure*,void*),void *ctx)
10969b94acceSBarry Smith {
10973a40ed3dSBarry Smith   PetscFunctionBegin;
109877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1099d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1100a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1101d4bb536fSBarry Smith   }
11029b94acceSBarry Smith   snes->computejacobian = func;
11039b94acceSBarry Smith   snes->jacP            = ctx;
11049b94acceSBarry Smith   snes->jacobian        = A;
11059b94acceSBarry Smith   snes->jacobian_pre    = B;
11063a40ed3dSBarry Smith   PetscFunctionReturn(0);
11079b94acceSBarry Smith }
11089b94acceSBarry Smith 
11095615d1e5SSatish Balay #undef __FUNC__
1110d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian"
111162fef451SLois Curfman McInnes /*@
111262fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
111362fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
111462fef451SLois Curfman McInnes 
111562fef451SLois Curfman McInnes    Input Parameter:
111662fef451SLois Curfman McInnes .  snes - the nonlinear solver context
111762fef451SLois Curfman McInnes 
111862fef451SLois Curfman McInnes    Output Parameters:
111962fef451SLois Curfman McInnes .  A - location to stash Hessian matrix (or PETSC_NULL)
112062fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
112162fef451SLois Curfman McInnes .  ctx - location to stash Hessian ctx (or PETSC_NULL)
112262fef451SLois Curfman McInnes 
1123fee21e36SBarry Smith    Not Collective, but Mat object is parallel if SNES object is parallel
1124fee21e36SBarry Smith 
112562fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
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 
1147*b2002411SLois Curfman McInnes    Input Parameters:
11489b94acceSBarry Smith .  snes - the SNES context
11498ddd3da0SLois Curfman McInnes .  x - the solution vector
11509b94acceSBarry Smith 
1151fee21e36SBarry Smith    Collective on SNES
1152fee21e36SBarry Smith 
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 
12509b94acceSBarry Smith    Input Parameter:
12519b94acceSBarry Smith .  snes - the SNES context
12529b94acceSBarry Smith 
1253fee21e36SBarry Smith    Collective on SNES
1254fee21e36SBarry 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 
12859b94acceSBarry Smith    Input Parameters:
12869b94acceSBarry Smith .  snes - the SNES context
128733174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
128833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
128933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
129033174efeSLois Curfman McInnes            of the change in the solution between steps
129133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
129233174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
12939b94acceSBarry Smith 
1294fee21e36SBarry Smith    Collective on SNES
1295fee21e36SBarry Smith 
129633174efeSLois Curfman McInnes    Options Database Keys:
129733174efeSLois Curfman McInnes $    -snes_atol <atol>
129833174efeSLois Curfman McInnes $    -snes_rtol <rtol>
129933174efeSLois Curfman McInnes $    -snes_stol <stol>
130033174efeSLois Curfman McInnes $    -snes_max_it <maxit>
130133174efeSLois Curfman McInnes $    -snes_max_funcs <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 
132833174efeSLois Curfman McInnes    Input Parameters:
132933174efeSLois Curfman McInnes .  snes - the SNES context
133033174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
133133174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
133233174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
133333174efeSLois Curfman McInnes            of the change in the solution between steps
133433174efeSLois Curfman McInnes .  maxit - maximum number of iterations
133533174efeSLois Curfman McInnes .  maxf - maximum number of function evaluations
133633174efeSLois Curfman McInnes 
1337fee21e36SBarry Smith    Not Collective
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 
13639b94acceSBarry Smith    Input Parameters:
13649b94acceSBarry Smith .  snes - the SNES context
13659b94acceSBarry Smith .  tol - tolerance
13669b94acceSBarry Smith 
1367fee21e36SBarry Smith    Collective on SNES
1368fee21e36SBarry Smith 
13699b94acceSBarry Smith    Options Database Key:
1370d7a720efSLois Curfman McInnes $    -snes_trtol <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 
13909b94acceSBarry Smith    Input Parameters:
13919b94acceSBarry Smith .  snes - the SNES context
13929b94acceSBarry Smith .  ftol - minimum function tolerance
13939b94acceSBarry Smith 
1394fee21e36SBarry Smith    Collective on SNES
1395fee21e36SBarry Smith 
13969b94acceSBarry Smith    Options Database Key:
1397d7a720efSLois Curfman McInnes $    -snes_fmin <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 
14249b94acceSBarry Smith    Input Parameters:
14259b94acceSBarry Smith .  snes - the SNES context
14269b94acceSBarry Smith .  func - monitoring routine
1427044dda88SLois Curfman McInnes .  mctx - [optional] user-defined context for private data for the
1428044dda88SLois Curfman McInnes           monitor routine (may be PETSC_NULL)
14299b94acceSBarry Smith 
1430fee21e36SBarry Smith    Collective on SNES
1431fee21e36SBarry Smith 
14329b94acceSBarry Smith    Calling sequence of func:
1433682d7d0cSBarry Smith    int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx)
14349b94acceSBarry Smith 
14359b94acceSBarry Smith $    snes - the SNES context
14369b94acceSBarry Smith $    its - iteration number
1437044dda88SLois Curfman McInnes $    mctx - [optional] monitoring context
14389b94acceSBarry Smith $
14399b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
14409b94acceSBarry Smith $    norm - 2-norm function value (may be estimated)
14419b94acceSBarry Smith $
14429b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
14439b94acceSBarry Smith $    norm - 2-norm gradient value (may be estimated)
14449b94acceSBarry Smith 
14459665c990SLois Curfman McInnes    Options Database Keys:
14469665c990SLois Curfman McInnes $    -snes_monitor        : sets SNESDefaultMonitor()
14479665c990SLois Curfman McInnes $    -snes_xmonitor       : sets line graph monitor,
14489665c990SLois Curfman McInnes $                           uses SNESLGMonitorCreate()
14499665c990SLois Curfman McInnes $    -snes_cancelmonitors : cancels all monitors that have
14509665c990SLois Curfman McInnes $                           been hardwired into a code by
14519665c990SLois Curfman McInnes $                           calls to SNESSetMonitor(), but
14529665c990SLois Curfman McInnes $                           does not cancel those set via
14539665c990SLois Curfman McInnes $                           the options database.
14549665c990SLois Curfman McInnes 
14559665c990SLois Curfman McInnes 
1456639f9d9dSBarry Smith    Notes:
14576bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
14586bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
14596bc08f3fSLois Curfman McInnes    order in which they were set.
1460639f9d9dSBarry Smith 
14619b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
14629b94acceSBarry Smith 
14635cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
14649b94acceSBarry Smith @*/
146574679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx )
14669b94acceSBarry Smith {
14673a40ed3dSBarry Smith   PetscFunctionBegin;
1468639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1469a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
1470639f9d9dSBarry Smith   }
1471639f9d9dSBarry Smith 
1472639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1473639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
14743a40ed3dSBarry Smith   PetscFunctionReturn(0);
14759b94acceSBarry Smith }
14769b94acceSBarry Smith 
14775615d1e5SSatish Balay #undef __FUNC__
14785cd90555SBarry Smith #define __FUNC__ "SNESClearMonitor"
14795cd90555SBarry Smith /*@C
14805cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
14815cd90555SBarry Smith 
14825cd90555SBarry Smith    Input Parameters:
14835cd90555SBarry Smith .  snes - the SNES context
14845cd90555SBarry Smith 
1485fee21e36SBarry Smith    Collective on SNES
1486fee21e36SBarry Smith 
14875cd90555SBarry Smith    Options Database:
14885cd90555SBarry Smith $    -snes_cancelmonitors : cancels all monitors that have
14895cd90555SBarry Smith $                           been hardwired into a code by
14905cd90555SBarry Smith $                           calls to SNESSetMonitor(), but
14915cd90555SBarry Smith $                           does not cancel those set via
14925cd90555SBarry Smith $                           the options database.
14935cd90555SBarry Smith 
14945cd90555SBarry Smith    Notes:
14955cd90555SBarry Smith      There is no way to clear one specific monitor from a SNES object.
14965cd90555SBarry Smith 
14975cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
14985cd90555SBarry Smith 
14995cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
15005cd90555SBarry Smith @*/
15015cd90555SBarry Smith int SNESClearMonitor( SNES snes )
15025cd90555SBarry Smith {
15035cd90555SBarry Smith   PetscFunctionBegin;
15045cd90555SBarry Smith   snes->numbermonitors = 0;
15055cd90555SBarry Smith   PetscFunctionReturn(0);
15065cd90555SBarry Smith }
15075cd90555SBarry Smith 
15085cd90555SBarry Smith #undef __FUNC__
1509d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest"
15109b94acceSBarry Smith /*@C
15119b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
15129b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
15139b94acceSBarry Smith 
15149b94acceSBarry Smith    Input Parameters:
15159b94acceSBarry Smith .  snes - the SNES context
15169b94acceSBarry Smith .  func - routine to test for convergence
1517044dda88SLois Curfman McInnes .  cctx - [optional] context for private data for the convergence routine
1518044dda88SLois Curfman McInnes           (may be PETSC_NULL)
15199b94acceSBarry Smith 
1520fee21e36SBarry Smith    Collective on SNES
1521fee21e36SBarry Smith 
15229b94acceSBarry Smith    Calling sequence of func:
15239b94acceSBarry Smith    int func (SNES snes,double xnorm,double gnorm,
15249b94acceSBarry Smith              double f,void *cctx)
15259b94acceSBarry Smith 
15269b94acceSBarry Smith $    snes - the SNES context
1527044dda88SLois Curfman McInnes $    cctx - [optional] convergence context
15289b94acceSBarry Smith $    xnorm - 2-norm of current iterate
15299b94acceSBarry Smith $
15309b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods:
15319b94acceSBarry Smith $    gnorm - 2-norm of current step
15329b94acceSBarry Smith $    f - 2-norm of function
15339b94acceSBarry Smith $
15349b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods:
15359b94acceSBarry Smith $    gnorm - 2-norm of current gradient
15369b94acceSBarry Smith $    f - function value
15379b94acceSBarry Smith 
15389b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
15399b94acceSBarry Smith 
154040191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
154140191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
15429b94acceSBarry Smith @*/
154374679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx)
15449b94acceSBarry Smith {
15453a40ed3dSBarry Smith   PetscFunctionBegin;
15469b94acceSBarry Smith   (snes)->converged = func;
15479b94acceSBarry Smith   (snes)->cnvP      = cctx;
15483a40ed3dSBarry Smith   PetscFunctionReturn(0);
15499b94acceSBarry Smith }
15509b94acceSBarry Smith 
15515615d1e5SSatish Balay #undef __FUNC__
1552d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory"
1553c9005455SLois Curfman McInnes /*@
1554c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1555c9005455SLois Curfman McInnes 
1556c9005455SLois Curfman McInnes    Input Parameters:
1557c9005455SLois Curfman McInnes .  snes - iterative context obtained from SNESCreate()
1558c9005455SLois Curfman McInnes .  a   - array to hold history
1559c9005455SLois Curfman McInnes .  na  - size of a
1560c9005455SLois Curfman McInnes 
1561fee21e36SBarry Smith    Collective on SNES
1562fee21e36SBarry Smith 
1563c9005455SLois Curfman McInnes    Notes:
1564c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1565c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1566c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1567c9005455SLois Curfman McInnes    at each step.
1568c9005455SLois Curfman McInnes 
1569c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1570c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1571c9005455SLois Curfman McInnes    during the section of code that is being timed.
1572c9005455SLois Curfman McInnes 
1573c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1574c9005455SLois Curfman McInnes @*/
1575c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na)
1576c9005455SLois Curfman McInnes {
15773a40ed3dSBarry Smith   PetscFunctionBegin;
1578c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1579c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1580c9005455SLois Curfman McInnes   snes->conv_hist      = a;
1581c9005455SLois Curfman McInnes   snes->conv_hist_size = na;
15823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1583c9005455SLois Curfman McInnes }
1584c9005455SLois Curfman McInnes 
1585c9005455SLois Curfman McInnes #undef __FUNC__
15865615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private"
15879b94acceSBarry Smith /*
15889b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
15899b94acceSBarry Smith    positive parameter delta.
15909b94acceSBarry Smith 
15919b94acceSBarry Smith     Input Parameters:
15929b94acceSBarry Smith .   snes - the SNES context
15939b94acceSBarry Smith .   y - approximate solution of linear system
15949b94acceSBarry Smith .   fnorm - 2-norm of current function
15959b94acceSBarry Smith .   delta - trust region size
15969b94acceSBarry Smith 
15979b94acceSBarry Smith     Output Parameters:
15989b94acceSBarry Smith .   gpnorm - predicted function norm at the new point, assuming local
15999b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
16009b94acceSBarry Smith     region, and exceeds zero otherwise.
16019b94acceSBarry Smith .   ynorm - 2-norm of the step
16029b94acceSBarry Smith 
16039b94acceSBarry Smith     Note:
160440191667SLois Curfman McInnes     For non-trust region methods such as SNES_EQ_LS, the parameter delta
16059b94acceSBarry Smith     is set to be the maximum allowable step size.
16069b94acceSBarry Smith 
16079b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
16089b94acceSBarry Smith */
16099b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
16109b94acceSBarry Smith                           double *gpnorm,double *ynorm)
16119b94acceSBarry Smith {
16129b94acceSBarry Smith   double norm;
16139b94acceSBarry Smith   Scalar cnorm;
16143a40ed3dSBarry Smith   int    ierr;
16153a40ed3dSBarry Smith 
16163a40ed3dSBarry Smith   PetscFunctionBegin;
16173a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
16189b94acceSBarry Smith   if (norm > *delta) {
16199b94acceSBarry Smith      norm = *delta/norm;
16209b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
16219b94acceSBarry Smith      cnorm = norm;
16229b94acceSBarry Smith      VecScale( &cnorm, y );
16239b94acceSBarry Smith      *ynorm = *delta;
16249b94acceSBarry Smith   } else {
16259b94acceSBarry Smith      *gpnorm = 0.0;
16269b94acceSBarry Smith      *ynorm = norm;
16279b94acceSBarry Smith   }
16283a40ed3dSBarry Smith   PetscFunctionReturn(0);
16299b94acceSBarry Smith }
16309b94acceSBarry Smith 
16315615d1e5SSatish Balay #undef __FUNC__
16325615d1e5SSatish Balay #define __FUNC__ "SNESSolve"
16339b94acceSBarry Smith /*@
16349b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
163563a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
16369b94acceSBarry Smith 
1637*b2002411SLois Curfman McInnes    Input Parameters:
16389b94acceSBarry Smith .  snes - the SNES context
16398ddd3da0SLois Curfman McInnes .  x - the solution vector
16409b94acceSBarry Smith 
16419b94acceSBarry Smith    Output Parameter:
1642*b2002411SLois Curfman McInnes .  its - number of iterations until termination
16439b94acceSBarry Smith 
1644fee21e36SBarry Smith    Collective on SNES
1645fee21e36SBarry Smith 
1646*b2002411SLois Curfman McInnes    Notes:
16478ddd3da0SLois Curfman McInnes    The user should initialize the vector, x, with the initial guess
16488ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
16498ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
16508ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
16518ddd3da0SLois Curfman McInnes 
16529b94acceSBarry Smith .keywords: SNES, nonlinear, solve
16539b94acceSBarry Smith 
165463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
16559b94acceSBarry Smith @*/
16568ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
16579b94acceSBarry Smith {
16583c7409f5SSatish Balay   int ierr, flg;
1659052efed2SBarry Smith 
16603a40ed3dSBarry Smith   PetscFunctionBegin;
166177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
166274679c65SBarry Smith   PetscValidIntPointer(its);
166382bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);}
1664c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
16659b94acceSBarry Smith   PLogEventBegin(SNES_Solve,snes,0,0,0);
1666c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
16679b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr);
16689b94acceSBarry Smith   PLogEventEnd(SNES_Solve,snes,0,0,0);
16693c7409f5SSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr);
16706d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); }
16713a40ed3dSBarry Smith   PetscFunctionReturn(0);
16729b94acceSBarry Smith }
16739b94acceSBarry Smith 
16749b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
16759b94acceSBarry Smith 
16765615d1e5SSatish Balay #undef __FUNC__
16775615d1e5SSatish Balay #define __FUNC__ "SNESSetType"
167882bf6240SBarry Smith /*@C
16794b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
16809b94acceSBarry Smith 
16819b94acceSBarry Smith    Input Parameters:
16829b94acceSBarry Smith .  snes - the SNES context
16839b94acceSBarry Smith .  method - a known method
16849b94acceSBarry Smith 
1685fee21e36SBarry Smith    Collective on SNES
1686fee21e36SBarry Smith 
1687ae12b187SLois Curfman McInnes   Options Database Command:
1688ae12b187SLois Curfman McInnes $ -snes_type  <method>
1689ae12b187SLois Curfman McInnes $    Use -help for a list of available methods
1690ae12b187SLois Curfman McInnes $    (for instance, ls or tr)
1691ae12b187SLois Curfman McInnes 
16929b94acceSBarry Smith    Notes:
16939b94acceSBarry Smith    See "petsc/include/snes.h" for available methods (for instance)
16949b94acceSBarry Smith $  Systems of nonlinear equations:
169540191667SLois Curfman McInnes $    SNES_EQ_LS - Newton's method with line search
169640191667SLois Curfman McInnes $    SNES_EQ_TR - Newton's method with trust region
16979b94acceSBarry Smith $  Unconstrained minimization:
169840191667SLois Curfman McInnes $    SNES_UM_TR - Newton's method with trust region
169940191667SLois Curfman McInnes $    SNES_UM_LS - Newton's method with line search
17009b94acceSBarry Smith 
1701ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1702ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1703ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1704ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1705ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1706ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1707ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1708ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1709ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1710ae12b187SLois Curfman McInnes   appropriate method.  In other words, this routine is for the advanced user.
1711a703fe33SLois Curfman McInnes 
1712f536c699SSatish Balay .keywords: SNES, set, method
17139b94acceSBarry Smith @*/
17144b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method)
17159b94acceSBarry Smith {
171684cb2905SBarry Smith   int ierr;
17179b94acceSBarry Smith   int (*r)(SNES);
17183a40ed3dSBarry Smith 
17193a40ed3dSBarry Smith   PetscFunctionBegin;
172077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
172182bf6240SBarry Smith 
172282bf6240SBarry Smith   if (!PetscStrcmp(snes->type_name,method)) PetscFunctionReturn(0);
172382bf6240SBarry Smith 
172482bf6240SBarry Smith   if (snes->setupcalled) {
1725e1311b90SBarry Smith     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
172682bf6240SBarry Smith     snes->data = 0;
172782bf6240SBarry Smith   }
17287d1a2b2bSBarry Smith 
17299b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
173082bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);}
173182bf6240SBarry Smith 
1732ecf371e4SBarry Smith   ierr =  DLRegisterFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr);
173382bf6240SBarry Smith 
17340452661fSBarry Smith   if (snes->data) PetscFree(snes->data);
173582bf6240SBarry Smith   snes->data = 0;
17363a40ed3dSBarry Smith   ierr = (*r)(snes); CHKERRQ(ierr);
173782bf6240SBarry Smith 
173882bf6240SBarry Smith   if (snes->type_name) PetscFree(snes->type_name);
173982bf6240SBarry Smith   snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name);
174082bf6240SBarry Smith   PetscStrcpy(snes->type_name,method);
174182bf6240SBarry Smith   snes->set_method_called = 1;
174282bf6240SBarry Smith 
17433a40ed3dSBarry Smith   PetscFunctionReturn(0);
17449b94acceSBarry Smith }
17459b94acceSBarry Smith 
1746a847f771SSatish Balay 
17479b94acceSBarry Smith /* --------------------------------------------------------------------- */
17485615d1e5SSatish Balay #undef __FUNC__
1749d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy"
17509b94acceSBarry Smith /*@C
17519b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
17529b94acceSBarry Smith    registered by SNESRegister().
17539b94acceSBarry Smith 
1754fee21e36SBarry Smith    Not Collective
1755fee21e36SBarry Smith 
17569b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
17579b94acceSBarry Smith 
17589b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
17599b94acceSBarry Smith @*/
1760cf256101SBarry Smith int SNESRegisterDestroy(void)
17619b94acceSBarry Smith {
176282bf6240SBarry Smith   int ierr;
176382bf6240SBarry Smith 
17643a40ed3dSBarry Smith   PetscFunctionBegin;
176582bf6240SBarry Smith   if (SNESList) {
176682bf6240SBarry Smith     ierr = DLRegisterDestroy( SNESList );CHKERRQ(ierr);
176782bf6240SBarry Smith     SNESList = 0;
17689b94acceSBarry Smith   }
176984cb2905SBarry Smith   SNESRegisterAllCalled = 0;
17703a40ed3dSBarry Smith   PetscFunctionReturn(0);
17719b94acceSBarry Smith }
17729b94acceSBarry Smith 
17735615d1e5SSatish Balay #undef __FUNC__
1774d4bb536fSBarry Smith #define __FUNC__ "SNESGetType"
17759b94acceSBarry Smith /*@C
17769a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
17779b94acceSBarry Smith 
17789b94acceSBarry Smith    Input Parameter:
17794b0e389bSBarry Smith .  snes - nonlinear solver context
17809b94acceSBarry Smith 
17819b94acceSBarry Smith    Output Parameter:
178282bf6240SBarry Smith .  method - SNES method (a charactor string)
17839b94acceSBarry Smith 
1784fee21e36SBarry Smith    Not Collective
1785fee21e36SBarry Smith 
17869b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name
17879b94acceSBarry Smith @*/
178882bf6240SBarry Smith int SNESGetType(SNES snes, SNESType *method)
17899b94acceSBarry Smith {
17903a40ed3dSBarry Smith   PetscFunctionBegin;
179182bf6240SBarry Smith   *method = snes->type_name;
17923a40ed3dSBarry Smith   PetscFunctionReturn(0);
17939b94acceSBarry Smith }
17949b94acceSBarry Smith 
17955615d1e5SSatish Balay #undef __FUNC__
1796d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution"
17979b94acceSBarry Smith /*@C
17989b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
17999b94acceSBarry Smith    stored.
18009b94acceSBarry Smith 
18019b94acceSBarry Smith    Input Parameter:
18029b94acceSBarry Smith .  snes - the SNES context
18039b94acceSBarry Smith 
18049b94acceSBarry Smith    Output Parameter:
18059b94acceSBarry Smith .  x - the solution
18069b94acceSBarry Smith 
1807fee21e36SBarry Smith    Not Collective, but Vec is parallel if SNES is parallel
1808fee21e36SBarry Smith 
18099b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
18109b94acceSBarry Smith 
1811a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
18129b94acceSBarry Smith @*/
18139b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
18149b94acceSBarry Smith {
18153a40ed3dSBarry Smith   PetscFunctionBegin;
181677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18179b94acceSBarry Smith   *x = snes->vec_sol_always;
18183a40ed3dSBarry Smith   PetscFunctionReturn(0);
18199b94acceSBarry Smith }
18209b94acceSBarry Smith 
18215615d1e5SSatish Balay #undef __FUNC__
1822d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate"
18239b94acceSBarry Smith /*@C
18249b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
18259b94acceSBarry Smith    stored.
18269b94acceSBarry Smith 
18279b94acceSBarry Smith    Input Parameter:
18289b94acceSBarry Smith .  snes - the SNES context
18299b94acceSBarry Smith 
18309b94acceSBarry Smith    Output Parameter:
18319b94acceSBarry Smith .  x - the solution update
18329b94acceSBarry Smith 
1833fee21e36SBarry Smith    Not Collective, but Vec is parallel if SNES is parallel
18349b94acceSBarry Smith 
18359b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
18369b94acceSBarry Smith 
18379b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
18389b94acceSBarry Smith @*/
18399b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
18409b94acceSBarry Smith {
18413a40ed3dSBarry Smith   PetscFunctionBegin;
184277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18439b94acceSBarry Smith   *x = snes->vec_sol_update_always;
18443a40ed3dSBarry Smith   PetscFunctionReturn(0);
18459b94acceSBarry Smith }
18469b94acceSBarry Smith 
18475615d1e5SSatish Balay #undef __FUNC__
1848d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction"
18499b94acceSBarry Smith /*@C
18503638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
18519b94acceSBarry Smith 
18529b94acceSBarry Smith    Input Parameter:
18539b94acceSBarry Smith .  snes - the SNES context
18549b94acceSBarry Smith 
18559b94acceSBarry Smith    Output Parameter:
18563638b69dSLois Curfman McInnes .  r - the function
18579b94acceSBarry Smith 
1858fee21e36SBarry Smith    Not Collective, but Vec is parallel if SNES is parallel
1859fee21e36SBarry Smith 
18609b94acceSBarry Smith    Notes:
18619b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
18629b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
18639b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
18649b94acceSBarry Smith 
1865a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
18669b94acceSBarry Smith 
186731615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
186831615d04SLois Curfman McInnes           SNESGetGradient()
18699b94acceSBarry Smith @*/
18709b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r)
18719b94acceSBarry Smith {
18723a40ed3dSBarry Smith   PetscFunctionBegin;
187377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1874a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1875a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1876a8c6a408SBarry Smith   }
18779b94acceSBarry Smith   *r = snes->vec_func_always;
18783a40ed3dSBarry Smith   PetscFunctionReturn(0);
18799b94acceSBarry Smith }
18809b94acceSBarry Smith 
18815615d1e5SSatish Balay #undef __FUNC__
1882d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient"
18839b94acceSBarry Smith /*@C
18843638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
18859b94acceSBarry Smith 
18869b94acceSBarry Smith    Input Parameter:
18879b94acceSBarry Smith .  snes - the SNES context
18889b94acceSBarry Smith 
18899b94acceSBarry Smith    Output Parameter:
18909b94acceSBarry Smith .  r - the gradient
18919b94acceSBarry Smith 
1892fee21e36SBarry Smith    Not Collective, but Vec is parallel if SNES is parallel
1893fee21e36SBarry Smith 
18949b94acceSBarry Smith    Notes:
18959b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
18969b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
18979b94acceSBarry Smith    SNESGetFunction().
18989b94acceSBarry Smith 
18999b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
19009b94acceSBarry Smith 
190131615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction()
19029b94acceSBarry Smith @*/
19039b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r)
19049b94acceSBarry Smith {
19053a40ed3dSBarry Smith   PetscFunctionBegin;
190677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
19073a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1908a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
19093a40ed3dSBarry Smith   }
19109b94acceSBarry Smith   *r = snes->vec_func_always;
19113a40ed3dSBarry Smith   PetscFunctionReturn(0);
19129b94acceSBarry Smith }
19139b94acceSBarry Smith 
19145615d1e5SSatish Balay #undef __FUNC__
1915d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction"
19169b94acceSBarry Smith /*@
19179b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
19189b94acceSBarry Smith    unconstrained minimization problems.
19199b94acceSBarry Smith 
19209b94acceSBarry Smith    Input Parameter:
19219b94acceSBarry Smith .  snes - the SNES context
19229b94acceSBarry Smith 
19239b94acceSBarry Smith    Output Parameter:
19249b94acceSBarry Smith .  r - the function
19259b94acceSBarry Smith 
1926fee21e36SBarry Smith    Not Collective, but Vec is parallel if SNES is parallel
1927fee21e36SBarry Smith 
19289b94acceSBarry Smith    Notes:
19299b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
19309b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
19319b94acceSBarry Smith    SNESGetFunction().
19329b94acceSBarry Smith 
19339b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
19349b94acceSBarry Smith 
193531615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction()
19369b94acceSBarry Smith @*/
19379b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r)
19389b94acceSBarry Smith {
19393a40ed3dSBarry Smith   PetscFunctionBegin;
194077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
194174679c65SBarry Smith   PetscValidScalarPointer(r);
19423a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1943a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
19443a40ed3dSBarry Smith   }
19459b94acceSBarry Smith   *r = snes->fc;
19463a40ed3dSBarry Smith   PetscFunctionReturn(0);
19479b94acceSBarry Smith }
19489b94acceSBarry Smith 
19495615d1e5SSatish Balay #undef __FUNC__
1950d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix"
19513c7409f5SSatish Balay /*@C
19523c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1953d850072dSLois Curfman McInnes    SNES options in the database.
19543c7409f5SSatish Balay 
19553c7409f5SSatish Balay    Input Parameter:
19563c7409f5SSatish Balay .  snes - the SNES context
19573c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
19583c7409f5SSatish Balay 
1959fee21e36SBarry Smith    Collective on SNES
1960fee21e36SBarry Smith 
1961d850072dSLois Curfman McInnes    Notes:
1962a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1963a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1964a83b1b31SSatish Balay    hyphen.
1965d850072dSLois Curfman McInnes 
19663c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1967a86d99e1SLois Curfman McInnes 
1968a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
19693c7409f5SSatish Balay @*/
19703c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
19713c7409f5SSatish Balay {
19723c7409f5SSatish Balay   int ierr;
19733c7409f5SSatish Balay 
19743a40ed3dSBarry Smith   PetscFunctionBegin;
197577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1976639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
19773c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
19783a40ed3dSBarry Smith   PetscFunctionReturn(0);
19793c7409f5SSatish Balay }
19803c7409f5SSatish Balay 
19815615d1e5SSatish Balay #undef __FUNC__
1982d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix"
19833c7409f5SSatish Balay /*@C
1984f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1985d850072dSLois Curfman McInnes    SNES options in the database.
19863c7409f5SSatish Balay 
19873c7409f5SSatish Balay    Input Parameter:
19883c7409f5SSatish Balay .  snes - the SNES context
19893c7409f5SSatish Balay .  prefix - the prefix to prepend to all option names
19903c7409f5SSatish Balay 
1991fee21e36SBarry Smith    Collective on SNES
1992fee21e36SBarry Smith 
1993d850072dSLois Curfman McInnes    Notes:
1994a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1995a83b1b31SSatish Balay    The first character of all runtime options is AUTOMATICALLY the
1996a83b1b31SSatish Balay    hyphen.
1997d850072dSLois Curfman McInnes 
19983c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
1999a86d99e1SLois Curfman McInnes 
2000a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20013c7409f5SSatish Balay @*/
20023c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
20033c7409f5SSatish Balay {
20043c7409f5SSatish Balay   int ierr;
20053c7409f5SSatish Balay 
20063a40ed3dSBarry Smith   PetscFunctionBegin;
200777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2008639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20093c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
20103a40ed3dSBarry Smith   PetscFunctionReturn(0);
20113c7409f5SSatish Balay }
20123c7409f5SSatish Balay 
20135615d1e5SSatish Balay #undef __FUNC__
2014d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix"
2015c83590e2SSatish Balay /*@
20163c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20173c7409f5SSatish Balay    SNES options in the database.
20183c7409f5SSatish Balay 
20193c7409f5SSatish Balay    Input Parameter:
20203c7409f5SSatish Balay .  snes - the SNES context
20213c7409f5SSatish Balay 
20223c7409f5SSatish Balay    Output Parameter:
20233c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20243c7409f5SSatish Balay 
2025fee21e36SBarry Smith    Not Collective
2026fee21e36SBarry Smith 
20273c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2028a86d99e1SLois Curfman McInnes 
2029a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20303c7409f5SSatish Balay @*/
20313c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
20323c7409f5SSatish Balay {
20333c7409f5SSatish Balay   int ierr;
20343c7409f5SSatish Balay 
20353a40ed3dSBarry Smith   PetscFunctionBegin;
203677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2037639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr);
20383a40ed3dSBarry Smith   PetscFunctionReturn(0);
20393c7409f5SSatish Balay }
20403c7409f5SSatish Balay 
2041ca161407SBarry Smith #undef __FUNC__
2042ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp"
2043ca161407SBarry Smith /*@
2044ca161407SBarry Smith    SNESPrintHelp - Prints all options for the SNES component.
2045ca161407SBarry Smith 
2046ca161407SBarry Smith    Input Parameter:
2047ca161407SBarry Smith .  snes - the SNES context
2048ca161407SBarry Smith 
2049fee21e36SBarry Smith    Collective on SNES
2050fee21e36SBarry Smith 
2051ca161407SBarry Smith    Options Database Keys:
2052ca161407SBarry Smith $  -help, -h
2053ca161407SBarry Smith 
2054ca161407SBarry Smith .keywords: SNES, nonlinear, help
2055ca161407SBarry Smith 
2056ca161407SBarry Smith .seealso: SNESSetFromOptions()
2057ca161407SBarry Smith @*/
2058ca161407SBarry Smith int SNESPrintHelp(SNES snes)
2059ca161407SBarry Smith {
2060ca161407SBarry Smith   char                p[64];
2061ca161407SBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
2062ca161407SBarry Smith   int                 ierr;
2063ca161407SBarry Smith 
2064ca161407SBarry Smith   PetscFunctionBegin;
2065ca161407SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2066ca161407SBarry Smith 
2067ca161407SBarry Smith   PetscStrcpy(p,"-");
2068ca161407SBarry Smith   if (snes->prefix) PetscStrcat(p, snes->prefix);
2069ca161407SBarry Smith 
2070ca161407SBarry Smith   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2071ca161407SBarry Smith 
207276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");
207382bf6240SBarry Smith   ierr = DLRegisterPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr);
207476be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);
207576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);
207676be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);
207776be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);
207876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);
207976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);
208076be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);
208176be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");
208276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);
208376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
2084ca161407SBarry Smith     residual norm at each iteration.\n",p);
208576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
2086ca161407SBarry Smith     residual norm for small residual norms. This is useful to conceal\n\
2087ca161407SBarry Smith     meaningless digits that may be different on different machines.\n",p);
208876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);
2089ca161407SBarry Smith   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
209076be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2091ca161407SBarry Smith      " Options for solving systems of nonlinear equations only:\n");
209276be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);
209376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);
209476be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);
209576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);
209676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);
209776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2098ca161407SBarry Smith      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);
209976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2100ca161407SBarry Smith      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);
210176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2102ca161407SBarry Smith      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);
210376be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2104ca161407SBarry Smith      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);
210576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2106ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);
210776be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2108ca161407SBarry Smith      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);
210976be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2110ca161407SBarry Smith      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);
2111ca161407SBarry Smith   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
211276be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,
2113ca161407SBarry Smith      " Options for solving unconstrained minimization problems only:\n");
211476be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);
211576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);
211676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);
2117ca161407SBarry Smith   }
211876be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p);
211976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"a particular method\n");
2120ca161407SBarry Smith   if (snes->printhelp) {
2121ca161407SBarry Smith     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2122ca161407SBarry Smith   }
2123ca161407SBarry Smith   PetscFunctionReturn(0);
2124ca161407SBarry Smith }
21253c7409f5SSatish Balay 
2126*b2002411SLois Curfman McInnes /*M
2127*b2002411SLois Curfman McInnes    SNESRegister - Adds a method to the nonlinear solver package.
21289b94acceSBarry Smith 
2129*b2002411SLois Curfman McInnes    Synopsis:
2130*b2002411SLois Curfman McInnes    SNESRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
21319b94acceSBarry Smith 
2132*b2002411SLois Curfman McInnes    Input Parameters:
2133*b2002411SLois Curfman McInnes .  name_solver - name of a new user-defined solver
2134*b2002411SLois Curfman McInnes .  path - path (either absolute or relative) the library containing this solver
2135*b2002411SLois Curfman McInnes .  name_create - name of routine to create method context
2136*b2002411SLois Curfman McInnes .  routine_create - routine to create method context
21379b94acceSBarry Smith 
2138*b2002411SLois Curfman McInnes    Notes:
2139*b2002411SLois Curfman McInnes    SNESRegister() may be called multiple times to add several user-defined solvers.
21403a40ed3dSBarry Smith 
2141*b2002411SLois Curfman McInnes    If dynamic libraries are used, then the fourth input argument (routine_create)
2142*b2002411SLois Curfman McInnes    is ignored.
2143*b2002411SLois Curfman McInnes 
2144*b2002411SLois Curfman McInnes    Sample usage:
2145*b2002411SLois Curfman McInnes    SNESRegister("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2146*b2002411SLois Curfman McInnes                 "MySolverCreate",MySolverCreate);
2147*b2002411SLois Curfman McInnes 
2148*b2002411SLois Curfman McInnes    Then, your solver can be chosen with the procedural interface via
2149*b2002411SLois Curfman McInnes $     SNESSetType(snes,"my_solver")
2150*b2002411SLois Curfman McInnes    or at runtime via the option
2151*b2002411SLois Curfman McInnes $     -snes_type my_solver
2152*b2002411SLois Curfman McInnes 
2153*b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register
2154*b2002411SLois Curfman McInnes 
2155*b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2156*b2002411SLois Curfman McInnes M*/
2157*b2002411SLois Curfman McInnes 
2158*b2002411SLois Curfman McInnes #undef __FUNC__
2159*b2002411SLois Curfman McInnes #define __FUNC__ "SNESRegister_Private"
2160*b2002411SLois Curfman McInnes int SNESRegister_Private(char *sname,char *path,char *name,int (*function)(SNES))
2161*b2002411SLois Curfman McInnes {
2162*b2002411SLois Curfman McInnes   char fullname[256];
2163*b2002411SLois Curfman McInnes   int  ierr;
2164*b2002411SLois Curfman McInnes 
2165*b2002411SLois Curfman McInnes   PetscFunctionBegin;
2166*b2002411SLois Curfman McInnes   PetscStrcpy(fullname,path); PetscStrcat(fullname,":");PetscStrcat(fullname,name);
2167*b2002411SLois Curfman McInnes   ierr = DLRegister_Private(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr);
2168*b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2169*b2002411SLois Curfman McInnes }
2170