14b9ad928SBarry Smith 2c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscksp.h" I*/ 3f3f935beSJed Brown /*I "petscpcmg.h" I*/ 44b9ad928SBarry Smith 54b9ad928SBarry Smith #undef __FUNCT__ 697177400SBarry Smith #define __FUNCT__ "PCMGDefaultResidual" 71f6cc5b2SSatish Balay /*@C 897177400SBarry Smith PCMGDefaultResidual - 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 2497177400SBarry Smith .seealso: PCMGSetResidual() 254b9ad928SBarry Smith @*/ 267087cfbeSBarry Smith PetscErrorCode PCMGDefaultResidual(Mat mat,Vec b,Vec x,Vec r) 274b9ad928SBarry Smith { 28dfbe8321SBarry Smith PetscErrorCode ierr; 294b9ad928SBarry Smith 304b9ad928SBarry Smith PetscFunctionBegin; 314b9ad928SBarry Smith ierr = MatMult(mat,x,r);CHKERRQ(ierr); 32efb30889SBarry Smith ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 334b9ad928SBarry Smith PetscFunctionReturn(0); 344b9ad928SBarry Smith } 354b9ad928SBarry Smith 364b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/ 374b9ad928SBarry Smith 384b9ad928SBarry Smith #undef __FUNCT__ 39d6337fedSHong Zhang #define __FUNCT__ "PCMGGetCoarseSolve" 40f39d8e23SSatish Balay /*@ 4197177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 424b9ad928SBarry Smith 434b9ad928SBarry Smith Not Collective 444b9ad928SBarry Smith 454b9ad928SBarry Smith Input Parameter: 464b9ad928SBarry Smith . pc - the multigrid context 474b9ad928SBarry Smith 484b9ad928SBarry Smith Output Parameter: 494b9ad928SBarry Smith . ksp - the coarse grid solver context 504b9ad928SBarry Smith 514b9ad928SBarry Smith Level: advanced 524b9ad928SBarry Smith 534b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid 544b9ad928SBarry Smith @*/ 557087cfbeSBarry Smith PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 564b9ad928SBarry Smith { 57f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 58f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 594b9ad928SBarry Smith 604b9ad928SBarry Smith PetscFunctionBegin; 61c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 62f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 634b9ad928SBarry Smith PetscFunctionReturn(0); 644b9ad928SBarry Smith } 654b9ad928SBarry Smith 664b9ad928SBarry Smith #undef __FUNCT__ 67d6337fedSHong Zhang #define __FUNCT__ "PCMGSetResidual" 684b9ad928SBarry Smith /*@C 6997177400SBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual 704b9ad928SBarry Smith on the lth level. 714b9ad928SBarry Smith 72ad4df100SBarry Smith Logically 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 77157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 78157726a2SBarry Smith previous one were provided then PCMGDefaultResidual() is used 794b9ad928SBarry Smith - mat - matrix associated with residual 804b9ad928SBarry Smith 814b9ad928SBarry Smith Level: advanced 824b9ad928SBarry Smith 834b9ad928SBarry Smith .keywords: MG, set, multigrid, residual, level 844b9ad928SBarry Smith 8597177400SBarry Smith .seealso: PCMGDefaultResidual() 864b9ad928SBarry Smith @*/ 877087cfbeSBarry Smith PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 884b9ad928SBarry Smith { 89f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 90f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 91298cc208SBarry Smith PetscErrorCode ierr; 924b9ad928SBarry Smith 934b9ad928SBarry Smith PetscFunctionBegin; 94c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 95e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 96157726a2SBarry Smith if (residual) { 97f3fbd535SBarry Smith mglevels[l]->residual = residual; 98157726a2SBarry Smith } if (!mglevels[l]->residual) { 99157726a2SBarry Smith mglevels[l]->residual = PCMGDefaultResidual; 100157726a2SBarry Smith } 101f3ae41bdSBarry Smith if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 102298cc208SBarry Smith ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 103f3fbd535SBarry Smith mglevels[l]->A = mat; 1044b9ad928SBarry Smith PetscFunctionReturn(0); 1054b9ad928SBarry Smith } 1064b9ad928SBarry Smith 1074b9ad928SBarry Smith #undef __FUNCT__ 108d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation" 1094b9ad928SBarry Smith /*@ 110aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 111bf5b2e24SBarry Smith interpolation from l-1 to the lth level 1124b9ad928SBarry Smith 113ad4df100SBarry Smith Logically Collective on PC and Mat 1144b9ad928SBarry Smith 1154b9ad928SBarry Smith Input Parameters: 1164b9ad928SBarry Smith + pc - the multigrid context 1174b9ad928SBarry Smith . mat - the interpolation operator 118bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 1194b9ad928SBarry Smith 1204b9ad928SBarry Smith Level: advanced 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith Notes: 1234b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 1244b9ad928SBarry Smith for the same level. 1254b9ad928SBarry Smith 1264b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1274b9ad928SBarry Smith out from the matrix size which one it is. 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith .keywords: multigrid, set, interpolate, level 1304b9ad928SBarry Smith 13197177400SBarry Smith .seealso: PCMGSetRestriction() 1324b9ad928SBarry Smith @*/ 1337087cfbeSBarry Smith PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 1344b9ad928SBarry Smith { 135f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 136f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 137fccaa45eSBarry Smith PetscErrorCode ierr; 1384b9ad928SBarry Smith 1394b9ad928SBarry Smith PetscFunctionBegin; 140c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 141e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 142e7e72b3dSBarry Smith if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 143c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1446bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr); 145f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 1464b9ad928SBarry Smith PetscFunctionReturn(0); 1474b9ad928SBarry Smith } 1484b9ad928SBarry Smith 1494b9ad928SBarry Smith #undef __FUNCT__ 150c88c5224SJed Brown #define __FUNCT__ "PCMGGetInterpolation" 151c88c5224SJed Brown /*@ 152c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 153c88c5224SJed Brown interpolation from l-1 to the lth level 154c88c5224SJed Brown 155c88c5224SJed Brown Logically Collective on PC 156c88c5224SJed Brown 157c88c5224SJed Brown Input Parameters: 158c88c5224SJed Brown + pc - the multigrid context 159c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 160c88c5224SJed Brown 161c88c5224SJed Brown Output Parameter: 162c88c5224SJed Brown . mat - the interpolation matrix 163c88c5224SJed Brown 164c88c5224SJed Brown Level: advanced 165c88c5224SJed Brown 166c88c5224SJed Brown .keywords: MG, get, multigrid, interpolation, level 167c88c5224SJed Brown 168c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 169c88c5224SJed Brown @*/ 170c88c5224SJed Brown PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 171c88c5224SJed Brown { 172c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 173c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 174c88c5224SJed Brown PetscErrorCode ierr; 175c88c5224SJed Brown 176c88c5224SJed Brown PetscFunctionBegin; 177c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 178c88c5224SJed Brown PetscValidPointer(mat,3); 179c88c5224SJed Brown if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 180c88c5224SJed Brown if (l <= 0 || mg->nlevels <= l) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 181c88c5224SJed Brown if (!mglevels[l]->interpolate) { 182c88c5224SJed Brown if (!mglevels[l]->restrct) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetInterpolation()"); 183c88c5224SJed Brown ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 184c88c5224SJed Brown } 185c88c5224SJed Brown *mat = mglevels[l]->interpolate; 186c88c5224SJed Brown PetscFunctionReturn(0); 187c88c5224SJed Brown } 188c88c5224SJed Brown 189c88c5224SJed Brown #undef __FUNCT__ 190d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction" 1914b9ad928SBarry Smith /*@ 19297177400SBarry Smith PCMGSetRestriction - Sets the function to be used to restrict vector 1934b9ad928SBarry Smith from level l to l-1. 1944b9ad928SBarry Smith 195ad4df100SBarry Smith Logically Collective on PC and Mat 1964b9ad928SBarry Smith 1974b9ad928SBarry Smith Input Parameters: 1984b9ad928SBarry Smith + pc - the multigrid context 199c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 200c88c5224SJed Brown - mat - the restriction matrix 2014b9ad928SBarry Smith 2024b9ad928SBarry Smith Level: advanced 2034b9ad928SBarry Smith 2044b9ad928SBarry Smith Notes: 2054b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 2064b9ad928SBarry Smith for the same level. 2074b9ad928SBarry Smith 2084b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 2094b9ad928SBarry Smith out from the matrix size which one it is. 2104b9ad928SBarry Smith 211aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 212fccaa45eSBarry Smith is used. 213fccaa45eSBarry Smith 2144b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level 2154b9ad928SBarry Smith 216aea2a34eSBarry Smith .seealso: PCMGSetInterpolation() 2174b9ad928SBarry Smith @*/ 2187087cfbeSBarry Smith PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 2194b9ad928SBarry Smith { 220fccaa45eSBarry Smith PetscErrorCode ierr; 221f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 222f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith PetscFunctionBegin; 225c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 226c88c5224SJed Brown PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 227e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 228e7e72b3dSBarry Smith if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 229c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 2306bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 231f3fbd535SBarry Smith mglevels[l]->restrct = mat; 2324b9ad928SBarry Smith PetscFunctionReturn(0); 2334b9ad928SBarry Smith } 2344b9ad928SBarry Smith 2354b9ad928SBarry Smith #undef __FUNCT__ 236c88c5224SJed Brown #define __FUNCT__ "PCMGGetRestriction" 237c88c5224SJed Brown /*@ 238c88c5224SJed Brown PCMGGetRestriction - Gets the function to be used to restrict vector 239c88c5224SJed Brown from level l to l-1. 240c88c5224SJed Brown 241c88c5224SJed Brown Logically Collective on PC and Mat 242c88c5224SJed Brown 243c88c5224SJed Brown Input Parameters: 244c88c5224SJed Brown + pc - the multigrid context 245c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 246c88c5224SJed Brown 247c88c5224SJed Brown Output Parameter: 248c88c5224SJed Brown . mat - the restriction matrix 249c88c5224SJed Brown 250c88c5224SJed Brown Level: advanced 251c88c5224SJed Brown 252c88c5224SJed Brown .keywords: MG, get, multigrid, restriction, level 253c88c5224SJed Brown 254c88c5224SJed Brown .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale() 255c88c5224SJed Brown @*/ 256c88c5224SJed Brown PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 257c88c5224SJed Brown { 258c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 259c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 260c88c5224SJed Brown PetscErrorCode ierr; 261c88c5224SJed Brown 262c88c5224SJed Brown PetscFunctionBegin; 263c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 264c88c5224SJed Brown PetscValidPointer(mat,3); 265c88c5224SJed Brown if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 266c88c5224SJed Brown if (l <= 0 || mg->nlevels <= l) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 267c88c5224SJed Brown if (!mglevels[l]->restrct) { 268c88c5224SJed Brown if (!mglevels[l]->interpolate) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 269c88c5224SJed Brown ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 270c88c5224SJed Brown } 271c88c5224SJed Brown *mat = mglevels[l]->restrct; 272c88c5224SJed Brown PetscFunctionReturn(0); 273c88c5224SJed Brown } 274c88c5224SJed Brown 275c88c5224SJed Brown #undef __FUNCT__ 27673250ac0SBarry Smith #define __FUNCT__ "PCMGSetRScale" 27773250ac0SBarry Smith /*@ 27873250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 27973250ac0SBarry Smith 280c88c5224SJed Brown Logically Collective on PC and Vec 281c88c5224SJed Brown 282c88c5224SJed Brown Input Parameters: 283c88c5224SJed Brown + pc - the multigrid context 284c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 285c88c5224SJed Brown . rscale - the scaling 286c88c5224SJed Brown 287c88c5224SJed Brown Level: advanced 288c88c5224SJed Brown 289c88c5224SJed Brown Notes: 290c88c5224SJed Brown When evaluating a function on a coarse level one does not want to do F( R * x) one does F( rscale * R * x) where rscale is 1 over the row sums of R. 291c88c5224SJed Brown 292c88c5224SJed Brown .keywords: MG, set, multigrid, restriction, level 293c88c5224SJed Brown 294c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale() 295c88c5224SJed Brown @*/ 296c88c5224SJed Brown PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 297c88c5224SJed Brown { 298c88c5224SJed Brown PetscErrorCode ierr; 299c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 300c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 301c88c5224SJed Brown 302c88c5224SJed Brown PetscFunctionBegin; 303c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 304c88c5224SJed Brown if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 305c88c5224SJed Brown if (l <= 0 || mg->nlevels <= l) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 306c88c5224SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 307c88c5224SJed Brown ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 308c88c5224SJed Brown mglevels[l]->rscale = rscale; 309c88c5224SJed Brown PetscFunctionReturn(0); 310c88c5224SJed Brown } 311c88c5224SJed Brown 312c88c5224SJed Brown #undef __FUNCT__ 313c88c5224SJed Brown #define __FUNCT__ "PCMGGetRScale" 314c88c5224SJed Brown /*@ 315c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 316c88c5224SJed Brown 317c88c5224SJed Brown Collective on PC 31873250ac0SBarry Smith 31973250ac0SBarry Smith Input Parameters: 32073250ac0SBarry Smith + pc - the multigrid context 32173250ac0SBarry Smith . rscale - the scaling 32273250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 32373250ac0SBarry Smith 32473250ac0SBarry Smith Level: advanced 32573250ac0SBarry Smith 32673250ac0SBarry Smith Notes: 32773250ac0SBarry Smith When evaluating a function on a coarse level one does not want to do F( R * x) one does F( rscale * R * x) where rscale is 1 over the row sums of R. 32873250ac0SBarry Smith 32973250ac0SBarry Smith .keywords: MG, set, multigrid, restriction, level 33073250ac0SBarry Smith 331c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGGetRestriction() 33273250ac0SBarry Smith @*/ 333c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 33473250ac0SBarry Smith { 33573250ac0SBarry Smith PetscErrorCode ierr; 33673250ac0SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 33773250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 33873250ac0SBarry Smith 33973250ac0SBarry Smith PetscFunctionBegin; 34073250ac0SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 34173250ac0SBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 342c88c5224SJed Brown if (l <= 0 || mg->nlevels <= l) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 343c88c5224SJed Brown if (!mglevels[l]->rscale) { 344c88c5224SJed Brown Mat R; 345c88c5224SJed Brown Vec X,Y,coarse,fine; 346c88c5224SJed Brown PetscInt M,N; 347c88c5224SJed Brown ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 348c88c5224SJed Brown ierr = MatGetVecs(R,&X,&Y);CHKERRQ(ierr); 349c88c5224SJed Brown ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 350c88c5224SJed Brown if (M < N) {fine = X; coarse = Y;} 351c88c5224SJed Brown else if (N < M) {fine = Y; coarse = X;} 352c88c5224SJed Brown else SETERRQ(((PetscObject)R)->comm,PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 353c88c5224SJed Brown ierr = VecSet(fine,1.);CHKERRQ(ierr); 354c88c5224SJed Brown ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 355c88c5224SJed Brown ierr = VecDestroy(&fine);CHKERRQ(ierr); 356c88c5224SJed Brown ierr = VecReciprocal(coarse);CHKERRQ(ierr); 357c88c5224SJed Brown mglevels[l]->rscale = coarse; 358c88c5224SJed Brown } 359c88c5224SJed Brown *rscale = mglevels[l]->rscale; 36073250ac0SBarry Smith PetscFunctionReturn(0); 36173250ac0SBarry Smith } 36273250ac0SBarry Smith 36373250ac0SBarry Smith #undef __FUNCT__ 364d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother" 365f39d8e23SSatish Balay /*@ 36697177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 36797177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 36897177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 3694b9ad928SBarry Smith post-smoothing. 3704b9ad928SBarry Smith 3714b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 3724b9ad928SBarry Smith 3734b9ad928SBarry Smith Input Parameters: 3744b9ad928SBarry Smith + pc - the multigrid context 3754b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 3764b9ad928SBarry Smith 3774b9ad928SBarry Smith Ouput Parameters: 3784b9ad928SBarry Smith . ksp - the smoother 3794b9ad928SBarry Smith 3804b9ad928SBarry Smith Level: advanced 3814b9ad928SBarry Smith 3824b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 3834b9ad928SBarry Smith 38497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 3854b9ad928SBarry Smith @*/ 3867087cfbeSBarry Smith PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 3874b9ad928SBarry Smith { 388f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 389f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith PetscFunctionBegin; 392c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 393f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 3944b9ad928SBarry Smith PetscFunctionReturn(0); 3954b9ad928SBarry Smith } 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith #undef __FUNCT__ 398d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp" 399f39d8e23SSatish Balay /*@ 40097177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 4014b9ad928SBarry Smith coarse grid correction (post-smoother). 4024b9ad928SBarry Smith 4034b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4044b9ad928SBarry Smith 4054b9ad928SBarry Smith Input Parameters: 4064b9ad928SBarry Smith + pc - the multigrid context 4074b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4084b9ad928SBarry Smith 4094b9ad928SBarry Smith Ouput Parameters: 4104b9ad928SBarry Smith . ksp - the smoother 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith Level: advanced 4134b9ad928SBarry Smith 41489cce641SBarry Smith Notes: calling this will result in a different pre and post smoother so you may need to 41589cce641SBarry Smith set options on the pre smoother also 41689cce641SBarry Smith 4174b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level 4184b9ad928SBarry Smith 41997177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 4204b9ad928SBarry Smith @*/ 4217087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 4224b9ad928SBarry Smith { 423f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 424f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 425dfbe8321SBarry Smith PetscErrorCode ierr; 426f69a0ea3SMatthew Knepley const char *prefix; 4274b9ad928SBarry Smith MPI_Comm comm; 4284b9ad928SBarry Smith 4294b9ad928SBarry Smith PetscFunctionBegin; 430c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4314b9ad928SBarry Smith /* 4324b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 4334b9ad928SBarry Smith Thus we check if a different one has already been allocated, 4344b9ad928SBarry Smith if not we allocate it. 4354b9ad928SBarry Smith */ 436e7e72b3dSBarry Smith if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 437f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 43819fd82e9SBarry Smith KSPType ksptype; 43919fd82e9SBarry Smith PCType pctype; 440336babb1SJed Brown PC ipc; 441336babb1SJed Brown PetscReal rtol,abstol,dtol; 442336babb1SJed Brown PetscInt maxits; 443336babb1SJed Brown KSPNormType normtype; 444f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 445f3fbd535SBarry Smith ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 446336babb1SJed Brown ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 447*3bf036e2SBarry Smith ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 448336babb1SJed Brown ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 449336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 450336babb1SJed Brown ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 451336babb1SJed Brown 452f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 453f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 454f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 455336babb1SJed Brown ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 456336babb1SJed Brown ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 457336babb1SJed Brown ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 458336babb1SJed Brown ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 459336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 460336babb1SJed Brown ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 461f3fbd535SBarry Smith ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr); 4624b9ad928SBarry Smith } 463f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 4644b9ad928SBarry Smith PetscFunctionReturn(0); 4654b9ad928SBarry Smith } 4664b9ad928SBarry Smith 4674b9ad928SBarry Smith #undef __FUNCT__ 4689dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown" 469f39d8e23SSatish Balay /*@ 47097177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 4714b9ad928SBarry Smith coarse grid correction (pre-smoother). 4724b9ad928SBarry Smith 4734b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4744b9ad928SBarry Smith 4754b9ad928SBarry Smith Input Parameters: 4764b9ad928SBarry Smith + pc - the multigrid context 4774b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4784b9ad928SBarry Smith 4794b9ad928SBarry Smith Ouput Parameters: 4804b9ad928SBarry Smith . ksp - the smoother 4814b9ad928SBarry Smith 4824b9ad928SBarry Smith Level: advanced 4834b9ad928SBarry Smith 48489cce641SBarry Smith Notes: calling this will result in a different pre and post smoother so you may need to 48589cce641SBarry Smith set options on the post smoother also 48689cce641SBarry Smith 4874b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 4884b9ad928SBarry Smith 48997177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 4904b9ad928SBarry Smith @*/ 4917087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 4924b9ad928SBarry Smith { 493dfbe8321SBarry Smith PetscErrorCode ierr; 494f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 495f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 4964b9ad928SBarry Smith 4974b9ad928SBarry Smith PetscFunctionBegin; 498c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4994b9ad928SBarry Smith /* make sure smoother up and down are different */ 500c5eb9154SBarry Smith if (l) { 50197177400SBarry Smith ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr); 502d8148a5aSMatthew Knepley } 503f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5044b9ad928SBarry Smith PetscFunctionReturn(0); 5054b9ad928SBarry Smith } 5064b9ad928SBarry Smith 5074b9ad928SBarry Smith #undef __FUNCT__ 5089dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel" 5094b9ad928SBarry Smith /*@ 51097177400SBarry Smith PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level. 5114b9ad928SBarry Smith 512ad4df100SBarry Smith Logically Collective on PC 5134b9ad928SBarry Smith 5144b9ad928SBarry Smith Input Parameters: 5154b9ad928SBarry Smith + pc - the multigrid context 5164b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 5174b9ad928SBarry Smith - n - the number of cycles 5184b9ad928SBarry Smith 5194b9ad928SBarry Smith Level: advanced 5204b9ad928SBarry Smith 5214b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 5224b9ad928SBarry Smith 52397177400SBarry Smith .seealso: PCMGSetCycles() 5244b9ad928SBarry Smith @*/ 5257087cfbeSBarry Smith PetscErrorCode PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c) 5264b9ad928SBarry Smith { 527f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 528f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5294b9ad928SBarry Smith 5304b9ad928SBarry Smith PetscFunctionBegin; 531c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 532e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 533c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,l,2); 534c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,c,3); 535f3fbd535SBarry Smith mglevels[l]->cycles = c; 5364b9ad928SBarry Smith PetscFunctionReturn(0); 5374b9ad928SBarry Smith } 5384b9ad928SBarry Smith 5394b9ad928SBarry Smith #undef __FUNCT__ 540d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs" 5414b9ad928SBarry Smith /*@ 54297177400SBarry Smith PCMGSetRhs - Sets the vector space to be used to store the right-hand side 543fccaa45eSBarry Smith on a particular level. 5444b9ad928SBarry Smith 545ad4df100SBarry Smith Logically Collective on PC and Vec 5464b9ad928SBarry Smith 5474b9ad928SBarry Smith Input Parameters: 5484b9ad928SBarry Smith + pc - the multigrid context 5494b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 5504b9ad928SBarry Smith - c - the space 5514b9ad928SBarry Smith 5524b9ad928SBarry Smith Level: advanced 5534b9ad928SBarry Smith 554fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 555fccaa45eSBarry Smith 556fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 557fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 558fccaa45eSBarry Smith 5594b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level 5604b9ad928SBarry Smith 56197177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 5624b9ad928SBarry Smith @*/ 5637087cfbeSBarry Smith PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 5644b9ad928SBarry Smith { 565fccaa45eSBarry Smith PetscErrorCode ierr; 566f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 567f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5684b9ad928SBarry Smith 5694b9ad928SBarry Smith PetscFunctionBegin; 570c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 571e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 572e7e72b3dSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 573c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 5746bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 575f3fbd535SBarry Smith mglevels[l]->b = c; 5764b9ad928SBarry Smith PetscFunctionReturn(0); 5774b9ad928SBarry Smith } 5784b9ad928SBarry Smith 5794b9ad928SBarry Smith #undef __FUNCT__ 580d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX" 5814b9ad928SBarry Smith /*@ 58297177400SBarry Smith PCMGSetX - Sets the vector space to be used to store the solution on a 583fccaa45eSBarry Smith particular level. 5844b9ad928SBarry Smith 585ad4df100SBarry Smith Logically Collective on PC and Vec 5864b9ad928SBarry Smith 5874b9ad928SBarry Smith Input Parameters: 5884b9ad928SBarry Smith + pc - the multigrid context 589251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 5904b9ad928SBarry Smith - c - the space 5914b9ad928SBarry Smith 5924b9ad928SBarry Smith Level: advanced 5934b9ad928SBarry Smith 594fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 595fccaa45eSBarry Smith 596fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 597fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 598fccaa45eSBarry Smith 5994b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level 6004b9ad928SBarry Smith 60197177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 6024b9ad928SBarry Smith @*/ 6037087cfbeSBarry Smith PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 6044b9ad928SBarry Smith { 605fccaa45eSBarry Smith PetscErrorCode ierr; 606f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 607f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6084b9ad928SBarry Smith 6094b9ad928SBarry Smith PetscFunctionBegin; 610c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 611e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 612e7e72b3dSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 613c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6146bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 615f3fbd535SBarry Smith mglevels[l]->x = c; 6164b9ad928SBarry Smith PetscFunctionReturn(0); 6174b9ad928SBarry Smith } 6184b9ad928SBarry Smith 6194b9ad928SBarry Smith #undef __FUNCT__ 620d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR" 6214b9ad928SBarry Smith /*@ 62297177400SBarry Smith PCMGSetR - Sets the vector space to be used to store the residual on a 623fccaa45eSBarry Smith particular level. 6244b9ad928SBarry Smith 625ad4df100SBarry Smith Logically Collective on PC and Vec 6264b9ad928SBarry Smith 6274b9ad928SBarry Smith Input Parameters: 6284b9ad928SBarry Smith + pc - the multigrid context 6294b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 6304b9ad928SBarry Smith - c - the space 6314b9ad928SBarry Smith 6324b9ad928SBarry Smith Level: advanced 6334b9ad928SBarry Smith 634fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 635fccaa45eSBarry Smith 636fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 637fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 638fccaa45eSBarry Smith 6394b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level 6404b9ad928SBarry Smith @*/ 6417087cfbeSBarry Smith PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 6424b9ad928SBarry Smith { 643fccaa45eSBarry Smith PetscErrorCode ierr; 644f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 645f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6464b9ad928SBarry Smith 6474b9ad928SBarry Smith PetscFunctionBegin; 648c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 649e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 650e7e72b3dSBarry Smith if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 651c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6526bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 653f3fbd535SBarry Smith mglevels[l]->r = c; 6544b9ad928SBarry Smith PetscFunctionReturn(0); 6554b9ad928SBarry Smith } 656