14b9ad928SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 34b9ad928SBarry Smith 4f9426fe0SMark Adams /* ---------------------------------------------------------------------------*/ 51f6cc5b2SSatish Balay /*@C 654b2cd4bSJed Brown PCMGResidualDefault - Default routine to calculate the residual. 74b9ad928SBarry Smith 84b9ad928SBarry Smith Collective on Mat and Vec 94b9ad928SBarry Smith 104b9ad928SBarry Smith Input Parameters: 114b9ad928SBarry Smith + mat - the matrix 124b9ad928SBarry Smith . b - the right-hand-side 134b9ad928SBarry Smith - x - the approximate solution 144b9ad928SBarry Smith 154b9ad928SBarry Smith Output Parameter: 164b9ad928SBarry Smith . r - location to store the residual 174b9ad928SBarry Smith 18d0e4de75SBarry Smith Level: developer 194b9ad928SBarry Smith 204b9ad928SBarry Smith .keywords: MG, default, multigrid, residual 214b9ad928SBarry Smith 2297177400SBarry Smith .seealso: PCMGSetResidual() 234b9ad928SBarry Smith @*/ 2454b2cd4bSJed Brown PetscErrorCode PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r) 254b9ad928SBarry Smith { 26dfbe8321SBarry Smith PetscErrorCode ierr; 274b9ad928SBarry Smith 284b9ad928SBarry Smith PetscFunctionBegin; 29f9426fe0SMark Adams ierr = MatResidual(mat,b,x,r);CHKERRQ(ierr); 304b9ad928SBarry Smith PetscFunctionReturn(0); 314b9ad928SBarry Smith } 324b9ad928SBarry Smith 33f39d8e23SSatish Balay /*@ 3497177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 354b9ad928SBarry Smith 364b9ad928SBarry Smith Not Collective 374b9ad928SBarry Smith 384b9ad928SBarry Smith Input Parameter: 394b9ad928SBarry Smith . pc - the multigrid context 404b9ad928SBarry Smith 414b9ad928SBarry Smith Output Parameter: 424b9ad928SBarry Smith . ksp - the coarse grid solver context 434b9ad928SBarry Smith 444b9ad928SBarry Smith Level: advanced 454b9ad928SBarry Smith 464b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid 474b9ad928SBarry Smith @*/ 487087cfbeSBarry Smith PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 494b9ad928SBarry Smith { 50f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 51f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 524b9ad928SBarry Smith 534b9ad928SBarry Smith PetscFunctionBegin; 54c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 55f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 564b9ad928SBarry Smith PetscFunctionReturn(0); 574b9ad928SBarry Smith } 584b9ad928SBarry Smith 594b9ad928SBarry Smith /*@C 6097177400SBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual 614b9ad928SBarry Smith on the lth level. 624b9ad928SBarry Smith 63ad4df100SBarry Smith Logically Collective on PC and Mat 644b9ad928SBarry Smith 654b9ad928SBarry Smith Input Parameters: 664b9ad928SBarry Smith + pc - the multigrid context 674b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 68157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 69d0e4de75SBarry Smith previous one were provided then a default is used 704b9ad928SBarry Smith - mat - matrix associated with residual 714b9ad928SBarry Smith 724b9ad928SBarry Smith Level: advanced 734b9ad928SBarry Smith 744b9ad928SBarry Smith .keywords: MG, set, multigrid, residual, level 754b9ad928SBarry Smith 7654b2cd4bSJed Brown .seealso: PCMGResidualDefault() 774b9ad928SBarry Smith @*/ 787087cfbeSBarry Smith PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 794b9ad928SBarry Smith { 80f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 81f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 82298cc208SBarry Smith PetscErrorCode ierr; 834b9ad928SBarry Smith 844b9ad928SBarry Smith PetscFunctionBegin; 85c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 86ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 872fa5cd67SKarl Rupp if (residual) mglevels[l]->residual = residual; 8854b2cd4bSJed Brown if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault; 89f3ae41bdSBarry Smith if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 90298cc208SBarry Smith ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 91f3fbd535SBarry Smith mglevels[l]->A = mat; 924b9ad928SBarry Smith PetscFunctionReturn(0); 934b9ad928SBarry Smith } 944b9ad928SBarry Smith 954b9ad928SBarry Smith /*@ 96aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 97bf5b2e24SBarry Smith interpolation from l-1 to the lth level 984b9ad928SBarry Smith 99ad4df100SBarry Smith Logically Collective on PC and Mat 1004b9ad928SBarry Smith 1014b9ad928SBarry Smith Input Parameters: 1024b9ad928SBarry Smith + pc - the multigrid context 1034b9ad928SBarry Smith . mat - the interpolation operator 104bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 1054b9ad928SBarry Smith 1064b9ad928SBarry Smith Level: advanced 1074b9ad928SBarry Smith 1084b9ad928SBarry Smith Notes: 1094b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 1104b9ad928SBarry Smith for the same level. 1114b9ad928SBarry Smith 1124b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1134b9ad928SBarry Smith out from the matrix size which one it is. 1144b9ad928SBarry Smith 1154b9ad928SBarry Smith .keywords: multigrid, set, interpolate, level 1164b9ad928SBarry Smith 11797177400SBarry Smith .seealso: PCMGSetRestriction() 1184b9ad928SBarry Smith @*/ 1197087cfbeSBarry Smith PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 1204b9ad928SBarry Smith { 121f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 122f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 123fccaa45eSBarry Smith PetscErrorCode ierr; 1244b9ad928SBarry Smith 1254b9ad928SBarry Smith PetscFunctionBegin; 126c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 127ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 128ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 129c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1306bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr); 1312fa5cd67SKarl Rupp 132f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 1334b9ad928SBarry Smith PetscFunctionReturn(0); 1344b9ad928SBarry Smith } 1354b9ad928SBarry Smith 136c88c5224SJed Brown /*@ 137c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 138c88c5224SJed Brown interpolation from l-1 to the lth level 139c88c5224SJed Brown 140c88c5224SJed Brown Logically Collective on PC 141c88c5224SJed Brown 142c88c5224SJed Brown Input Parameters: 143c88c5224SJed Brown + pc - the multigrid context 144c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 145c88c5224SJed Brown 146c88c5224SJed Brown Output Parameter: 1473ad4599aSBarry Smith . mat - the interpolation matrix, can be NULL 148c88c5224SJed Brown 149c88c5224SJed Brown Level: advanced 150c88c5224SJed Brown 151c88c5224SJed Brown .keywords: MG, get, multigrid, interpolation, level 152c88c5224SJed Brown 153c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 154c88c5224SJed Brown @*/ 155c88c5224SJed Brown PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 156c88c5224SJed Brown { 157c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 158c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 159c88c5224SJed Brown PetscErrorCode ierr; 160c88c5224SJed Brown 161c88c5224SJed Brown PetscFunctionBegin; 162c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1633ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 164ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 165ce94432eSBarry 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); 166c88c5224SJed Brown if (!mglevels[l]->interpolate) { 1675aa31b60SBarry Smith if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()"); 168c88c5224SJed Brown ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 169c88c5224SJed Brown } 1703ad4599aSBarry Smith if (mat) *mat = mglevels[l]->interpolate; 171c88c5224SJed Brown PetscFunctionReturn(0); 172c88c5224SJed Brown } 173c88c5224SJed Brown 1744b9ad928SBarry Smith /*@ 17597177400SBarry Smith PCMGSetRestriction - Sets the function to be used to restrict vector 1764b9ad928SBarry Smith from level l to l-1. 1774b9ad928SBarry Smith 178ad4df100SBarry Smith Logically Collective on PC and Mat 1794b9ad928SBarry Smith 1804b9ad928SBarry Smith Input Parameters: 1814b9ad928SBarry Smith + pc - the multigrid context 182c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 183c88c5224SJed Brown - mat - the restriction matrix 1844b9ad928SBarry Smith 1854b9ad928SBarry Smith Level: advanced 1864b9ad928SBarry Smith 1874b9ad928SBarry Smith Notes: 1884b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 1894b9ad928SBarry Smith for the same level. 1904b9ad928SBarry Smith 1914b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 1924b9ad928SBarry Smith out from the matrix size which one it is. 1934b9ad928SBarry Smith 194aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 195fccaa45eSBarry Smith is used. 196fccaa45eSBarry Smith 1974b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level 1984b9ad928SBarry Smith 199aea2a34eSBarry Smith .seealso: PCMGSetInterpolation() 2004b9ad928SBarry Smith @*/ 2017087cfbeSBarry Smith PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 2024b9ad928SBarry Smith { 203fccaa45eSBarry Smith PetscErrorCode ierr; 204f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 205f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2064b9ad928SBarry Smith 2074b9ad928SBarry Smith PetscFunctionBegin; 208c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 209c88c5224SJed Brown PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 210ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 211ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 212c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 2136bf464f9SBarry Smith ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 2142fa5cd67SKarl Rupp 215f3fbd535SBarry Smith mglevels[l]->restrct = mat; 2164b9ad928SBarry Smith PetscFunctionReturn(0); 2174b9ad928SBarry Smith } 2184b9ad928SBarry Smith 219c88c5224SJed Brown /*@ 220c88c5224SJed Brown PCMGGetRestriction - Gets the function to be used to restrict vector 221c88c5224SJed Brown from level l to l-1. 222c88c5224SJed Brown 223c88c5224SJed Brown Logically Collective on PC and Mat 224c88c5224SJed Brown 225c88c5224SJed Brown Input Parameters: 226c88c5224SJed Brown + pc - the multigrid context 227c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 228c88c5224SJed Brown 229c88c5224SJed Brown Output Parameter: 230c88c5224SJed Brown . mat - the restriction matrix 231c88c5224SJed Brown 232c88c5224SJed Brown Level: advanced 233c88c5224SJed Brown 234c88c5224SJed Brown .keywords: MG, get, multigrid, restriction, level 235c88c5224SJed Brown 236c88c5224SJed Brown .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale() 237c88c5224SJed Brown @*/ 238c88c5224SJed Brown PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 239c88c5224SJed Brown { 240c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 241c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 242c88c5224SJed Brown PetscErrorCode ierr; 243c88c5224SJed Brown 244c88c5224SJed Brown PetscFunctionBegin; 245c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2463ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 247ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 248ce94432eSBarry 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); 249c88c5224SJed Brown if (!mglevels[l]->restrct) { 250ce94432eSBarry Smith if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 251c88c5224SJed Brown ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 252c88c5224SJed Brown } 2533ad4599aSBarry Smith if (mat) *mat = mglevels[l]->restrct; 254c88c5224SJed Brown PetscFunctionReturn(0); 255c88c5224SJed Brown } 256c88c5224SJed Brown 25773250ac0SBarry Smith /*@ 25873250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 25973250ac0SBarry Smith 260c88c5224SJed Brown Logically Collective on PC and Vec 261c88c5224SJed Brown 262c88c5224SJed Brown Input Parameters: 263c88c5224SJed Brown + pc - the multigrid context 264c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 265c88c5224SJed Brown . rscale - the scaling 266c88c5224SJed Brown 267c88c5224SJed Brown Level: advanced 268c88c5224SJed Brown 269c88c5224SJed Brown Notes: 270c88c5224SJed 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. 271c88c5224SJed Brown 272c88c5224SJed Brown .keywords: MG, set, multigrid, restriction, level 273c88c5224SJed Brown 274c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale() 275c88c5224SJed Brown @*/ 276c88c5224SJed Brown PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 277c88c5224SJed Brown { 278c88c5224SJed Brown PetscErrorCode ierr; 279c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 280c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 281c88c5224SJed Brown 282c88c5224SJed Brown PetscFunctionBegin; 283c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 284ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 285ce94432eSBarry 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); 286c88c5224SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 287c88c5224SJed Brown ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 2882fa5cd67SKarl Rupp 289c88c5224SJed Brown mglevels[l]->rscale = rscale; 290c88c5224SJed Brown PetscFunctionReturn(0); 291c88c5224SJed Brown } 292c88c5224SJed Brown 293c88c5224SJed Brown /*@ 294c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 295c88c5224SJed Brown 296c88c5224SJed Brown Collective on PC 29773250ac0SBarry Smith 29873250ac0SBarry Smith Input Parameters: 29973250ac0SBarry Smith + pc - the multigrid context 30073250ac0SBarry Smith . rscale - the scaling 30173250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 30273250ac0SBarry Smith 30373250ac0SBarry Smith Level: advanced 30473250ac0SBarry Smith 30573250ac0SBarry Smith Notes: 30673250ac0SBarry 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. 30773250ac0SBarry Smith 30873250ac0SBarry Smith .keywords: MG, set, multigrid, restriction, level 30973250ac0SBarry Smith 310c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGGetRestriction() 31173250ac0SBarry Smith @*/ 312c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 31373250ac0SBarry Smith { 31473250ac0SBarry Smith PetscErrorCode ierr; 31573250ac0SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 31673250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 31773250ac0SBarry Smith 31873250ac0SBarry Smith PetscFunctionBegin; 31973250ac0SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 320ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 321ce94432eSBarry 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); 322c88c5224SJed Brown if (!mglevels[l]->rscale) { 323c88c5224SJed Brown Mat R; 324c88c5224SJed Brown Vec X,Y,coarse,fine; 325c88c5224SJed Brown PetscInt M,N; 326c88c5224SJed Brown ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 3272a7a6963SBarry Smith ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr); 328c88c5224SJed Brown ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 3292fa5cd67SKarl Rupp if (M < N) { 3302fa5cd67SKarl Rupp fine = X; 3312fa5cd67SKarl Rupp coarse = Y; 3322fa5cd67SKarl Rupp } else if (N < M) { 3332fa5cd67SKarl Rupp fine = Y; coarse = X; 334ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 335c88c5224SJed Brown ierr = VecSet(fine,1.);CHKERRQ(ierr); 336c88c5224SJed Brown ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 337c88c5224SJed Brown ierr = VecDestroy(&fine);CHKERRQ(ierr); 338c88c5224SJed Brown ierr = VecReciprocal(coarse);CHKERRQ(ierr); 339c88c5224SJed Brown mglevels[l]->rscale = coarse; 340c88c5224SJed Brown } 341c88c5224SJed Brown *rscale = mglevels[l]->rscale; 34273250ac0SBarry Smith PetscFunctionReturn(0); 34373250ac0SBarry Smith } 34473250ac0SBarry Smith 345f39d8e23SSatish Balay /*@ 34697177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 34797177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 34897177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 3494b9ad928SBarry Smith post-smoothing. 3504b9ad928SBarry Smith 3514b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith Input Parameters: 3544b9ad928SBarry Smith + pc - the multigrid context 3554b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 3564b9ad928SBarry Smith 3574b9ad928SBarry Smith Ouput Parameters: 3584b9ad928SBarry Smith . ksp - the smoother 3594b9ad928SBarry Smith 36057420d5bSBarry Smith Notes: 36157420d5bSBarry 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. 36257420d5bSBarry Smith You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc) 36357420d5bSBarry Smith and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing. 36457420d5bSBarry Smith 3654b9ad928SBarry Smith Level: advanced 3664b9ad928SBarry Smith 3674b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 3684b9ad928SBarry Smith 36997177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 3704b9ad928SBarry Smith @*/ 3717087cfbeSBarry Smith PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 3724b9ad928SBarry Smith { 373f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 374f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 3754b9ad928SBarry Smith 3764b9ad928SBarry Smith PetscFunctionBegin; 377c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 378f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 3794b9ad928SBarry Smith PetscFunctionReturn(0); 3804b9ad928SBarry Smith } 3814b9ad928SBarry Smith 382f39d8e23SSatish Balay /*@ 38397177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 3844b9ad928SBarry Smith coarse grid correction (post-smoother). 3854b9ad928SBarry Smith 3864b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith Input Parameters: 3894b9ad928SBarry Smith + pc - the multigrid context 3904b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 3914b9ad928SBarry Smith 3924b9ad928SBarry Smith Ouput Parameters: 3934b9ad928SBarry Smith . ksp - the smoother 3944b9ad928SBarry Smith 3954b9ad928SBarry Smith Level: advanced 3964b9ad928SBarry Smith 39789cce641SBarry Smith Notes: calling this will result in a different pre and post smoother so you may need to 39889cce641SBarry Smith set options on the pre smoother also 39989cce641SBarry Smith 4004b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level 4014b9ad928SBarry Smith 40297177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 4034b9ad928SBarry Smith @*/ 4047087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 4054b9ad928SBarry Smith { 406f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 407f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 408dfbe8321SBarry Smith PetscErrorCode ierr; 409f69a0ea3SMatthew Knepley const char *prefix; 4104b9ad928SBarry Smith MPI_Comm comm; 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith PetscFunctionBegin; 413c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4144b9ad928SBarry Smith /* 4154b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 4164b9ad928SBarry Smith Thus we check if a different one has already been allocated, 4174b9ad928SBarry Smith if not we allocate it. 4184b9ad928SBarry Smith */ 419ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 420f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 42119fd82e9SBarry Smith KSPType ksptype; 42219fd82e9SBarry Smith PCType pctype; 423336babb1SJed Brown PC ipc; 424336babb1SJed Brown PetscReal rtol,abstol,dtol; 425336babb1SJed Brown PetscInt maxits; 426336babb1SJed Brown KSPNormType normtype; 427f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 428f3fbd535SBarry Smith ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 429336babb1SJed Brown ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 4303bf036e2SBarry Smith ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 431336babb1SJed Brown ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 432336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 433336babb1SJed Brown ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 434336babb1SJed Brown 435f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 436422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr); 437f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 438f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 439336babb1SJed Brown ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 440336babb1SJed Brown ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 441336babb1SJed Brown ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 4420059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 443336babb1SJed Brown ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 444336babb1SJed Brown ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 4453bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr); 446ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr); 4474b9ad928SBarry Smith } 448f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 4494b9ad928SBarry Smith PetscFunctionReturn(0); 4504b9ad928SBarry Smith } 4514b9ad928SBarry Smith 452f39d8e23SSatish Balay /*@ 45397177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 4544b9ad928SBarry Smith coarse grid correction (pre-smoother). 4554b9ad928SBarry Smith 4564b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 4574b9ad928SBarry Smith 4584b9ad928SBarry Smith Input Parameters: 4594b9ad928SBarry Smith + pc - the multigrid context 4604b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 4614b9ad928SBarry Smith 4624b9ad928SBarry Smith Ouput Parameters: 4634b9ad928SBarry Smith . ksp - the smoother 4644b9ad928SBarry Smith 4654b9ad928SBarry Smith Level: advanced 4664b9ad928SBarry Smith 46789cce641SBarry Smith Notes: calling this will result in a different pre and post smoother so you may need to 46889cce641SBarry Smith set options on the post smoother also 46989cce641SBarry Smith 4704b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 4714b9ad928SBarry Smith 47297177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 4734b9ad928SBarry Smith @*/ 4747087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 4754b9ad928SBarry Smith { 476dfbe8321SBarry Smith PetscErrorCode ierr; 477f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 478f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 4794b9ad928SBarry Smith 4804b9ad928SBarry Smith PetscFunctionBegin; 481c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4824b9ad928SBarry Smith /* make sure smoother up and down are different */ 483c5eb9154SBarry Smith if (l) { 4840298fd71SBarry Smith ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr); 485d8148a5aSMatthew Knepley } 486f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 4874b9ad928SBarry Smith PetscFunctionReturn(0); 4884b9ad928SBarry Smith } 4894b9ad928SBarry Smith 4904b9ad928SBarry Smith /*@ 491*cab31ae5SJed Brown PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 4924b9ad928SBarry Smith 493ad4df100SBarry Smith Logically Collective on PC 4944b9ad928SBarry Smith 4954b9ad928SBarry Smith Input Parameters: 4964b9ad928SBarry Smith + pc - the multigrid context 497c1cbb1deSBarry Smith . l - the level (0 is coarsest) 498c1cbb1deSBarry Smith - c - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 4994b9ad928SBarry Smith 5004b9ad928SBarry Smith Level: advanced 5014b9ad928SBarry Smith 5024b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 5034b9ad928SBarry Smith 504c1cbb1deSBarry Smith .seealso: PCMGSetCycleType() 5054b9ad928SBarry Smith @*/ 506c1cbb1deSBarry Smith PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c) 5074b9ad928SBarry Smith { 508f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 509f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5104b9ad928SBarry Smith 5114b9ad928SBarry Smith PetscFunctionBegin; 512c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 513ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 514c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,l,2); 515c51679f6SSatish Balay PetscValidLogicalCollectiveEnum(pc,c,3); 516f3fbd535SBarry Smith mglevels[l]->cycles = c; 5174b9ad928SBarry Smith PetscFunctionReturn(0); 5184b9ad928SBarry Smith } 5194b9ad928SBarry Smith 5204b9ad928SBarry Smith /*@ 52197177400SBarry Smith PCMGSetRhs - Sets the vector space to be used to store the right-hand side 522fccaa45eSBarry Smith on a particular level. 5234b9ad928SBarry Smith 524ad4df100SBarry Smith Logically Collective on PC and Vec 5254b9ad928SBarry Smith 5264b9ad928SBarry Smith Input Parameters: 5274b9ad928SBarry Smith + pc - the multigrid context 5284b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 5294b9ad928SBarry Smith - c - the space 5304b9ad928SBarry Smith 5314b9ad928SBarry Smith Level: advanced 5324b9ad928SBarry Smith 533fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 534fccaa45eSBarry Smith 535fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 536fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 537fccaa45eSBarry Smith 5384b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level 5394b9ad928SBarry Smith 54097177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 5414b9ad928SBarry Smith @*/ 5427087cfbeSBarry Smith PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 5434b9ad928SBarry Smith { 544fccaa45eSBarry Smith PetscErrorCode ierr; 545f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 546f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5474b9ad928SBarry Smith 5484b9ad928SBarry Smith PetscFunctionBegin; 549c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 550ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 551ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 552c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 5536bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 5542fa5cd67SKarl Rupp 555f3fbd535SBarry Smith mglevels[l]->b = c; 5564b9ad928SBarry Smith PetscFunctionReturn(0); 5574b9ad928SBarry Smith } 5584b9ad928SBarry Smith 5594b9ad928SBarry Smith /*@ 56097177400SBarry Smith PCMGSetX - Sets the vector space to be used to store the solution on a 561fccaa45eSBarry Smith particular level. 5624b9ad928SBarry Smith 563ad4df100SBarry Smith Logically Collective on PC and Vec 5644b9ad928SBarry Smith 5654b9ad928SBarry Smith Input Parameters: 5664b9ad928SBarry Smith + pc - the multigrid context 567251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 5684b9ad928SBarry Smith - c - the space 5694b9ad928SBarry Smith 5704b9ad928SBarry Smith Level: advanced 5714b9ad928SBarry Smith 572fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 573fccaa45eSBarry Smith 574fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 575fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 576fccaa45eSBarry Smith 5774b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level 5784b9ad928SBarry Smith 57997177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 5804b9ad928SBarry Smith @*/ 5817087cfbeSBarry Smith PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 5824b9ad928SBarry Smith { 583fccaa45eSBarry Smith PetscErrorCode ierr; 584f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 585f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5864b9ad928SBarry Smith 5874b9ad928SBarry Smith PetscFunctionBegin; 588c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 589ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 590ce94432eSBarry Smith if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 591c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 5926bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 5932fa5cd67SKarl Rupp 594f3fbd535SBarry Smith mglevels[l]->x = c; 5954b9ad928SBarry Smith PetscFunctionReturn(0); 5964b9ad928SBarry Smith } 5974b9ad928SBarry Smith 5984b9ad928SBarry Smith /*@ 59997177400SBarry Smith PCMGSetR - Sets the vector space to be used to store the residual on a 600fccaa45eSBarry Smith particular level. 6014b9ad928SBarry Smith 602ad4df100SBarry Smith Logically Collective on PC and Vec 6034b9ad928SBarry Smith 6044b9ad928SBarry Smith Input Parameters: 6054b9ad928SBarry Smith + pc - the multigrid context 6064b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 6074b9ad928SBarry Smith - c - the space 6084b9ad928SBarry Smith 6094b9ad928SBarry Smith Level: advanced 6104b9ad928SBarry Smith 611fccaa45eSBarry Smith Notes: If this is not provided PETSc will automatically generate one. 612fccaa45eSBarry Smith 613fccaa45eSBarry Smith You do not need to keep a reference to this vector if you do 614fccaa45eSBarry Smith not need it PCDestroy() will properly free it. 615fccaa45eSBarry Smith 6164b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level 6174b9ad928SBarry Smith @*/ 6187087cfbeSBarry Smith PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 6194b9ad928SBarry Smith { 620fccaa45eSBarry Smith PetscErrorCode ierr; 621f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 622f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6234b9ad928SBarry Smith 6244b9ad928SBarry Smith PetscFunctionBegin; 625c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 626ce94432eSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 627ce94432eSBarry Smith if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 628c3122656SLisandro Dalcin ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 6296bf464f9SBarry Smith ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 6302fa5cd67SKarl Rupp 631f3fbd535SBarry Smith mglevels[l]->r = c; 6324b9ad928SBarry Smith PetscFunctionReturn(0); 6334b9ad928SBarry Smith } 634