1dba47a55SKris Buschelman #define PETSCKSP_DLL 24b9ad928SBarry Smith 34b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscksp.h" I*/ 44b9ad928SBarry Smith /*I "petscmg.h" I*/ 54b9ad928SBarry Smith 64b9ad928SBarry Smith #undef __FUNCT__ 74b9ad928SBarry Smith #define __FUNCT__ "MGDefaultResidual" 84b9ad928SBarry Smith /*@C 94b9ad928SBarry Smith MGDefaultResidual - Default routine to calculate the residual. 104b9ad928SBarry Smith 114b9ad928SBarry Smith Collective on Mat and Vec 124b9ad928SBarry Smith 134b9ad928SBarry Smith Input Parameters: 144b9ad928SBarry Smith + mat - the matrix 154b9ad928SBarry Smith . b - the right-hand-side 164b9ad928SBarry Smith - x - the approximate solution 174b9ad928SBarry Smith 184b9ad928SBarry Smith Output Parameter: 194b9ad928SBarry Smith . r - location to store the residual 204b9ad928SBarry Smith 214b9ad928SBarry Smith Level: advanced 224b9ad928SBarry Smith 234b9ad928SBarry Smith .keywords: MG, default, multigrid, residual 244b9ad928SBarry Smith 254b9ad928SBarry Smith .seealso: MGSetResidual() 264b9ad928SBarry Smith @*/ 27dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGDefaultResidual(Mat mat,Vec b,Vec x,Vec r) 284b9ad928SBarry Smith { 29dfbe8321SBarry Smith PetscErrorCode ierr; 304b9ad928SBarry Smith PetscScalar mone = -1.0; 314b9ad928SBarry Smith 324b9ad928SBarry Smith PetscFunctionBegin; 334b9ad928SBarry Smith ierr = MatMult(mat,x,r);CHKERRQ(ierr); 344b9ad928SBarry Smith ierr = VecAYPX(&mone,b,r);CHKERRQ(ierr); 354b9ad928SBarry Smith PetscFunctionReturn(0); 364b9ad928SBarry Smith } 374b9ad928SBarry Smith 384b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/ 394b9ad928SBarry Smith 404b9ad928SBarry Smith #undef __FUNCT__ 414b9ad928SBarry Smith #define __FUNCT__ "MGGetCoarseSolve" 424b9ad928SBarry Smith /*@C 434b9ad928SBarry Smith MGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 444b9ad928SBarry Smith 454b9ad928SBarry Smith Not Collective 464b9ad928SBarry Smith 474b9ad928SBarry Smith Input Parameter: 484b9ad928SBarry Smith . pc - the multigrid context 494b9ad928SBarry Smith 504b9ad928SBarry Smith Output Parameter: 514b9ad928SBarry Smith . ksp - the coarse grid solver context 524b9ad928SBarry Smith 534b9ad928SBarry Smith Level: advanced 544b9ad928SBarry Smith 554b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid 564b9ad928SBarry Smith @*/ 57dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGGetCoarseSolve(PC pc,KSP *ksp) 584b9ad928SBarry Smith { 594b9ad928SBarry Smith MG *mg = (MG*)pc->data; 604b9ad928SBarry Smith 614b9ad928SBarry Smith PetscFunctionBegin; 624b9ad928SBarry Smith *ksp = mg[0]->smoothd; 634b9ad928SBarry Smith PetscFunctionReturn(0); 644b9ad928SBarry Smith } 654b9ad928SBarry Smith 664b9ad928SBarry Smith #undef __FUNCT__ 674b9ad928SBarry Smith #define __FUNCT__ "MGSetResidual" 684b9ad928SBarry Smith /*@C 694b9ad928SBarry Smith MGSetResidual - Sets the function to be used to calculate the residual 704b9ad928SBarry Smith on the lth level. 714b9ad928SBarry Smith 724b9ad928SBarry Smith Collective on PC and Mat 734b9ad928SBarry Smith 744b9ad928SBarry Smith Input Parameters: 754b9ad928SBarry Smith + pc - the multigrid context 764b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 774b9ad928SBarry Smith . residual - function used to form residual (usually MGDefaultResidual) 784b9ad928SBarry Smith - mat - matrix associated with residual 794b9ad928SBarry Smith 804b9ad928SBarry Smith Level: advanced 814b9ad928SBarry Smith 824b9ad928SBarry Smith .keywords: MG, set, multigrid, residual, level 834b9ad928SBarry Smith 844b9ad928SBarry Smith .seealso: MGDefaultResidual() 854b9ad928SBarry Smith @*/ 86dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 874b9ad928SBarry Smith { 884b9ad928SBarry Smith MG *mg = (MG*)pc->data; 894b9ad928SBarry Smith 904b9ad928SBarry Smith PetscFunctionBegin; 914b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 924b9ad928SBarry Smith 934b9ad928SBarry Smith mg[l]->residual = residual; 944b9ad928SBarry Smith mg[l]->A = mat; 954b9ad928SBarry Smith PetscFunctionReturn(0); 964b9ad928SBarry Smith } 974b9ad928SBarry Smith 984b9ad928SBarry Smith #undef __FUNCT__ 994b9ad928SBarry Smith #define __FUNCT__ "MGSetInterpolate" 1004b9ad928SBarry Smith /*@ 1014b9ad928SBarry Smith MGSetInterpolate - Sets the function to be used to calculate the 1024b9ad928SBarry Smith interpolation on the lth level. 1034b9ad928SBarry Smith 1044b9ad928SBarry Smith Collective on PC and Mat 1054b9ad928SBarry Smith 1064b9ad928SBarry Smith Input Parameters: 1074b9ad928SBarry Smith + pc - the multigrid context 1084b9ad928SBarry Smith . mat - the interpolation operator 1094b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 1104b9ad928SBarry Smith 1114b9ad928SBarry Smith Level: advanced 1124b9ad928SBarry Smith 1134b9ad928SBarry Smith Notes: 1144b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 1154b9ad928SBarry Smith for the same level. 1164b9ad928SBarry Smith 1174b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1184b9ad928SBarry Smith out from the matrix size which one it is. 1194b9ad928SBarry Smith 120fccaa45eSBarry Smith If you do not set this, the transpose of the Mat set with MGSetRestriction() 121fccaa45eSBarry Smith is used. 122fccaa45eSBarry Smith 1234b9ad928SBarry Smith .keywords: multigrid, set, interpolate, level 1244b9ad928SBarry Smith 1254b9ad928SBarry Smith .seealso: MGSetRestriction() 1264b9ad928SBarry Smith @*/ 127dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetInterpolate(PC pc,PetscInt l,Mat mat) 1284b9ad928SBarry Smith { 1294b9ad928SBarry Smith MG *mg = (MG*)pc->data; 130fccaa45eSBarry Smith PetscErrorCode ierr; 1314b9ad928SBarry Smith 1324b9ad928SBarry Smith PetscFunctionBegin; 1334b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 134*630ba2eeSBarry Smith if (!l) SETTERQ(PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 1354b9ad928SBarry Smith mg[l]->interpolate = mat; 136fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1374b9ad928SBarry Smith PetscFunctionReturn(0); 1384b9ad928SBarry Smith } 1394b9ad928SBarry Smith 1404b9ad928SBarry Smith #undef __FUNCT__ 1414b9ad928SBarry Smith #define __FUNCT__ "MGSetRestriction" 1424b9ad928SBarry Smith /*@ 1434b9ad928SBarry Smith MGSetRestriction - Sets the function to be used to restrict vector 1444b9ad928SBarry Smith from level l to l-1. 1454b9ad928SBarry Smith 1464b9ad928SBarry Smith Collective on PC and Mat 1474b9ad928SBarry Smith 1484b9ad928SBarry Smith Input Parameters: 1494b9ad928SBarry Smith + pc - the multigrid context 1504b9ad928SBarry Smith . mat - the restriction matrix 1514b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 1524b9ad928SBarry Smith 1534b9ad928SBarry Smith Level: advanced 1544b9ad928SBarry Smith 1554b9ad928SBarry Smith Notes: 1564b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 1574b9ad928SBarry Smith for the same level. 1584b9ad928SBarry Smith 1594b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1604b9ad928SBarry Smith out from the matrix size which one it is. 1614b9ad928SBarry Smith 162fccaa45eSBarry Smith If you do not set this, the transpose of the Mat set with MGSetInterpolate() 163fccaa45eSBarry Smith is used. 164fccaa45eSBarry Smith 1654b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level 1664b9ad928SBarry Smith 1674b9ad928SBarry Smith .seealso: MGSetInterpolate() 1684b9ad928SBarry Smith @*/ 169dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetRestriction(PC pc,PetscInt l,Mat mat) 1704b9ad928SBarry Smith { 171fccaa45eSBarry Smith PetscErrorCode ierr; 1724b9ad928SBarry Smith MG *mg = (MG*)pc->data; 1734b9ad928SBarry Smith 1744b9ad928SBarry Smith PetscFunctionBegin; 1754b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 176*630ba2eeSBarry Smith if (!l) SETTERQ(PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 1774b9ad928SBarry Smith mg[l]->restrct = mat; 178fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1794b9ad928SBarry Smith PetscFunctionReturn(0); 1804b9ad928SBarry Smith } 1814b9ad928SBarry Smith 1824b9ad928SBarry Smith #undef __FUNCT__ 1834b9ad928SBarry Smith #define __FUNCT__ "MGGetSmoother" 1844b9ad928SBarry Smith /*@C 1854b9ad928SBarry Smith MGGetSmoother - Gets the KSP context to be used as smoother for 1864b9ad928SBarry Smith both pre- and post-smoothing. Call both MGGetSmootherUp() and 1874b9ad928SBarry Smith MGGetSmootherDown() to use different functions for pre- and 1884b9ad928SBarry Smith post-smoothing. 1894b9ad928SBarry Smith 1904b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 1914b9ad928SBarry Smith 1924b9ad928SBarry Smith Input Parameters: 1934b9ad928SBarry Smith + pc - the multigrid context 1944b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 1954b9ad928SBarry Smith 1964b9ad928SBarry Smith Ouput Parameters: 1974b9ad928SBarry Smith . ksp - the smoother 1984b9ad928SBarry Smith 1994b9ad928SBarry Smith Level: advanced 2004b9ad928SBarry Smith 2014b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 2024b9ad928SBarry Smith 2034b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmootherDown() 2044b9ad928SBarry Smith @*/ 205dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGGetSmoother(PC pc,PetscInt l,KSP *ksp) 2064b9ad928SBarry Smith { 2074b9ad928SBarry Smith MG *mg = (MG*)pc->data; 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith PetscFunctionBegin; 2104b9ad928SBarry Smith *ksp = mg[l]->smoothd; 2114b9ad928SBarry Smith PetscFunctionReturn(0); 2124b9ad928SBarry Smith } 2134b9ad928SBarry Smith 2144b9ad928SBarry Smith #undef __FUNCT__ 2154b9ad928SBarry Smith #define __FUNCT__ "MGGetSmootherUp" 2164b9ad928SBarry Smith /*@C 2174b9ad928SBarry Smith MGGetSmootherUp - Gets the KSP context to be used as smoother after 2184b9ad928SBarry Smith coarse grid correction (post-smoother). 2194b9ad928SBarry Smith 2204b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 2214b9ad928SBarry Smith 2224b9ad928SBarry Smith Input Parameters: 2234b9ad928SBarry Smith + pc - the multigrid context 2244b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 2254b9ad928SBarry Smith 2264b9ad928SBarry Smith Ouput Parameters: 2274b9ad928SBarry Smith . ksp - the smoother 2284b9ad928SBarry Smith 2294b9ad928SBarry Smith Level: advanced 2304b9ad928SBarry Smith 2314b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level 2324b9ad928SBarry Smith 2334b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmootherDown() 2344b9ad928SBarry Smith @*/ 235dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 2364b9ad928SBarry Smith { 2374b9ad928SBarry Smith MG *mg = (MG*)pc->data; 238dfbe8321SBarry Smith PetscErrorCode ierr; 2394b9ad928SBarry Smith char *prefix; 2404b9ad928SBarry Smith MPI_Comm comm; 2414b9ad928SBarry Smith 2424b9ad928SBarry Smith PetscFunctionBegin; 2434b9ad928SBarry Smith /* 2444b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 2454b9ad928SBarry Smith Thus we check if a different one has already been allocated, 2464b9ad928SBarry Smith if not we allocate it. 2474b9ad928SBarry Smith */ 248b03c7568SBarry Smith if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 2494b9ad928SBarry Smith if (mg[l]->smoothu == mg[l]->smoothd) { 2504b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[l]->smoothd,&comm);CHKERRQ(ierr); 2512e70eadcSBarry Smith ierr = KSPGetOptionsPrefix(mg[l]->smoothd,&prefix);CHKERRQ(ierr); 2524b9ad928SBarry Smith ierr = KSPCreate(comm,&mg[l]->smoothu);CHKERRQ(ierr); 253d7d82daaSBarry Smith ierr = KSPSetTolerances(mg[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 2544b9ad928SBarry Smith ierr = KSPSetOptionsPrefix(mg[l]->smoothu,prefix);CHKERRQ(ierr); 25552e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,mg[l]->smoothu);CHKERRQ(ierr); 2564b9ad928SBarry Smith } 2574b9ad928SBarry Smith if (ksp) *ksp = mg[l]->smoothu; 2584b9ad928SBarry Smith PetscFunctionReturn(0); 2594b9ad928SBarry Smith } 2604b9ad928SBarry Smith 2614b9ad928SBarry Smith #undef __FUNCT__ 2624b9ad928SBarry Smith #define __FUNCT__ "MGGetSmootherDown" 2634b9ad928SBarry Smith /*@C 2644b9ad928SBarry Smith MGGetSmootherDown - Gets the KSP context to be used as smoother before 2654b9ad928SBarry Smith coarse grid correction (pre-smoother). 2664b9ad928SBarry Smith 2674b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 2684b9ad928SBarry Smith 2694b9ad928SBarry Smith Input Parameters: 2704b9ad928SBarry Smith + pc - the multigrid context 2714b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 2724b9ad928SBarry Smith 2734b9ad928SBarry Smith Ouput Parameters: 2744b9ad928SBarry Smith . ksp - the smoother 2754b9ad928SBarry Smith 2764b9ad928SBarry Smith Level: advanced 2774b9ad928SBarry Smith 2784b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 2794b9ad928SBarry Smith 2804b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmoother() 2814b9ad928SBarry Smith @*/ 282dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 2834b9ad928SBarry Smith { 284dfbe8321SBarry Smith PetscErrorCode ierr; 2854b9ad928SBarry Smith MG *mg = (MG*)pc->data; 2864b9ad928SBarry Smith 2874b9ad928SBarry Smith PetscFunctionBegin; 2884b9ad928SBarry Smith /* make sure smoother up and down are different */ 2894b9ad928SBarry Smith ierr = MGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr); 2904b9ad928SBarry Smith *ksp = mg[l]->smoothd; 2914b9ad928SBarry Smith PetscFunctionReturn(0); 2924b9ad928SBarry Smith } 2934b9ad928SBarry Smith 2944b9ad928SBarry Smith #undef __FUNCT__ 2954b9ad928SBarry Smith #define __FUNCT__ "MGSetCyclesOnLevel" 2964b9ad928SBarry Smith /*@ 2974b9ad928SBarry Smith MGSetCyclesOnLevel - Sets the number of cycles to run on this level. 2984b9ad928SBarry Smith 2994b9ad928SBarry Smith Collective on PC 3004b9ad928SBarry Smith 3014b9ad928SBarry Smith Input Parameters: 3024b9ad928SBarry Smith + pc - the multigrid context 3034b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3044b9ad928SBarry Smith - n - the number of cycles 3054b9ad928SBarry Smith 3064b9ad928SBarry Smith Level: advanced 3074b9ad928SBarry Smith 3084b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 3094b9ad928SBarry Smith 3104b9ad928SBarry Smith .seealso: MGSetCycles() 3114b9ad928SBarry Smith @*/ 312dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c) 3134b9ad928SBarry Smith { 3144b9ad928SBarry Smith MG *mg = (MG*)pc->data; 3154b9ad928SBarry Smith 3164b9ad928SBarry Smith PetscFunctionBegin; 3174b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 3184b9ad928SBarry Smith mg[l]->cycles = c; 3194b9ad928SBarry Smith PetscFunctionReturn(0); 3204b9ad928SBarry Smith } 3214b9ad928SBarry Smith 3224b9ad928SBarry Smith #undef __FUNCT__ 3234b9ad928SBarry Smith #define __FUNCT__ "MGSetRhs" 3244b9ad928SBarry Smith /*@ 3254b9ad928SBarry Smith MGSetRhs - Sets the vector space to be used to store the right-hand side 326fccaa45eSBarry Smith on a particular level. 3274b9ad928SBarry Smith 3284b9ad928SBarry Smith Collective on PC and Vec 3294b9ad928SBarry Smith 3304b9ad928SBarry Smith Input Parameters: 3314b9ad928SBarry Smith + pc - the multigrid context 3324b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3334b9ad928SBarry Smith - c - the space 3344b9ad928SBarry Smith 3354b9ad928SBarry Smith Level: advanced 3364b9ad928SBarry Smith 337fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 338fccaa45eSBarry Smith 339fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 340fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 341fccaa45eSBarry Smith 3424b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level 3434b9ad928SBarry Smith 3444b9ad928SBarry Smith .seealso: MGSetX(), MGSetR() 3454b9ad928SBarry Smith @*/ 346dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetRhs(PC pc,PetscInt l,Vec c) 3474b9ad928SBarry Smith { 348fccaa45eSBarry Smith PetscErrorCode ierr; 3494b9ad928SBarry Smith MG *mg = (MG*)pc->data; 3504b9ad928SBarry Smith 3514b9ad928SBarry Smith PetscFunctionBegin; 3524b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 353fccaa45eSBarry Smith if (l == mg[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 3544b9ad928SBarry Smith mg[l]->b = c; 355fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 3564b9ad928SBarry Smith PetscFunctionReturn(0); 3574b9ad928SBarry Smith } 3584b9ad928SBarry Smith 3594b9ad928SBarry Smith #undef __FUNCT__ 3604b9ad928SBarry Smith #define __FUNCT__ "MGSetX" 3614b9ad928SBarry Smith /*@ 3624b9ad928SBarry Smith MGSetX - Sets the vector space to be used to store the solution on a 363fccaa45eSBarry Smith particular level. 3644b9ad928SBarry Smith 3654b9ad928SBarry Smith Collective on PC and Vec 3664b9ad928SBarry Smith 3674b9ad928SBarry Smith Input Parameters: 3684b9ad928SBarry Smith + pc - the multigrid context 3694b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3704b9ad928SBarry Smith - c - the space 3714b9ad928SBarry Smith 3724b9ad928SBarry Smith Level: advanced 3734b9ad928SBarry Smith 374fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 375fccaa45eSBarry Smith 376fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 377fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 378fccaa45eSBarry Smith 3794b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level 3804b9ad928SBarry Smith 3814b9ad928SBarry Smith .seealso: MGSetRhs(), MGSetR() 3824b9ad928SBarry Smith @*/ 383dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetX(PC pc,PetscInt l,Vec c) 3844b9ad928SBarry Smith { 385fccaa45eSBarry Smith PetscErrorCode ierr; 3864b9ad928SBarry Smith MG *mg = (MG*)pc->data; 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith PetscFunctionBegin; 3894b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 390fccaa45eSBarry Smith if (l == mg[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 3914b9ad928SBarry Smith mg[l]->x = c; 392fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 3934b9ad928SBarry Smith PetscFunctionReturn(0); 3944b9ad928SBarry Smith } 3954b9ad928SBarry Smith 3964b9ad928SBarry Smith #undef __FUNCT__ 3974b9ad928SBarry Smith #define __FUNCT__ "MGSetR" 3984b9ad928SBarry Smith /*@ 3994b9ad928SBarry Smith MGSetR - Sets the vector space to be used to store the residual on a 400fccaa45eSBarry Smith particular level. 4014b9ad928SBarry Smith 4024b9ad928SBarry Smith Collective on PC and Vec 4034b9ad928SBarry Smith 4044b9ad928SBarry Smith Input Parameters: 4054b9ad928SBarry Smith + pc - the multigrid context 4064b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 4074b9ad928SBarry Smith - c - the space 4084b9ad928SBarry Smith 4094b9ad928SBarry Smith Level: advanced 4104b9ad928SBarry Smith 411fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 412fccaa45eSBarry Smith 413fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 414fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 415fccaa45eSBarry Smith 4164b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level 4174b9ad928SBarry Smith @*/ 418dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetR(PC pc,PetscInt l,Vec c) 4194b9ad928SBarry Smith { 420fccaa45eSBarry Smith PetscErrorCode ierr; 4214b9ad928SBarry Smith MG *mg = (MG*)pc->data; 4224b9ad928SBarry Smith 4234b9ad928SBarry Smith PetscFunctionBegin; 4244b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 425*630ba2eeSBarry Smith if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 4264b9ad928SBarry Smith mg[l]->r = c; 427fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 4284b9ad928SBarry Smith PetscFunctionReturn(0); 4294b9ad928SBarry Smith } 4304b9ad928SBarry Smith 4314b9ad928SBarry Smith 4324b9ad928SBarry Smith 4334b9ad928SBarry Smith 4344b9ad928SBarry Smith 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith 4374b9ad928SBarry Smith 438