xref: /petsc/src/snes/mf/snesmfj.c (revision c481317fb29dd7d10f52c4080c7dd13f6ffd59ed)
181e6777dSBarry Smith 
281e6777dSBarry Smith #ifndef lint
3*c481317fSBarry Smith static char vcid[] = "$Id: snesmfj.c,v 1.49 1997/04/10 00:06:10 bsmith Exp bsmith $";
481e6777dSBarry Smith #endif
581e6777dSBarry Smith 
670f55243SBarry Smith #include "src/snes/snesimpl.h"   /*I  "snes.h"   I*/
7eb9086c3SLois Curfman McInnes #include "pinclude/pviewer.h"
8*c481317fSBarry Smith #include <math.h>
981e6777dSBarry Smith 
1050361f65SLois Curfman McInnes typedef struct {  /* default context for matrix-free SNES */
11eb9086c3SLois Curfman McInnes   SNES        snes;      /* SNES context */
12eb9086c3SLois Curfman McInnes   Vec         w;         /* work vector */
13eb9086c3SLois Curfman McInnes   PCNullSpace sp;        /* null space context */
14eb9086c3SLois Curfman McInnes   double      error_rel; /* square root of relative error in computing function */
15639f9d9dSBarry Smith   double      umin;      /* minimum allowable u'a value relative to |u|_1 */
1639e2f89bSBarry Smith } MFCtx_Private;
1739e2f89bSBarry Smith 
185615d1e5SSatish Balay #undef __FUNC__
195eea60f9SBarry Smith #define __FUNC__ "SNESMatrixFreeDestroy_Private" /* ADIC Ignore */
20fae171e0SBarry Smith int SNESMatrixFreeDestroy_Private(PetscObject obj)
21b9fa9cd0SBarry Smith {
22b9fa9cd0SBarry Smith   int           ierr;
23fae171e0SBarry Smith   Mat           mat = (Mat) obj;
24fae171e0SBarry Smith   MFCtx_Private *ctx;
25fae171e0SBarry Smith 
26fae171e0SBarry Smith   ierr = MatShellGetContext(mat,(void **)&ctx);
27b9fa9cd0SBarry Smith   ierr = VecDestroy(ctx->w); CHKERRQ(ierr);
28b4fd4287SBarry Smith   if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp); CHKERRQ(ierr);}
29fae171e0SBarry Smith   PetscFree(ctx);
30b9fa9cd0SBarry Smith   return 0;
31b9fa9cd0SBarry Smith }
3250361f65SLois Curfman McInnes 
335615d1e5SSatish Balay #undef __FUNC__
345eea60f9SBarry Smith #define __FUNC__ "SNESMatrixFreeView_Private" /* ADIC Ignore */
3539e2f89bSBarry Smith /*
361d1bbde2SLois Curfman McInnes    SNESMatrixFreeView_Private - Views matrix-free parameters.
3739e2f89bSBarry Smith  */
38eb9086c3SLois Curfman McInnes int SNESMatrixFreeView_Private(Mat J,Viewer viewer)
39eb9086c3SLois Curfman McInnes {
40eb9086c3SLois Curfman McInnes   int           ierr;
41eb9086c3SLois Curfman McInnes   MFCtx_Private *ctx;
42eb9086c3SLois Curfman McInnes   MPI_Comm      comm;
43eb9086c3SLois Curfman McInnes   FILE          *fd;
44eb9086c3SLois Curfman McInnes   ViewerType    vtype;
45eb9086c3SLois Curfman McInnes 
46eb9086c3SLois Curfman McInnes   PetscObjectGetComm((PetscObject)J,&comm);
47eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
48eb9086c3SLois Curfman McInnes   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
49eb9086c3SLois Curfman McInnes   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
50eb9086c3SLois Curfman McInnes   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
51eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");
52eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
53eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"    umin=%g (minimum iterate parameter)\n",ctx->umin);
54eb9086c3SLois Curfman McInnes   }
55eb9086c3SLois Curfman McInnes   return 0;
56eb9086c3SLois Curfman McInnes }
57eb9086c3SLois Curfman McInnes 
58*c481317fSBarry Smith extern int VecDot_Seq(Vec,Vec,double *);
59*c481317fSBarry Smith extern int VecNorm_Seq(Vec,NormType,double *);
60*c481317fSBarry Smith 
615615d1e5SSatish Balay #undef __FUNC__
625615d1e5SSatish Balay #define __FUNC__ "SNESMatrixFreeMult_Private"
63eb9086c3SLois Curfman McInnes /*
64eb9086c3SLois Curfman McInnes   SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector
65eb9086c3SLois Curfman McInnes   product, y = F'(u)*a:
66eb9086c3SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
67eb9086c3SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
68eb9086c3SLois Curfman McInnes         u = current iterate
69eb9086c3SLois Curfman McInnes         h = difference interval
70eb9086c3SLois Curfman McInnes */
71eb9086c3SLois Curfman McInnes int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y)
7239e2f89bSBarry Smith {
73fae171e0SBarry Smith   MFCtx_Private *ctx;
74fae171e0SBarry Smith   SNES          snes;
75*c481317fSBarry Smith   double        ovalues[3],values[3],norm, sum, umin;
765334005bSBarry Smith   Scalar        h, dot, mone = -1.0;
77fae171e0SBarry Smith   Vec           w,U,F;
7883e56afdSLois Curfman McInnes   int           ierr, (*eval_fct)(SNES,Vec,Vec);
79*c481317fSBarry Smith   MPI_Comm      comm;
8039e2f89bSBarry Smith 
81*c481317fSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
826e38a044SBarry Smith   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
83fae171e0SBarry Smith   snes = ctx->snes;
84fae171e0SBarry Smith   w    = ctx->w;
85eb9086c3SLois Curfman McInnes   umin = ctx->umin;
86fae171e0SBarry Smith 
8757c37382SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
8857c37382SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
8957c37382SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
9057c37382SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
9157c37382SLois Curfman McInnes 
9278b31e54SBarry Smith   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
9383e56afdSLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
9483e56afdSLois Curfman McInnes     eval_fct = SNESComputeFunction;
9578b31e54SBarry Smith     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
9683e56afdSLois Curfman McInnes   }
9783e56afdSLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
9883e56afdSLois Curfman McInnes     eval_fct = SNESComputeGradient;
9983e56afdSLois Curfman McInnes     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
10083e56afdSLois Curfman McInnes   }
101e3372554SBarry Smith   else SETERRQ(1,0,"Invalid method class");
10250361f65SLois Curfman McInnes 
103eb9086c3SLois Curfman McInnes   /* Determine a "good" step size, h */
104*c481317fSBarry Smith 
105*c481317fSBarry Smith   /*
106eb9086c3SLois Curfman McInnes     ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
107eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
108eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
109*c481317fSBarry Smith   */
110*c481317fSBarry Smith 
111*c481317fSBarry Smith   /*
112*c481317fSBarry Smith      Call the Seq Vector routines and then do
113*c481317fSBarry Smith     a single reduction to reduce the number of communications
114*c481317fSBarry Smith      required
115*c481317fSBarry Smith   */
116*c481317fSBarry Smith 
117*c481317fSBarry Smith   PLogEventBegin(VEC_Dot,U,a,0,0);
118*c481317fSBarry Smith   ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr);
119*c481317fSBarry Smith   PLogEventEnd(VEC_Dot,U,a,0,0);
120*c481317fSBarry Smith   PLogEventBegin(VEC_Norm,a,0,0,0);
121*c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
122*c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
123*c481317fSBarry Smith   ovalues[2] = ovalues[2]*ovalues[2];
124*c481317fSBarry Smith   MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );
125*c481317fSBarry Smith   dot = values[0]; sum = values[1]; norm = sqrt(values[2]);
126*c481317fSBarry Smith   PLogEventEnd(VEC_Norm,a,0,0,0);
127*c481317fSBarry Smith 
128*c481317fSBarry Smith   /*
129*c481317fSBarry Smith      the plogeventbegin() below should really be above,
130*c481317fSBarry Smith      but they cannot be nested so it excludes the time
131*c481317fSBarry Smith      to compute h
132*c481317fSBarry Smith   */
133*c481317fSBarry Smith   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
134eb9086c3SLois Curfman McInnes 
135eb9086c3SLois Curfman McInnes   /* Safeguard for step sizes too small */
136edd2f0e1SBarry Smith   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
13719a167f6SBarry Smith #if defined(PETSC_COMPLEX)
138eb9086c3SLois Curfman McInnes   else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum;
139eb9086c3SLois Curfman McInnes   else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum;
14019a167f6SBarry Smith #else
141eb9086c3SLois Curfman McInnes   else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
142eb9086c3SLois Curfman McInnes   else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
14319a167f6SBarry Smith #endif
144eb9086c3SLois Curfman McInnes   h = ctx->error_rel*dot/(norm*norm);
14539e2f89bSBarry Smith 
146eb9086c3SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
147eb9086c3SLois Curfman McInnes   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
14883e56afdSLois Curfman McInnes   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
149195ee392SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
1505334005bSBarry Smith   h = 1.0/h;
151195ee392SLois Curfman McInnes   ierr = VecScale(&h,y); CHKERRQ(ierr);
152b4fd4287SBarry Smith   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
15357c37382SLois Curfman McInnes 
154eb9086c3SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
15539e2f89bSBarry Smith   return 0;
15639e2f89bSBarry Smith }
15783e56afdSLois Curfman McInnes 
1585615d1e5SSatish Balay #undef __FUNC__
1595615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatCreate"
1604b828684SBarry Smith /*@C
1615392566eSBarry Smith    SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix
16250361f65SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
16350361f65SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
1645392566eSBarry Smith 
1655392566eSBarry Smith    Input Parameters:
166eb9086c3SLois Curfman McInnes .  snes - the SNES context
1675392566eSBarry Smith .  x - vector where SNES solution is to be stored.
1685392566eSBarry Smith 
169eb9086c3SLois Curfman McInnes    Output Parameter:
1705392566eSBarry Smith .  J - the matrix-free matrix
1715392566eSBarry Smith 
17265afa06eSLois Curfman McInnes    Notes:
17365afa06eSLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
17465afa06eSLois Curfman McInnes    and work space for performing finite difference approximations of
1753d5db8e1SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
1763d5db8e1SLois Curfman McInnes 
1773d5db8e1SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
1783d5db8e1SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
1793d5db8e1SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
1803d5db8e1SLois Curfman McInnes $   where
1813d5db8e1SLois Curfman McInnes $        error_rel = square root of relative error in
1823d5db8e1SLois Curfman McInnes $                    function evaluation
1833d5db8e1SLois Curfman McInnes $        umin = minimum iterate parameter
1843d5db8e1SLois Curfman McInnes 
1853d5db8e1SLois Curfman McInnes    The user can set these parameters via SNESSetMatrixFreeParameters().
1863d5db8e1SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
18765afa06eSLois Curfman McInnes 
18865afa06eSLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
18965afa06eSLois Curfman McInnes    matrix context.
19065afa06eSLois Curfman McInnes 
1913d5db8e1SLois Curfman McInnes    Options Database Keys:
1923d5db8e1SLois Curfman McInnes $  -snes_mf_err <error_rel>
1933d5db8e1SLois Curfman McInnes $  -snes_mf_unim <umin>
1943d5db8e1SLois Curfman McInnes 
19565afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
19665afa06eSLois Curfman McInnes 
1973d5db8e1SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters()
1985392566eSBarry Smith @*/
1995392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J)
2005392566eSBarry Smith {
2015392566eSBarry Smith   MPI_Comm      comm;
2025392566eSBarry Smith   MFCtx_Private *mfctx;
203eb9086c3SLois Curfman McInnes   int           n, nloc, ierr, flg;
2043e966953SLois Curfman McInnes   char          p[64];
2055392566eSBarry Smith 
2060452661fSBarry Smith   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
207464493b3SBarry Smith   PLogObjectMemory(snes,sizeof(MFCtx_Private));
208b4fd4287SBarry Smith   mfctx->sp   = 0;
2095392566eSBarry Smith   mfctx->snes = snes;
210eb9086c3SLois Curfman McInnes   mfctx->error_rel = 1.e-8; /* assumes double precision */
211639f9d9dSBarry Smith   mfctx->umin      = 1.e-6;
212eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
213eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
214639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
2153e966953SLois Curfman McInnes   PetscStrcpy(p,"-");
2163e966953SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
217639f9d9dSBarry Smith   if (flg) {
2183e966953SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
2193e966953SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin);
220639f9d9dSBarry Smith   }
2215392566eSBarry Smith   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
222195ee392SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
223195ee392SLois Curfman McInnes   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
2247ddc982cSLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
225029af93fSBarry Smith   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
2261c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr);
2271c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr);
2281c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr);
229b9fa9cd0SBarry Smith   PLogObjectParent(*J,mfctx->w);
230d370d78aSBarry Smith   PLogObjectParent(snes,*J);
23181e6777dSBarry Smith   return 0;
23281e6777dSBarry Smith }
23381e6777dSBarry Smith 
2345615d1e5SSatish Balay #undef __FUNC__
2355615d1e5SSatish Balay #define __FUNC__ "SNESSetMatrixFreeParameters"
236b4fd4287SBarry Smith /*@
237eb9086c3SLois Curfman McInnes    SNESSetMatrixFreeParameters - Sets the parameters for the approximation of
238eb9086c3SLois Curfman McInnes    matrix-vector products using finite differences.
239eb9086c3SLois Curfman McInnes 
240639f9d9dSBarry Smith $       J(u)*a = [J(u+h*a) - J(u)]/h where
241639f9d9dSBarry Smith $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
242639f9d9dSBarry Smith $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
243639f9d9dSBarry Smith $
244eb9086c3SLois Curfman McInnes    Input Parameters:
245eb9086c3SLois Curfman McInnes .  snes - the SNES context
246ccb6204bSLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
247ccb6204bSLois Curfman McInnes      the relative error in the function evaluations)
248eb9086c3SLois Curfman McInnes .  umin - minimum allowable u-value
249eb9086c3SLois Curfman McInnes 
250eb9086c3SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
2513d5db8e1SLois Curfman McInnes 
2523d5db8e1SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate()
253eb9086c3SLois Curfman McInnes @*/
254eb9086c3SLois Curfman McInnes int SNESSetMatrixFreeParameters(SNES snes,double error,double umin)
255eb9086c3SLois Curfman McInnes {
256eb9086c3SLois Curfman McInnes   MFCtx_Private *ctx;
257eb9086c3SLois Curfman McInnes   int           ierr;
258eb9086c3SLois Curfman McInnes   Mat           mat;
259eb9086c3SLois Curfman McInnes 
260eb9086c3SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
261eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
262eb9086c3SLois Curfman McInnes   if (ctx) {
263eb9086c3SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
264eb9086c3SLois Curfman McInnes     if (umin != PETSC_DEFAULT)  ctx->umin = umin;
265eb9086c3SLois Curfman McInnes   }
266eb9086c3SLois Curfman McInnes   return 0;
267eb9086c3SLois Curfman McInnes }
268eb9086c3SLois Curfman McInnes 
2695615d1e5SSatish Balay #undef __FUNC__
2705615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace"
271eb9086c3SLois Curfman McInnes /*@
27221a45821SLois Curfman McInnes    SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that
273f525115eSLois Curfman McInnes    an operator is supposed to have.  Since roundoff will create a
274b4fd4287SBarry Smith    small component in the null space, if you know the null space
275b4fd4287SBarry Smith    you may have it automatically removed.
276b4fd4287SBarry Smith 
277b4fd4287SBarry Smith    Input Parameters:
27821a45821SLois Curfman McInnes .  J - the matrix-free matrix context
27921a45821SLois Curfman McInnes .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
280b4fd4287SBarry Smith .  n - number of vectors (excluding constant vector) in null space
28121a45821SLois Curfman McInnes .  vecs - the vectors that span the null space (excluding the constant vector);
282b4fd4287SBarry Smith .         these vectors must be orthonormal
283b4fd4287SBarry Smith 
284b4fd4287SBarry Smith .keywords: SNES, matrix-free, null space
285b4fd4287SBarry Smith @*/
286b4fd4287SBarry Smith int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
287b4fd4287SBarry Smith {
288b4fd4287SBarry Smith   int           ierr;
289b4fd4287SBarry Smith   MFCtx_Private *ctx;
290b4fd4287SBarry Smith   MPI_Comm      comm;
291b4fd4287SBarry Smith 
292b4fd4287SBarry Smith   PetscObjectGetComm((PetscObject)J,&comm);
293b4fd4287SBarry Smith 
294b4fd4287SBarry Smith   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
295b4fd4287SBarry Smith   /* no context indicates that it is not the "matrix free" matrix type */
296b4fd4287SBarry Smith   if (!ctx) return 0;
297b4fd4287SBarry Smith   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
298b4fd4287SBarry Smith   return 0;
299b4fd4287SBarry Smith }
300b4fd4287SBarry Smith 
3015334005bSBarry Smith 
3025334005bSBarry Smith 
3035334005bSBarry Smith 
304