14b9ad928SBarry Smith 2c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscksp.h" I*/ 3f3f935beSJed Brown /*I "petscpcmg.h" I*/ 44b9ad928SBarry Smith 54b9ad928SBarry Smith #undef __FUNCT__ 697177400SBarry Smith #define __FUNCT__ "PCMGDefaultResidual" 71f6cc5b2SSatish Balay /*@C 897177400SBarry Smith PCMGDefaultResidual - Default routine to calculate the residual. 94b9ad928SBarry Smith 104b9ad928SBarry Smith Collective on Mat and Vec 114b9ad928SBarry Smith 124b9ad928SBarry Smith Input Parameters: 134b9ad928SBarry Smith + mat - the matrix 144b9ad928SBarry Smith . b - the right-hand-side 154b9ad928SBarry Smith - x - the approximate solution 164b9ad928SBarry Smith 174b9ad928SBarry Smith Output Parameter: 184b9ad928SBarry Smith . r - location to store the residual 194b9ad928SBarry Smith 204b9ad928SBarry Smith Level: advanced 214b9ad928SBarry Smith 224b9ad928SBarry Smith .keywords: MG, default, multigrid, residual 234b9ad928SBarry Smith 2497177400SBarry Smith .seealso: PCMGSetResidual() 254b9ad928SBarry Smith @*/ 267087cfbeSBarry Smith PetscErrorCode PCMGDefaultResidual(Mat mat,Vec b,Vec x,Vec r) 274b9ad928SBarry Smith { 28dfbe8321SBarry Smith PetscErrorCode ierr; 294b9ad928SBarry Smith 304b9ad928SBarry Smith PetscFunctionBegin; 314b9ad928SBarry Smith ierr = MatMult(mat,x,r);CHKERRQ(ierr); 32efb30889SBarry Smith ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 334b9ad928SBarry Smith PetscFunctionReturn(0); 344b9ad928SBarry Smith } 354b9ad928SBarry Smith 364b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/ 374b9ad928SBarry Smith 384b9ad928SBarry Smith #undef __FUNCT__ 39d6337fedSHong Zhang #define __FUNCT__ "PCMGGetCoarseSolve" 40f39d8e23SSatish Balay /*@ 4197177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 424b9ad928SBarry Smith 434b9ad928SBarry Smith Not Collective 444b9ad928SBarry Smith 454b9ad928SBarry Smith Input Parameter: 464b9ad928SBarry Smith . pc - the multigrid context 474b9ad928SBarry Smith 484b9ad928SBarry Smith Output Parameter: 494b9ad928SBarry Smith . ksp - the coarse grid solver context 504b9ad928SBarry Smith 514b9ad928SBarry Smith Level: advanced 524b9ad928SBarry Smith 534b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid 544b9ad928SBarry Smith @*/ 557087cfbeSBarry Smith PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 564b9ad928SBarry Smith { 57f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 58f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 594b9ad928SBarry Smith 604b9ad928SBarry Smith PetscFunctionBegin; 61c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 62f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 634b9ad928SBarry Smith PetscFunctionReturn(0); 644b9ad928SBarry Smith } 654b9ad928SBarry Smith 664b9ad928SBarry Smith #undef __FUNCT__ 67d6337fedSHong Zhang #define __FUNCT__ "PCMGSetResidual" 684b9ad928SBarry Smith /*@C 6997177400SBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual 704b9ad928SBarry Smith on the lth level. 714b9ad928SBarry Smith 72ad4df100SBarry Smith Logically Collective on PC and Mat 734b9ad928SBarry Smith 744b9ad928SBarry Smith Input Parameters: 754b9ad928SBarry Smith + pc - the multigrid context 764b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 77157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 78157726a2SBarry Smith previous one were provided then PCMGDefaultResidual() is used 794b9ad928SBarry Smith - mat - matrix associated with residual 804b9ad928SBarry Smith 814b9ad928SBarry Smith Level: advanced 824b9ad928SBarry Smith 834b9ad928SBarry Smith .keywords: MG, set, multigrid, residual, level 844b9ad928SBarry Smith 8597177400SBarry Smith .seealso: PCMGDefaultResidual() 864b9ad928SBarry Smith @*/ 877087cfbeSBarry Smith PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 884b9ad928SBarry Smith { 89f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 90f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 91298cc208SBarry Smith PetscErrorCode ierr; 924b9ad928SBarry Smith 934b9ad928SBarry Smith PetscFunctionBegin; 94c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 95*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 962fa5cd67SKarl Rupp if (residual) mglevels[l]->residual = residual; 972fa5cd67SKarl Rupp if (!mglevels[l]->residual) mglevels[l]->residual = PCMGDefaultResidual; 98f3ae41bdSBarry Smith if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 99298cc208SBarry Smith ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 1002fa5cd67SKarl Rupp 101f3fbd535SBarry Smith mglevels[l]->A = mat; 1024b9ad928SBarry Smith PetscFunctionReturn(0); 1034b9ad928SBarry Smith } 1044b9ad928SBarry Smith 1054b9ad928SBarry Smith #undef __FUNCT__ 106d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation" 1074b9ad928SBarry Smith /*@ 108aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 109bf5b2e24SBarry Smith interpolation from l-1 to the lth level 1104b9ad928SBarry Smith 111ad4df100SBarry Smith Logically Collective on PC and Mat 1124b9ad928SBarry Smith 1134b9ad928SBarry Smith Input Parameters: 1144b9ad928SBarry Smith + pc - the multigrid context 1154b9ad928SBarry Smith . mat - the interpolation operator 116bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 1174b9ad928SBarry Smith 1184b9ad928SBarry Smith Level: advanced 1194b9ad928SBarry Smith 1204b9ad928SBarry Smith Notes: 1214b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 1224b9ad928SBarry Smith for the same level. 1234b9ad928SBarry Smith 1244b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1254b9ad928SBarry Smith out from the matrix size which one it is. 1264b9ad928SBarry Smith 1274b9ad928SBarry Smith .keywords: multigrid, set, interpolate, level 1284b9ad928SBarry Smith 12997177400SBarry Smith .seealso: PCMGSetRestriction() 1304b9ad928SBarry Smith @*/ 1317087cfbeSBarry Smith PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 1324b9ad928SBarry Smith { 133f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 134f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 135fccaa45eSBarry Smith PetscErrorCode ierr; 1364b9ad928SBarry Smith 1374b9ad928SBarry Smith PetscFunctionBegin; 138c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 139*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 140*ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 141c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1426bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr); 1432fa5cd67SKarl Rupp 144f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 1454b9ad928SBarry Smith PetscFunctionReturn(0); 1464b9ad928SBarry Smith } 1474b9ad928SBarry Smith 1484b9ad928SBarry Smith #undef __FUNCT__ 149c88c5224SJed Brown #define __FUNCT__ "PCMGGetInterpolation" 150c88c5224SJed Brown /*@ 151c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 152c88c5224SJed Brown interpolation from l-1 to the lth level 153c88c5224SJed Brown 154c88c5224SJed Brown Logically Collective on PC 155c88c5224SJed Brown 156c88c5224SJed Brown Input Parameters: 157c88c5224SJed Brown + pc - the multigrid context 158c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 159c88c5224SJed Brown 160c88c5224SJed Brown Output Parameter: 161c88c5224SJed Brown . mat - the interpolation matrix 162c88c5224SJed Brown 163c88c5224SJed Brown Level: advanced 164c88c5224SJed Brown 165c88c5224SJed Brown .keywords: MG, get, multigrid, interpolation, level 166c88c5224SJed Brown 167c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 168c88c5224SJed Brown @*/ 169c88c5224SJed Brown PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 170c88c5224SJed Brown { 171c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 172c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 173c88c5224SJed Brown PetscErrorCode ierr; 174c88c5224SJed Brown 175c88c5224SJed Brown PetscFunctionBegin; 176c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 177c88c5224SJed Brown PetscValidPointer(mat,3); 178*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 179*ce94432eSBarry 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); 180c88c5224SJed Brown if (!mglevels[l]->interpolate) { 181*ce94432eSBarry Smith if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetInterpolation()"); 182c88c5224SJed Brown ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 183c88c5224SJed Brown } 184c88c5224SJed Brown *mat = mglevels[l]->interpolate; 185c88c5224SJed Brown PetscFunctionReturn(0); 186c88c5224SJed Brown } 187c88c5224SJed Brown 188c88c5224SJed Brown #undef __FUNCT__ 189d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction" 1904b9ad928SBarry Smith /*@ 19197177400SBarry Smith PCMGSetRestriction - Sets the function to be used to restrict vector 1924b9ad928SBarry Smith from level l to l-1. 1934b9ad928SBarry Smith 194ad4df100SBarry Smith Logically Collective on PC and Mat 1954b9ad928SBarry Smith 1964b9ad928SBarry Smith Input Parameters: 1974b9ad928SBarry Smith + pc - the multigrid context 198c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 199c88c5224SJed Brown - mat - the restriction matrix 2004b9ad928SBarry Smith 2014b9ad928SBarry Smith Level: advanced 2024b9ad928SBarry Smith 2034b9ad928SBarry Smith Notes: 2044b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 2054b9ad928SBarry Smith for the same level. 2064b9ad928SBarry Smith 2074b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 2084b9ad928SBarry Smith out from the matrix size which one it is. 2094b9ad928SBarry Smith 210aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 211fccaa45eSBarry Smith is used. 212fccaa45eSBarry Smith 2134b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level 2144b9ad928SBarry Smith 215aea2a34eSBarry Smith .seealso: PCMGSetInterpolation() 2164b9ad928SBarry Smith @*/ 2177087cfbeSBarry Smith PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 2184b9ad928SBarry Smith { 219fccaa45eSBarry Smith PetscErrorCode ierr; 220f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 221f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2224b9ad928SBarry Smith 2234b9ad928SBarry Smith PetscFunctionBegin; 224c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 225c88c5224SJed Brown PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 226*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 227*ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 228c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 2296bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 2302fa5cd67SKarl Rupp 231f3fbd535SBarry Smith mglevels[l]->restrct = mat; 2324b9ad928SBarry Smith PetscFunctionReturn(0); 2334b9ad928SBarry Smith } 2344b9ad928SBarry Smith 2354b9ad928SBarry Smith #undef __FUNCT__ 236c88c5224SJed Brown #define __FUNCT__ "PCMGGetRestriction" 237c88c5224SJed Brown /*@ 238c88c5224SJed Brown PCMGGetRestriction - Gets the function to be used to restrict vector 239c88c5224SJed Brown from level l to l-1. 240c88c5224SJed Brown 241c88c5224SJed Brown Logically Collective on PC and Mat 242c88c5224SJed Brown 243c88c5224SJed Brown Input Parameters: 244c88c5224SJed Brown + pc - the multigrid context 245c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 246c88c5224SJed Brown 247c88c5224SJed Brown Output Parameter: 248c88c5224SJed Brown . mat - the restriction matrix 249c88c5224SJed Brown 250c88c5224SJed Brown Level: advanced 251c88c5224SJed Brown 252c88c5224SJed Brown .keywords: MG, get, multigrid, restriction, level 253c88c5224SJed Brown 254c88c5224SJed Brown .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale() 255c88c5224SJed Brown @*/ 256c88c5224SJed Brown PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 257c88c5224SJed Brown { 258c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 259c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 260c88c5224SJed Brown PetscErrorCode ierr; 261c88c5224SJed Brown 262c88c5224SJed Brown PetscFunctionBegin; 263c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 264c88c5224SJed Brown PetscValidPointer(mat,3); 265*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 266*ce94432eSBarry 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); 267c88c5224SJed Brown if (!mglevels[l]->restrct) { 268*ce94432eSBarry Smith if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 269c88c5224SJed Brown ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 270c88c5224SJed Brown } 271c88c5224SJed Brown *mat = mglevels[l]->restrct; 272c88c5224SJed Brown PetscFunctionReturn(0); 273c88c5224SJed Brown } 274c88c5224SJed Brown 275c88c5224SJed Brown #undef __FUNCT__ 27673250ac0SBarry Smith #define __FUNCT__ "PCMGSetRScale" 27773250ac0SBarry Smith /*@ 27873250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 27973250ac0SBarry Smith 280c88c5224SJed Brown Logically Collective on PC and Vec 281c88c5224SJed Brown 282c88c5224SJed Brown Input Parameters: 283c88c5224SJed Brown + pc - the multigrid context 284c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 285c88c5224SJed Brown . rscale - the scaling 286c88c5224SJed Brown 287c88c5224SJed Brown Level: advanced 288c88c5224SJed Brown 289c88c5224SJed Brown Notes: 290c88c5224SJed Brown When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R. 291c88c5224SJed Brown 292c88c5224SJed Brown .keywords: MG, set, multigrid, restriction, level 293c88c5224SJed Brown 294c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale() 295c88c5224SJed Brown @*/ 296c88c5224SJed Brown PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 297c88c5224SJed Brown { 298c88c5224SJed Brown PetscErrorCode ierr; 299c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 300c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 301c88c5224SJed Brown 302c88c5224SJed Brown PetscFunctionBegin; 303c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 304*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 305*ce94432eSBarry 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); 306c88c5224SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 307c88c5224SJed Brown ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 3082fa5cd67SKarl Rupp 309c88c5224SJed Brown mglevels[l]->rscale = rscale; 310c88c5224SJed Brown PetscFunctionReturn(0); 311c88c5224SJed Brown } 312c88c5224SJed Brown 313c88c5224SJed Brown #undef __FUNCT__ 314c88c5224SJed Brown #define __FUNCT__ "PCMGGetRScale" 315c88c5224SJed Brown /*@ 316c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 317c88c5224SJed Brown 318c88c5224SJed Brown Collective on PC 31973250ac0SBarry Smith 32073250ac0SBarry Smith Input Parameters: 32173250ac0SBarry Smith + pc - the multigrid context 32273250ac0SBarry Smith . rscale - the scaling 32373250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 32473250ac0SBarry Smith 32573250ac0SBarry Smith Level: advanced 32673250ac0SBarry Smith 32773250ac0SBarry Smith Notes: 32873250ac0SBarry 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. 32973250ac0SBarry Smith 33073250ac0SBarry Smith .keywords: MG, set, multigrid, restriction, level 33173250ac0SBarry Smith 332c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGGetRestriction() 33373250ac0SBarry Smith @*/ 334c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 33573250ac0SBarry Smith { 33673250ac0SBarry Smith PetscErrorCode ierr; 33773250ac0SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 33873250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 33973250ac0SBarry Smith 34073250ac0SBarry Smith PetscFunctionBegin; 34173250ac0SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 342*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 343*ce94432eSBarry 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); 344c88c5224SJed Brown if (!mglevels[l]->rscale) { 345c88c5224SJed Brown Mat R; 346c88c5224SJed Brown Vec X,Y,coarse,fine; 347c88c5224SJed Brown PetscInt M,N; 348c88c5224SJed Brown ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 349c88c5224SJed Brown ierr = MatGetVecs(R,&X,&Y);CHKERRQ(ierr); 350c88c5224SJed Brown ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 3512fa5cd67SKarl Rupp if (M < N) { 3522fa5cd67SKarl Rupp fine = X; 3532fa5cd67SKarl Rupp coarse = Y; 3542fa5cd67SKarl Rupp } else if (N < M) { 3552fa5cd67SKarl Rupp fine = Y; coarse = X; 356*ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 357c88c5224SJed Brown ierr = VecSet(fine,1.);CHKERRQ(ierr); 358c88c5224SJed Brown ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 359c88c5224SJed Brown ierr = VecDestroy(&fine);CHKERRQ(ierr); 360c88c5224SJed Brown ierr = VecReciprocal(coarse);CHKERRQ(ierr); 361c88c5224SJed Brown mglevels[l]->rscale = coarse; 362c88c5224SJed Brown } 363c88c5224SJed Brown *rscale = mglevels[l]->rscale; 36473250ac0SBarry Smith PetscFunctionReturn(0); 36573250ac0SBarry Smith } 36673250ac0SBarry Smith 36773250ac0SBarry Smith #undef __FUNCT__ 368d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother" 369f39d8e23SSatish Balay /*@ 37097177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 37197177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 37297177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 3734b9ad928SBarry Smith post-smoothing. 3744b9ad928SBarry Smith 3754b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 3764b9ad928SBarry Smith 3774b9ad928SBarry Smith Input Parameters: 3784b9ad928SBarry Smith + pc - the multigrid context 3794b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 3804b9ad928SBarry Smith 3814b9ad928SBarry Smith Ouput Parameters: 3824b9ad928SBarry Smith . ksp - the smoother 3834b9ad928SBarry Smith 3844b9ad928SBarry Smith Level: advanced 3854b9ad928SBarry Smith 3864b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 3874b9ad928SBarry Smith 38897177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 3894b9ad928SBarry Smith @*/ 3907087cfbeSBarry Smith PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 3914b9ad928SBarry Smith { 392f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 393f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 3944b9ad928SBarry Smith 3954b9ad928SBarry Smith PetscFunctionBegin; 396c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 397f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 3984b9ad928SBarry Smith PetscFunctionReturn(0); 3994b9ad928SBarry Smith } 4004b9ad928SBarry Smith 4014b9ad928SBarry Smith #undef __FUNCT__ 402d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp" 403f39d8e23SSatish Balay /*@ 40497177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 4054b9ad928SBarry Smith coarse grid correction (post-smoother). 4064b9ad928SBarry Smith 4074b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4084b9ad928SBarry Smith 4094b9ad928SBarry Smith Input Parameters: 4104b9ad928SBarry Smith + pc - the multigrid context 4114b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4124b9ad928SBarry Smith 4134b9ad928SBarry Smith Ouput Parameters: 4144b9ad928SBarry Smith . ksp - the smoother 4154b9ad928SBarry Smith 4164b9ad928SBarry Smith Level: advanced 4174b9ad928SBarry Smith 41889cce641SBarry Smith Notes: calling this will result in a different pre and post smoother so you may need to 41989cce641SBarry Smith set options on the pre smoother also 42089cce641SBarry Smith 4214b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level 4224b9ad928SBarry Smith 42397177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 4244b9ad928SBarry Smith @*/ 4257087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 4264b9ad928SBarry Smith { 427f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 428f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 429dfbe8321SBarry Smith PetscErrorCode ierr; 430f69a0ea3SMatthew Knepley const char *prefix; 4314b9ad928SBarry Smith MPI_Comm comm; 4324b9ad928SBarry Smith 4334b9ad928SBarry Smith PetscFunctionBegin; 434c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4354b9ad928SBarry Smith /* 4364b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 4374b9ad928SBarry Smith Thus we check if a different one has already been allocated, 4384b9ad928SBarry Smith if not we allocate it. 4394b9ad928SBarry Smith */ 440*ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 441f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 44219fd82e9SBarry Smith KSPType ksptype; 44319fd82e9SBarry Smith PCType pctype; 444336babb1SJed Brown PC ipc; 445336babb1SJed Brown PetscReal rtol,abstol,dtol; 446336babb1SJed Brown PetscInt maxits; 447336babb1SJed Brown KSPNormType normtype; 448f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 449f3fbd535SBarry Smith ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 450336babb1SJed Brown ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 4513bf036e2SBarry Smith ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 452336babb1SJed Brown ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 453336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 454336babb1SJed Brown ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 455336babb1SJed Brown 456f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 457f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 458f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 459336babb1SJed Brown ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 460336babb1SJed Brown ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 461336babb1SJed Brown ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 4620298fd71SBarry Smith ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPSkipConverged,NULL,NULL);CHKERRQ(ierr); 463336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 464336babb1SJed Brown ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 465f3fbd535SBarry Smith ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr); 4664b9ad928SBarry Smith } 467f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 4684b9ad928SBarry Smith PetscFunctionReturn(0); 4694b9ad928SBarry Smith } 4704b9ad928SBarry Smith 4714b9ad928SBarry Smith #undef __FUNCT__ 4729dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown" 473f39d8e23SSatish Balay /*@ 47497177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 4754b9ad928SBarry Smith coarse grid correction (pre-smoother). 4764b9ad928SBarry Smith 4774b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4784b9ad928SBarry Smith 4794b9ad928SBarry Smith Input Parameters: 4804b9ad928SBarry Smith + pc - the multigrid context 4814b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4824b9ad928SBarry Smith 4834b9ad928SBarry Smith Ouput Parameters: 4844b9ad928SBarry Smith . ksp - the smoother 4854b9ad928SBarry Smith 4864b9ad928SBarry Smith Level: advanced 4874b9ad928SBarry Smith 48889cce641SBarry Smith Notes: calling this will result in a different pre and post smoother so you may need to 48989cce641SBarry Smith set options on the post smoother also 49089cce641SBarry Smith 4914b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 4924b9ad928SBarry Smith 49397177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 4944b9ad928SBarry Smith @*/ 4957087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 4964b9ad928SBarry Smith { 497dfbe8321SBarry Smith PetscErrorCode ierr; 498f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 499f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5004b9ad928SBarry Smith 5014b9ad928SBarry Smith PetscFunctionBegin; 502c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5034b9ad928SBarry Smith /* make sure smoother up and down are different */ 504c5eb9154SBarry Smith if (l) { 5050298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr); 506d8148a5aSMatthew Knepley } 507f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5084b9ad928SBarry Smith PetscFunctionReturn(0); 5094b9ad928SBarry Smith } 5104b9ad928SBarry Smith 5114b9ad928SBarry Smith #undef __FUNCT__ 5129dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel" 5134b9ad928SBarry Smith /*@ 51497177400SBarry Smith PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level. 5154b9ad928SBarry Smith 516ad4df100SBarry Smith Logically Collective on PC 5174b9ad928SBarry Smith 5184b9ad928SBarry Smith Input Parameters: 5194b9ad928SBarry Smith + pc - the multigrid context 5204b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 5214b9ad928SBarry Smith - n - the number of cycles 5224b9ad928SBarry Smith 5234b9ad928SBarry Smith Level: advanced 5244b9ad928SBarry Smith 5254b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 5264b9ad928SBarry Smith 52797177400SBarry Smith .seealso: PCMGSetCycles() 5284b9ad928SBarry Smith @*/ 5297087cfbeSBarry Smith PetscErrorCode PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c) 5304b9ad928SBarry Smith { 531f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 532f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5334b9ad928SBarry Smith 5344b9ad928SBarry Smith PetscFunctionBegin; 535c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 536*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 537c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,l,2); 538c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,c,3); 539f3fbd535SBarry Smith mglevels[l]->cycles = c; 5404b9ad928SBarry Smith PetscFunctionReturn(0); 5414b9ad928SBarry Smith } 5424b9ad928SBarry Smith 5434b9ad928SBarry Smith #undef __FUNCT__ 544d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs" 5454b9ad928SBarry Smith /*@ 54697177400SBarry Smith PCMGSetRhs - Sets the vector space to be used to store the right-hand side 547fccaa45eSBarry Smith on a particular level. 5484b9ad928SBarry Smith 549ad4df100SBarry Smith Logically Collective on PC and Vec 5504b9ad928SBarry Smith 5514b9ad928SBarry Smith Input Parameters: 5524b9ad928SBarry Smith + pc - the multigrid context 5534b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 5544b9ad928SBarry Smith - c - the space 5554b9ad928SBarry Smith 5564b9ad928SBarry Smith Level: advanced 5574b9ad928SBarry Smith 558fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 559fccaa45eSBarry Smith 560fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 561fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 562fccaa45eSBarry Smith 5634b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level 5644b9ad928SBarry Smith 56597177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 5664b9ad928SBarry Smith @*/ 5677087cfbeSBarry Smith PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 5684b9ad928SBarry Smith { 569fccaa45eSBarry Smith PetscErrorCode ierr; 570f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 571f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5724b9ad928SBarry Smith 5734b9ad928SBarry Smith PetscFunctionBegin; 574c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 575*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 576*ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 577c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 5786bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 5792fa5cd67SKarl Rupp 580f3fbd535SBarry Smith mglevels[l]->b = c; 5814b9ad928SBarry Smith PetscFunctionReturn(0); 5824b9ad928SBarry Smith } 5834b9ad928SBarry Smith 5844b9ad928SBarry Smith #undef __FUNCT__ 585d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX" 5864b9ad928SBarry Smith /*@ 58797177400SBarry Smith PCMGSetX - Sets the vector space to be used to store the solution on a 588fccaa45eSBarry Smith particular level. 5894b9ad928SBarry Smith 590ad4df100SBarry Smith Logically Collective on PC and Vec 5914b9ad928SBarry Smith 5924b9ad928SBarry Smith Input Parameters: 5934b9ad928SBarry Smith + pc - the multigrid context 594251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 5954b9ad928SBarry Smith - c - the space 5964b9ad928SBarry Smith 5974b9ad928SBarry Smith Level: advanced 5984b9ad928SBarry Smith 599fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 600fccaa45eSBarry Smith 601fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 602fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 603fccaa45eSBarry Smith 6044b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level 6054b9ad928SBarry Smith 60697177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 6074b9ad928SBarry Smith @*/ 6087087cfbeSBarry Smith PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 6094b9ad928SBarry Smith { 610fccaa45eSBarry Smith PetscErrorCode ierr; 611f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 612f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6134b9ad928SBarry Smith 6144b9ad928SBarry Smith PetscFunctionBegin; 615c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 616*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 617*ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 618c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6196bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 6202fa5cd67SKarl Rupp 621f3fbd535SBarry Smith mglevels[l]->x = c; 6224b9ad928SBarry Smith PetscFunctionReturn(0); 6234b9ad928SBarry Smith } 6244b9ad928SBarry Smith 6254b9ad928SBarry Smith #undef __FUNCT__ 626d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR" 6274b9ad928SBarry Smith /*@ 62897177400SBarry Smith PCMGSetR - Sets the vector space to be used to store the residual on a 629fccaa45eSBarry Smith particular level. 6304b9ad928SBarry Smith 631ad4df100SBarry Smith Logically Collective on PC and Vec 6324b9ad928SBarry Smith 6334b9ad928SBarry Smith Input Parameters: 6344b9ad928SBarry Smith + pc - the multigrid context 6354b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 6364b9ad928SBarry Smith - c - the space 6374b9ad928SBarry Smith 6384b9ad928SBarry Smith Level: advanced 6394b9ad928SBarry Smith 640fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 641fccaa45eSBarry Smith 642fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 643fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 644fccaa45eSBarry Smith 6454b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level 6464b9ad928SBarry Smith @*/ 6477087cfbeSBarry Smith PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 6484b9ad928SBarry Smith { 649fccaa45eSBarry Smith PetscErrorCode ierr; 650f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 651f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6524b9ad928SBarry Smith 6534b9ad928SBarry Smith PetscFunctionBegin; 654c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 655*ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 656*ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 657c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6586bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 6592fa5cd67SKarl Rupp 660f3fbd535SBarry Smith mglevels[l]->r = c; 6614b9ad928SBarry Smith PetscFunctionReturn(0); 6624b9ad928SBarry Smith } 663