14b9ad928SBarry Smith 24b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscksp.h" I*/ 34b9ad928SBarry Smith /*I "petscmg.h" I*/ 44b9ad928SBarry Smith 54b9ad928SBarry Smith #undef __FUNCT__ 64b9ad928SBarry Smith #define __FUNCT__ "MGDefaultResidual" 74b9ad928SBarry Smith /*@C 84b9ad928SBarry Smith MGDefaultResidual - Default routine to calculate the residual. 94b9ad928SBarry Smith 104b9ad928SBarry Smith Collective on Mat and Vec 114b9ad928SBarry Smith 124b9ad928SBarry Smith Input Parameters: 134b9ad928SBarry Smith + mat - the matrix 144b9ad928SBarry Smith . b - the right-hand-side 154b9ad928SBarry Smith - x - the approximate solution 164b9ad928SBarry Smith 174b9ad928SBarry Smith Output Parameter: 184b9ad928SBarry Smith . r - location to store the residual 194b9ad928SBarry Smith 204b9ad928SBarry Smith Level: advanced 214b9ad928SBarry Smith 224b9ad928SBarry Smith .keywords: MG, default, multigrid, residual 234b9ad928SBarry Smith 244b9ad928SBarry Smith .seealso: MGSetResidual() 254b9ad928SBarry Smith @*/ 26dfbe8321SBarry Smith PetscErrorCode MGDefaultResidual(Mat mat,Vec b,Vec x,Vec r) 274b9ad928SBarry Smith { 28dfbe8321SBarry Smith PetscErrorCode ierr; 294b9ad928SBarry Smith PetscScalar mone = -1.0; 304b9ad928SBarry Smith 314b9ad928SBarry Smith PetscFunctionBegin; 324b9ad928SBarry Smith ierr = MatMult(mat,x,r);CHKERRQ(ierr); 334b9ad928SBarry Smith ierr = VecAYPX(&mone,b,r);CHKERRQ(ierr); 344b9ad928SBarry Smith PetscFunctionReturn(0); 354b9ad928SBarry Smith } 364b9ad928SBarry Smith 374b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/ 384b9ad928SBarry Smith 394b9ad928SBarry Smith #undef __FUNCT__ 404b9ad928SBarry Smith #define __FUNCT__ "MGGetCoarseSolve" 414b9ad928SBarry Smith /*@C 424b9ad928SBarry Smith MGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 434b9ad928SBarry Smith 444b9ad928SBarry Smith Not Collective 454b9ad928SBarry Smith 464b9ad928SBarry Smith Input Parameter: 474b9ad928SBarry Smith . pc - the multigrid context 484b9ad928SBarry Smith 494b9ad928SBarry Smith Output Parameter: 504b9ad928SBarry Smith . ksp - the coarse grid solver context 514b9ad928SBarry Smith 524b9ad928SBarry Smith Level: advanced 534b9ad928SBarry Smith 544b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid 554b9ad928SBarry Smith @*/ 56dfbe8321SBarry Smith PetscErrorCode MGGetCoarseSolve(PC pc,KSP *ksp) 574b9ad928SBarry Smith { 584b9ad928SBarry Smith MG *mg = (MG*)pc->data; 594b9ad928SBarry Smith 604b9ad928SBarry Smith PetscFunctionBegin; 614b9ad928SBarry Smith *ksp = mg[0]->smoothd; 624b9ad928SBarry Smith PetscFunctionReturn(0); 634b9ad928SBarry Smith } 644b9ad928SBarry Smith 654b9ad928SBarry Smith #undef __FUNCT__ 664b9ad928SBarry Smith #define __FUNCT__ "MGSetResidual" 674b9ad928SBarry Smith /*@C 684b9ad928SBarry Smith MGSetResidual - Sets the function to be used to calculate the residual 694b9ad928SBarry Smith on the lth level. 704b9ad928SBarry Smith 714b9ad928SBarry Smith Collective on PC and Mat 724b9ad928SBarry Smith 734b9ad928SBarry Smith Input Parameters: 744b9ad928SBarry Smith + pc - the multigrid context 754b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 764b9ad928SBarry Smith . residual - function used to form residual (usually MGDefaultResidual) 774b9ad928SBarry Smith - mat - matrix associated with residual 784b9ad928SBarry Smith 794b9ad928SBarry Smith Level: advanced 804b9ad928SBarry Smith 814b9ad928SBarry Smith .keywords: MG, set, multigrid, residual, level 824b9ad928SBarry Smith 834b9ad928SBarry Smith .seealso: MGDefaultResidual() 844b9ad928SBarry Smith @*/ 8579416396SBarry Smith PetscErrorCode MGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 864b9ad928SBarry Smith { 874b9ad928SBarry Smith MG *mg = (MG*)pc->data; 884b9ad928SBarry Smith 894b9ad928SBarry Smith PetscFunctionBegin; 904b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 914b9ad928SBarry Smith 924b9ad928SBarry Smith mg[l]->residual = residual; 934b9ad928SBarry Smith mg[l]->A = mat; 944b9ad928SBarry Smith PetscFunctionReturn(0); 954b9ad928SBarry Smith } 964b9ad928SBarry Smith 974b9ad928SBarry Smith #undef __FUNCT__ 984b9ad928SBarry Smith #define __FUNCT__ "MGSetInterpolate" 994b9ad928SBarry Smith /*@ 1004b9ad928SBarry Smith MGSetInterpolate - Sets the function to be used to calculate the 1014b9ad928SBarry Smith interpolation on the lth level. 1024b9ad928SBarry Smith 1034b9ad928SBarry Smith Collective on PC and Mat 1044b9ad928SBarry Smith 1054b9ad928SBarry Smith Input Parameters: 1064b9ad928SBarry Smith + pc - the multigrid context 1074b9ad928SBarry Smith . mat - the interpolation operator 1084b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 1094b9ad928SBarry Smith 1104b9ad928SBarry Smith Level: advanced 1114b9ad928SBarry Smith 1124b9ad928SBarry Smith Notes: 1134b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 1144b9ad928SBarry Smith for the same level. 1154b9ad928SBarry Smith 1164b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1174b9ad928SBarry Smith out from the matrix size which one it is. 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith .keywords: multigrid, set, interpolate, level 1204b9ad928SBarry Smith 1214b9ad928SBarry Smith .seealso: MGSetRestriction() 1224b9ad928SBarry Smith @*/ 12379416396SBarry Smith PetscErrorCode MGSetInterpolate(PC pc,PetscInt l,Mat mat) 1244b9ad928SBarry Smith { 1254b9ad928SBarry Smith MG *mg = (MG*)pc->data; 1264b9ad928SBarry Smith 1274b9ad928SBarry Smith PetscFunctionBegin; 1284b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1294b9ad928SBarry Smith mg[l]->interpolate = mat; 1304b9ad928SBarry Smith PetscFunctionReturn(0); 1314b9ad928SBarry Smith } 1324b9ad928SBarry Smith 1334b9ad928SBarry Smith #undef __FUNCT__ 1344b9ad928SBarry Smith #define __FUNCT__ "MGSetRestriction" 1354b9ad928SBarry Smith /*@ 1364b9ad928SBarry Smith MGSetRestriction - Sets the function to be used to restrict vector 1374b9ad928SBarry Smith from level l to l-1. 1384b9ad928SBarry Smith 1394b9ad928SBarry Smith Collective on PC and Mat 1404b9ad928SBarry Smith 1414b9ad928SBarry Smith Input Parameters: 1424b9ad928SBarry Smith + pc - the multigrid context 1434b9ad928SBarry Smith . mat - the restriction matrix 1444b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 1454b9ad928SBarry Smith 1464b9ad928SBarry Smith Level: advanced 1474b9ad928SBarry Smith 1484b9ad928SBarry Smith Notes: 1494b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 1504b9ad928SBarry Smith for the same level. 1514b9ad928SBarry Smith 1524b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1534b9ad928SBarry Smith out from the matrix size which one it is. 1544b9ad928SBarry Smith 1554b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level 1564b9ad928SBarry Smith 1574b9ad928SBarry Smith .seealso: MGSetInterpolate() 1584b9ad928SBarry Smith @*/ 15979416396SBarry Smith PetscErrorCode MGSetRestriction(PC pc,PetscInt l,Mat mat) 1604b9ad928SBarry Smith { 1614b9ad928SBarry Smith MG *mg = (MG*)pc->data; 1624b9ad928SBarry Smith 1634b9ad928SBarry Smith PetscFunctionBegin; 1644b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1654b9ad928SBarry Smith mg[l]->restrct = mat; 1664b9ad928SBarry Smith PetscFunctionReturn(0); 1674b9ad928SBarry Smith } 1684b9ad928SBarry Smith 1694b9ad928SBarry Smith #undef __FUNCT__ 1704b9ad928SBarry Smith #define __FUNCT__ "MGGetSmoother" 1714b9ad928SBarry Smith /*@C 1724b9ad928SBarry Smith MGGetSmoother - Gets the KSP context to be used as smoother for 1734b9ad928SBarry Smith both pre- and post-smoothing. Call both MGGetSmootherUp() and 1744b9ad928SBarry Smith MGGetSmootherDown() to use different functions for pre- and 1754b9ad928SBarry Smith post-smoothing. 1764b9ad928SBarry Smith 1774b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 1784b9ad928SBarry Smith 1794b9ad928SBarry Smith Input Parameters: 1804b9ad928SBarry Smith + pc - the multigrid context 1814b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 1824b9ad928SBarry Smith 1834b9ad928SBarry Smith Ouput Parameters: 1844b9ad928SBarry Smith . ksp - the smoother 1854b9ad928SBarry Smith 1864b9ad928SBarry Smith Level: advanced 1874b9ad928SBarry Smith 1884b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 1894b9ad928SBarry Smith 1904b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmootherDown() 1914b9ad928SBarry Smith @*/ 19279416396SBarry Smith PetscErrorCode MGGetSmoother(PC pc,PetscInt l,KSP *ksp) 1934b9ad928SBarry Smith { 1944b9ad928SBarry Smith MG *mg = (MG*)pc->data; 1954b9ad928SBarry Smith 1964b9ad928SBarry Smith PetscFunctionBegin; 1974b9ad928SBarry Smith *ksp = mg[l]->smoothd; 1984b9ad928SBarry Smith PetscFunctionReturn(0); 1994b9ad928SBarry Smith } 2004b9ad928SBarry Smith 2014b9ad928SBarry Smith #undef __FUNCT__ 2024b9ad928SBarry Smith #define __FUNCT__ "MGGetSmootherUp" 2034b9ad928SBarry Smith /*@C 2044b9ad928SBarry Smith MGGetSmootherUp - Gets the KSP context to be used as smoother after 2054b9ad928SBarry Smith coarse grid correction (post-smoother). 2064b9ad928SBarry Smith 2074b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith Input Parameters: 2104b9ad928SBarry Smith + pc - the multigrid context 2114b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 2124b9ad928SBarry Smith 2134b9ad928SBarry Smith Ouput Parameters: 2144b9ad928SBarry Smith . ksp - the smoother 2154b9ad928SBarry Smith 2164b9ad928SBarry Smith Level: advanced 2174b9ad928SBarry Smith 2184b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level 2194b9ad928SBarry Smith 2204b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmootherDown() 2214b9ad928SBarry Smith @*/ 22279416396SBarry Smith PetscErrorCode MGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 2234b9ad928SBarry Smith { 2244b9ad928SBarry Smith MG *mg = (MG*)pc->data; 225dfbe8321SBarry Smith PetscErrorCode ierr; 2264b9ad928SBarry Smith char *prefix; 2274b9ad928SBarry Smith MPI_Comm comm; 2284b9ad928SBarry Smith 2294b9ad928SBarry Smith PetscFunctionBegin; 2304b9ad928SBarry Smith /* 2314b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 2324b9ad928SBarry Smith Thus we check if a different one has already been allocated, 2334b9ad928SBarry Smith if not we allocate it. 2344b9ad928SBarry Smith */ 235b03c7568SBarry Smith if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 2364b9ad928SBarry Smith if (mg[l]->smoothu == mg[l]->smoothd) { 2374b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[l]->smoothd,&comm);CHKERRQ(ierr); 2382e70eadcSBarry Smith ierr = KSPGetOptionsPrefix(mg[l]->smoothd,&prefix);CHKERRQ(ierr); 2394b9ad928SBarry Smith ierr = KSPCreate(comm,&mg[l]->smoothu);CHKERRQ(ierr); 240d7d82daaSBarry Smith ierr = KSPSetTolerances(mg[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 2414b9ad928SBarry Smith ierr = KSPSetOptionsPrefix(mg[l]->smoothu,prefix);CHKERRQ(ierr); 242*52e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,mg[l]->smoothu);CHKERRQ(ierr); 2434b9ad928SBarry Smith } 2444b9ad928SBarry Smith if (ksp) *ksp = mg[l]->smoothu; 2454b9ad928SBarry Smith PetscFunctionReturn(0); 2464b9ad928SBarry Smith } 2474b9ad928SBarry Smith 2484b9ad928SBarry Smith #undef __FUNCT__ 2494b9ad928SBarry Smith #define __FUNCT__ "MGGetSmootherDown" 2504b9ad928SBarry Smith /*@C 2514b9ad928SBarry Smith MGGetSmootherDown - Gets the KSP context to be used as smoother before 2524b9ad928SBarry Smith coarse grid correction (pre-smoother). 2534b9ad928SBarry Smith 2544b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 2554b9ad928SBarry Smith 2564b9ad928SBarry Smith Input Parameters: 2574b9ad928SBarry Smith + pc - the multigrid context 2584b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 2594b9ad928SBarry Smith 2604b9ad928SBarry Smith Ouput Parameters: 2614b9ad928SBarry Smith . ksp - the smoother 2624b9ad928SBarry Smith 2634b9ad928SBarry Smith Level: advanced 2644b9ad928SBarry Smith 2654b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 2664b9ad928SBarry Smith 2674b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmoother() 2684b9ad928SBarry Smith @*/ 26979416396SBarry Smith PetscErrorCode MGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 2704b9ad928SBarry Smith { 271dfbe8321SBarry Smith PetscErrorCode ierr; 2724b9ad928SBarry Smith MG *mg = (MG*)pc->data; 2734b9ad928SBarry Smith 2744b9ad928SBarry Smith PetscFunctionBegin; 2754b9ad928SBarry Smith /* make sure smoother up and down are different */ 2764b9ad928SBarry Smith ierr = MGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr); 2774b9ad928SBarry Smith *ksp = mg[l]->smoothd; 2784b9ad928SBarry Smith PetscFunctionReturn(0); 2794b9ad928SBarry Smith } 2804b9ad928SBarry Smith 2814b9ad928SBarry Smith #undef __FUNCT__ 2824b9ad928SBarry Smith #define __FUNCT__ "MGSetCyclesOnLevel" 2834b9ad928SBarry Smith /*@ 2844b9ad928SBarry Smith MGSetCyclesOnLevel - Sets the number of cycles to run on this level. 2854b9ad928SBarry Smith 2864b9ad928SBarry Smith Collective on PC 2874b9ad928SBarry Smith 2884b9ad928SBarry Smith Input Parameters: 2894b9ad928SBarry Smith + pc - the multigrid context 2904b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 2914b9ad928SBarry Smith - n - the number of cycles 2924b9ad928SBarry Smith 2934b9ad928SBarry Smith Level: advanced 2944b9ad928SBarry Smith 2954b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 2964b9ad928SBarry Smith 2974b9ad928SBarry Smith .seealso: MGSetCycles() 2984b9ad928SBarry Smith @*/ 29979416396SBarry Smith PetscErrorCode MGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c) 3004b9ad928SBarry Smith { 3014b9ad928SBarry Smith MG *mg = (MG*)pc->data; 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith PetscFunctionBegin; 3044b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 3054b9ad928SBarry Smith mg[l]->cycles = c; 3064b9ad928SBarry Smith PetscFunctionReturn(0); 3074b9ad928SBarry Smith } 3084b9ad928SBarry Smith 3094b9ad928SBarry Smith #undef __FUNCT__ 3104b9ad928SBarry Smith #define __FUNCT__ "MGSetRhs" 3114b9ad928SBarry Smith /*@ 3124b9ad928SBarry Smith MGSetRhs - Sets the vector space to be used to store the right-hand side 3134b9ad928SBarry Smith on a particular level. The user should free this space at the conclusion 3144b9ad928SBarry Smith of multigrid use. 3154b9ad928SBarry Smith 3164b9ad928SBarry Smith Collective on PC and Vec 3174b9ad928SBarry Smith 3184b9ad928SBarry Smith Input Parameters: 3194b9ad928SBarry Smith + pc - the multigrid context 3204b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3214b9ad928SBarry Smith - c - the space 3224b9ad928SBarry Smith 3234b9ad928SBarry Smith Level: advanced 3244b9ad928SBarry Smith 3254b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level 3264b9ad928SBarry Smith 3274b9ad928SBarry Smith .seealso: MGSetX(), MGSetR() 3284b9ad928SBarry Smith @*/ 32979416396SBarry Smith PetscErrorCode MGSetRhs(PC pc,PetscInt l,Vec c) 3304b9ad928SBarry Smith { 3314b9ad928SBarry Smith MG *mg = (MG*)pc->data; 3324b9ad928SBarry Smith 3334b9ad928SBarry Smith PetscFunctionBegin; 3344b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 3354b9ad928SBarry Smith mg[l]->b = c; 3364b9ad928SBarry Smith PetscFunctionReturn(0); 3374b9ad928SBarry Smith } 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith #undef __FUNCT__ 3404b9ad928SBarry Smith #define __FUNCT__ "MGSetX" 3414b9ad928SBarry Smith /*@ 3424b9ad928SBarry Smith MGSetX - Sets the vector space to be used to store the solution on a 3434b9ad928SBarry Smith particular level. The user should free this space at the conclusion 3444b9ad928SBarry Smith of multigrid use. 3454b9ad928SBarry Smith 3464b9ad928SBarry Smith Collective on PC and Vec 3474b9ad928SBarry Smith 3484b9ad928SBarry Smith Input Parameters: 3494b9ad928SBarry Smith + pc - the multigrid context 3504b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3514b9ad928SBarry Smith - c - the space 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith Level: advanced 3544b9ad928SBarry Smith 3554b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level 3564b9ad928SBarry Smith 3574b9ad928SBarry Smith .seealso: MGSetRhs(), MGSetR() 3584b9ad928SBarry Smith @*/ 35979416396SBarry Smith PetscErrorCode MGSetX(PC pc,PetscInt l,Vec c) 3604b9ad928SBarry Smith { 3614b9ad928SBarry Smith MG *mg = (MG*)pc->data; 3624b9ad928SBarry Smith 3634b9ad928SBarry Smith PetscFunctionBegin; 3644b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 3654b9ad928SBarry Smith mg[l]->x = c; 3664b9ad928SBarry Smith PetscFunctionReturn(0); 3674b9ad928SBarry Smith } 3684b9ad928SBarry Smith 3694b9ad928SBarry Smith #undef __FUNCT__ 3704b9ad928SBarry Smith #define __FUNCT__ "MGSetR" 3714b9ad928SBarry Smith /*@ 3724b9ad928SBarry Smith MGSetR - Sets the vector space to be used to store the residual on a 3734b9ad928SBarry Smith particular level. The user should free this space at the conclusion of 3744b9ad928SBarry Smith multigrid use. 3754b9ad928SBarry Smith 3764b9ad928SBarry Smith Collective on PC and Vec 3774b9ad928SBarry Smith 3784b9ad928SBarry Smith Input Parameters: 3794b9ad928SBarry Smith + pc - the multigrid context 3804b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3814b9ad928SBarry Smith - c - the space 3824b9ad928SBarry Smith 3834b9ad928SBarry Smith Level: advanced 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level 3864b9ad928SBarry Smith @*/ 38779416396SBarry Smith PetscErrorCode MGSetR(PC pc,PetscInt l,Vec c) 3884b9ad928SBarry Smith { 3894b9ad928SBarry Smith MG *mg = (MG*)pc->data; 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith PetscFunctionBegin; 3924b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 3934b9ad928SBarry Smith mg[l]->r = c; 3944b9ad928SBarry Smith PetscFunctionReturn(0); 3954b9ad928SBarry Smith } 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith 3984b9ad928SBarry Smith 3994b9ad928SBarry Smith 4004b9ad928SBarry Smith 4014b9ad928SBarry Smith 4024b9ad928SBarry Smith 4034b9ad928SBarry Smith 404