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 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 61d083f849SBarry Smith Logically Collective on PC 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 95d083f849SBarry Smith Logically Collective on PC 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 1308a2c336bSFande Kong /*@ 13195750439SFande Kong PCMGSetOperators - Sets operator and preconditioning matrix for lth level 1328a2c336bSFande Kong 133d083f849SBarry Smith Logically Collective on PC 1348a2c336bSFande Kong 1358a2c336bSFande Kong Input Parameters: 1368a2c336bSFande Kong + pc - the multigrid context 1378a2c336bSFande Kong . Amat - the operator 1388a2c336bSFande Kong . pmat - the preconditioning operator 13995750439SFande Kong - l - the level (0 is the coarsest) to supply 1408a2c336bSFande Kong 1418a2c336bSFande Kong Level: advanced 1428a2c336bSFande Kong 1438a2c336bSFande Kong .keywords: multigrid, set, interpolate, level 1448a2c336bSFande Kong 1458a2c336bSFande Kong .seealso: PCMGSetRestriction(), PCMGSetInterpolation() 1468a2c336bSFande Kong @*/ 1478a2c336bSFande Kong PetscErrorCode PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat) 148360ee056SFande Kong { 149360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 150360ee056SFande Kong PC_MG_Levels **mglevels = mg->levels; 151360ee056SFande Kong PetscErrorCode ierr; 152360ee056SFande Kong 153360ee056SFande Kong PetscFunctionBegin; 154360ee056SFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1558a2c336bSFande Kong PetscValidHeaderSpecific(Amat,MAT_CLASSID,3); 1568a2c336bSFande Kong PetscValidHeaderSpecific(Pmat,MAT_CLASSID,4); 157360ee056SFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1588a2c336bSFande Kong ierr = KSPSetOperators(mglevels[l]->smoothd,Amat,Pmat);CHKERRQ(ierr); 159360ee056SFande Kong PetscFunctionReturn(0); 160360ee056SFande Kong } 161360ee056SFande Kong 162c88c5224SJed Brown /*@ 163c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 164c88c5224SJed Brown interpolation from l-1 to the lth level 165c88c5224SJed Brown 166c88c5224SJed Brown Logically Collective on PC 167c88c5224SJed Brown 168c88c5224SJed Brown Input Parameters: 169c88c5224SJed Brown + pc - the multigrid context 170c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 171c88c5224SJed Brown 172c88c5224SJed Brown Output Parameter: 1733ad4599aSBarry Smith . mat - the interpolation matrix, can be NULL 174c88c5224SJed Brown 175c88c5224SJed Brown Level: advanced 176c88c5224SJed Brown 177c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 178c88c5224SJed Brown @*/ 179c88c5224SJed Brown PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 180c88c5224SJed Brown { 181c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 182c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 183c88c5224SJed Brown PetscErrorCode ierr; 184c88c5224SJed Brown 185c88c5224SJed Brown PetscFunctionBegin; 186c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1873ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 188ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 189ce94432eSBarry 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); 190c88c5224SJed Brown if (!mglevels[l]->interpolate) { 1915aa31b60SBarry Smith if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()"); 192c88c5224SJed Brown ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 193c88c5224SJed Brown } 1943ad4599aSBarry Smith if (mat) *mat = mglevels[l]->interpolate; 195c88c5224SJed Brown PetscFunctionReturn(0); 196c88c5224SJed Brown } 197c88c5224SJed Brown 1988a2c336bSFande Kong /*@ 199eab52d2bSLawrence Mitchell PCMGSetRestriction - Sets the function to be used to restrict dual vectors 2004b9ad928SBarry Smith from level l to l-1. 2014b9ad928SBarry Smith 202d083f849SBarry Smith Logically Collective on PC 2034b9ad928SBarry Smith 2044b9ad928SBarry Smith Input Parameters: 2054b9ad928SBarry Smith + pc - the multigrid context 206c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 207c88c5224SJed Brown - mat - the restriction matrix 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith Level: advanced 2104b9ad928SBarry Smith 2114b9ad928SBarry Smith Notes: 2124b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 2134b9ad928SBarry Smith for the same level. 2144b9ad928SBarry Smith 2154b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 2164b9ad928SBarry Smith out from the matrix size which one it is. 2174b9ad928SBarry Smith 218aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 219fccaa45eSBarry Smith is used. 220fccaa45eSBarry Smith 221*a4355517SPierre Jolivet .seealso: PCMGSetInterpolation() 2224b9ad928SBarry Smith @*/ 2237087cfbeSBarry Smith PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 2244b9ad928SBarry Smith { 225fccaa45eSBarry Smith PetscErrorCode ierr; 226f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 227f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2284b9ad928SBarry Smith 2294b9ad928SBarry Smith PetscFunctionBegin; 230c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 231c88c5224SJed Brown PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 232ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 233ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 234c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 2356bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 2362fa5cd67SKarl Rupp 237f3fbd535SBarry Smith mglevels[l]->restrct = mat; 2384b9ad928SBarry Smith PetscFunctionReturn(0); 2394b9ad928SBarry Smith } 2404b9ad928SBarry Smith 241c88c5224SJed Brown /*@ 242eab52d2bSLawrence Mitchell PCMGGetRestriction - Gets the function to be used to restrict dual vectors 243c88c5224SJed Brown from level l to l-1. 244c88c5224SJed Brown 245d083f849SBarry Smith Logically Collective on PC 246c88c5224SJed Brown 247c88c5224SJed Brown Input Parameters: 248c88c5224SJed Brown + pc - the multigrid context 249c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 250c88c5224SJed Brown 251c88c5224SJed Brown Output Parameter: 252c88c5224SJed Brown . mat - the restriction matrix 253c88c5224SJed Brown 254c88c5224SJed Brown Level: advanced 255c88c5224SJed Brown 256eab52d2bSLawrence Mitchell .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection() 257c88c5224SJed Brown @*/ 258c88c5224SJed Brown PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 259c88c5224SJed Brown { 260c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 261c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 262c88c5224SJed Brown PetscErrorCode ierr; 263c88c5224SJed Brown 264c88c5224SJed Brown PetscFunctionBegin; 265c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2663ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 267ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 268ce94432eSBarry 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); 269c88c5224SJed Brown if (!mglevels[l]->restrct) { 270ce94432eSBarry Smith if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 271c88c5224SJed Brown ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 272c88c5224SJed Brown } 2733ad4599aSBarry Smith if (mat) *mat = mglevels[l]->restrct; 274c88c5224SJed Brown PetscFunctionReturn(0); 275c88c5224SJed Brown } 276c88c5224SJed Brown 27773250ac0SBarry Smith /*@ 27873250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 27973250ac0SBarry Smith 280d083f849SBarry Smith Logically Collective on PC 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: 290eab52d2bSLawrence 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. 291c88c5224SJed Brown 292eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection() 293c88c5224SJed Brown @*/ 294c88c5224SJed Brown PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 295c88c5224SJed Brown { 296c88c5224SJed Brown PetscErrorCode ierr; 297c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 298c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 299c88c5224SJed Brown 300c88c5224SJed Brown PetscFunctionBegin; 301c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 302ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 303ce94432eSBarry 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); 304c88c5224SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 305c88c5224SJed Brown ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 3062fa5cd67SKarl Rupp 307c88c5224SJed Brown mglevels[l]->rscale = rscale; 308c88c5224SJed Brown PetscFunctionReturn(0); 309c88c5224SJed Brown } 310c88c5224SJed Brown 311c88c5224SJed Brown /*@ 312c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 313c88c5224SJed Brown 314c88c5224SJed Brown Collective on PC 31573250ac0SBarry Smith 31673250ac0SBarry Smith Input Parameters: 31773250ac0SBarry Smith + pc - the multigrid context 31873250ac0SBarry Smith . rscale - the scaling 31973250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 32073250ac0SBarry Smith 32173250ac0SBarry Smith Level: advanced 32273250ac0SBarry Smith 32373250ac0SBarry Smith Notes: 324eab52d2bSLawrence 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. 32573250ac0SBarry Smith 326eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection() 32773250ac0SBarry Smith @*/ 328c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 32973250ac0SBarry Smith { 33073250ac0SBarry Smith PetscErrorCode ierr; 33173250ac0SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 33273250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 33373250ac0SBarry Smith 33473250ac0SBarry Smith PetscFunctionBegin; 33573250ac0SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 336ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 337ce94432eSBarry 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); 338c88c5224SJed Brown if (!mglevels[l]->rscale) { 339c88c5224SJed Brown Mat R; 340c88c5224SJed Brown Vec X,Y,coarse,fine; 341c88c5224SJed Brown PetscInt M,N; 342c88c5224SJed Brown ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 3432a7a6963SBarry Smith ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr); 344c88c5224SJed Brown ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 3452fa5cd67SKarl Rupp if (M < N) { 3462fa5cd67SKarl Rupp fine = X; 3472fa5cd67SKarl Rupp coarse = Y; 3482fa5cd67SKarl Rupp } else if (N < M) { 3492fa5cd67SKarl Rupp fine = Y; coarse = X; 350ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 351c88c5224SJed Brown ierr = VecSet(fine,1.);CHKERRQ(ierr); 352c88c5224SJed Brown ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 353c88c5224SJed Brown ierr = VecDestroy(&fine);CHKERRQ(ierr); 354c88c5224SJed Brown ierr = VecReciprocal(coarse);CHKERRQ(ierr); 355c88c5224SJed Brown mglevels[l]->rscale = coarse; 356c88c5224SJed Brown } 357c88c5224SJed Brown *rscale = mglevels[l]->rscale; 35873250ac0SBarry Smith PetscFunctionReturn(0); 35973250ac0SBarry Smith } 36073250ac0SBarry Smith 361f39d8e23SSatish Balay /*@ 362eab52d2bSLawrence Mitchell PCMGSetInjection - Sets the function to be used to inject primal vectors 363eab52d2bSLawrence Mitchell from level l to l-1. 364eab52d2bSLawrence Mitchell 365d083f849SBarry Smith Logically Collective on PC 366eab52d2bSLawrence Mitchell 367eab52d2bSLawrence Mitchell Input Parameters: 368eab52d2bSLawrence Mitchell + pc - the multigrid context 369eab52d2bSLawrence Mitchell . l - the level (0 is coarsest) to supply [Do not supply 0] 370eab52d2bSLawrence Mitchell - mat - the injection matrix 371eab52d2bSLawrence Mitchell 372eab52d2bSLawrence Mitchell Level: advanced 373eab52d2bSLawrence Mitchell 374eab52d2bSLawrence Mitchell .seealso: PCMGSetRestriction() 375eab52d2bSLawrence Mitchell @*/ 376eab52d2bSLawrence Mitchell PetscErrorCode PCMGSetInjection(PC pc,PetscInt l,Mat mat) 377eab52d2bSLawrence Mitchell { 378eab52d2bSLawrence Mitchell PetscErrorCode ierr; 379eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 380eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 381eab52d2bSLawrence Mitchell 382eab52d2bSLawrence Mitchell PetscFunctionBegin; 383eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 384eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 385eab52d2bSLawrence Mitchell if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 386eab52d2bSLawrence Mitchell if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 387eab52d2bSLawrence Mitchell ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 388eab52d2bSLawrence Mitchell ierr = MatDestroy(&mglevels[l]->inject);CHKERRQ(ierr); 389eab52d2bSLawrence Mitchell 390eab52d2bSLawrence Mitchell mglevels[l]->inject = mat; 391eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 392eab52d2bSLawrence Mitchell } 393eab52d2bSLawrence Mitchell 394eab52d2bSLawrence Mitchell /*@ 395eab52d2bSLawrence Mitchell PCMGGetInjection - Gets the function to be used to inject primal vectors 396eab52d2bSLawrence Mitchell from level l to l-1. 397eab52d2bSLawrence Mitchell 398d083f849SBarry Smith Logically Collective on PC 399eab52d2bSLawrence Mitchell 400eab52d2bSLawrence Mitchell Input Parameters: 401eab52d2bSLawrence Mitchell + pc - the multigrid context 402eab52d2bSLawrence Mitchell - l - the level (0 is coarsest) to supply [Do not supply 0] 403eab52d2bSLawrence Mitchell 404eab52d2bSLawrence Mitchell Output Parameter: 40599a38656SLawrence Mitchell . mat - the restriction matrix (may be NULL if no injection is available). 406eab52d2bSLawrence Mitchell 407eab52d2bSLawrence Mitchell Level: advanced 408eab52d2bSLawrence Mitchell 409eab52d2bSLawrence Mitchell .seealso: PCMGSetInjection(), PCMGetGetRestriction() 410eab52d2bSLawrence Mitchell @*/ 411eab52d2bSLawrence Mitchell PetscErrorCode PCMGGetInjection(PC pc,PetscInt l,Mat *mat) 412eab52d2bSLawrence Mitchell { 413eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 414eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 415eab52d2bSLawrence Mitchell 416eab52d2bSLawrence Mitchell PetscFunctionBegin; 417eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 418eab52d2bSLawrence Mitchell if (mat) PetscValidPointer(mat,3); 419eab52d2bSLawrence Mitchell if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 420eab52d2bSLawrence 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); 421eab52d2bSLawrence Mitchell if (mat) *mat = mglevels[l]->inject; 422eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 423eab52d2bSLawrence Mitchell } 424eab52d2bSLawrence Mitchell 425eab52d2bSLawrence Mitchell /*@ 42697177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 42797177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 42897177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 4294b9ad928SBarry Smith post-smoothing. 4304b9ad928SBarry Smith 4314b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4324b9ad928SBarry Smith 4334b9ad928SBarry Smith Input Parameters: 4344b9ad928SBarry Smith + pc - the multigrid context 4354b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4364b9ad928SBarry Smith 4374b9ad928SBarry Smith Ouput Parameters: 4384b9ad928SBarry Smith . ksp - the smoother 4394b9ad928SBarry Smith 44057420d5bSBarry Smith Notes: 44157420d5bSBarry 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. 44257420d5bSBarry Smith You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc) 44357420d5bSBarry Smith and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing. 44457420d5bSBarry Smith 4454b9ad928SBarry Smith Level: advanced 4464b9ad928SBarry Smith 44728668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve() 4484b9ad928SBarry Smith @*/ 4497087cfbeSBarry Smith PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 4504b9ad928SBarry Smith { 451f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 452f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 4534b9ad928SBarry Smith 4544b9ad928SBarry Smith PetscFunctionBegin; 455c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 456f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 4574b9ad928SBarry Smith PetscFunctionReturn(0); 4584b9ad928SBarry Smith } 4594b9ad928SBarry Smith 460f39d8e23SSatish Balay /*@ 46197177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 4624b9ad928SBarry Smith coarse grid correction (post-smoother). 4634b9ad928SBarry Smith 4644b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4654b9ad928SBarry Smith 4664b9ad928SBarry Smith Input Parameters: 4674b9ad928SBarry Smith + pc - the multigrid context 4684b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4694b9ad928SBarry Smith 4704b9ad928SBarry Smith Ouput Parameters: 4714b9ad928SBarry Smith . ksp - the smoother 4724b9ad928SBarry Smith 4734b9ad928SBarry Smith Level: advanced 4744b9ad928SBarry Smith 47595452b02SPatrick Sanan Notes: 47695452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 47789cce641SBarry Smith set options on the pre smoother also 47889cce641SBarry Smith 47997177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 4804b9ad928SBarry Smith @*/ 4817087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 4824b9ad928SBarry Smith { 483f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 484f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 485dfbe8321SBarry Smith PetscErrorCode ierr; 486f69a0ea3SMatthew Knepley const char *prefix; 4874b9ad928SBarry Smith MPI_Comm comm; 4884b9ad928SBarry Smith 4894b9ad928SBarry Smith PetscFunctionBegin; 490c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4914b9ad928SBarry Smith /* 4924b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 4934b9ad928SBarry Smith Thus we check if a different one has already been allocated, 4944b9ad928SBarry Smith if not we allocate it. 4954b9ad928SBarry Smith */ 496ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 497f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 49819fd82e9SBarry Smith KSPType ksptype; 49919fd82e9SBarry Smith PCType pctype; 500336babb1SJed Brown PC ipc; 501336babb1SJed Brown PetscReal rtol,abstol,dtol; 502336babb1SJed Brown PetscInt maxits; 503336babb1SJed Brown KSPNormType normtype; 504f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 505f3fbd535SBarry Smith ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 506336babb1SJed Brown ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 5073bf036e2SBarry Smith ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 508336babb1SJed Brown ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 509336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 510336babb1SJed Brown ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 511336babb1SJed Brown 512f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 513422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr); 514f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 515f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 516336babb1SJed Brown ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 517336babb1SJed Brown ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 518336babb1SJed Brown ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 5190059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 520336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 521336babb1SJed Brown ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 5223bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr); 523ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr); 5244b9ad928SBarry Smith } 525f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 5264b9ad928SBarry Smith PetscFunctionReturn(0); 5274b9ad928SBarry Smith } 5284b9ad928SBarry Smith 529f39d8e23SSatish Balay /*@ 53097177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 5314b9ad928SBarry Smith coarse grid correction (pre-smoother). 5324b9ad928SBarry Smith 5334b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 5344b9ad928SBarry Smith 5354b9ad928SBarry Smith Input Parameters: 5364b9ad928SBarry Smith + pc - the multigrid context 5374b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5384b9ad928SBarry Smith 5394b9ad928SBarry Smith Ouput Parameters: 5404b9ad928SBarry Smith . ksp - the smoother 5414b9ad928SBarry Smith 5424b9ad928SBarry Smith Level: advanced 5434b9ad928SBarry Smith 54495452b02SPatrick Sanan Notes: 54595452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 54689cce641SBarry Smith set options on the post smoother also 54789cce641SBarry Smith 54897177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 5494b9ad928SBarry Smith @*/ 5507087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 5514b9ad928SBarry Smith { 552dfbe8321SBarry Smith PetscErrorCode ierr; 553f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 554f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5554b9ad928SBarry Smith 5564b9ad928SBarry Smith PetscFunctionBegin; 557c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5584b9ad928SBarry Smith /* make sure smoother up and down are different */ 559c5eb9154SBarry Smith if (l) { 5600298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr); 561d8148a5aSMatthew Knepley } 562f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5634b9ad928SBarry Smith PetscFunctionReturn(0); 5644b9ad928SBarry Smith } 5654b9ad928SBarry Smith 5664b9ad928SBarry Smith /*@ 567cab31ae5SJed Brown PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 5684b9ad928SBarry Smith 569ad4df100SBarry Smith Logically Collective on PC 5704b9ad928SBarry Smith 5714b9ad928SBarry Smith Input Parameters: 5724b9ad928SBarry Smith + pc - the multigrid context 573c1cbb1deSBarry Smith . l - the level (0 is coarsest) 574c1cbb1deSBarry Smith - c - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 5754b9ad928SBarry Smith 5764b9ad928SBarry Smith Level: advanced 5774b9ad928SBarry Smith 578c1cbb1deSBarry Smith .seealso: PCMGSetCycleType() 5794b9ad928SBarry Smith @*/ 580c1cbb1deSBarry Smith PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c) 5814b9ad928SBarry Smith { 582f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 583f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5844b9ad928SBarry Smith 5854b9ad928SBarry Smith PetscFunctionBegin; 586c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 587ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 588c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,l,2); 589c51679f6SSatish Balay PetscValidLogicalCollectiveEnum(pc,c,3); 590f3fbd535SBarry Smith mglevels[l]->cycles = c; 5914b9ad928SBarry Smith PetscFunctionReturn(0); 5924b9ad928SBarry Smith } 5934b9ad928SBarry Smith 5944b9ad928SBarry Smith /*@ 59597177400SBarry Smith PCMGSetRhs - Sets the vector space to be used to store the right-hand side 596fccaa45eSBarry Smith on a particular level. 5974b9ad928SBarry Smith 598d083f849SBarry Smith Logically Collective on PC 5994b9ad928SBarry Smith 6004b9ad928SBarry Smith Input Parameters: 6014b9ad928SBarry Smith + pc - the multigrid context 6024b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 6034b9ad928SBarry Smith - c - the space 6044b9ad928SBarry Smith 6054b9ad928SBarry Smith Level: advanced 6064b9ad928SBarry Smith 60795452b02SPatrick Sanan Notes: 60895452b02SPatrick Sanan If this is not provided PETSc will automatically generate one. 609fccaa45eSBarry Smith 610fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 611fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 612fccaa45eSBarry Smith 61397177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 6144b9ad928SBarry Smith @*/ 6157087cfbeSBarry Smith PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 6164b9ad928SBarry Smith { 617fccaa45eSBarry Smith PetscErrorCode ierr; 618f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 619f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6204b9ad928SBarry Smith 6214b9ad928SBarry Smith PetscFunctionBegin; 622c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 623ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 624ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 625c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6266bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 6272fa5cd67SKarl Rupp 628f3fbd535SBarry Smith mglevels[l]->b = c; 6294b9ad928SBarry Smith PetscFunctionReturn(0); 6304b9ad928SBarry Smith } 6314b9ad928SBarry Smith 6324b9ad928SBarry Smith /*@ 63397177400SBarry Smith PCMGSetX - Sets the vector space to be used to store the solution on a 634fccaa45eSBarry Smith particular level. 6354b9ad928SBarry Smith 636d083f849SBarry Smith Logically Collective on PC 6374b9ad928SBarry Smith 6384b9ad928SBarry Smith Input Parameters: 6394b9ad928SBarry Smith + pc - the multigrid context 640251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 6414b9ad928SBarry Smith - c - the space 6424b9ad928SBarry Smith 6434b9ad928SBarry Smith Level: advanced 6444b9ad928SBarry Smith 64595452b02SPatrick Sanan Notes: 64695452b02SPatrick Sanan If this is not provided PETSc will automatically generate one. 647fccaa45eSBarry Smith 648fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 649fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 650fccaa45eSBarry Smith 65197177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 6524b9ad928SBarry Smith @*/ 6537087cfbeSBarry Smith PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 6544b9ad928SBarry Smith { 655fccaa45eSBarry Smith PetscErrorCode ierr; 656f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 657f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6584b9ad928SBarry Smith 6594b9ad928SBarry Smith PetscFunctionBegin; 660c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 661ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 662ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 663c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6646bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 6652fa5cd67SKarl Rupp 666f3fbd535SBarry Smith mglevels[l]->x = c; 6674b9ad928SBarry Smith PetscFunctionReturn(0); 6684b9ad928SBarry Smith } 6694b9ad928SBarry Smith 6704b9ad928SBarry Smith /*@ 67197177400SBarry Smith PCMGSetR - Sets the vector space to be used to store the residual on a 672fccaa45eSBarry Smith particular level. 6734b9ad928SBarry Smith 674d083f849SBarry Smith Logically Collective on PC 6754b9ad928SBarry Smith 6764b9ad928SBarry Smith Input Parameters: 6774b9ad928SBarry Smith + pc - the multigrid context 6784b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 6794b9ad928SBarry Smith - c - the space 6804b9ad928SBarry Smith 6814b9ad928SBarry Smith Level: advanced 6824b9ad928SBarry Smith 68395452b02SPatrick Sanan Notes: 68495452b02SPatrick Sanan If this is not provided PETSc will automatically generate one. 685fccaa45eSBarry Smith 686fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 687fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 688fccaa45eSBarry Smith 6894b9ad928SBarry Smith @*/ 6907087cfbeSBarry Smith PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 6914b9ad928SBarry Smith { 692fccaa45eSBarry Smith PetscErrorCode ierr; 693f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 694f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6954b9ad928SBarry Smith 6964b9ad928SBarry Smith PetscFunctionBegin; 697c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 698ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 699ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 700c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 7016bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 7022fa5cd67SKarl Rupp 703f3fbd535SBarry Smith mglevels[l]->r = c; 7044b9ad928SBarry Smith PetscFunctionReturn(0); 7054b9ad928SBarry Smith } 706