14b9ad928SBarry Smith 2c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscksp.h" I*/ 34b9ad928SBarry Smith 44b9ad928SBarry Smith #undef __FUNCT__ 5d0e4de75SBarry Smith #define __FUNCT__ "PCMGResidual_Default" 61f6cc5b2SSatish Balay /*@C 7d0e4de75SBarry Smith PCMGResidual_Default - Default routine to calculate the residual. 84b9ad928SBarry Smith 94b9ad928SBarry Smith Collective on Mat and Vec 104b9ad928SBarry Smith 114b9ad928SBarry Smith Input Parameters: 124b9ad928SBarry Smith + mat - the matrix 134b9ad928SBarry Smith . b - the right-hand-side 144b9ad928SBarry Smith - x - the approximate solution 154b9ad928SBarry Smith 164b9ad928SBarry Smith Output Parameter: 174b9ad928SBarry Smith . r - location to store the residual 184b9ad928SBarry Smith 19d0e4de75SBarry Smith Level: developer 204b9ad928SBarry Smith 214b9ad928SBarry Smith .keywords: MG, default, multigrid, residual 224b9ad928SBarry Smith 2397177400SBarry Smith .seealso: PCMGSetResidual() 244b9ad928SBarry Smith @*/ 25d0e4de75SBarry Smith PetscErrorCode PCMGResidual_Default(Mat mat,Vec b,Vec x,Vec r) 264b9ad928SBarry Smith { 27dfbe8321SBarry Smith PetscErrorCode ierr; 284b9ad928SBarry Smith 294b9ad928SBarry Smith PetscFunctionBegin; 304b9ad928SBarry Smith ierr = MatMult(mat,x,r);CHKERRQ(ierr); 31efb30889SBarry Smith ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 324b9ad928SBarry Smith PetscFunctionReturn(0); 334b9ad928SBarry Smith } 344b9ad928SBarry Smith 354b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/ 364b9ad928SBarry Smith 374b9ad928SBarry Smith #undef __FUNCT__ 38d6337fedSHong Zhang #define __FUNCT__ "PCMGGetCoarseSolve" 39f39d8e23SSatish Balay /*@ 4097177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 414b9ad928SBarry Smith 424b9ad928SBarry Smith Not Collective 434b9ad928SBarry Smith 444b9ad928SBarry Smith Input Parameter: 454b9ad928SBarry Smith . pc - the multigrid context 464b9ad928SBarry Smith 474b9ad928SBarry Smith Output Parameter: 484b9ad928SBarry Smith . ksp - the coarse grid solver context 494b9ad928SBarry Smith 504b9ad928SBarry Smith Level: advanced 514b9ad928SBarry Smith 524b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid 534b9ad928SBarry Smith @*/ 547087cfbeSBarry Smith PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 554b9ad928SBarry Smith { 56f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 57f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 584b9ad928SBarry Smith 594b9ad928SBarry Smith PetscFunctionBegin; 60c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 61f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 624b9ad928SBarry Smith PetscFunctionReturn(0); 634b9ad928SBarry Smith } 644b9ad928SBarry Smith 654b9ad928SBarry Smith #undef __FUNCT__ 66d6337fedSHong Zhang #define __FUNCT__ "PCMGSetResidual" 674b9ad928SBarry Smith /*@C 6897177400SBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual 694b9ad928SBarry Smith on the lth level. 704b9ad928SBarry Smith 71ad4df100SBarry Smith Logically Collective on PC and Mat 724b9ad928SBarry Smith 734b9ad928SBarry Smith Input Parameters: 744b9ad928SBarry Smith + pc - the multigrid context 754b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 76157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 77d0e4de75SBarry Smith previous one were provided then a default is used 784b9ad928SBarry Smith - mat - matrix associated with residual 794b9ad928SBarry Smith 804b9ad928SBarry Smith Level: advanced 814b9ad928SBarry Smith 824b9ad928SBarry Smith .keywords: MG, set, multigrid, residual, level 834b9ad928SBarry Smith 84d0e4de75SBarry Smith .seealso: PCMGResidual_Default() 854b9ad928SBarry Smith @*/ 867087cfbeSBarry Smith PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 874b9ad928SBarry Smith { 88f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 89f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 90298cc208SBarry Smith PetscErrorCode ierr; 914b9ad928SBarry Smith 924b9ad928SBarry Smith PetscFunctionBegin; 93c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 94ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 952fa5cd67SKarl Rupp if (residual) mglevels[l]->residual = residual; 96d0e4de75SBarry Smith if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidual_Default; 97f3ae41bdSBarry Smith if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 98298cc208SBarry Smith ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 992fa5cd67SKarl Rupp 100f3fbd535SBarry Smith mglevels[l]->A = mat; 1014b9ad928SBarry Smith PetscFunctionReturn(0); 1024b9ad928SBarry Smith } 1034b9ad928SBarry Smith 1044b9ad928SBarry Smith #undef __FUNCT__ 105d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation" 1064b9ad928SBarry Smith /*@ 107aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 108bf5b2e24SBarry Smith interpolation from l-1 to the lth level 1094b9ad928SBarry Smith 110ad4df100SBarry Smith Logically Collective on PC and Mat 1114b9ad928SBarry Smith 1124b9ad928SBarry Smith Input Parameters: 1134b9ad928SBarry Smith + pc - the multigrid context 1144b9ad928SBarry Smith . mat - the interpolation operator 115bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 1164b9ad928SBarry Smith 1174b9ad928SBarry Smith Level: advanced 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith Notes: 1204b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 1214b9ad928SBarry Smith for the same level. 1224b9ad928SBarry Smith 1234b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1244b9ad928SBarry Smith out from the matrix size which one it is. 1254b9ad928SBarry Smith 1264b9ad928SBarry Smith .keywords: multigrid, set, interpolate, level 1274b9ad928SBarry Smith 12897177400SBarry Smith .seealso: PCMGSetRestriction() 1294b9ad928SBarry Smith @*/ 1307087cfbeSBarry Smith PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 1314b9ad928SBarry Smith { 132f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 133f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 134fccaa45eSBarry Smith PetscErrorCode ierr; 1354b9ad928SBarry Smith 1364b9ad928SBarry Smith PetscFunctionBegin; 137c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 138ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 139ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 140c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1416bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr); 1422fa5cd67SKarl Rupp 143f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 1444b9ad928SBarry Smith PetscFunctionReturn(0); 1454b9ad928SBarry Smith } 1464b9ad928SBarry Smith 1474b9ad928SBarry Smith #undef __FUNCT__ 148c88c5224SJed Brown #define __FUNCT__ "PCMGGetInterpolation" 149c88c5224SJed Brown /*@ 150c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 151c88c5224SJed Brown interpolation from l-1 to the lth level 152c88c5224SJed Brown 153c88c5224SJed Brown Logically Collective on PC 154c88c5224SJed Brown 155c88c5224SJed Brown Input Parameters: 156c88c5224SJed Brown + pc - the multigrid context 157c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 158c88c5224SJed Brown 159c88c5224SJed Brown Output Parameter: 160c88c5224SJed Brown . mat - the interpolation matrix 161c88c5224SJed Brown 162c88c5224SJed Brown Level: advanced 163c88c5224SJed Brown 164c88c5224SJed Brown .keywords: MG, get, multigrid, interpolation, level 165c88c5224SJed Brown 166c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 167c88c5224SJed Brown @*/ 168c88c5224SJed Brown PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 169c88c5224SJed Brown { 170c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 171c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 172c88c5224SJed Brown PetscErrorCode ierr; 173c88c5224SJed Brown 174c88c5224SJed Brown PetscFunctionBegin; 175c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 176c88c5224SJed Brown PetscValidPointer(mat,3); 177ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 178ce94432eSBarry 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); 179c88c5224SJed Brown if (!mglevels[l]->interpolate) { 180ce94432eSBarry Smith if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetInterpolation()"); 181c88c5224SJed Brown ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 182c88c5224SJed Brown } 183c88c5224SJed Brown *mat = mglevels[l]->interpolate; 184c88c5224SJed Brown PetscFunctionReturn(0); 185c88c5224SJed Brown } 186c88c5224SJed Brown 187c88c5224SJed Brown #undef __FUNCT__ 188d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction" 1894b9ad928SBarry Smith /*@ 19097177400SBarry Smith PCMGSetRestriction - Sets the function to be used to restrict vector 1914b9ad928SBarry Smith from level l to l-1. 1924b9ad928SBarry Smith 193ad4df100SBarry Smith Logically Collective on PC and Mat 1944b9ad928SBarry Smith 1954b9ad928SBarry Smith Input Parameters: 1964b9ad928SBarry Smith + pc - the multigrid context 197c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 198c88c5224SJed Brown - mat - the restriction matrix 1994b9ad928SBarry Smith 2004b9ad928SBarry Smith Level: advanced 2014b9ad928SBarry Smith 2024b9ad928SBarry Smith Notes: 2034b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 2044b9ad928SBarry Smith for the same level. 2054b9ad928SBarry Smith 2064b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 2074b9ad928SBarry Smith out from the matrix size which one it is. 2084b9ad928SBarry Smith 209aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 210fccaa45eSBarry Smith is used. 211fccaa45eSBarry Smith 2124b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level 2134b9ad928SBarry Smith 214aea2a34eSBarry Smith .seealso: PCMGSetInterpolation() 2154b9ad928SBarry Smith @*/ 2167087cfbeSBarry Smith PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 2174b9ad928SBarry Smith { 218fccaa45eSBarry Smith PetscErrorCode ierr; 219f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 220f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2214b9ad928SBarry Smith 2224b9ad928SBarry Smith PetscFunctionBegin; 223c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 224c88c5224SJed Brown PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 225ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 226ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 227c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 2286bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 2292fa5cd67SKarl Rupp 230f3fbd535SBarry Smith mglevels[l]->restrct = mat; 2314b9ad928SBarry Smith PetscFunctionReturn(0); 2324b9ad928SBarry Smith } 2334b9ad928SBarry Smith 2344b9ad928SBarry Smith #undef __FUNCT__ 235c88c5224SJed Brown #define __FUNCT__ "PCMGGetRestriction" 236c88c5224SJed Brown /*@ 237c88c5224SJed Brown PCMGGetRestriction - Gets the function to be used to restrict vector 238c88c5224SJed Brown from level l to l-1. 239c88c5224SJed Brown 240c88c5224SJed Brown Logically Collective on PC and Mat 241c88c5224SJed Brown 242c88c5224SJed Brown Input Parameters: 243c88c5224SJed Brown + pc - the multigrid context 244c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 245c88c5224SJed Brown 246c88c5224SJed Brown Output Parameter: 247c88c5224SJed Brown . mat - the restriction matrix 248c88c5224SJed Brown 249c88c5224SJed Brown Level: advanced 250c88c5224SJed Brown 251c88c5224SJed Brown .keywords: MG, get, multigrid, restriction, level 252c88c5224SJed Brown 253c88c5224SJed Brown .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale() 254c88c5224SJed Brown @*/ 255c88c5224SJed Brown PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 256c88c5224SJed Brown { 257c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 258c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 259c88c5224SJed Brown PetscErrorCode ierr; 260c88c5224SJed Brown 261c88c5224SJed Brown PetscFunctionBegin; 262c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 263c88c5224SJed Brown PetscValidPointer(mat,3); 264ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 265ce94432eSBarry 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); 266c88c5224SJed Brown if (!mglevels[l]->restrct) { 267ce94432eSBarry Smith if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 268c88c5224SJed Brown ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 269c88c5224SJed Brown } 270c88c5224SJed Brown *mat = mglevels[l]->restrct; 271c88c5224SJed Brown PetscFunctionReturn(0); 272c88c5224SJed Brown } 273c88c5224SJed Brown 274c88c5224SJed Brown #undef __FUNCT__ 27573250ac0SBarry Smith #define __FUNCT__ "PCMGSetRScale" 27673250ac0SBarry Smith /*@ 27773250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 27873250ac0SBarry Smith 279c88c5224SJed Brown Logically Collective on PC and Vec 280c88c5224SJed Brown 281c88c5224SJed Brown Input Parameters: 282c88c5224SJed Brown + pc - the multigrid context 283c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 284c88c5224SJed Brown . rscale - the scaling 285c88c5224SJed Brown 286c88c5224SJed Brown Level: advanced 287c88c5224SJed Brown 288c88c5224SJed Brown Notes: 289c88c5224SJed Brown When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R. 290c88c5224SJed Brown 291c88c5224SJed Brown .keywords: MG, set, multigrid, restriction, level 292c88c5224SJed Brown 293c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale() 294c88c5224SJed Brown @*/ 295c88c5224SJed Brown PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 296c88c5224SJed Brown { 297c88c5224SJed Brown PetscErrorCode ierr; 298c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 299c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 300c88c5224SJed Brown 301c88c5224SJed Brown PetscFunctionBegin; 302c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 303ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 304ce94432eSBarry 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); 305c88c5224SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 306c88c5224SJed Brown ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 3072fa5cd67SKarl Rupp 308c88c5224SJed Brown mglevels[l]->rscale = rscale; 309c88c5224SJed Brown PetscFunctionReturn(0); 310c88c5224SJed Brown } 311c88c5224SJed Brown 312c88c5224SJed Brown #undef __FUNCT__ 313c88c5224SJed Brown #define __FUNCT__ "PCMGGetRScale" 314c88c5224SJed Brown /*@ 315c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 316c88c5224SJed Brown 317c88c5224SJed Brown Collective on PC 31873250ac0SBarry Smith 31973250ac0SBarry Smith Input Parameters: 32073250ac0SBarry Smith + pc - the multigrid context 32173250ac0SBarry Smith . rscale - the scaling 32273250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 32373250ac0SBarry Smith 32473250ac0SBarry Smith Level: advanced 32573250ac0SBarry Smith 32673250ac0SBarry Smith Notes: 32773250ac0SBarry Smith When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R. 32873250ac0SBarry Smith 32973250ac0SBarry Smith .keywords: MG, set, multigrid, restriction, level 33073250ac0SBarry Smith 331c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGGetRestriction() 33273250ac0SBarry Smith @*/ 333c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 33473250ac0SBarry Smith { 33573250ac0SBarry Smith PetscErrorCode ierr; 33673250ac0SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 33773250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 33873250ac0SBarry Smith 33973250ac0SBarry Smith PetscFunctionBegin; 34073250ac0SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 341ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 342ce94432eSBarry 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); 343c88c5224SJed Brown if (!mglevels[l]->rscale) { 344c88c5224SJed Brown Mat R; 345c88c5224SJed Brown Vec X,Y,coarse,fine; 346c88c5224SJed Brown PetscInt M,N; 347c88c5224SJed Brown ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 348c88c5224SJed Brown ierr = MatGetVecs(R,&X,&Y);CHKERRQ(ierr); 349c88c5224SJed Brown ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 3502fa5cd67SKarl Rupp if (M < N) { 3512fa5cd67SKarl Rupp fine = X; 3522fa5cd67SKarl Rupp coarse = Y; 3532fa5cd67SKarl Rupp } else if (N < M) { 3542fa5cd67SKarl Rupp fine = Y; coarse = X; 355ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 356c88c5224SJed Brown ierr = VecSet(fine,1.);CHKERRQ(ierr); 357c88c5224SJed Brown ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 358c88c5224SJed Brown ierr = VecDestroy(&fine);CHKERRQ(ierr); 359c88c5224SJed Brown ierr = VecReciprocal(coarse);CHKERRQ(ierr); 360c88c5224SJed Brown mglevels[l]->rscale = coarse; 361c88c5224SJed Brown } 362c88c5224SJed Brown *rscale = mglevels[l]->rscale; 36373250ac0SBarry Smith PetscFunctionReturn(0); 36473250ac0SBarry Smith } 36573250ac0SBarry Smith 36673250ac0SBarry Smith #undef __FUNCT__ 367d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother" 368f39d8e23SSatish Balay /*@ 36997177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 37097177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 37197177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 3724b9ad928SBarry Smith post-smoothing. 3734b9ad928SBarry Smith 3744b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 3754b9ad928SBarry Smith 3764b9ad928SBarry Smith Input Parameters: 3774b9ad928SBarry Smith + pc - the multigrid context 3784b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 3794b9ad928SBarry Smith 3804b9ad928SBarry Smith Ouput Parameters: 3814b9ad928SBarry Smith . ksp - the smoother 3824b9ad928SBarry Smith 3834b9ad928SBarry Smith Level: advanced 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 3864b9ad928SBarry Smith 38797177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 3884b9ad928SBarry Smith @*/ 3897087cfbeSBarry Smith PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 3904b9ad928SBarry Smith { 391f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 392f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 3934b9ad928SBarry Smith 3944b9ad928SBarry Smith PetscFunctionBegin; 395c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 396f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 3974b9ad928SBarry Smith PetscFunctionReturn(0); 3984b9ad928SBarry Smith } 3994b9ad928SBarry Smith 4004b9ad928SBarry Smith #undef __FUNCT__ 401d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp" 402f39d8e23SSatish Balay /*@ 40397177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 4044b9ad928SBarry Smith coarse grid correction (post-smoother). 4054b9ad928SBarry Smith 4064b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4074b9ad928SBarry Smith 4084b9ad928SBarry Smith Input Parameters: 4094b9ad928SBarry Smith + pc - the multigrid context 4104b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith Ouput Parameters: 4134b9ad928SBarry Smith . ksp - the smoother 4144b9ad928SBarry Smith 4154b9ad928SBarry Smith Level: advanced 4164b9ad928SBarry Smith 41789cce641SBarry Smith Notes: calling this will result in a different pre and post smoother so you may need to 41889cce641SBarry Smith set options on the pre smoother also 41989cce641SBarry Smith 4204b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level 4214b9ad928SBarry Smith 42297177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 4234b9ad928SBarry Smith @*/ 4247087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 4254b9ad928SBarry Smith { 426f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 427f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 428dfbe8321SBarry Smith PetscErrorCode ierr; 429f69a0ea3SMatthew Knepley const char *prefix; 4304b9ad928SBarry Smith MPI_Comm comm; 4314b9ad928SBarry Smith 4324b9ad928SBarry Smith PetscFunctionBegin; 433c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4344b9ad928SBarry Smith /* 4354b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 4364b9ad928SBarry Smith Thus we check if a different one has already been allocated, 4374b9ad928SBarry Smith if not we allocate it. 4384b9ad928SBarry Smith */ 439ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 440f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 44119fd82e9SBarry Smith KSPType ksptype; 44219fd82e9SBarry Smith PCType pctype; 443336babb1SJed Brown PC ipc; 444336babb1SJed Brown PetscReal rtol,abstol,dtol; 445336babb1SJed Brown PetscInt maxits; 446336babb1SJed Brown KSPNormType normtype; 447f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 448f3fbd535SBarry Smith ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 449336babb1SJed Brown ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 4503bf036e2SBarry Smith ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 451336babb1SJed Brown ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 452336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 453336babb1SJed Brown ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 454336babb1SJed Brown 455f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 456f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 457f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 458336babb1SJed Brown ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 459336babb1SJed Brown ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 460336babb1SJed Brown ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 4610298fd71SBarry Smith ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPSkipConverged,NULL,NULL);CHKERRQ(ierr); 462336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 463336babb1SJed Brown ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 464*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr); 4654b9ad928SBarry Smith } 466f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 4674b9ad928SBarry Smith PetscFunctionReturn(0); 4684b9ad928SBarry Smith } 4694b9ad928SBarry Smith 4704b9ad928SBarry Smith #undef __FUNCT__ 4719dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown" 472f39d8e23SSatish Balay /*@ 47397177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 4744b9ad928SBarry Smith coarse grid correction (pre-smoother). 4754b9ad928SBarry Smith 4764b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4774b9ad928SBarry Smith 4784b9ad928SBarry Smith Input Parameters: 4794b9ad928SBarry Smith + pc - the multigrid context 4804b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4814b9ad928SBarry Smith 4824b9ad928SBarry Smith Ouput Parameters: 4834b9ad928SBarry Smith . ksp - the smoother 4844b9ad928SBarry Smith 4854b9ad928SBarry Smith Level: advanced 4864b9ad928SBarry Smith 48789cce641SBarry Smith Notes: calling this will result in a different pre and post smoother so you may need to 48889cce641SBarry Smith set options on the post smoother also 48989cce641SBarry Smith 4904b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 4914b9ad928SBarry Smith 49297177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 4934b9ad928SBarry Smith @*/ 4947087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 4954b9ad928SBarry Smith { 496dfbe8321SBarry Smith PetscErrorCode ierr; 497f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 498f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 4994b9ad928SBarry Smith 5004b9ad928SBarry Smith PetscFunctionBegin; 501c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5024b9ad928SBarry Smith /* make sure smoother up and down are different */ 503c5eb9154SBarry Smith if (l) { 5040298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr); 505d8148a5aSMatthew Knepley } 506f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5074b9ad928SBarry Smith PetscFunctionReturn(0); 5084b9ad928SBarry Smith } 5094b9ad928SBarry Smith 5104b9ad928SBarry Smith #undef __FUNCT__ 5119dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel" 5124b9ad928SBarry Smith /*@ 51397177400SBarry Smith PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level. 5144b9ad928SBarry Smith 515ad4df100SBarry Smith Logically Collective on PC 5164b9ad928SBarry Smith 5174b9ad928SBarry Smith Input Parameters: 5184b9ad928SBarry Smith + pc - the multigrid context 5194b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 5204b9ad928SBarry Smith - n - the number of cycles 5214b9ad928SBarry Smith 5224b9ad928SBarry Smith Level: advanced 5234b9ad928SBarry Smith 5244b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 5254b9ad928SBarry Smith 52697177400SBarry Smith .seealso: PCMGSetCycles() 5274b9ad928SBarry Smith @*/ 5287087cfbeSBarry Smith PetscErrorCode PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c) 5294b9ad928SBarry Smith { 530f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 531f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5324b9ad928SBarry Smith 5334b9ad928SBarry Smith PetscFunctionBegin; 534c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 535ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 536c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,l,2); 537c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,c,3); 538f3fbd535SBarry Smith mglevels[l]->cycles = c; 5394b9ad928SBarry Smith PetscFunctionReturn(0); 5404b9ad928SBarry Smith } 5414b9ad928SBarry Smith 5424b9ad928SBarry Smith #undef __FUNCT__ 543d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs" 5444b9ad928SBarry Smith /*@ 54597177400SBarry Smith PCMGSetRhs - Sets the vector space to be used to store the right-hand side 546fccaa45eSBarry Smith on a particular level. 5474b9ad928SBarry Smith 548ad4df100SBarry Smith Logically Collective on PC and Vec 5494b9ad928SBarry Smith 5504b9ad928SBarry Smith Input Parameters: 5514b9ad928SBarry Smith + pc - the multigrid context 5524b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 5534b9ad928SBarry Smith - c - the space 5544b9ad928SBarry Smith 5554b9ad928SBarry Smith Level: advanced 5564b9ad928SBarry Smith 557fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 558fccaa45eSBarry Smith 559fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 560fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 561fccaa45eSBarry Smith 5624b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level 5634b9ad928SBarry Smith 56497177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 5654b9ad928SBarry Smith @*/ 5667087cfbeSBarry Smith PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 5674b9ad928SBarry Smith { 568fccaa45eSBarry Smith PetscErrorCode ierr; 569f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 570f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5714b9ad928SBarry Smith 5724b9ad928SBarry Smith PetscFunctionBegin; 573c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 574ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 575ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 576c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 5776bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 5782fa5cd67SKarl Rupp 579f3fbd535SBarry Smith mglevels[l]->b = c; 5804b9ad928SBarry Smith PetscFunctionReturn(0); 5814b9ad928SBarry Smith } 5824b9ad928SBarry Smith 5834b9ad928SBarry Smith #undef __FUNCT__ 584d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX" 5854b9ad928SBarry Smith /*@ 58697177400SBarry Smith PCMGSetX - Sets the vector space to be used to store the solution on a 587fccaa45eSBarry Smith particular level. 5884b9ad928SBarry Smith 589ad4df100SBarry Smith Logically Collective on PC and Vec 5904b9ad928SBarry Smith 5914b9ad928SBarry Smith Input Parameters: 5924b9ad928SBarry Smith + pc - the multigrid context 593251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 5944b9ad928SBarry Smith - c - the space 5954b9ad928SBarry Smith 5964b9ad928SBarry Smith Level: advanced 5974b9ad928SBarry Smith 598fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 599fccaa45eSBarry Smith 600fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 601fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 602fccaa45eSBarry Smith 6034b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level 6044b9ad928SBarry Smith 60597177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 6064b9ad928SBarry Smith @*/ 6077087cfbeSBarry Smith PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 6084b9ad928SBarry Smith { 609fccaa45eSBarry Smith PetscErrorCode ierr; 610f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 611f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6124b9ad928SBarry Smith 6134b9ad928SBarry Smith PetscFunctionBegin; 614c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 615ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 616ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 617c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6186bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 6192fa5cd67SKarl Rupp 620f3fbd535SBarry Smith mglevels[l]->x = c; 6214b9ad928SBarry Smith PetscFunctionReturn(0); 6224b9ad928SBarry Smith } 6234b9ad928SBarry Smith 6244b9ad928SBarry Smith #undef __FUNCT__ 625d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR" 6264b9ad928SBarry Smith /*@ 62797177400SBarry Smith PCMGSetR - Sets the vector space to be used to store the residual on a 628fccaa45eSBarry Smith particular level. 6294b9ad928SBarry Smith 630ad4df100SBarry Smith Logically Collective on PC and Vec 6314b9ad928SBarry Smith 6324b9ad928SBarry Smith Input Parameters: 6334b9ad928SBarry Smith + pc - the multigrid context 6344b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 6354b9ad928SBarry Smith - c - the space 6364b9ad928SBarry Smith 6374b9ad928SBarry Smith Level: advanced 6384b9ad928SBarry Smith 639fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 640fccaa45eSBarry Smith 641fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 642fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 643fccaa45eSBarry Smith 6444b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level 6454b9ad928SBarry Smith @*/ 6467087cfbeSBarry Smith PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 6474b9ad928SBarry Smith { 648fccaa45eSBarry Smith PetscErrorCode ierr; 649f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 650f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6514b9ad928SBarry Smith 6524b9ad928SBarry Smith PetscFunctionBegin; 653c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 654ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 655ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 656c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6576bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 6582fa5cd67SKarl Rupp 659f3fbd535SBarry Smith mglevels[l]->r = c; 6604b9ad928SBarry Smith PetscFunctionReturn(0); 6614b9ad928SBarry Smith } 662