14b9ad928SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 34b9ad928SBarry Smith 4f9426fe0SMark Adams /* ---------------------------------------------------------------------------*/ 51f6cc5b2SSatish Balay /*@C 654b2cd4bSJed Brown PCMGResidualDefault - Default routine to calculate the residual. 74b9ad928SBarry Smith 8d083f849SBarry Smith Collective on Mat 94b9ad928SBarry Smith 104b9ad928SBarry Smith Input Parameters: 114b9ad928SBarry Smith + mat - the matrix 124b9ad928SBarry Smith . b - the right-hand-side 134b9ad928SBarry Smith - x - the approximate solution 144b9ad928SBarry Smith 154b9ad928SBarry Smith Output Parameter: 164b9ad928SBarry Smith . r - location to store the residual 174b9ad928SBarry Smith 18d0e4de75SBarry Smith Level: developer 194b9ad928SBarry Smith 2097177400SBarry Smith .seealso: PCMGSetResidual() 214b9ad928SBarry Smith @*/ 2254b2cd4bSJed Brown PetscErrorCode PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r) 234b9ad928SBarry Smith { 24dfbe8321SBarry Smith PetscErrorCode ierr; 254b9ad928SBarry Smith 264b9ad928SBarry Smith PetscFunctionBegin; 27f9426fe0SMark Adams ierr = MatResidual(mat,b,x,r);CHKERRQ(ierr); 284b9ad928SBarry Smith PetscFunctionReturn(0); 294b9ad928SBarry Smith } 304b9ad928SBarry Smith 31*fcb023d4SJed Brown /*@C 32*fcb023d4SJed Brown PCMGResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system 33*fcb023d4SJed Brown 34*fcb023d4SJed Brown Collective on Mat 35*fcb023d4SJed Brown 36*fcb023d4SJed Brown Input Parameters: 37*fcb023d4SJed Brown + mat - the matrix 38*fcb023d4SJed Brown . b - the right-hand-side 39*fcb023d4SJed Brown - x - the approximate solution 40*fcb023d4SJed Brown 41*fcb023d4SJed Brown Output Parameter: 42*fcb023d4SJed Brown . r - location to store the residual 43*fcb023d4SJed Brown 44*fcb023d4SJed Brown Level: developer 45*fcb023d4SJed Brown 46*fcb023d4SJed Brown .seealso: PCMGSetResidualTranspose() 47*fcb023d4SJed Brown @*/ 48*fcb023d4SJed Brown PetscErrorCode PCMGResidualTransposeDefault(Mat mat,Vec b,Vec x,Vec r) 49*fcb023d4SJed Brown { 50*fcb023d4SJed Brown PetscErrorCode ierr; 51*fcb023d4SJed Brown 52*fcb023d4SJed Brown PetscFunctionBegin; 53*fcb023d4SJed Brown ierr = MatMultTranspose(mat,x,r);CHKERRQ(ierr); 54*fcb023d4SJed Brown ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 55*fcb023d4SJed Brown PetscFunctionReturn(0); 56*fcb023d4SJed Brown } 57*fcb023d4SJed Brown 58f39d8e23SSatish Balay /*@ 5997177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 604b9ad928SBarry Smith 614b9ad928SBarry Smith Not Collective 624b9ad928SBarry Smith 634b9ad928SBarry Smith Input Parameter: 644b9ad928SBarry Smith . pc - the multigrid context 654b9ad928SBarry Smith 664b9ad928SBarry Smith Output Parameter: 674b9ad928SBarry Smith . ksp - the coarse grid solver context 684b9ad928SBarry Smith 694b9ad928SBarry Smith Level: advanced 704b9ad928SBarry Smith 7128668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetSmoother() 724b9ad928SBarry Smith @*/ 737087cfbeSBarry Smith PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 744b9ad928SBarry Smith { 75f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 76f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 774b9ad928SBarry Smith 784b9ad928SBarry Smith PetscFunctionBegin; 79c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 80f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 814b9ad928SBarry Smith PetscFunctionReturn(0); 824b9ad928SBarry Smith } 834b9ad928SBarry Smith 844b9ad928SBarry Smith /*@C 8597177400SBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual 864b9ad928SBarry Smith on the lth level. 874b9ad928SBarry Smith 88d083f849SBarry Smith Logically Collective on PC 894b9ad928SBarry Smith 904b9ad928SBarry Smith Input Parameters: 914b9ad928SBarry Smith + pc - the multigrid context 924b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 93157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 94d0e4de75SBarry Smith previous one were provided then a default is used 954b9ad928SBarry Smith - mat - matrix associated with residual 964b9ad928SBarry Smith 974b9ad928SBarry Smith Level: advanced 984b9ad928SBarry Smith 9954b2cd4bSJed Brown .seealso: PCMGResidualDefault() 1004b9ad928SBarry Smith @*/ 1017087cfbeSBarry Smith PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 1024b9ad928SBarry Smith { 103f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 104f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 105298cc208SBarry Smith PetscErrorCode ierr; 1064b9ad928SBarry Smith 1074b9ad928SBarry Smith PetscFunctionBegin; 108c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 109ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1102fa5cd67SKarl Rupp if (residual) mglevels[l]->residual = residual; 11154b2cd4bSJed Brown if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault; 112f3ae41bdSBarry Smith if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 113298cc208SBarry Smith ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 114f3fbd535SBarry Smith mglevels[l]->A = mat; 1154b9ad928SBarry Smith PetscFunctionReturn(0); 1164b9ad928SBarry Smith } 1174b9ad928SBarry Smith 118*fcb023d4SJed Brown /*@C 119*fcb023d4SJed Brown PCMGSetResidualTranspose - Sets the function to be used to calculate the residual of the transposed linear system 120*fcb023d4SJed Brown on the lth level. 121*fcb023d4SJed Brown 122*fcb023d4SJed Brown Logically Collective on PC 123*fcb023d4SJed Brown 124*fcb023d4SJed Brown Input Parameters: 125*fcb023d4SJed Brown + pc - the multigrid context 126*fcb023d4SJed Brown . l - the level (0 is coarsest) to supply 127*fcb023d4SJed Brown . residualt - function used to form transpose of residual, if none is provided the previously provide one is used, if no 128*fcb023d4SJed Brown previous one were provided then a default is used 129*fcb023d4SJed Brown - mat - matrix associated with residual 130*fcb023d4SJed Brown 131*fcb023d4SJed Brown Level: advanced 132*fcb023d4SJed Brown 133*fcb023d4SJed Brown .seealso: PCMGResidualTransposeDefault() 134*fcb023d4SJed Brown @*/ 135*fcb023d4SJed Brown PetscErrorCode PCMGSetResidualTranspose(PC pc,PetscInt l,PetscErrorCode (*residualt)(Mat,Vec,Vec,Vec),Mat mat) 136*fcb023d4SJed Brown { 137*fcb023d4SJed Brown PC_MG *mg = (PC_MG*)pc->data; 138*fcb023d4SJed Brown PC_MG_Levels **mglevels = mg->levels; 139*fcb023d4SJed Brown PetscErrorCode ierr; 140*fcb023d4SJed Brown 141*fcb023d4SJed Brown PetscFunctionBegin; 142*fcb023d4SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 143*fcb023d4SJed Brown if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 144*fcb023d4SJed Brown if (residualt) mglevels[l]->residualtranspose = residualt; 145*fcb023d4SJed Brown if (!mglevels[l]->residualtranspose) mglevels[l]->residualtranspose = PCMGResidualTransposeDefault; 146*fcb023d4SJed Brown if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 147*fcb023d4SJed Brown ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 148*fcb023d4SJed Brown mglevels[l]->A = mat; 149*fcb023d4SJed Brown PetscFunctionReturn(0); 150*fcb023d4SJed Brown } 151*fcb023d4SJed Brown 1524b9ad928SBarry Smith /*@ 153aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 154bf5b2e24SBarry Smith interpolation from l-1 to the lth level 1554b9ad928SBarry Smith 156d083f849SBarry Smith Logically Collective on PC 1574b9ad928SBarry Smith 1584b9ad928SBarry Smith Input Parameters: 1594b9ad928SBarry Smith + pc - the multigrid context 1604b9ad928SBarry Smith . mat - the interpolation operator 161bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 1624b9ad928SBarry Smith 1634b9ad928SBarry Smith Level: advanced 1644b9ad928SBarry Smith 1654b9ad928SBarry Smith Notes: 1664b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 1674b9ad928SBarry Smith for the same level. 1684b9ad928SBarry Smith 1694b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1704b9ad928SBarry Smith out from the matrix size which one it is. 1714b9ad928SBarry Smith 17297177400SBarry Smith .seealso: PCMGSetRestriction() 1734b9ad928SBarry Smith @*/ 1747087cfbeSBarry Smith PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 1754b9ad928SBarry Smith { 176f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 177f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 178fccaa45eSBarry Smith PetscErrorCode ierr; 1794b9ad928SBarry Smith 1804b9ad928SBarry Smith PetscFunctionBegin; 181c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 182ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 183ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 184c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1856bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr); 1862fa5cd67SKarl Rupp 187f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 1884b9ad928SBarry Smith PetscFunctionReturn(0); 1894b9ad928SBarry Smith } 1904b9ad928SBarry Smith 1918a2c336bSFande Kong /*@ 19295750439SFande Kong PCMGSetOperators - Sets operator and preconditioning matrix for lth level 1938a2c336bSFande Kong 194d083f849SBarry Smith Logically Collective on PC 1958a2c336bSFande Kong 1968a2c336bSFande Kong Input Parameters: 1978a2c336bSFande Kong + pc - the multigrid context 1988a2c336bSFande Kong . Amat - the operator 1998a2c336bSFande Kong . pmat - the preconditioning operator 20095750439SFande Kong - l - the level (0 is the coarsest) to supply 2018a2c336bSFande Kong 2028a2c336bSFande Kong Level: advanced 2038a2c336bSFande Kong 2048a2c336bSFande Kong .keywords: multigrid, set, interpolate, level 2058a2c336bSFande Kong 2068a2c336bSFande Kong .seealso: PCMGSetRestriction(), PCMGSetInterpolation() 2078a2c336bSFande Kong @*/ 2088a2c336bSFande Kong PetscErrorCode PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat) 209360ee056SFande Kong { 210360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 211360ee056SFande Kong PC_MG_Levels **mglevels = mg->levels; 212360ee056SFande Kong PetscErrorCode ierr; 213360ee056SFande Kong 214360ee056SFande Kong PetscFunctionBegin; 215360ee056SFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2168a2c336bSFande Kong PetscValidHeaderSpecific(Amat,MAT_CLASSID,3); 2178a2c336bSFande Kong PetscValidHeaderSpecific(Pmat,MAT_CLASSID,4); 218360ee056SFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 2198a2c336bSFande Kong ierr = KSPSetOperators(mglevels[l]->smoothd,Amat,Pmat);CHKERRQ(ierr); 220360ee056SFande Kong PetscFunctionReturn(0); 221360ee056SFande Kong } 222360ee056SFande Kong 223c88c5224SJed Brown /*@ 224c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 225c88c5224SJed Brown interpolation from l-1 to the lth level 226c88c5224SJed Brown 227c88c5224SJed Brown Logically Collective on PC 228c88c5224SJed Brown 229c88c5224SJed Brown Input Parameters: 230c88c5224SJed Brown + pc - the multigrid context 231c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 232c88c5224SJed Brown 233c88c5224SJed Brown Output Parameter: 2343ad4599aSBarry Smith . mat - the interpolation matrix, can be NULL 235c88c5224SJed Brown 236c88c5224SJed Brown Level: advanced 237c88c5224SJed Brown 238c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 239c88c5224SJed Brown @*/ 240c88c5224SJed Brown PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 241c88c5224SJed Brown { 242c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 243c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 244c88c5224SJed Brown PetscErrorCode ierr; 245c88c5224SJed Brown 246c88c5224SJed Brown PetscFunctionBegin; 247c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2483ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 249ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 250ce94432eSBarry Smith if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 251c88c5224SJed Brown if (!mglevels[l]->interpolate) { 2525aa31b60SBarry Smith if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()"); 253c88c5224SJed Brown ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 254c88c5224SJed Brown } 2553ad4599aSBarry Smith if (mat) *mat = mglevels[l]->interpolate; 256c88c5224SJed Brown PetscFunctionReturn(0); 257c88c5224SJed Brown } 258c88c5224SJed Brown 2598a2c336bSFande Kong /*@ 260eab52d2bSLawrence Mitchell PCMGSetRestriction - Sets the function to be used to restrict dual vectors 2614b9ad928SBarry Smith from level l to l-1. 2624b9ad928SBarry Smith 263d083f849SBarry Smith Logically Collective on PC 2644b9ad928SBarry Smith 2654b9ad928SBarry Smith Input Parameters: 2664b9ad928SBarry Smith + pc - the multigrid context 267c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 268c88c5224SJed Brown - mat - the restriction matrix 2694b9ad928SBarry Smith 2704b9ad928SBarry Smith Level: advanced 2714b9ad928SBarry Smith 2724b9ad928SBarry Smith Notes: 2734b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 2744b9ad928SBarry Smith for the same level. 2754b9ad928SBarry Smith 2764b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 2774b9ad928SBarry Smith out from the matrix size which one it is. 2784b9ad928SBarry Smith 279aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 280fccaa45eSBarry Smith is used. 281fccaa45eSBarry Smith 282a4355517SPierre Jolivet .seealso: PCMGSetInterpolation() 2834b9ad928SBarry Smith @*/ 2847087cfbeSBarry Smith PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 2854b9ad928SBarry Smith { 286fccaa45eSBarry Smith PetscErrorCode ierr; 287f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 288f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2894b9ad928SBarry Smith 2904b9ad928SBarry Smith PetscFunctionBegin; 291c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 292c88c5224SJed Brown PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 293ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 294ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 295c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 2966bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 2972fa5cd67SKarl Rupp 298f3fbd535SBarry Smith mglevels[l]->restrct = mat; 2994b9ad928SBarry Smith PetscFunctionReturn(0); 3004b9ad928SBarry Smith } 3014b9ad928SBarry Smith 302c88c5224SJed Brown /*@ 303eab52d2bSLawrence Mitchell PCMGGetRestriction - Gets the function to be used to restrict dual vectors 304c88c5224SJed Brown from level l to l-1. 305c88c5224SJed Brown 306d083f849SBarry Smith Logically Collective on PC 307c88c5224SJed Brown 308c88c5224SJed Brown Input Parameters: 309c88c5224SJed Brown + pc - the multigrid context 310c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 311c88c5224SJed Brown 312c88c5224SJed Brown Output Parameter: 313c88c5224SJed Brown . mat - the restriction matrix 314c88c5224SJed Brown 315c88c5224SJed Brown Level: advanced 316c88c5224SJed Brown 317eab52d2bSLawrence Mitchell .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection() 318c88c5224SJed Brown @*/ 319c88c5224SJed Brown PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 320c88c5224SJed Brown { 321c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 322c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 323c88c5224SJed Brown PetscErrorCode ierr; 324c88c5224SJed Brown 325c88c5224SJed Brown PetscFunctionBegin; 326c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3273ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 328ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 329ce94432eSBarry Smith if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 330c88c5224SJed Brown if (!mglevels[l]->restrct) { 331ce94432eSBarry Smith if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 332c88c5224SJed Brown ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 333c88c5224SJed Brown } 3343ad4599aSBarry Smith if (mat) *mat = mglevels[l]->restrct; 335c88c5224SJed Brown PetscFunctionReturn(0); 336c88c5224SJed Brown } 337c88c5224SJed Brown 33873250ac0SBarry Smith /*@ 33973250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 34073250ac0SBarry Smith 341d083f849SBarry Smith Logically Collective on PC 342c88c5224SJed Brown 343c88c5224SJed Brown Input Parameters: 344c88c5224SJed Brown + pc - the multigrid context 345c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 346c88c5224SJed Brown . rscale - the scaling 347c88c5224SJed Brown 348c88c5224SJed Brown Level: advanced 349c88c5224SJed Brown 350c88c5224SJed Brown Notes: 351eab52d2bSLawrence Mitchell 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. It is preferable to use PCMGSetInjection() to control moving primal vectors. 352c88c5224SJed Brown 353eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection() 354c88c5224SJed Brown @*/ 355c88c5224SJed Brown PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 356c88c5224SJed Brown { 357c88c5224SJed Brown PetscErrorCode ierr; 358c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 359c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 360c88c5224SJed Brown 361c88c5224SJed Brown PetscFunctionBegin; 362c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 363ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 364ce94432eSBarry Smith if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 365c88c5224SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 366c88c5224SJed Brown ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 3672fa5cd67SKarl Rupp 368c88c5224SJed Brown mglevels[l]->rscale = rscale; 369c88c5224SJed Brown PetscFunctionReturn(0); 370c88c5224SJed Brown } 371c88c5224SJed Brown 372c88c5224SJed Brown /*@ 373c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 374c88c5224SJed Brown 375c88c5224SJed Brown Collective on PC 37673250ac0SBarry Smith 37773250ac0SBarry Smith Input Parameters: 37873250ac0SBarry Smith + pc - the multigrid context 37973250ac0SBarry Smith . rscale - the scaling 38073250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 38173250ac0SBarry Smith 38273250ac0SBarry Smith Level: advanced 38373250ac0SBarry Smith 38473250ac0SBarry Smith Notes: 385eab52d2bSLawrence Mitchell 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. It is preferable to use PCMGGetInjection() to control moving primal vectors. 38673250ac0SBarry Smith 387eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection() 38873250ac0SBarry Smith @*/ 389c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 39073250ac0SBarry Smith { 39173250ac0SBarry Smith PetscErrorCode ierr; 39273250ac0SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 39373250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 39473250ac0SBarry Smith 39573250ac0SBarry Smith PetscFunctionBegin; 39673250ac0SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 397ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 398ce94432eSBarry Smith if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 399c88c5224SJed Brown if (!mglevels[l]->rscale) { 400c88c5224SJed Brown Mat R; 401c88c5224SJed Brown Vec X,Y,coarse,fine; 402c88c5224SJed Brown PetscInt M,N; 403c88c5224SJed Brown ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 4042a7a6963SBarry Smith ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr); 405c88c5224SJed Brown ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 4062fa5cd67SKarl Rupp if (M < N) { 4072fa5cd67SKarl Rupp fine = X; 4082fa5cd67SKarl Rupp coarse = Y; 4092fa5cd67SKarl Rupp } else if (N < M) { 4102fa5cd67SKarl Rupp fine = Y; coarse = X; 411ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 412c88c5224SJed Brown ierr = VecSet(fine,1.);CHKERRQ(ierr); 413c88c5224SJed Brown ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 414c88c5224SJed Brown ierr = VecDestroy(&fine);CHKERRQ(ierr); 415c88c5224SJed Brown ierr = VecReciprocal(coarse);CHKERRQ(ierr); 416c88c5224SJed Brown mglevels[l]->rscale = coarse; 417c88c5224SJed Brown } 418c88c5224SJed Brown *rscale = mglevels[l]->rscale; 41973250ac0SBarry Smith PetscFunctionReturn(0); 42073250ac0SBarry Smith } 42173250ac0SBarry Smith 422f39d8e23SSatish Balay /*@ 423eab52d2bSLawrence Mitchell PCMGSetInjection - Sets the function to be used to inject primal vectors 424eab52d2bSLawrence Mitchell from level l to l-1. 425eab52d2bSLawrence Mitchell 426d083f849SBarry Smith Logically Collective on PC 427eab52d2bSLawrence Mitchell 428eab52d2bSLawrence Mitchell Input Parameters: 429eab52d2bSLawrence Mitchell + pc - the multigrid context 430eab52d2bSLawrence Mitchell . l - the level (0 is coarsest) to supply [Do not supply 0] 431eab52d2bSLawrence Mitchell - mat - the injection matrix 432eab52d2bSLawrence Mitchell 433eab52d2bSLawrence Mitchell Level: advanced 434eab52d2bSLawrence Mitchell 435eab52d2bSLawrence Mitchell .seealso: PCMGSetRestriction() 436eab52d2bSLawrence Mitchell @*/ 437eab52d2bSLawrence Mitchell PetscErrorCode PCMGSetInjection(PC pc,PetscInt l,Mat mat) 438eab52d2bSLawrence Mitchell { 439eab52d2bSLawrence Mitchell PetscErrorCode ierr; 440eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 441eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 442eab52d2bSLawrence Mitchell 443eab52d2bSLawrence Mitchell PetscFunctionBegin; 444eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 445eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 446eab52d2bSLawrence Mitchell if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 447eab52d2bSLawrence Mitchell if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 448eab52d2bSLawrence Mitchell ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 449eab52d2bSLawrence Mitchell ierr = MatDestroy(&mglevels[l]->inject);CHKERRQ(ierr); 450eab52d2bSLawrence Mitchell 451eab52d2bSLawrence Mitchell mglevels[l]->inject = mat; 452eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 453eab52d2bSLawrence Mitchell } 454eab52d2bSLawrence Mitchell 455eab52d2bSLawrence Mitchell /*@ 456eab52d2bSLawrence Mitchell PCMGGetInjection - Gets the function to be used to inject primal vectors 457eab52d2bSLawrence Mitchell from level l to l-1. 458eab52d2bSLawrence Mitchell 459d083f849SBarry Smith Logically Collective on PC 460eab52d2bSLawrence Mitchell 461eab52d2bSLawrence Mitchell Input Parameters: 462eab52d2bSLawrence Mitchell + pc - the multigrid context 463eab52d2bSLawrence Mitchell - l - the level (0 is coarsest) to supply [Do not supply 0] 464eab52d2bSLawrence Mitchell 465eab52d2bSLawrence Mitchell Output Parameter: 46699a38656SLawrence Mitchell . mat - the restriction matrix (may be NULL if no injection is available). 467eab52d2bSLawrence Mitchell 468eab52d2bSLawrence Mitchell Level: advanced 469eab52d2bSLawrence Mitchell 470eab52d2bSLawrence Mitchell .seealso: PCMGSetInjection(), PCMGetGetRestriction() 471eab52d2bSLawrence Mitchell @*/ 472eab52d2bSLawrence Mitchell PetscErrorCode PCMGGetInjection(PC pc,PetscInt l,Mat *mat) 473eab52d2bSLawrence Mitchell { 474eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 475eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 476eab52d2bSLawrence Mitchell 477eab52d2bSLawrence Mitchell PetscFunctionBegin; 478eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 479eab52d2bSLawrence Mitchell if (mat) PetscValidPointer(mat,3); 480eab52d2bSLawrence Mitchell if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 481eab52d2bSLawrence Mitchell if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 482eab52d2bSLawrence Mitchell if (mat) *mat = mglevels[l]->inject; 483eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 484eab52d2bSLawrence Mitchell } 485eab52d2bSLawrence Mitchell 486eab52d2bSLawrence Mitchell /*@ 48797177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 48897177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 48997177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 4904b9ad928SBarry Smith post-smoothing. 4914b9ad928SBarry Smith 4924b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4934b9ad928SBarry Smith 4944b9ad928SBarry Smith Input Parameters: 4954b9ad928SBarry Smith + pc - the multigrid context 4964b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4974b9ad928SBarry Smith 4984b9ad928SBarry Smith Ouput Parameters: 4994b9ad928SBarry Smith . ksp - the smoother 5004b9ad928SBarry Smith 50157420d5bSBarry Smith Notes: 50257420d5bSBarry Smith Once you have called this routine, you can call KSPSetOperators(ksp,...) on the resulting ksp to provide the operators for the smoother for this level. 50357420d5bSBarry Smith You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc) 50457420d5bSBarry Smith and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing. 50557420d5bSBarry Smith 5064b9ad928SBarry Smith Level: advanced 5074b9ad928SBarry Smith 50828668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve() 5094b9ad928SBarry Smith @*/ 5107087cfbeSBarry Smith PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 5114b9ad928SBarry Smith { 512f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 513f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5144b9ad928SBarry Smith 5154b9ad928SBarry Smith PetscFunctionBegin; 516c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 517f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5184b9ad928SBarry Smith PetscFunctionReturn(0); 5194b9ad928SBarry Smith } 5204b9ad928SBarry Smith 521f39d8e23SSatish Balay /*@ 52297177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 5234b9ad928SBarry Smith coarse grid correction (post-smoother). 5244b9ad928SBarry Smith 5254b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 5264b9ad928SBarry Smith 5274b9ad928SBarry Smith Input Parameters: 5284b9ad928SBarry Smith + pc - the multigrid context 5294b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5304b9ad928SBarry Smith 5314b9ad928SBarry Smith Ouput Parameters: 5324b9ad928SBarry Smith . ksp - the smoother 5334b9ad928SBarry Smith 5344b9ad928SBarry Smith Level: advanced 5354b9ad928SBarry Smith 53695452b02SPatrick Sanan Notes: 53795452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 53889cce641SBarry Smith set options on the pre smoother also 53989cce641SBarry Smith 54097177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 5414b9ad928SBarry Smith @*/ 5427087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 5434b9ad928SBarry Smith { 544f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 545f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 546dfbe8321SBarry Smith PetscErrorCode ierr; 547f69a0ea3SMatthew Knepley const char *prefix; 5484b9ad928SBarry Smith MPI_Comm comm; 5494b9ad928SBarry Smith 5504b9ad928SBarry Smith PetscFunctionBegin; 551c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5524b9ad928SBarry Smith /* 5534b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 5544b9ad928SBarry Smith Thus we check if a different one has already been allocated, 5554b9ad928SBarry Smith if not we allocate it. 5564b9ad928SBarry Smith */ 557ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 558f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 55919fd82e9SBarry Smith KSPType ksptype; 56019fd82e9SBarry Smith PCType pctype; 561336babb1SJed Brown PC ipc; 562336babb1SJed Brown PetscReal rtol,abstol,dtol; 563336babb1SJed Brown PetscInt maxits; 564336babb1SJed Brown KSPNormType normtype; 565f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 566f3fbd535SBarry Smith ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 567336babb1SJed Brown ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 5683bf036e2SBarry Smith ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 569336babb1SJed Brown ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 570336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 571336babb1SJed Brown ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 572336babb1SJed Brown 573f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 574422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr); 575f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 576f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 577336babb1SJed Brown ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 578336babb1SJed Brown ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 579336babb1SJed Brown ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 5800059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 581336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 582336babb1SJed Brown ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 5833bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr); 584ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr); 5854b9ad928SBarry Smith } 586f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 5874b9ad928SBarry Smith PetscFunctionReturn(0); 5884b9ad928SBarry Smith } 5894b9ad928SBarry Smith 590f39d8e23SSatish Balay /*@ 59197177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 5924b9ad928SBarry Smith coarse grid correction (pre-smoother). 5934b9ad928SBarry Smith 5944b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 5954b9ad928SBarry Smith 5964b9ad928SBarry Smith Input Parameters: 5974b9ad928SBarry Smith + pc - the multigrid context 5984b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5994b9ad928SBarry Smith 6004b9ad928SBarry Smith Ouput Parameters: 6014b9ad928SBarry Smith . ksp - the smoother 6024b9ad928SBarry Smith 6034b9ad928SBarry Smith Level: advanced 6044b9ad928SBarry Smith 60595452b02SPatrick Sanan Notes: 60695452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 60789cce641SBarry Smith set options on the post smoother also 60889cce641SBarry Smith 60997177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 6104b9ad928SBarry Smith @*/ 6117087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 6124b9ad928SBarry Smith { 613dfbe8321SBarry Smith PetscErrorCode ierr; 614f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 615f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6164b9ad928SBarry Smith 6174b9ad928SBarry Smith PetscFunctionBegin; 618c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 6194b9ad928SBarry Smith /* make sure smoother up and down are different */ 620c5eb9154SBarry Smith if (l) { 6210298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr); 622d8148a5aSMatthew Knepley } 623f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 6244b9ad928SBarry Smith PetscFunctionReturn(0); 6254b9ad928SBarry Smith } 6264b9ad928SBarry Smith 6274b9ad928SBarry Smith /*@ 628cab31ae5SJed Brown PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 6294b9ad928SBarry Smith 630ad4df100SBarry Smith Logically Collective on PC 6314b9ad928SBarry Smith 6324b9ad928SBarry Smith Input Parameters: 6334b9ad928SBarry Smith + pc - the multigrid context 634c1cbb1deSBarry Smith . l - the level (0 is coarsest) 635c1cbb1deSBarry Smith - c - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 6364b9ad928SBarry Smith 6374b9ad928SBarry Smith Level: advanced 6384b9ad928SBarry Smith 639c1cbb1deSBarry Smith .seealso: PCMGSetCycleType() 6404b9ad928SBarry Smith @*/ 641c1cbb1deSBarry Smith PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c) 6424b9ad928SBarry Smith { 643f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 644f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6454b9ad928SBarry Smith 6464b9ad928SBarry Smith PetscFunctionBegin; 647c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 648ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 649c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,l,2); 650c51679f6SSatish Balay PetscValidLogicalCollectiveEnum(pc,c,3); 651f3fbd535SBarry Smith mglevels[l]->cycles = c; 6524b9ad928SBarry Smith PetscFunctionReturn(0); 6534b9ad928SBarry Smith } 6544b9ad928SBarry Smith 6554b9ad928SBarry Smith /*@ 656f3b08a26SMatthew G. Knepley PCMGSetRhs - Sets the vector to be used to store the right-hand side on a particular level. 6574b9ad928SBarry Smith 658d083f849SBarry Smith Logically Collective on PC 6594b9ad928SBarry Smith 6604b9ad928SBarry Smith Input Parameters: 6614b9ad928SBarry Smith + pc - the multigrid context 6624b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 663f3b08a26SMatthew G. Knepley - c - the Vec 6644b9ad928SBarry Smith 6654b9ad928SBarry Smith Level: advanced 6664b9ad928SBarry Smith 66795452b02SPatrick Sanan Notes: 668f3b08a26SMatthew G. Knepley If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it. 669fccaa45eSBarry Smith 670f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, right-hand-side, rhs, level 67197177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 6724b9ad928SBarry Smith @*/ 6737087cfbeSBarry Smith PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 6744b9ad928SBarry Smith { 675fccaa45eSBarry Smith PetscErrorCode ierr; 676f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 677f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6784b9ad928SBarry Smith 6794b9ad928SBarry Smith PetscFunctionBegin; 680c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 681ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 682ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 683c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6846bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 6852fa5cd67SKarl Rupp 686f3fbd535SBarry Smith mglevels[l]->b = c; 6874b9ad928SBarry Smith PetscFunctionReturn(0); 6884b9ad928SBarry Smith } 6894b9ad928SBarry Smith 6904b9ad928SBarry Smith /*@ 691f3b08a26SMatthew G. Knepley PCMGSetX - Sets the vector to be used to store the solution on a particular level. 6924b9ad928SBarry Smith 693d083f849SBarry Smith Logically Collective on PC 6944b9ad928SBarry Smith 6954b9ad928SBarry Smith Input Parameters: 6964b9ad928SBarry Smith + pc - the multigrid context 697251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 698f3b08a26SMatthew G. Knepley - c - the Vec 6994b9ad928SBarry Smith 7004b9ad928SBarry Smith Level: advanced 7014b9ad928SBarry Smith 70295452b02SPatrick Sanan Notes: 703f3b08a26SMatthew G. Knepley If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it. 704fccaa45eSBarry Smith 705f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, solution, level 70697177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 7074b9ad928SBarry Smith @*/ 7087087cfbeSBarry Smith PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 7094b9ad928SBarry Smith { 710fccaa45eSBarry Smith PetscErrorCode ierr; 711f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 712f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7134b9ad928SBarry Smith 7144b9ad928SBarry Smith PetscFunctionBegin; 715c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 716ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 717ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 718c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 7196bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 7202fa5cd67SKarl Rupp 721f3fbd535SBarry Smith mglevels[l]->x = c; 7224b9ad928SBarry Smith PetscFunctionReturn(0); 7234b9ad928SBarry Smith } 7244b9ad928SBarry Smith 7254b9ad928SBarry Smith /*@ 726f3b08a26SMatthew G. Knepley PCMGSetR - Sets the vector to be used to store the residual on a particular level. 7274b9ad928SBarry Smith 728d083f849SBarry Smith Logically Collective on PC 7294b9ad928SBarry Smith 7304b9ad928SBarry Smith Input Parameters: 7314b9ad928SBarry Smith + pc - the multigrid context 7324b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 733f3b08a26SMatthew G. Knepley - c - the Vec 7344b9ad928SBarry Smith 7354b9ad928SBarry Smith Level: advanced 7364b9ad928SBarry Smith 73795452b02SPatrick Sanan Notes: 738f3b08a26SMatthew G. Knepley If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it. 739fccaa45eSBarry Smith 740f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, residual, level 741f3b08a26SMatthew G. Knepley .seealso: PCMGSetRhs(), PCMGSetX() 7424b9ad928SBarry Smith @*/ 7437087cfbeSBarry Smith PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 7444b9ad928SBarry Smith { 745fccaa45eSBarry Smith PetscErrorCode ierr; 746f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 747f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7484b9ad928SBarry Smith 7494b9ad928SBarry Smith PetscFunctionBegin; 750c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 751ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 752ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 753c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 7546bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 7552fa5cd67SKarl Rupp 756f3fbd535SBarry Smith mglevels[l]->r = c; 7574b9ad928SBarry Smith PetscFunctionReturn(0); 7584b9ad928SBarry Smith } 759