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__ 797177400SBarry Smith #define __FUNCT__ "PCMGDefaultResidual" 8*f39d8e23SSatish Balay /*@ 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 PetscScalar mone = -1.0; 314b9ad928SBarry Smith 324b9ad928SBarry Smith PetscFunctionBegin; 334b9ad928SBarry Smith ierr = MatMult(mat,x,r);CHKERRQ(ierr); 342dcb1b2aSMatthew Knepley ierr = VecAYPX(r,mone,b);CHKERRQ(ierr); 354b9ad928SBarry Smith PetscFunctionReturn(0); 364b9ad928SBarry Smith } 374b9ad928SBarry Smith 384b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/ 394b9ad928SBarry Smith 404b9ad928SBarry Smith #undef __FUNCT__ 414b9ad928SBarry Smith #define __FUNCT__ "MGGetCoarseSolve" 42*f39d8e23SSatish Balay /*@ 4397177400SBarry Smith PCMGGetCoarseSolve - 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 @*/ 5797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetCoarseSolve(PC pc,KSP *ksp) 584b9ad928SBarry Smith { 599dcbbd2bSBarry Smith PC_MG **mg = (PC_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 6997177400SBarry Smith PCMGSetResidual - 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 7797177400SBarry Smith . residual - function used to form residual (usually PCMGDefaultResidual) 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 8497177400SBarry Smith .seealso: PCMGDefaultResidual() 854b9ad928SBarry Smith @*/ 8697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 874b9ad928SBarry Smith { 889dcbbd2bSBarry Smith PC_MG **mg = (PC_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 /*@ 10197177400SBarry Smith PCMGSetInterpolate - 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 12097177400SBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetRestriction() 121fccaa45eSBarry Smith is used. 122fccaa45eSBarry Smith 1234b9ad928SBarry Smith .keywords: multigrid, set, interpolate, level 1244b9ad928SBarry Smith 12597177400SBarry Smith .seealso: PCMGSetRestriction() 1264b9ad928SBarry Smith @*/ 12797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetInterpolate(PC pc,PetscInt l,Mat mat) 1284b9ad928SBarry Smith { 1299dcbbd2bSBarry Smith PC_MG **mg = (PC_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"); 1346d9142f8SRichard Katz if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 135c132512eSBarry Smith if (mg[l]->interpolate) {ierr = MatDestroy(mg[l]->interpolate);CHKERRQ(ierr);} 1364b9ad928SBarry Smith mg[l]->interpolate = mat; 137fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1384b9ad928SBarry Smith PetscFunctionReturn(0); 1394b9ad928SBarry Smith } 1404b9ad928SBarry Smith 1414b9ad928SBarry Smith #undef __FUNCT__ 1424b9ad928SBarry Smith #define __FUNCT__ "MGSetRestriction" 1434b9ad928SBarry Smith /*@ 14497177400SBarry Smith PCMGSetRestriction - Sets the function to be used to restrict vector 1454b9ad928SBarry Smith from level l to l-1. 1464b9ad928SBarry Smith 1474b9ad928SBarry Smith Collective on PC and Mat 1484b9ad928SBarry Smith 1494b9ad928SBarry Smith Input Parameters: 1504b9ad928SBarry Smith + pc - the multigrid context 1514b9ad928SBarry Smith . mat - the restriction matrix 1524b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 1534b9ad928SBarry Smith 1544b9ad928SBarry Smith Level: advanced 1554b9ad928SBarry Smith 1564b9ad928SBarry Smith Notes: 1574b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 1584b9ad928SBarry Smith for the same level. 1594b9ad928SBarry Smith 1604b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1614b9ad928SBarry Smith out from the matrix size which one it is. 1624b9ad928SBarry Smith 16397177400SBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolate() 164fccaa45eSBarry Smith is used. 165fccaa45eSBarry Smith 1664b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level 1674b9ad928SBarry Smith 16897177400SBarry Smith .seealso: PCMGSetInterpolate() 1694b9ad928SBarry Smith @*/ 17097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 1714b9ad928SBarry Smith { 172fccaa45eSBarry Smith PetscErrorCode ierr; 1739dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 1744b9ad928SBarry Smith 1754b9ad928SBarry Smith PetscFunctionBegin; 1764b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1776d9142f8SRichard Katz if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 178c132512eSBarry Smith if (mg[l]->restrct) {ierr = MatDestroy(mg[l]->restrct);CHKERRQ(ierr);} 1794b9ad928SBarry Smith mg[l]->restrct = mat; 180fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1814b9ad928SBarry Smith PetscFunctionReturn(0); 1824b9ad928SBarry Smith } 1834b9ad928SBarry Smith 1844b9ad928SBarry Smith #undef __FUNCT__ 1854b9ad928SBarry Smith #define __FUNCT__ "MGGetSmoother" 186*f39d8e23SSatish Balay /*@ 18797177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 18897177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 18997177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 1904b9ad928SBarry Smith post-smoothing. 1914b9ad928SBarry Smith 1924b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 1934b9ad928SBarry Smith 1944b9ad928SBarry Smith Input Parameters: 1954b9ad928SBarry Smith + pc - the multigrid context 1964b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 1974b9ad928SBarry Smith 1984b9ad928SBarry Smith Ouput Parameters: 1994b9ad928SBarry Smith . ksp - the smoother 2004b9ad928SBarry Smith 2014b9ad928SBarry Smith Level: advanced 2024b9ad928SBarry Smith 2034b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 2044b9ad928SBarry Smith 20597177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 2064b9ad928SBarry Smith @*/ 20797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 2084b9ad928SBarry Smith { 2099dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 2104b9ad928SBarry Smith 2114b9ad928SBarry Smith PetscFunctionBegin; 2124b9ad928SBarry Smith *ksp = mg[l]->smoothd; 2134b9ad928SBarry Smith PetscFunctionReturn(0); 2144b9ad928SBarry Smith } 2154b9ad928SBarry Smith 2164b9ad928SBarry Smith #undef __FUNCT__ 2174b9ad928SBarry Smith #define __FUNCT__ "MGGetSmootherUp" 218*f39d8e23SSatish Balay /*@ 21997177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 2204b9ad928SBarry Smith coarse grid correction (post-smoother). 2214b9ad928SBarry Smith 2224b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith Input Parameters: 2254b9ad928SBarry Smith + pc - the multigrid context 2264b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 2274b9ad928SBarry Smith 2284b9ad928SBarry Smith Ouput Parameters: 2294b9ad928SBarry Smith . ksp - the smoother 2304b9ad928SBarry Smith 2314b9ad928SBarry Smith Level: advanced 2324b9ad928SBarry Smith 2334b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level 2344b9ad928SBarry Smith 23597177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 2364b9ad928SBarry Smith @*/ 23797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 2384b9ad928SBarry Smith { 2399dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 240dfbe8321SBarry Smith PetscErrorCode ierr; 241f69a0ea3SMatthew Knepley const char *prefix; 2424b9ad928SBarry Smith MPI_Comm comm; 2434b9ad928SBarry Smith 2444b9ad928SBarry Smith PetscFunctionBegin; 2454b9ad928SBarry Smith /* 2464b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 2474b9ad928SBarry Smith Thus we check if a different one has already been allocated, 2484b9ad928SBarry Smith if not we allocate it. 2494b9ad928SBarry Smith */ 250b03c7568SBarry Smith if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 2514b9ad928SBarry Smith if (mg[l]->smoothu == mg[l]->smoothd) { 2524b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[l]->smoothd,&comm);CHKERRQ(ierr); 2532e70eadcSBarry Smith ierr = KSPGetOptionsPrefix(mg[l]->smoothd,&prefix);CHKERRQ(ierr); 2544b9ad928SBarry Smith ierr = KSPCreate(comm,&mg[l]->smoothu);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" 265*f39d8e23SSatish 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 */ 29197177400SBarry Smith ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr); 2924b9ad928SBarry Smith *ksp = mg[l]->smoothd; 2934b9ad928SBarry Smith PetscFunctionReturn(0); 2944b9ad928SBarry Smith } 2954b9ad928SBarry Smith 2964b9ad928SBarry Smith #undef __FUNCT__ 2979dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel" 2984b9ad928SBarry Smith /*@ 29997177400SBarry Smith PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level. 3004b9ad928SBarry Smith 3014b9ad928SBarry Smith Collective on PC 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith Input Parameters: 3044b9ad928SBarry Smith + pc - the multigrid context 3054b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3064b9ad928SBarry Smith - n - the number of cycles 3074b9ad928SBarry Smith 3084b9ad928SBarry Smith Level: advanced 3094b9ad928SBarry Smith 3104b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 3114b9ad928SBarry Smith 31297177400SBarry Smith .seealso: PCMGSetCycles() 3134b9ad928SBarry Smith @*/ 31497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c) 3154b9ad928SBarry Smith { 3169dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 3174b9ad928SBarry Smith 3184b9ad928SBarry Smith PetscFunctionBegin; 3194b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 3204b9ad928SBarry Smith mg[l]->cycles = c; 3214b9ad928SBarry Smith PetscFunctionReturn(0); 3224b9ad928SBarry Smith } 3234b9ad928SBarry Smith 3244b9ad928SBarry Smith #undef __FUNCT__ 3254b9ad928SBarry Smith #define __FUNCT__ "MGSetRhs" 3264b9ad928SBarry Smith /*@ 32797177400SBarry Smith PCMGSetRhs - Sets the vector space to be used to store the right-hand side 328fccaa45eSBarry Smith on a particular level. 3294b9ad928SBarry Smith 3304b9ad928SBarry Smith Collective on PC and Vec 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith Input Parameters: 3334b9ad928SBarry Smith + pc - the multigrid context 3344b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3354b9ad928SBarry Smith - c - the space 3364b9ad928SBarry Smith 3374b9ad928SBarry Smith Level: advanced 3384b9ad928SBarry Smith 339fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 340fccaa45eSBarry Smith 341fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 342fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 343fccaa45eSBarry Smith 3444b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level 3454b9ad928SBarry Smith 34697177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 3474b9ad928SBarry Smith @*/ 34897177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRhs(PC pc,PetscInt l,Vec c) 3494b9ad928SBarry Smith { 350fccaa45eSBarry Smith PetscErrorCode ierr; 3519dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith PetscFunctionBegin; 3544b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 355fccaa45eSBarry Smith if (l == mg[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 356c132512eSBarry Smith if (mg[l]->b) {ierr = VecDestroy(mg[l]->b);CHKERRQ(ierr);} 3574b9ad928SBarry Smith mg[l]->b = c; 358fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 3594b9ad928SBarry Smith PetscFunctionReturn(0); 3604b9ad928SBarry Smith } 3614b9ad928SBarry Smith 3624b9ad928SBarry Smith #undef __FUNCT__ 3634b9ad928SBarry Smith #define __FUNCT__ "MGSetX" 3644b9ad928SBarry Smith /*@ 36597177400SBarry Smith PCMGSetX - Sets the vector space to be used to store the solution on a 366fccaa45eSBarry Smith particular level. 3674b9ad928SBarry Smith 3684b9ad928SBarry Smith Collective on PC and Vec 3694b9ad928SBarry Smith 3704b9ad928SBarry Smith Input Parameters: 3714b9ad928SBarry Smith + pc - the multigrid context 3724b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 3734b9ad928SBarry Smith - c - the space 3744b9ad928SBarry Smith 3754b9ad928SBarry Smith Level: advanced 3764b9ad928SBarry Smith 377fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 378fccaa45eSBarry Smith 379fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 380fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 381fccaa45eSBarry Smith 3824b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level 3834b9ad928SBarry Smith 38497177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 3854b9ad928SBarry Smith @*/ 38697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetX(PC pc,PetscInt l,Vec c) 3874b9ad928SBarry Smith { 388fccaa45eSBarry Smith PetscErrorCode ierr; 3899dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith PetscFunctionBegin; 3924b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 393fccaa45eSBarry Smith if (l == mg[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 394c132512eSBarry Smith if (mg[l]->x) {ierr = VecDestroy(mg[l]->x);CHKERRQ(ierr);} 3954b9ad928SBarry Smith mg[l]->x = c; 396fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 3974b9ad928SBarry Smith PetscFunctionReturn(0); 3984b9ad928SBarry Smith } 3994b9ad928SBarry Smith 4004b9ad928SBarry Smith #undef __FUNCT__ 4014b9ad928SBarry Smith #define __FUNCT__ "MGSetR" 4024b9ad928SBarry Smith /*@ 40397177400SBarry Smith PCMGSetR - Sets the vector space to be used to store the residual on a 404fccaa45eSBarry Smith particular level. 4054b9ad928SBarry Smith 4064b9ad928SBarry Smith Collective on PC and Vec 4074b9ad928SBarry Smith 4084b9ad928SBarry Smith Input Parameters: 4094b9ad928SBarry Smith + pc - the multigrid context 4104b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 4114b9ad928SBarry Smith - c - the space 4124b9ad928SBarry Smith 4134b9ad928SBarry Smith Level: advanced 4144b9ad928SBarry Smith 415fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 416fccaa45eSBarry Smith 417fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 418fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 419fccaa45eSBarry Smith 4204b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level 4214b9ad928SBarry Smith @*/ 42297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetR(PC pc,PetscInt l,Vec c) 4234b9ad928SBarry Smith { 424fccaa45eSBarry Smith PetscErrorCode ierr; 4259dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 4264b9ad928SBarry Smith 4274b9ad928SBarry Smith PetscFunctionBegin; 4284b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 429630ba2eeSBarry Smith if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 430c132512eSBarry Smith if (mg[l]->r) {ierr = VecDestroy(mg[l]->r);CHKERRQ(ierr);} 4314b9ad928SBarry Smith mg[l]->r = c; 432fccaa45eSBarry Smith ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 4334b9ad928SBarry Smith PetscFunctionReturn(0); 4344b9ad928SBarry Smith } 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith 4374b9ad928SBarry Smith 4384b9ad928SBarry Smith 4394b9ad928SBarry Smith 4404b9ad928SBarry Smith 4414b9ad928SBarry Smith 4424b9ad928SBarry Smith 443