1dba47a55SKris Buschelman #define PETSCKSP_DLL 24b9ad928SBarry Smith 3*7c4f633dSBarry 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__ 797177400SBarry Smith #define __FUNCT__ "PCMGDefaultResidual" 81f6cc5b2SSatish Balay /*@C 997177400SBarry Smith PCMGDefaultResidual - 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 2597177400SBarry Smith .seealso: PCMGSetResidual() 264b9ad928SBarry Smith @*/ 2797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGDefaultResidual(Mat mat,Vec b,Vec x,Vec r) 284b9ad928SBarry Smith { 29dfbe8321SBarry Smith PetscErrorCode ierr; 304b9ad928SBarry Smith 314b9ad928SBarry Smith PetscFunctionBegin; 324b9ad928SBarry Smith ierr = MatMult(mat,x,r);CHKERRQ(ierr); 33efb30889SBarry Smith ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 344b9ad928SBarry Smith PetscFunctionReturn(0); 354b9ad928SBarry Smith } 364b9ad928SBarry Smith 374b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/ 384b9ad928SBarry Smith 394b9ad928SBarry Smith #undef __FUNCT__ 40d6337fedSHong Zhang #define __FUNCT__ "PCMGGetCoarseSolve" 41f39d8e23SSatish Balay /*@ 4297177400SBarry Smith PCMGGetCoarseSolve - 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 @*/ 5697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetCoarseSolve(PC pc,KSP *ksp) 574b9ad928SBarry Smith { 589dcbbd2bSBarry Smith PC_MG **mg = (PC_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__ 66d6337fedSHong Zhang #define __FUNCT__ "PCMGSetResidual" 674b9ad928SBarry Smith /*@C 6897177400SBarry Smith PCMGSetResidual - 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 7697177400SBarry Smith . residual - function used to form residual (usually PCMGDefaultResidual) 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 8397177400SBarry Smith .seealso: PCMGDefaultResidual() 844b9ad928SBarry Smith @*/ 8597177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 864b9ad928SBarry Smith { 879dcbbd2bSBarry Smith PC_MG **mg = (PC_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__ 98d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation" 994b9ad928SBarry Smith /*@ 100aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 101bf5b2e24SBarry Smith interpolation from l-1 to 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 108bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 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 11997177400SBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetRestriction() 120fccaa45eSBarry Smith is used. 121fccaa45eSBarry Smith 1224b9ad928SBarry Smith .keywords: multigrid, set, interpolate, level 1234b9ad928SBarry Smith 12497177400SBarry Smith .seealso: PCMGSetRestriction() 1254b9ad928SBarry Smith @*/ 126aea2a34eSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 1274b9ad928SBarry Smith { 1289dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 129fccaa45eSBarry Smith PetscErrorCode ierr; 1304b9ad928SBarry Smith 1314b9ad928SBarry Smith PetscFunctionBegin; 1324b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1336d9142f8SRichard Katz if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 134c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 135c132512eSBarry Smith if (mg[l]->interpolate) {ierr = MatDestroy(mg[l]->interpolate);CHKERRQ(ierr);} 1364b9ad928SBarry Smith mg[l]->interpolate = mat; 1374b9ad928SBarry Smith PetscFunctionReturn(0); 1384b9ad928SBarry Smith } 1394b9ad928SBarry Smith 1404b9ad928SBarry Smith #undef __FUNCT__ 141d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction" 1424b9ad928SBarry Smith /*@ 14397177400SBarry Smith PCMGSetRestriction - 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 151bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 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 162aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 163fccaa45eSBarry Smith is used. 164fccaa45eSBarry Smith 1654b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level 1664b9ad928SBarry Smith 167aea2a34eSBarry Smith .seealso: PCMGSetInterpolation() 1684b9ad928SBarry Smith @*/ 16997177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 1704b9ad928SBarry Smith { 171fccaa45eSBarry Smith PetscErrorCode ierr; 1729dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 1734b9ad928SBarry Smith 1744b9ad928SBarry Smith PetscFunctionBegin; 1754b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1766d9142f8SRichard Katz if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 177c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 178c132512eSBarry Smith if (mg[l]->restrct) {ierr = MatDestroy(mg[l]->restrct);CHKERRQ(ierr);} 1794b9ad928SBarry Smith mg[l]->restrct = mat; 1804b9ad928SBarry Smith PetscFunctionReturn(0); 1814b9ad928SBarry Smith } 1824b9ad928SBarry Smith 1834b9ad928SBarry Smith #undef __FUNCT__ 184d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother" 185f39d8e23SSatish Balay /*@ 18697177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 18797177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 18897177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 1894b9ad928SBarry Smith post-smoothing. 1904b9ad928SBarry Smith 1914b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 1924b9ad928SBarry Smith 1934b9ad928SBarry Smith Input Parameters: 1944b9ad928SBarry Smith + pc - the multigrid context 1954b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 1964b9ad928SBarry Smith 1974b9ad928SBarry Smith Ouput Parameters: 1984b9ad928SBarry Smith . ksp - the smoother 1994b9ad928SBarry Smith 2004b9ad928SBarry Smith Level: advanced 2014b9ad928SBarry Smith 2024b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 2034b9ad928SBarry Smith 20497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 2054b9ad928SBarry Smith @*/ 20697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 2074b9ad928SBarry Smith { 2089dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 2094b9ad928SBarry Smith 2104b9ad928SBarry Smith PetscFunctionBegin; 2114b9ad928SBarry Smith *ksp = mg[l]->smoothd; 2124b9ad928SBarry Smith PetscFunctionReturn(0); 2134b9ad928SBarry Smith } 2144b9ad928SBarry Smith 2154b9ad928SBarry Smith #undef __FUNCT__ 216d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp" 217f39d8e23SSatish Balay /*@ 21897177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 2194b9ad928SBarry Smith coarse grid correction (post-smoother). 2204b9ad928SBarry Smith 2214b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 2224b9ad928SBarry Smith 2234b9ad928SBarry Smith Input Parameters: 2244b9ad928SBarry Smith + pc - the multigrid context 2254b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 2264b9ad928SBarry Smith 2274b9ad928SBarry Smith Ouput Parameters: 2284b9ad928SBarry Smith . ksp - the smoother 2294b9ad928SBarry Smith 2304b9ad928SBarry Smith Level: advanced 2314b9ad928SBarry Smith 2324b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level 2334b9ad928SBarry Smith 23497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 2354b9ad928SBarry Smith @*/ 23697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 2374b9ad928SBarry Smith { 2389dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 239dfbe8321SBarry Smith PetscErrorCode ierr; 240f69a0ea3SMatthew Knepley const char *prefix; 2414b9ad928SBarry Smith MPI_Comm comm; 2424b9ad928SBarry Smith 2434b9ad928SBarry Smith PetscFunctionBegin; 2444b9ad928SBarry Smith /* 2454b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 2464b9ad928SBarry Smith Thus we check if a different one has already been allocated, 2474b9ad928SBarry Smith if not we allocate it. 2484b9ad928SBarry Smith */ 249b03c7568SBarry Smith if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 2504b9ad928SBarry Smith if (mg[l]->smoothu == mg[l]->smoothd) { 2514b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[l]->smoothd,&comm);CHKERRQ(ierr); 2522e70eadcSBarry Smith ierr = KSPGetOptionsPrefix(mg[l]->smoothd,&prefix);CHKERRQ(ierr); 2534b9ad928SBarry Smith ierr = KSPCreate(comm,&mg[l]->smoothu);CHKERRQ(ierr); 2541cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mg[l]->smoothu,(PetscObject)pc,mg[0]->levels-l);CHKERRQ(ierr); 255d7d82daaSBarry Smith ierr = KSPSetTolerances(mg[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 2564b9ad928SBarry Smith ierr = KSPSetOptionsPrefix(mg[l]->smoothu,prefix);CHKERRQ(ierr); 25752e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,mg[l]->smoothu);CHKERRQ(ierr); 2584b9ad928SBarry Smith } 2594b9ad928SBarry Smith if (ksp) *ksp = mg[l]->smoothu; 2604b9ad928SBarry Smith PetscFunctionReturn(0); 2614b9ad928SBarry Smith } 2624b9ad928SBarry Smith 2634b9ad928SBarry Smith #undef __FUNCT__ 2649dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown" 265f39d8e23SSatish Balay /*@ 26697177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 2674b9ad928SBarry Smith coarse grid correction (pre-smoother). 2684b9ad928SBarry Smith 2694b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 2704b9ad928SBarry Smith 2714b9ad928SBarry Smith Input Parameters: 2724b9ad928SBarry Smith + pc - the multigrid context 2734b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 2744b9ad928SBarry Smith 2754b9ad928SBarry Smith Ouput Parameters: 2764b9ad928SBarry Smith . ksp - the smoother 2774b9ad928SBarry Smith 2784b9ad928SBarry Smith Level: advanced 2794b9ad928SBarry Smith 2804b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 2814b9ad928SBarry Smith 28297177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 2834b9ad928SBarry Smith @*/ 28497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 2854b9ad928SBarry Smith { 286dfbe8321SBarry Smith PetscErrorCode ierr; 2879dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 2884b9ad928SBarry Smith 2894b9ad928SBarry Smith PetscFunctionBegin; 2904b9ad928SBarry Smith /* make sure smoother up and down are different */ 291d8148a5aSMatthew Knepley if (l != 0) { 29297177400SBarry Smith ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr); 293d8148a5aSMatthew Knepley } 2944b9ad928SBarry Smith *ksp = mg[l]->smoothd; 2954b9ad928SBarry Smith PetscFunctionReturn(0); 2964b9ad928SBarry Smith } 2974b9ad928SBarry Smith 2984b9ad928SBarry Smith #undef __FUNCT__ 2999dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel" 3004b9ad928SBarry Smith /*@ 30197177400SBarry Smith PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level. 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith Collective on PC 3044b9ad928SBarry Smith 3054b9ad928SBarry Smith Input Parameters: 3064b9ad928SBarry Smith + pc - the multigrid context 3074b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3084b9ad928SBarry Smith - n - the number of cycles 3094b9ad928SBarry Smith 3104b9ad928SBarry Smith Level: advanced 3114b9ad928SBarry Smith 3124b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 3134b9ad928SBarry Smith 31497177400SBarry Smith .seealso: PCMGSetCycles() 3154b9ad928SBarry Smith @*/ 31697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c) 3174b9ad928SBarry Smith { 3189dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 3194b9ad928SBarry Smith 3204b9ad928SBarry Smith PetscFunctionBegin; 3214b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 3224b9ad928SBarry Smith mg[l]->cycles = c; 3234b9ad928SBarry Smith PetscFunctionReturn(0); 3244b9ad928SBarry Smith } 3254b9ad928SBarry Smith 3264b9ad928SBarry Smith #undef __FUNCT__ 327d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs" 3284b9ad928SBarry Smith /*@ 32997177400SBarry Smith PCMGSetRhs - Sets the vector space to be used to store the right-hand side 330fccaa45eSBarry Smith on a particular level. 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith Collective on PC and Vec 3334b9ad928SBarry Smith 3344b9ad928SBarry Smith Input Parameters: 3354b9ad928SBarry Smith + pc - the multigrid context 3364b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3374b9ad928SBarry Smith - c - the space 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith Level: advanced 3404b9ad928SBarry Smith 341fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 342fccaa45eSBarry Smith 343fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 344fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 345fccaa45eSBarry Smith 3464b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level 3474b9ad928SBarry Smith 34897177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 3494b9ad928SBarry Smith @*/ 35097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRhs(PC pc,PetscInt l,Vec c) 3514b9ad928SBarry Smith { 352fccaa45eSBarry Smith PetscErrorCode ierr; 3539dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 3544b9ad928SBarry Smith 3554b9ad928SBarry Smith PetscFunctionBegin; 3564b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 357fccaa45eSBarry Smith if (l == mg[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 358c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 359c132512eSBarry Smith if (mg[l]->b) {ierr = VecDestroy(mg[l]->b);CHKERRQ(ierr);} 3604b9ad928SBarry Smith mg[l]->b = c; 3614b9ad928SBarry Smith PetscFunctionReturn(0); 3624b9ad928SBarry Smith } 3634b9ad928SBarry Smith 3644b9ad928SBarry Smith #undef __FUNCT__ 365d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX" 3664b9ad928SBarry Smith /*@ 36797177400SBarry Smith PCMGSetX - Sets the vector space to be used to store the solution on a 368fccaa45eSBarry Smith particular level. 3694b9ad928SBarry Smith 3704b9ad928SBarry Smith Collective on PC and Vec 3714b9ad928SBarry Smith 3724b9ad928SBarry Smith Input Parameters: 3734b9ad928SBarry Smith + pc - the multigrid context 3744b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3754b9ad928SBarry Smith - c - the space 3764b9ad928SBarry Smith 3774b9ad928SBarry Smith Level: advanced 3784b9ad928SBarry Smith 379fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 380fccaa45eSBarry Smith 381fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 382fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 383fccaa45eSBarry Smith 3844b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level 3854b9ad928SBarry Smith 38697177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 3874b9ad928SBarry Smith @*/ 38897177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetX(PC pc,PetscInt l,Vec c) 3894b9ad928SBarry Smith { 390fccaa45eSBarry Smith PetscErrorCode ierr; 3919dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 3924b9ad928SBarry Smith 3934b9ad928SBarry Smith PetscFunctionBegin; 3944b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 395fccaa45eSBarry Smith if (l == mg[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 396c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 397c132512eSBarry Smith if (mg[l]->x) {ierr = VecDestroy(mg[l]->x);CHKERRQ(ierr);} 3984b9ad928SBarry Smith mg[l]->x = c; 3994b9ad928SBarry Smith PetscFunctionReturn(0); 4004b9ad928SBarry Smith } 4014b9ad928SBarry Smith 4024b9ad928SBarry Smith #undef __FUNCT__ 403d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR" 4044b9ad928SBarry Smith /*@ 40597177400SBarry Smith PCMGSetR - Sets the vector space to be used to store the residual on a 406fccaa45eSBarry Smith particular level. 4074b9ad928SBarry Smith 4084b9ad928SBarry Smith Collective on PC and Vec 4094b9ad928SBarry Smith 4104b9ad928SBarry Smith Input Parameters: 4114b9ad928SBarry Smith + pc - the multigrid context 4124b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 4134b9ad928SBarry Smith - c - the space 4144b9ad928SBarry Smith 4154b9ad928SBarry Smith Level: advanced 4164b9ad928SBarry Smith 417fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 418fccaa45eSBarry Smith 419fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 420fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 421fccaa45eSBarry Smith 4224b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level 4234b9ad928SBarry Smith @*/ 42497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetR(PC pc,PetscInt l,Vec c) 4254b9ad928SBarry Smith { 426fccaa45eSBarry Smith PetscErrorCode ierr; 4279dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 4284b9ad928SBarry Smith 4294b9ad928SBarry Smith PetscFunctionBegin; 4304b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 431630ba2eeSBarry Smith if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 432c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 433c132512eSBarry Smith if (mg[l]->r) {ierr = VecDestroy(mg[l]->r);CHKERRQ(ierr);} 4344b9ad928SBarry Smith mg[l]->r = c; 4354b9ad928SBarry Smith PetscFunctionReturn(0); 4364b9ad928SBarry Smith } 437