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 84b9ad928SBarry Smith Collective on Mat and Vec 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 31f39d8e23SSatish Balay /*@ 3297177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 334b9ad928SBarry Smith 344b9ad928SBarry Smith Not Collective 354b9ad928SBarry Smith 364b9ad928SBarry Smith Input Parameter: 374b9ad928SBarry Smith . pc - the multigrid context 384b9ad928SBarry Smith 394b9ad928SBarry Smith Output Parameter: 404b9ad928SBarry Smith . ksp - the coarse grid solver context 414b9ad928SBarry Smith 424b9ad928SBarry Smith Level: advanced 434b9ad928SBarry Smith 4428668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetSmoother() 454b9ad928SBarry Smith @*/ 467087cfbeSBarry Smith PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 474b9ad928SBarry Smith { 48f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 49f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 504b9ad928SBarry Smith 514b9ad928SBarry Smith PetscFunctionBegin; 52c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 53f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 544b9ad928SBarry Smith PetscFunctionReturn(0); 554b9ad928SBarry Smith } 564b9ad928SBarry Smith 574b9ad928SBarry Smith /*@C 5897177400SBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual 594b9ad928SBarry Smith on the lth level. 604b9ad928SBarry Smith 61ad4df100SBarry Smith Logically Collective on PC and Mat 624b9ad928SBarry Smith 634b9ad928SBarry Smith Input Parameters: 644b9ad928SBarry Smith + pc - the multigrid context 654b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 66157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 67d0e4de75SBarry Smith previous one were provided then a default is used 684b9ad928SBarry Smith - mat - matrix associated with residual 694b9ad928SBarry Smith 704b9ad928SBarry Smith Level: advanced 714b9ad928SBarry Smith 7254b2cd4bSJed Brown .seealso: PCMGResidualDefault() 734b9ad928SBarry Smith @*/ 747087cfbeSBarry Smith PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 754b9ad928SBarry Smith { 76f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 77f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 78298cc208SBarry Smith PetscErrorCode ierr; 794b9ad928SBarry Smith 804b9ad928SBarry Smith PetscFunctionBegin; 81c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 82ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 832fa5cd67SKarl Rupp if (residual) mglevels[l]->residual = residual; 8454b2cd4bSJed Brown if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault; 85f3ae41bdSBarry Smith if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 86298cc208SBarry Smith ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 87f3fbd535SBarry Smith mglevels[l]->A = mat; 884b9ad928SBarry Smith PetscFunctionReturn(0); 894b9ad928SBarry Smith } 904b9ad928SBarry Smith 914b9ad928SBarry Smith /*@ 92aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 93bf5b2e24SBarry Smith interpolation from l-1 to the lth level 944b9ad928SBarry Smith 95ad4df100SBarry Smith Logically Collective on PC and Mat 964b9ad928SBarry Smith 974b9ad928SBarry Smith Input Parameters: 984b9ad928SBarry Smith + pc - the multigrid context 994b9ad928SBarry Smith . mat - the interpolation operator 100bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 1014b9ad928SBarry Smith 1024b9ad928SBarry Smith Level: advanced 1034b9ad928SBarry Smith 1044b9ad928SBarry Smith Notes: 1054b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 1064b9ad928SBarry Smith for the same level. 1074b9ad928SBarry Smith 1084b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1094b9ad928SBarry Smith out from the matrix size which one it is. 1104b9ad928SBarry Smith 11197177400SBarry Smith .seealso: PCMGSetRestriction() 1124b9ad928SBarry Smith @*/ 1137087cfbeSBarry Smith PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 1144b9ad928SBarry Smith { 115f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 116f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 117fccaa45eSBarry Smith PetscErrorCode ierr; 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith PetscFunctionBegin; 120c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 121ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 122ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 123c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1246bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr); 1252fa5cd67SKarl Rupp 126f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 1274b9ad928SBarry Smith PetscFunctionReturn(0); 1284b9ad928SBarry Smith } 1294b9ad928SBarry Smith 130*360ee056SFande Kong PetscErrorCode PCMGSetOperators(PC pc,PetscInt l,Mat mat) 131*360ee056SFande Kong { 132*360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 133*360ee056SFande Kong PC_MG_Levels **mglevels = mg->levels; 134*360ee056SFande Kong PetscErrorCode ierr; 135*360ee056SFande Kong 136*360ee056SFande Kong PetscFunctionBegin; 137*360ee056SFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 138*360ee056SFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 139*360ee056SFande Kong ierr = KSPSetOperators(mglevels[l]->smoothd,mat,mat);CHKERRQ(ierr); 140*360ee056SFande Kong PetscFunctionReturn(0); 141*360ee056SFande Kong } 142*360ee056SFande Kong 143c88c5224SJed Brown /*@ 144c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 145c88c5224SJed Brown interpolation from l-1 to the lth level 146c88c5224SJed Brown 147c88c5224SJed Brown Logically Collective on PC 148c88c5224SJed Brown 149c88c5224SJed Brown Input Parameters: 150c88c5224SJed Brown + pc - the multigrid context 151c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 152c88c5224SJed Brown 153c88c5224SJed Brown Output Parameter: 1543ad4599aSBarry Smith . mat - the interpolation matrix, can be NULL 155c88c5224SJed Brown 156c88c5224SJed Brown Level: advanced 157c88c5224SJed Brown 158c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 159c88c5224SJed Brown @*/ 160c88c5224SJed Brown PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 161c88c5224SJed Brown { 162c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 163c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 164c88c5224SJed Brown PetscErrorCode ierr; 165c88c5224SJed Brown 166c88c5224SJed Brown PetscFunctionBegin; 167c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1683ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 169ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 170ce94432eSBarry 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); 171c88c5224SJed Brown if (!mglevels[l]->interpolate) { 1725aa31b60SBarry Smith if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()"); 173c88c5224SJed Brown ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 174c88c5224SJed Brown } 1753ad4599aSBarry Smith if (mat) *mat = mglevels[l]->interpolate; 176c88c5224SJed Brown PetscFunctionReturn(0); 177c88c5224SJed Brown } 178c88c5224SJed Brown 1794b9ad928SBarry Smith /*@ 180eab52d2bSLawrence Mitchell PCMGSetRestriction - Sets the function to be used to restrict dual vectors 1814b9ad928SBarry Smith from level l to l-1. 1824b9ad928SBarry Smith 183ad4df100SBarry Smith Logically Collective on PC and Mat 1844b9ad928SBarry Smith 1854b9ad928SBarry Smith Input Parameters: 1864b9ad928SBarry Smith + pc - the multigrid context 187c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 188c88c5224SJed Brown - mat - the restriction matrix 1894b9ad928SBarry Smith 1904b9ad928SBarry Smith Level: advanced 1914b9ad928SBarry Smith 1924b9ad928SBarry Smith Notes: 1934b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 1944b9ad928SBarry Smith for the same level. 1954b9ad928SBarry Smith 1964b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1974b9ad928SBarry Smith out from the matrix size which one it is. 1984b9ad928SBarry Smith 199aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 200fccaa45eSBarry Smith is used. 201fccaa45eSBarry Smith 202eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGetSetInjection() 2034b9ad928SBarry Smith @*/ 2047087cfbeSBarry Smith PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 2054b9ad928SBarry Smith { 206fccaa45eSBarry Smith PetscErrorCode ierr; 207f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 208f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2094b9ad928SBarry Smith 2104b9ad928SBarry Smith PetscFunctionBegin; 211c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 212c88c5224SJed Brown PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 213ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 214ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 215c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 2166bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 2172fa5cd67SKarl Rupp 218f3fbd535SBarry Smith mglevels[l]->restrct = mat; 2194b9ad928SBarry Smith PetscFunctionReturn(0); 2204b9ad928SBarry Smith } 2214b9ad928SBarry Smith 222c88c5224SJed Brown /*@ 223eab52d2bSLawrence Mitchell PCMGGetRestriction - Gets the function to be used to restrict dual vectors 224c88c5224SJed Brown from level l to l-1. 225c88c5224SJed Brown 226c88c5224SJed Brown Logically Collective on PC and Mat 227c88c5224SJed Brown 228c88c5224SJed Brown Input Parameters: 229c88c5224SJed Brown + pc - the multigrid context 230c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 231c88c5224SJed Brown 232c88c5224SJed Brown Output Parameter: 233c88c5224SJed Brown . mat - the restriction matrix 234c88c5224SJed Brown 235c88c5224SJed Brown Level: advanced 236c88c5224SJed Brown 237eab52d2bSLawrence Mitchell .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection() 238c88c5224SJed Brown @*/ 239c88c5224SJed Brown PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 240c88c5224SJed Brown { 241c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 242c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 243c88c5224SJed Brown PetscErrorCode ierr; 244c88c5224SJed Brown 245c88c5224SJed Brown PetscFunctionBegin; 246c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2473ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 248ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 249ce94432eSBarry 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); 250c88c5224SJed Brown if (!mglevels[l]->restrct) { 251ce94432eSBarry Smith if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 252c88c5224SJed Brown ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 253c88c5224SJed Brown } 2543ad4599aSBarry Smith if (mat) *mat = mglevels[l]->restrct; 255c88c5224SJed Brown PetscFunctionReturn(0); 256c88c5224SJed Brown } 257c88c5224SJed Brown 25873250ac0SBarry Smith /*@ 25973250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 26073250ac0SBarry Smith 261c88c5224SJed Brown Logically Collective on PC and Vec 262c88c5224SJed Brown 263c88c5224SJed Brown Input Parameters: 264c88c5224SJed Brown + pc - the multigrid context 265c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 266c88c5224SJed Brown . rscale - the scaling 267c88c5224SJed Brown 268c88c5224SJed Brown Level: advanced 269c88c5224SJed Brown 270c88c5224SJed Brown Notes: 271eab52d2bSLawrence 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. 272c88c5224SJed Brown 273eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection() 274c88c5224SJed Brown @*/ 275c88c5224SJed Brown PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 276c88c5224SJed Brown { 277c88c5224SJed Brown PetscErrorCode ierr; 278c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 279c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 280c88c5224SJed Brown 281c88c5224SJed Brown PetscFunctionBegin; 282c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 283ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 284ce94432eSBarry 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); 285c88c5224SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 286c88c5224SJed Brown ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 2872fa5cd67SKarl Rupp 288c88c5224SJed Brown mglevels[l]->rscale = rscale; 289c88c5224SJed Brown PetscFunctionReturn(0); 290c88c5224SJed Brown } 291c88c5224SJed Brown 292c88c5224SJed Brown /*@ 293c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 294c88c5224SJed Brown 295c88c5224SJed Brown Collective on PC 29673250ac0SBarry Smith 29773250ac0SBarry Smith Input Parameters: 29873250ac0SBarry Smith + pc - the multigrid context 29973250ac0SBarry Smith . rscale - the scaling 30073250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 30173250ac0SBarry Smith 30273250ac0SBarry Smith Level: advanced 30373250ac0SBarry Smith 30473250ac0SBarry Smith Notes: 305eab52d2bSLawrence 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. 30673250ac0SBarry Smith 307eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection() 30873250ac0SBarry Smith @*/ 309c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 31073250ac0SBarry Smith { 31173250ac0SBarry Smith PetscErrorCode ierr; 31273250ac0SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 31373250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 31473250ac0SBarry Smith 31573250ac0SBarry Smith PetscFunctionBegin; 31673250ac0SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 317ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 318ce94432eSBarry 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); 319c88c5224SJed Brown if (!mglevels[l]->rscale) { 320c88c5224SJed Brown Mat R; 321c88c5224SJed Brown Vec X,Y,coarse,fine; 322c88c5224SJed Brown PetscInt M,N; 323c88c5224SJed Brown ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 3242a7a6963SBarry Smith ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr); 325c88c5224SJed Brown ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 3262fa5cd67SKarl Rupp if (M < N) { 3272fa5cd67SKarl Rupp fine = X; 3282fa5cd67SKarl Rupp coarse = Y; 3292fa5cd67SKarl Rupp } else if (N < M) { 3302fa5cd67SKarl Rupp fine = Y; coarse = X; 331ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 332c88c5224SJed Brown ierr = VecSet(fine,1.);CHKERRQ(ierr); 333c88c5224SJed Brown ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 334c88c5224SJed Brown ierr = VecDestroy(&fine);CHKERRQ(ierr); 335c88c5224SJed Brown ierr = VecReciprocal(coarse);CHKERRQ(ierr); 336c88c5224SJed Brown mglevels[l]->rscale = coarse; 337c88c5224SJed Brown } 338c88c5224SJed Brown *rscale = mglevels[l]->rscale; 33973250ac0SBarry Smith PetscFunctionReturn(0); 34073250ac0SBarry Smith } 34173250ac0SBarry Smith 342f39d8e23SSatish Balay /*@ 343eab52d2bSLawrence Mitchell PCMGSetInjection - Sets the function to be used to inject primal vectors 344eab52d2bSLawrence Mitchell from level l to l-1. 345eab52d2bSLawrence Mitchell 346eab52d2bSLawrence Mitchell Logically Collective on PC and Mat 347eab52d2bSLawrence Mitchell 348eab52d2bSLawrence Mitchell Input Parameters: 349eab52d2bSLawrence Mitchell + pc - the multigrid context 350eab52d2bSLawrence Mitchell . l - the level (0 is coarsest) to supply [Do not supply 0] 351eab52d2bSLawrence Mitchell - mat - the injection matrix 352eab52d2bSLawrence Mitchell 353eab52d2bSLawrence Mitchell Level: advanced 354eab52d2bSLawrence Mitchell 355eab52d2bSLawrence Mitchell .seealso: PCMGSetRestriction() 356eab52d2bSLawrence Mitchell @*/ 357eab52d2bSLawrence Mitchell PetscErrorCode PCMGSetInjection(PC pc,PetscInt l,Mat mat) 358eab52d2bSLawrence Mitchell { 359eab52d2bSLawrence Mitchell PetscErrorCode ierr; 360eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 361eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 362eab52d2bSLawrence Mitchell 363eab52d2bSLawrence Mitchell PetscFunctionBegin; 364eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 365eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 366eab52d2bSLawrence Mitchell if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 367eab52d2bSLawrence Mitchell if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 368eab52d2bSLawrence Mitchell ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 369eab52d2bSLawrence Mitchell ierr = MatDestroy(&mglevels[l]->inject);CHKERRQ(ierr); 370eab52d2bSLawrence Mitchell 371eab52d2bSLawrence Mitchell mglevels[l]->inject = mat; 372eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 373eab52d2bSLawrence Mitchell } 374eab52d2bSLawrence Mitchell 375eab52d2bSLawrence Mitchell /*@ 376eab52d2bSLawrence Mitchell PCMGGetInjection - Gets the function to be used to inject primal vectors 377eab52d2bSLawrence Mitchell from level l to l-1. 378eab52d2bSLawrence Mitchell 379eab52d2bSLawrence Mitchell Logically Collective on PC and Mat 380eab52d2bSLawrence Mitchell 381eab52d2bSLawrence Mitchell Input Parameters: 382eab52d2bSLawrence Mitchell + pc - the multigrid context 383eab52d2bSLawrence Mitchell - l - the level (0 is coarsest) to supply [Do not supply 0] 384eab52d2bSLawrence Mitchell 385eab52d2bSLawrence Mitchell Output Parameter: 38699a38656SLawrence Mitchell . mat - the restriction matrix (may be NULL if no injection is available). 387eab52d2bSLawrence Mitchell 388eab52d2bSLawrence Mitchell Level: advanced 389eab52d2bSLawrence Mitchell 390eab52d2bSLawrence Mitchell .seealso: PCMGSetInjection(), PCMGetGetRestriction() 391eab52d2bSLawrence Mitchell @*/ 392eab52d2bSLawrence Mitchell PetscErrorCode PCMGGetInjection(PC pc,PetscInt l,Mat *mat) 393eab52d2bSLawrence Mitchell { 394eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 395eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 396eab52d2bSLawrence Mitchell 397eab52d2bSLawrence Mitchell PetscFunctionBegin; 398eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 399eab52d2bSLawrence Mitchell if (mat) PetscValidPointer(mat,3); 400eab52d2bSLawrence Mitchell if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 401eab52d2bSLawrence 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); 402eab52d2bSLawrence Mitchell if (mat) *mat = mglevels[l]->inject; 403eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 404eab52d2bSLawrence Mitchell } 405eab52d2bSLawrence Mitchell 406eab52d2bSLawrence Mitchell /*@ 40797177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 40897177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 40997177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 4104b9ad928SBarry Smith post-smoothing. 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith Input Parameters: 4154b9ad928SBarry Smith + pc - the multigrid context 4164b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4174b9ad928SBarry Smith 4184b9ad928SBarry Smith Ouput Parameters: 4194b9ad928SBarry Smith . ksp - the smoother 4204b9ad928SBarry Smith 42157420d5bSBarry Smith Notes: 42257420d5bSBarry 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. 42357420d5bSBarry Smith You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc) 42457420d5bSBarry Smith and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing. 42557420d5bSBarry Smith 4264b9ad928SBarry Smith Level: advanced 4274b9ad928SBarry Smith 42828668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve() 4294b9ad928SBarry Smith @*/ 4307087cfbeSBarry Smith PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 4314b9ad928SBarry Smith { 432f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 433f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 4344b9ad928SBarry Smith 4354b9ad928SBarry Smith PetscFunctionBegin; 436c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 437f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 4384b9ad928SBarry Smith PetscFunctionReturn(0); 4394b9ad928SBarry Smith } 4404b9ad928SBarry Smith 441f39d8e23SSatish Balay /*@ 44297177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 4434b9ad928SBarry Smith coarse grid correction (post-smoother). 4444b9ad928SBarry Smith 4454b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4464b9ad928SBarry Smith 4474b9ad928SBarry Smith Input Parameters: 4484b9ad928SBarry Smith + pc - the multigrid context 4494b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4504b9ad928SBarry Smith 4514b9ad928SBarry Smith Ouput Parameters: 4524b9ad928SBarry Smith . ksp - the smoother 4534b9ad928SBarry Smith 4544b9ad928SBarry Smith Level: advanced 4554b9ad928SBarry Smith 45695452b02SPatrick Sanan Notes: 45795452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 45889cce641SBarry Smith set options on the pre smoother also 45989cce641SBarry Smith 46097177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 4614b9ad928SBarry Smith @*/ 4627087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 4634b9ad928SBarry Smith { 464f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 465f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 466dfbe8321SBarry Smith PetscErrorCode ierr; 467f69a0ea3SMatthew Knepley const char *prefix; 4684b9ad928SBarry Smith MPI_Comm comm; 4694b9ad928SBarry Smith 4704b9ad928SBarry Smith PetscFunctionBegin; 471c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4724b9ad928SBarry Smith /* 4734b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 4744b9ad928SBarry Smith Thus we check if a different one has already been allocated, 4754b9ad928SBarry Smith if not we allocate it. 4764b9ad928SBarry Smith */ 477ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 478f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 47919fd82e9SBarry Smith KSPType ksptype; 48019fd82e9SBarry Smith PCType pctype; 481336babb1SJed Brown PC ipc; 482336babb1SJed Brown PetscReal rtol,abstol,dtol; 483336babb1SJed Brown PetscInt maxits; 484336babb1SJed Brown KSPNormType normtype; 485f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 486f3fbd535SBarry Smith ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 487336babb1SJed Brown ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 4883bf036e2SBarry Smith ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 489336babb1SJed Brown ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 490336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 491336babb1SJed Brown ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 492336babb1SJed Brown 493f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 494422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr); 495f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 496f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 497336babb1SJed Brown ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 498336babb1SJed Brown ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 499336babb1SJed Brown ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 5000059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 501336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 502336babb1SJed Brown ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 5033bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr); 504ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr); 5054b9ad928SBarry Smith } 506f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 5074b9ad928SBarry Smith PetscFunctionReturn(0); 5084b9ad928SBarry Smith } 5094b9ad928SBarry Smith 510f39d8e23SSatish Balay /*@ 51197177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 5124b9ad928SBarry Smith coarse grid correction (pre-smoother). 5134b9ad928SBarry Smith 5144b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 5154b9ad928SBarry Smith 5164b9ad928SBarry Smith Input Parameters: 5174b9ad928SBarry Smith + pc - the multigrid context 5184b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5194b9ad928SBarry Smith 5204b9ad928SBarry Smith Ouput Parameters: 5214b9ad928SBarry Smith . ksp - the smoother 5224b9ad928SBarry Smith 5234b9ad928SBarry Smith Level: advanced 5244b9ad928SBarry Smith 52595452b02SPatrick Sanan Notes: 52695452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 52789cce641SBarry Smith set options on the post smoother also 52889cce641SBarry Smith 52997177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 5304b9ad928SBarry Smith @*/ 5317087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 5324b9ad928SBarry Smith { 533dfbe8321SBarry Smith PetscErrorCode ierr; 534f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 535f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5364b9ad928SBarry Smith 5374b9ad928SBarry Smith PetscFunctionBegin; 538c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5394b9ad928SBarry Smith /* make sure smoother up and down are different */ 540c5eb9154SBarry Smith if (l) { 5410298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr); 542d8148a5aSMatthew Knepley } 543f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5444b9ad928SBarry Smith PetscFunctionReturn(0); 5454b9ad928SBarry Smith } 5464b9ad928SBarry Smith 5474b9ad928SBarry Smith /*@ 548cab31ae5SJed Brown PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 5494b9ad928SBarry Smith 550ad4df100SBarry Smith Logically Collective on PC 5514b9ad928SBarry Smith 5524b9ad928SBarry Smith Input Parameters: 5534b9ad928SBarry Smith + pc - the multigrid context 554c1cbb1deSBarry Smith . l - the level (0 is coarsest) 555c1cbb1deSBarry Smith - c - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 5564b9ad928SBarry Smith 5574b9ad928SBarry Smith Level: advanced 5584b9ad928SBarry Smith 559c1cbb1deSBarry Smith .seealso: PCMGSetCycleType() 5604b9ad928SBarry Smith @*/ 561c1cbb1deSBarry Smith PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c) 5624b9ad928SBarry Smith { 563f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 564f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5654b9ad928SBarry Smith 5664b9ad928SBarry Smith PetscFunctionBegin; 567c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 568ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 569c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,l,2); 570c51679f6SSatish Balay PetscValidLogicalCollectiveEnum(pc,c,3); 571f3fbd535SBarry Smith mglevels[l]->cycles = c; 5724b9ad928SBarry Smith PetscFunctionReturn(0); 5734b9ad928SBarry Smith } 5744b9ad928SBarry Smith 5754b9ad928SBarry Smith /*@ 57697177400SBarry Smith PCMGSetRhs - Sets the vector space to be used to store the right-hand side 577fccaa45eSBarry Smith on a particular level. 5784b9ad928SBarry Smith 579ad4df100SBarry Smith Logically Collective on PC and Vec 5804b9ad928SBarry Smith 5814b9ad928SBarry Smith Input Parameters: 5824b9ad928SBarry Smith + pc - the multigrid context 5834b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 5844b9ad928SBarry Smith - c - the space 5854b9ad928SBarry Smith 5864b9ad928SBarry Smith Level: advanced 5874b9ad928SBarry Smith 58895452b02SPatrick Sanan Notes: 58995452b02SPatrick Sanan If this is not provided PETSc will automatically generate one. 590fccaa45eSBarry Smith 591fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 592fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 593fccaa45eSBarry Smith 59497177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 5954b9ad928SBarry Smith @*/ 5967087cfbeSBarry Smith PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 5974b9ad928SBarry Smith { 598fccaa45eSBarry Smith PetscErrorCode ierr; 599f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 600f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6014b9ad928SBarry Smith 6024b9ad928SBarry Smith PetscFunctionBegin; 603c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 604ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 605ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 606c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6076bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 6082fa5cd67SKarl Rupp 609f3fbd535SBarry Smith mglevels[l]->b = c; 6104b9ad928SBarry Smith PetscFunctionReturn(0); 6114b9ad928SBarry Smith } 6124b9ad928SBarry Smith 6134b9ad928SBarry Smith /*@ 61497177400SBarry Smith PCMGSetX - Sets the vector space to be used to store the solution on a 615fccaa45eSBarry Smith particular level. 6164b9ad928SBarry Smith 617ad4df100SBarry Smith Logically Collective on PC and Vec 6184b9ad928SBarry Smith 6194b9ad928SBarry Smith Input Parameters: 6204b9ad928SBarry Smith + pc - the multigrid context 621251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 6224b9ad928SBarry Smith - c - the space 6234b9ad928SBarry Smith 6244b9ad928SBarry Smith Level: advanced 6254b9ad928SBarry Smith 62695452b02SPatrick Sanan Notes: 62795452b02SPatrick Sanan If this is not provided PETSc will automatically generate one. 628fccaa45eSBarry Smith 629fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 630fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 631fccaa45eSBarry Smith 63297177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 6334b9ad928SBarry Smith @*/ 6347087cfbeSBarry Smith PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 6354b9ad928SBarry Smith { 636fccaa45eSBarry Smith PetscErrorCode ierr; 637f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 638f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6394b9ad928SBarry Smith 6404b9ad928SBarry Smith PetscFunctionBegin; 641c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 642ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 643ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 644c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6456bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 6462fa5cd67SKarl Rupp 647f3fbd535SBarry Smith mglevels[l]->x = c; 6484b9ad928SBarry Smith PetscFunctionReturn(0); 6494b9ad928SBarry Smith } 6504b9ad928SBarry Smith 6514b9ad928SBarry Smith /*@ 65297177400SBarry Smith PCMGSetR - Sets the vector space to be used to store the residual on a 653fccaa45eSBarry Smith particular level. 6544b9ad928SBarry Smith 655ad4df100SBarry Smith Logically Collective on PC and Vec 6564b9ad928SBarry Smith 6574b9ad928SBarry Smith Input Parameters: 6584b9ad928SBarry Smith + pc - the multigrid context 6594b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 6604b9ad928SBarry Smith - c - the space 6614b9ad928SBarry Smith 6624b9ad928SBarry Smith Level: advanced 6634b9ad928SBarry Smith 66495452b02SPatrick Sanan Notes: 66595452b02SPatrick Sanan If this is not provided PETSc will automatically generate one. 666fccaa45eSBarry Smith 667fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 668fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 669fccaa45eSBarry Smith 6704b9ad928SBarry Smith @*/ 6717087cfbeSBarry Smith PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 6724b9ad928SBarry Smith { 673fccaa45eSBarry Smith PetscErrorCode ierr; 674f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 675f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6764b9ad928SBarry Smith 6774b9ad928SBarry Smith PetscFunctionBegin; 678c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 679ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 680ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 681c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6826bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 6832fa5cd67SKarl Rupp 684f3fbd535SBarry Smith mglevels[l]->r = c; 6854b9ad928SBarry Smith PetscFunctionReturn(0); 6864b9ad928SBarry Smith } 687