14b9ad928SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 34b9ad928SBarry Smith 4f9426fe0SMark Adams /* ---------------------------------------------------------------------------*/ 51f6cc5b2SSatish Balay /*@C 654b2cd4bSJed Brown PCMGResidualDefault - Default routine to calculate the residual. 74b9ad928SBarry Smith 8d083f849SBarry Smith Collective on Mat 94b9ad928SBarry Smith 104b9ad928SBarry Smith Input Parameters: 114b9ad928SBarry Smith + mat - the matrix 124b9ad928SBarry Smith . b - the right-hand-side 134b9ad928SBarry Smith - x - the approximate solution 144b9ad928SBarry Smith 154b9ad928SBarry Smith Output Parameter: 164b9ad928SBarry Smith . r - location to store the residual 174b9ad928SBarry Smith 18d0e4de75SBarry Smith Level: developer 194b9ad928SBarry Smith 2097177400SBarry Smith .seealso: PCMGSetResidual() 214b9ad928SBarry Smith @*/ 2254b2cd4bSJed Brown PetscErrorCode PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r) 234b9ad928SBarry Smith { 244b9ad928SBarry Smith PetscFunctionBegin; 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatResidual(mat,b,x,r)); 264b9ad928SBarry Smith PetscFunctionReturn(0); 274b9ad928SBarry Smith } 284b9ad928SBarry Smith 29fcb023d4SJed Brown /*@C 30fcb023d4SJed Brown PCMGResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system 31fcb023d4SJed Brown 32fcb023d4SJed Brown Collective on Mat 33fcb023d4SJed Brown 34fcb023d4SJed Brown Input Parameters: 35fcb023d4SJed Brown + mat - the matrix 36fcb023d4SJed Brown . b - the right-hand-side 37fcb023d4SJed Brown - x - the approximate solution 38fcb023d4SJed Brown 39fcb023d4SJed Brown Output Parameter: 40fcb023d4SJed Brown . r - location to store the residual 41fcb023d4SJed Brown 42fcb023d4SJed Brown Level: developer 43fcb023d4SJed Brown 44fcb023d4SJed Brown .seealso: PCMGSetResidualTranspose() 45fcb023d4SJed Brown @*/ 46fcb023d4SJed Brown PetscErrorCode PCMGResidualTransposeDefault(Mat mat,Vec b,Vec x,Vec r) 47fcb023d4SJed Brown { 48fcb023d4SJed Brown PetscFunctionBegin; 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(mat,x,r)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(r,-1.0,b)); 51fcb023d4SJed Brown PetscFunctionReturn(0); 52fcb023d4SJed Brown } 53fcb023d4SJed Brown 5430b0564aSStefano Zampini /*@C 5530b0564aSStefano Zampini PCMGMatResidualDefault - Default routine to calculate the residual. 5630b0564aSStefano Zampini 5730b0564aSStefano Zampini Collective on Mat 5830b0564aSStefano Zampini 5930b0564aSStefano Zampini Input Parameters: 6030b0564aSStefano Zampini + mat - the matrix 6130b0564aSStefano Zampini . b - the right-hand-side 6230b0564aSStefano Zampini - x - the approximate solution 6330b0564aSStefano Zampini 6430b0564aSStefano Zampini Output Parameter: 6530b0564aSStefano Zampini . r - location to store the residual 6630b0564aSStefano Zampini 6730b0564aSStefano Zampini Level: developer 6830b0564aSStefano Zampini 6930b0564aSStefano Zampini .seealso: PCMGSetMatResidual() 7030b0564aSStefano Zampini @*/ 7130b0564aSStefano Zampini PetscErrorCode PCMGMatResidualDefault(Mat mat,Mat b,Mat x,Mat r) 7230b0564aSStefano Zampini { 7330b0564aSStefano Zampini PetscFunctionBegin; 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(mat,x,MAT_REUSE_MATRIX,PETSC_DEFAULT,&r)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatAYPX(r,-1.0,b,UNKNOWN_NONZERO_PATTERN)); 7630b0564aSStefano Zampini PetscFunctionReturn(0); 7730b0564aSStefano Zampini } 7830b0564aSStefano Zampini 7930b0564aSStefano Zampini /*@C 8030b0564aSStefano Zampini PCMGMatResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system 8130b0564aSStefano Zampini 8230b0564aSStefano Zampini Collective on Mat 8330b0564aSStefano Zampini 8430b0564aSStefano Zampini Input Parameters: 8530b0564aSStefano Zampini + mat - the matrix 8630b0564aSStefano Zampini . b - the right-hand-side 8730b0564aSStefano Zampini - x - the approximate solution 8830b0564aSStefano Zampini 8930b0564aSStefano Zampini Output Parameter: 9030b0564aSStefano Zampini . r - location to store the residual 9130b0564aSStefano Zampini 9230b0564aSStefano Zampini Level: developer 9330b0564aSStefano Zampini 9430b0564aSStefano Zampini .seealso: PCMGSetMatResidualTranspose() 9530b0564aSStefano Zampini @*/ 9630b0564aSStefano Zampini PetscErrorCode PCMGMatResidualTransposeDefault(Mat mat,Mat b,Mat x,Mat r) 9730b0564aSStefano Zampini { 9830b0564aSStefano Zampini PetscFunctionBegin; 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(mat,x,MAT_REUSE_MATRIX,PETSC_DEFAULT,&r)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatAYPX(r,-1.0,b,UNKNOWN_NONZERO_PATTERN)); 10130b0564aSStefano Zampini PetscFunctionReturn(0); 10230b0564aSStefano Zampini } 103f39d8e23SSatish Balay /*@ 10497177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 1054b9ad928SBarry Smith 1064b9ad928SBarry Smith Not Collective 1074b9ad928SBarry Smith 1084b9ad928SBarry Smith Input Parameter: 1094b9ad928SBarry Smith . pc - the multigrid context 1104b9ad928SBarry Smith 1114b9ad928SBarry Smith Output Parameter: 1124b9ad928SBarry Smith . ksp - the coarse grid solver context 1134b9ad928SBarry Smith 1144b9ad928SBarry Smith Level: advanced 1154b9ad928SBarry Smith 11628668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetSmoother() 1174b9ad928SBarry Smith @*/ 1187087cfbeSBarry Smith PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 1194b9ad928SBarry Smith { 120f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 121f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1224b9ad928SBarry Smith 1234b9ad928SBarry Smith PetscFunctionBegin; 124c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 125f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 1264b9ad928SBarry Smith PetscFunctionReturn(0); 1274b9ad928SBarry Smith } 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith /*@C 13097177400SBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual 1314b9ad928SBarry Smith on the lth level. 1324b9ad928SBarry Smith 133d083f849SBarry Smith Logically Collective on PC 1344b9ad928SBarry Smith 1354b9ad928SBarry Smith Input Parameters: 1364b9ad928SBarry Smith + pc - the multigrid context 1374b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 138157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 139d0e4de75SBarry Smith previous one were provided then a default is used 1404b9ad928SBarry Smith - mat - matrix associated with residual 1414b9ad928SBarry Smith 1424b9ad928SBarry Smith Level: advanced 1434b9ad928SBarry Smith 14454b2cd4bSJed Brown .seealso: PCMGResidualDefault() 1454b9ad928SBarry Smith @*/ 1467087cfbeSBarry Smith PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 1474b9ad928SBarry Smith { 148f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 149f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1504b9ad928SBarry Smith 1514b9ad928SBarry Smith PetscFunctionBegin; 152c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 153*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1542fa5cd67SKarl Rupp if (residual) mglevels[l]->residual = residual; 15554b2cd4bSJed Brown if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault; 15630b0564aSStefano Zampini mglevels[l]->matresidual = PCMGMatResidualDefault; 1575f80ce2aSJacob Faibussowitsch if (mat) CHKERRQ(PetscObjectReference((PetscObject)mat)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mglevels[l]->A)); 159f3fbd535SBarry Smith mglevels[l]->A = mat; 1604b9ad928SBarry Smith PetscFunctionReturn(0); 1614b9ad928SBarry Smith } 1624b9ad928SBarry Smith 163fcb023d4SJed Brown /*@C 164fcb023d4SJed Brown PCMGSetResidualTranspose - Sets the function to be used to calculate the residual of the transposed linear system 165fcb023d4SJed Brown on the lth level. 166fcb023d4SJed Brown 167fcb023d4SJed Brown Logically Collective on PC 168fcb023d4SJed Brown 169fcb023d4SJed Brown Input Parameters: 170fcb023d4SJed Brown + pc - the multigrid context 171fcb023d4SJed Brown . l - the level (0 is coarsest) to supply 172fcb023d4SJed Brown . residualt - function used to form transpose of residual, if none is provided the previously provide one is used, if no 173fcb023d4SJed Brown previous one were provided then a default is used 174fcb023d4SJed Brown - mat - matrix associated with residual 175fcb023d4SJed Brown 176fcb023d4SJed Brown Level: advanced 177fcb023d4SJed Brown 178fcb023d4SJed Brown .seealso: PCMGResidualTransposeDefault() 179fcb023d4SJed Brown @*/ 180fcb023d4SJed Brown PetscErrorCode PCMGSetResidualTranspose(PC pc,PetscInt l,PetscErrorCode (*residualt)(Mat,Vec,Vec,Vec),Mat mat) 181fcb023d4SJed Brown { 182fcb023d4SJed Brown PC_MG *mg = (PC_MG*)pc->data; 183fcb023d4SJed Brown PC_MG_Levels **mglevels = mg->levels; 184fcb023d4SJed Brown 185fcb023d4SJed Brown PetscFunctionBegin; 186fcb023d4SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 187*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 188fcb023d4SJed Brown if (residualt) mglevels[l]->residualtranspose = residualt; 189fcb023d4SJed Brown if (!mglevels[l]->residualtranspose) mglevels[l]->residualtranspose = PCMGResidualTransposeDefault; 19030b0564aSStefano Zampini mglevels[l]->matresidualtranspose = PCMGMatResidualTransposeDefault; 1915f80ce2aSJacob Faibussowitsch if (mat) CHKERRQ(PetscObjectReference((PetscObject)mat)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mglevels[l]->A)); 193fcb023d4SJed Brown mglevels[l]->A = mat; 194fcb023d4SJed Brown PetscFunctionReturn(0); 195fcb023d4SJed Brown } 196fcb023d4SJed Brown 1974b9ad928SBarry Smith /*@ 198aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 199bf5b2e24SBarry Smith interpolation from l-1 to the lth level 2004b9ad928SBarry Smith 201d083f849SBarry Smith Logically Collective on PC 2024b9ad928SBarry Smith 2034b9ad928SBarry Smith Input Parameters: 2044b9ad928SBarry Smith + pc - the multigrid context 2054b9ad928SBarry Smith . mat - the interpolation operator 206bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 2074b9ad928SBarry Smith 2084b9ad928SBarry Smith Level: advanced 2094b9ad928SBarry Smith 2104b9ad928SBarry Smith Notes: 2114b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 2124b9ad928SBarry Smith for the same level. 2134b9ad928SBarry Smith 2144b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 2154b9ad928SBarry Smith out from the matrix size which one it is. 2164b9ad928SBarry Smith 21797177400SBarry Smith .seealso: PCMGSetRestriction() 2184b9ad928SBarry Smith @*/ 2197087cfbeSBarry Smith PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 2204b9ad928SBarry Smith { 221f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 222f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith PetscFunctionBegin; 225c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 226*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 227*28b400f6SJacob Faibussowitsch PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)mat)); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mglevels[l]->interpolate)); 2302fa5cd67SKarl Rupp 231f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 2324b9ad928SBarry Smith PetscFunctionReturn(0); 2334b9ad928SBarry Smith } 2344b9ad928SBarry Smith 2358a2c336bSFande Kong /*@ 23695750439SFande Kong PCMGSetOperators - Sets operator and preconditioning matrix for lth level 2378a2c336bSFande Kong 238d083f849SBarry Smith Logically Collective on PC 2398a2c336bSFande Kong 2408a2c336bSFande Kong Input Parameters: 2418a2c336bSFande Kong + pc - the multigrid context 2428a2c336bSFande Kong . Amat - the operator 2438a2c336bSFande Kong . pmat - the preconditioning operator 24495750439SFande Kong - l - the level (0 is the coarsest) to supply 2458a2c336bSFande Kong 2468a2c336bSFande Kong Level: advanced 2478a2c336bSFande Kong 2488a2c336bSFande Kong .keywords: multigrid, set, interpolate, level 2498a2c336bSFande Kong 2508a2c336bSFande Kong .seealso: PCMGSetRestriction(), PCMGSetInterpolation() 2518a2c336bSFande Kong @*/ 2528a2c336bSFande Kong PetscErrorCode PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat) 253360ee056SFande Kong { 254360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 255360ee056SFande Kong PC_MG_Levels **mglevels = mg->levels; 256360ee056SFande Kong 257360ee056SFande Kong PetscFunctionBegin; 258360ee056SFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2598a2c336bSFande Kong PetscValidHeaderSpecific(Amat,MAT_CLASSID,3); 2608a2c336bSFande Kong PetscValidHeaderSpecific(Pmat,MAT_CLASSID,4); 261*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(mglevels[l]->smoothd,Amat,Pmat)); 263360ee056SFande Kong PetscFunctionReturn(0); 264360ee056SFande Kong } 265360ee056SFande Kong 266c88c5224SJed Brown /*@ 267c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 268c88c5224SJed Brown interpolation from l-1 to the lth level 269c88c5224SJed Brown 270c88c5224SJed Brown Logically Collective on PC 271c88c5224SJed Brown 272c88c5224SJed Brown Input Parameters: 273c88c5224SJed Brown + pc - the multigrid context 274c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 275c88c5224SJed Brown 276c88c5224SJed Brown Output Parameter: 2773ad4599aSBarry Smith . mat - the interpolation matrix, can be NULL 278c88c5224SJed Brown 279c88c5224SJed Brown Level: advanced 280c88c5224SJed Brown 281c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 282c88c5224SJed Brown @*/ 283c88c5224SJed Brown PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 284c88c5224SJed Brown { 285c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 286c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 287c88c5224SJed Brown 288c88c5224SJed Brown PetscFunctionBegin; 289c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2903ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 291*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 2922c71b3e2SJacob Faibussowitsch PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 293c88c5224SJed Brown if (!mglevels[l]->interpolate) { 2942c71b3e2SJacob Faibussowitsch PetscCheckFalse(!mglevels[l]->restrct,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()"); 2955f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetInterpolation(pc,l,mglevels[l]->restrct)); 296c88c5224SJed Brown } 2973ad4599aSBarry Smith if (mat) *mat = mglevels[l]->interpolate; 298c88c5224SJed Brown PetscFunctionReturn(0); 299c88c5224SJed Brown } 300c88c5224SJed Brown 3018a2c336bSFande Kong /*@ 302eab52d2bSLawrence Mitchell PCMGSetRestriction - Sets the function to be used to restrict dual vectors 3034b9ad928SBarry Smith from level l to l-1. 3044b9ad928SBarry Smith 305d083f849SBarry Smith Logically Collective on PC 3064b9ad928SBarry Smith 3074b9ad928SBarry Smith Input Parameters: 3084b9ad928SBarry Smith + pc - the multigrid context 309c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 310c88c5224SJed Brown - mat - the restriction matrix 3114b9ad928SBarry Smith 3124b9ad928SBarry Smith Level: advanced 3134b9ad928SBarry Smith 3144b9ad928SBarry Smith Notes: 3154b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 3164b9ad928SBarry Smith for the same level. 3174b9ad928SBarry Smith 3184b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 3194b9ad928SBarry Smith out from the matrix size which one it is. 3204b9ad928SBarry Smith 321aea2a34eSBarry Smith If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 322fccaa45eSBarry Smith is used. 323fccaa45eSBarry Smith 324a4355517SPierre Jolivet .seealso: PCMGSetInterpolation() 3254b9ad928SBarry Smith @*/ 3267087cfbeSBarry Smith PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 3274b9ad928SBarry Smith { 328f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 329f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 3304b9ad928SBarry Smith 3314b9ad928SBarry Smith PetscFunctionBegin; 332c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 333c88c5224SJed Brown PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 334*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 335*28b400f6SJacob Faibussowitsch PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 3365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)mat)); 3375f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mglevels[l]->restrct)); 3382fa5cd67SKarl Rupp 339f3fbd535SBarry Smith mglevels[l]->restrct = mat; 3404b9ad928SBarry Smith PetscFunctionReturn(0); 3414b9ad928SBarry Smith } 3424b9ad928SBarry Smith 343c88c5224SJed Brown /*@ 344eab52d2bSLawrence Mitchell PCMGGetRestriction - Gets the function to be used to restrict dual vectors 345c88c5224SJed Brown from level l to l-1. 346c88c5224SJed Brown 347d083f849SBarry Smith Logically Collective on PC 348c88c5224SJed Brown 349c88c5224SJed Brown Input Parameters: 350c88c5224SJed Brown + pc - the multigrid context 351c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 352c88c5224SJed Brown 353c88c5224SJed Brown Output Parameter: 354c88c5224SJed Brown . mat - the restriction matrix 355c88c5224SJed Brown 356c88c5224SJed Brown Level: advanced 357c88c5224SJed Brown 358eab52d2bSLawrence Mitchell .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection() 359c88c5224SJed Brown @*/ 360c88c5224SJed Brown PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 361c88c5224SJed Brown { 362c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 363c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 364c88c5224SJed Brown 365c88c5224SJed Brown PetscFunctionBegin; 366c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3673ad4599aSBarry Smith if (mat) PetscValidPointer(mat,3); 368*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 3692c71b3e2SJacob Faibussowitsch PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 370c88c5224SJed Brown if (!mglevels[l]->restrct) { 3712c71b3e2SJacob Faibussowitsch PetscCheckFalse(!mglevels[l]->interpolate,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 3725f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGSetRestriction(pc,l,mglevels[l]->interpolate)); 373c88c5224SJed Brown } 3743ad4599aSBarry Smith if (mat) *mat = mglevels[l]->restrct; 375c88c5224SJed Brown PetscFunctionReturn(0); 376c88c5224SJed Brown } 377c88c5224SJed Brown 37873250ac0SBarry Smith /*@ 37973250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 38073250ac0SBarry Smith 381d083f849SBarry Smith Logically Collective on PC 382c88c5224SJed Brown 383c88c5224SJed Brown Input Parameters: 384c88c5224SJed Brown + pc - the multigrid context 385c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 386c88c5224SJed Brown . rscale - the scaling 387c88c5224SJed Brown 388c88c5224SJed Brown Level: advanced 389c88c5224SJed Brown 390c88c5224SJed Brown Notes: 391eab52d2bSLawrence Mitchell When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R. It is preferable to use PCMGSetInjection() to control moving primal vectors. 392c88c5224SJed Brown 393eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection() 394c88c5224SJed Brown @*/ 395c88c5224SJed Brown PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 396c88c5224SJed Brown { 397c88c5224SJed Brown PC_MG *mg = (PC_MG*)pc->data; 398c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 399c88c5224SJed Brown 400c88c5224SJed Brown PetscFunctionBegin; 401c88c5224SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 402*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 4032c71b3e2SJacob Faibussowitsch PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 4045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)rscale)); 4055f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&mglevels[l]->rscale)); 4062fa5cd67SKarl Rupp 407c88c5224SJed Brown mglevels[l]->rscale = rscale; 408c88c5224SJed Brown PetscFunctionReturn(0); 409c88c5224SJed Brown } 410c88c5224SJed Brown 411c88c5224SJed Brown /*@ 412c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 413c88c5224SJed Brown 414c88c5224SJed Brown Collective on PC 41573250ac0SBarry Smith 41673250ac0SBarry Smith Input Parameters: 41773250ac0SBarry Smith + pc - the multigrid context 41873250ac0SBarry Smith . rscale - the scaling 41973250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 42073250ac0SBarry Smith 42173250ac0SBarry Smith Level: advanced 42273250ac0SBarry Smith 42373250ac0SBarry Smith Notes: 424eab52d2bSLawrence Mitchell When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R. It is preferable to use PCMGGetInjection() to control moving primal vectors. 42573250ac0SBarry Smith 426eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection() 42773250ac0SBarry Smith @*/ 428c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 42973250ac0SBarry Smith { 43073250ac0SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 43173250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 43273250ac0SBarry Smith 43373250ac0SBarry Smith PetscFunctionBegin; 43473250ac0SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 435*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 4362c71b3e2SJacob Faibussowitsch PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 437c88c5224SJed Brown if (!mglevels[l]->rscale) { 438c88c5224SJed Brown Mat R; 439c88c5224SJed Brown Vec X,Y,coarse,fine; 440c88c5224SJed Brown PetscInt M,N; 4415f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGGetRestriction(pc,l,&R)); 4425f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(R,&X,&Y)); 4435f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(R,&M,&N)); 4442fa5cd67SKarl Rupp if (M < N) { 4452fa5cd67SKarl Rupp fine = X; 4462fa5cd67SKarl Rupp coarse = Y; 4472fa5cd67SKarl Rupp } else if (N < M) { 4482fa5cd67SKarl Rupp fine = Y; coarse = X; 449ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 4505f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(fine,1.)); 4515f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestrict(R,fine,coarse)); 4525f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&fine)); 4535f80ce2aSJacob Faibussowitsch CHKERRQ(VecReciprocal(coarse)); 454c88c5224SJed Brown mglevels[l]->rscale = coarse; 455c88c5224SJed Brown } 456c88c5224SJed Brown *rscale = mglevels[l]->rscale; 45773250ac0SBarry Smith PetscFunctionReturn(0); 45873250ac0SBarry Smith } 45973250ac0SBarry Smith 460f39d8e23SSatish Balay /*@ 461eab52d2bSLawrence Mitchell PCMGSetInjection - Sets the function to be used to inject primal vectors 462eab52d2bSLawrence Mitchell from level l to l-1. 463eab52d2bSLawrence Mitchell 464d083f849SBarry Smith Logically Collective on PC 465eab52d2bSLawrence Mitchell 466eab52d2bSLawrence Mitchell Input Parameters: 467eab52d2bSLawrence Mitchell + pc - the multigrid context 468eab52d2bSLawrence Mitchell . l - the level (0 is coarsest) to supply [Do not supply 0] 469eab52d2bSLawrence Mitchell - mat - the injection matrix 470eab52d2bSLawrence Mitchell 471eab52d2bSLawrence Mitchell Level: advanced 472eab52d2bSLawrence Mitchell 473eab52d2bSLawrence Mitchell .seealso: PCMGSetRestriction() 474eab52d2bSLawrence Mitchell @*/ 475eab52d2bSLawrence Mitchell PetscErrorCode PCMGSetInjection(PC pc,PetscInt l,Mat mat) 476eab52d2bSLawrence Mitchell { 477eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 478eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 479eab52d2bSLawrence Mitchell 480eab52d2bSLawrence Mitchell PetscFunctionBegin; 481eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 482eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 483*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 484*28b400f6SJacob Faibussowitsch PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 4855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)mat)); 4865f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mglevels[l]->inject)); 487eab52d2bSLawrence Mitchell 488eab52d2bSLawrence Mitchell mglevels[l]->inject = mat; 489eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 490eab52d2bSLawrence Mitchell } 491eab52d2bSLawrence Mitchell 492eab52d2bSLawrence Mitchell /*@ 493eab52d2bSLawrence Mitchell PCMGGetInjection - Gets the function to be used to inject primal vectors 494eab52d2bSLawrence Mitchell from level l to l-1. 495eab52d2bSLawrence Mitchell 496d083f849SBarry Smith Logically Collective on PC 497eab52d2bSLawrence Mitchell 498eab52d2bSLawrence Mitchell Input Parameters: 499eab52d2bSLawrence Mitchell + pc - the multigrid context 500eab52d2bSLawrence Mitchell - l - the level (0 is coarsest) to supply [Do not supply 0] 501eab52d2bSLawrence Mitchell 502eab52d2bSLawrence Mitchell Output Parameter: 50399a38656SLawrence Mitchell . mat - the restriction matrix (may be NULL if no injection is available). 504eab52d2bSLawrence Mitchell 505eab52d2bSLawrence Mitchell Level: advanced 506eab52d2bSLawrence Mitchell 507eab52d2bSLawrence Mitchell .seealso: PCMGSetInjection(), PCMGetGetRestriction() 508eab52d2bSLawrence Mitchell @*/ 509eab52d2bSLawrence Mitchell PetscErrorCode PCMGGetInjection(PC pc,PetscInt l,Mat *mat) 510eab52d2bSLawrence Mitchell { 511eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG*)pc->data; 512eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 513eab52d2bSLawrence Mitchell 514eab52d2bSLawrence Mitchell PetscFunctionBegin; 515eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc,PC_CLASSID,1); 516eab52d2bSLawrence Mitchell if (mat) PetscValidPointer(mat,3); 517*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 5182c71b3e2SJacob Faibussowitsch PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 519eab52d2bSLawrence Mitchell if (mat) *mat = mglevels[l]->inject; 520eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 521eab52d2bSLawrence Mitchell } 522eab52d2bSLawrence Mitchell 523eab52d2bSLawrence Mitchell /*@ 52497177400SBarry Smith PCMGGetSmoother - Gets the KSP context to be used as smoother for 52597177400SBarry Smith both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 52697177400SBarry Smith PCMGGetSmootherDown() to use different functions for pre- and 5274b9ad928SBarry Smith post-smoothing. 5284b9ad928SBarry Smith 5294b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 5304b9ad928SBarry Smith 5314b9ad928SBarry Smith Input Parameters: 5324b9ad928SBarry Smith + pc - the multigrid context 5334b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5344b9ad928SBarry Smith 53501d2d390SJose E. Roman Output Parameter: 5364b9ad928SBarry Smith . ksp - the smoother 5374b9ad928SBarry Smith 53857420d5bSBarry Smith Notes: 53957420d5bSBarry 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. 54057420d5bSBarry Smith You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc) 54157420d5bSBarry Smith and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing. 54257420d5bSBarry Smith 5434b9ad928SBarry Smith Level: advanced 5444b9ad928SBarry Smith 54528668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve() 5464b9ad928SBarry Smith @*/ 5477087cfbeSBarry Smith PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 5484b9ad928SBarry Smith { 549f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 550f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5514b9ad928SBarry Smith 5524b9ad928SBarry Smith PetscFunctionBegin; 553c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 554f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5554b9ad928SBarry Smith PetscFunctionReturn(0); 5564b9ad928SBarry Smith } 5574b9ad928SBarry Smith 558f39d8e23SSatish Balay /*@ 55997177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 5604b9ad928SBarry Smith coarse grid correction (post-smoother). 5614b9ad928SBarry Smith 5624b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 5634b9ad928SBarry Smith 5644b9ad928SBarry Smith Input Parameters: 5654b9ad928SBarry Smith + pc - the multigrid context 5664b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5674b9ad928SBarry Smith 56801d2d390SJose E. Roman Output Parameter: 5694b9ad928SBarry Smith . ksp - the smoother 5704b9ad928SBarry Smith 5714b9ad928SBarry Smith Level: advanced 5724b9ad928SBarry Smith 57395452b02SPatrick Sanan Notes: 57495452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 57589cce641SBarry Smith set options on the pre smoother also 57689cce641SBarry Smith 57797177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 5784b9ad928SBarry Smith @*/ 5797087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 5804b9ad928SBarry Smith { 581f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 582f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 583f69a0ea3SMatthew Knepley const char *prefix; 5844b9ad928SBarry Smith MPI_Comm comm; 5854b9ad928SBarry Smith 5864b9ad928SBarry Smith PetscFunctionBegin; 587c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5884b9ad928SBarry Smith /* 5894b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 5904b9ad928SBarry Smith Thus we check if a different one has already been allocated, 5914b9ad928SBarry Smith if not we allocate it. 5924b9ad928SBarry Smith */ 593*28b400f6SJacob Faibussowitsch PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 594f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 59519fd82e9SBarry Smith KSPType ksptype; 59619fd82e9SBarry Smith PCType pctype; 597336babb1SJed Brown PC ipc; 598336babb1SJed Brown PetscReal rtol,abstol,dtol; 599336babb1SJed Brown PetscInt maxits; 600336babb1SJed Brown KSPNormType normtype; 6015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm)); 6025f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix)); 6035f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits)); 6045f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetType(mglevels[l]->smoothd,&ksptype)); 6055f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetNormType(mglevels[l]->smoothd,&normtype)); 6065f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(mglevels[l]->smoothd,&ipc)); 6075f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetType(ipc,&pctype)); 608336babb1SJed Brown 6095f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreate(comm,&mglevels[l]->smoothu)); 6105f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure)); 6115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l)); 6125f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix)); 6135f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits)); 6145f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(mglevels[l]->smoothu,ksptype)); 6155f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetNormType(mglevels[l]->smoothu,normtype)); 6165f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL)); 6175f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(mglevels[l]->smoothu,&ipc)); 6185f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(ipc,pctype)); 6195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu)); 6205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level)); 6214b9ad928SBarry Smith } 622f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 6234b9ad928SBarry Smith PetscFunctionReturn(0); 6244b9ad928SBarry Smith } 6254b9ad928SBarry Smith 626f39d8e23SSatish Balay /*@ 62797177400SBarry Smith PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 6284b9ad928SBarry Smith coarse grid correction (pre-smoother). 6294b9ad928SBarry Smith 6304b9ad928SBarry Smith Not Collective, KSP returned is parallel if PC is 6314b9ad928SBarry Smith 6324b9ad928SBarry Smith Input Parameters: 6334b9ad928SBarry Smith + pc - the multigrid context 6344b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 6354b9ad928SBarry Smith 63601d2d390SJose E. Roman Output Parameter: 6374b9ad928SBarry Smith . ksp - the smoother 6384b9ad928SBarry Smith 6394b9ad928SBarry Smith Level: advanced 6404b9ad928SBarry Smith 64195452b02SPatrick Sanan Notes: 64295452b02SPatrick Sanan calling this will result in a different pre and post smoother so you may need to 64389cce641SBarry Smith set options on the post smoother also 64489cce641SBarry Smith 64597177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 6464b9ad928SBarry Smith @*/ 6477087cfbeSBarry Smith PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 6484b9ad928SBarry Smith { 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); 6544b9ad928SBarry Smith /* make sure smoother up and down are different */ 655c5eb9154SBarry Smith if (l) { 6565f80ce2aSJacob Faibussowitsch CHKERRQ(PCMGGetSmootherUp(pc,l,NULL)); 657d8148a5aSMatthew Knepley } 658f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 6594b9ad928SBarry Smith PetscFunctionReturn(0); 6604b9ad928SBarry Smith } 6614b9ad928SBarry Smith 6624b9ad928SBarry Smith /*@ 663cab31ae5SJed Brown PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 6644b9ad928SBarry Smith 665ad4df100SBarry Smith Logically Collective on PC 6664b9ad928SBarry Smith 6674b9ad928SBarry Smith Input Parameters: 6684b9ad928SBarry Smith + pc - the multigrid context 669c1cbb1deSBarry Smith . l - the level (0 is coarsest) 670c1cbb1deSBarry Smith - c - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 6714b9ad928SBarry Smith 6724b9ad928SBarry Smith Level: advanced 6734b9ad928SBarry Smith 674c1cbb1deSBarry Smith .seealso: PCMGSetCycleType() 6754b9ad928SBarry Smith @*/ 676c1cbb1deSBarry Smith PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c) 6774b9ad928SBarry Smith { 678f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 679f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6804b9ad928SBarry Smith 6814b9ad928SBarry Smith PetscFunctionBegin; 682c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 683*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 684c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,l,2); 685c51679f6SSatish Balay PetscValidLogicalCollectiveEnum(pc,c,3); 686f3fbd535SBarry Smith mglevels[l]->cycles = c; 6874b9ad928SBarry Smith PetscFunctionReturn(0); 6884b9ad928SBarry Smith } 6894b9ad928SBarry Smith 6904b9ad928SBarry Smith /*@ 691f3b08a26SMatthew G. Knepley PCMGSetRhs - Sets the vector to be used to store the right-hand side on a particular level. 6924b9ad928SBarry Smith 693d083f849SBarry Smith Logically Collective on PC 6944b9ad928SBarry Smith 6954b9ad928SBarry Smith Input Parameters: 6964b9ad928SBarry Smith + pc - the multigrid context 6974b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 698f3b08a26SMatthew G. Knepley - c - the Vec 6994b9ad928SBarry Smith 7004b9ad928SBarry Smith Level: advanced 7014b9ad928SBarry Smith 70295452b02SPatrick Sanan Notes: 703f3b08a26SMatthew G. Knepley If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it. 704fccaa45eSBarry Smith 705f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, right-hand-side, rhs, level 70697177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR() 7074b9ad928SBarry Smith @*/ 7087087cfbeSBarry Smith PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 7094b9ad928SBarry Smith { 710f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 711f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7124b9ad928SBarry Smith 7134b9ad928SBarry Smith PetscFunctionBegin; 714c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 715*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 7162c71b3e2SJacob Faibussowitsch PetscCheckFalse(l == mglevels[0]->levels-1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 7175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)c)); 7185f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&mglevels[l]->b)); 7192fa5cd67SKarl Rupp 720f3fbd535SBarry Smith mglevels[l]->b = c; 7214b9ad928SBarry Smith PetscFunctionReturn(0); 7224b9ad928SBarry Smith } 7234b9ad928SBarry Smith 7244b9ad928SBarry Smith /*@ 725f3b08a26SMatthew G. Knepley PCMGSetX - Sets the vector to be used to store the solution on a particular level. 7264b9ad928SBarry Smith 727d083f849SBarry Smith Logically Collective on PC 7284b9ad928SBarry Smith 7294b9ad928SBarry Smith Input Parameters: 7304b9ad928SBarry Smith + pc - the multigrid context 731251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 732f3b08a26SMatthew G. Knepley - c - the Vec 7334b9ad928SBarry Smith 7344b9ad928SBarry Smith Level: advanced 7354b9ad928SBarry Smith 73695452b02SPatrick Sanan Notes: 737f3b08a26SMatthew G. Knepley If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it. 738fccaa45eSBarry Smith 739f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, solution, level 74097177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR() 7414b9ad928SBarry Smith @*/ 7427087cfbeSBarry Smith PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 7434b9ad928SBarry Smith { 744f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 745f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7464b9ad928SBarry Smith 7474b9ad928SBarry Smith PetscFunctionBegin; 748c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 749*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 7502c71b3e2SJacob Faibussowitsch PetscCheckFalse(l == mglevels[0]->levels-1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 7515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)c)); 7525f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&mglevels[l]->x)); 7532fa5cd67SKarl Rupp 754f3fbd535SBarry Smith mglevels[l]->x = c; 7554b9ad928SBarry Smith PetscFunctionReturn(0); 7564b9ad928SBarry Smith } 7574b9ad928SBarry Smith 7584b9ad928SBarry Smith /*@ 759f3b08a26SMatthew G. Knepley PCMGSetR - Sets the vector to be used to store the residual on a particular level. 7604b9ad928SBarry Smith 761d083f849SBarry Smith Logically Collective on PC 7624b9ad928SBarry Smith 7634b9ad928SBarry Smith Input Parameters: 7644b9ad928SBarry Smith + pc - the multigrid context 7654b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 766f3b08a26SMatthew G. Knepley - c - the Vec 7674b9ad928SBarry Smith 7684b9ad928SBarry Smith Level: advanced 7694b9ad928SBarry Smith 77095452b02SPatrick Sanan Notes: 771f3b08a26SMatthew G. Knepley If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it. 772fccaa45eSBarry Smith 773f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, residual, level 774f3b08a26SMatthew G. Knepley .seealso: PCMGSetRhs(), PCMGSetX() 7754b9ad928SBarry Smith @*/ 7767087cfbeSBarry Smith PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 7774b9ad928SBarry Smith { 778f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 779f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7804b9ad928SBarry Smith 7814b9ad928SBarry Smith PetscFunctionBegin; 782c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 783*28b400f6SJacob Faibussowitsch PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 784*28b400f6SJacob Faibussowitsch PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 7855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)c)); 7865f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&mglevels[l]->r)); 7872fa5cd67SKarl Rupp 788f3fbd535SBarry Smith mglevels[l]->r = c; 7894b9ad928SBarry Smith PetscFunctionReturn(0); 7904b9ad928SBarry Smith } 791