14b9ad928SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 34b9ad928SBarry Smith 41f6cc5b2SSatish Balay /*@C 554b2cd4bSJed Brown PCMGResidualDefault - Default routine to calculate the residual. 64b9ad928SBarry Smith 7*f1580f4eSBarry Smith Collective on mat 84b9ad928SBarry Smith 94b9ad928SBarry Smith Input Parameters: 104b9ad928SBarry Smith + mat - the matrix 114b9ad928SBarry Smith . b - the right-hand-side 124b9ad928SBarry Smith - x - the approximate solution 134b9ad928SBarry Smith 144b9ad928SBarry Smith Output Parameter: 154b9ad928SBarry Smith . r - location to store the residual 164b9ad928SBarry Smith 17d0e4de75SBarry Smith Level: developer 184b9ad928SBarry Smith 19*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetResidual()`, `PCMGSetMatResidual()` 204b9ad928SBarry Smith @*/ 219371c9d4SSatish Balay PetscErrorCode PCMGResidualDefault(Mat mat, Vec b, Vec x, Vec r) { 224b9ad928SBarry Smith PetscFunctionBegin; 239566063dSJacob Faibussowitsch PetscCall(MatResidual(mat, b, x, r)); 244b9ad928SBarry Smith PetscFunctionReturn(0); 254b9ad928SBarry Smith } 264b9ad928SBarry Smith 27fcb023d4SJed Brown /*@C 28fcb023d4SJed Brown PCMGResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system 29fcb023d4SJed Brown 30*f1580f4eSBarry Smith Collective on mat 31fcb023d4SJed Brown 32fcb023d4SJed Brown Input Parameters: 33fcb023d4SJed Brown + mat - the matrix 34fcb023d4SJed Brown . b - the right-hand-side 35fcb023d4SJed Brown - x - the approximate solution 36fcb023d4SJed Brown 37fcb023d4SJed Brown Output Parameter: 38fcb023d4SJed Brown . r - location to store the residual 39fcb023d4SJed Brown 40fcb023d4SJed Brown Level: developer 41fcb023d4SJed Brown 42*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetResidualTranspose()`, `PCMGMatResidualTransposeDefault()` 43fcb023d4SJed Brown @*/ 449371c9d4SSatish Balay PetscErrorCode PCMGResidualTransposeDefault(Mat mat, Vec b, Vec x, Vec r) { 45fcb023d4SJed Brown PetscFunctionBegin; 469566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat, x, r)); 479566063dSJacob Faibussowitsch PetscCall(VecAYPX(r, -1.0, b)); 48fcb023d4SJed Brown PetscFunctionReturn(0); 49fcb023d4SJed Brown } 50fcb023d4SJed Brown 5130b0564aSStefano Zampini /*@C 5230b0564aSStefano Zampini PCMGMatResidualDefault - Default routine to calculate the residual. 5330b0564aSStefano Zampini 54*f1580f4eSBarry Smith Collective on mat 5530b0564aSStefano Zampini 5630b0564aSStefano Zampini Input Parameters: 5730b0564aSStefano Zampini + mat - the matrix 5830b0564aSStefano Zampini . b - the right-hand-side 5930b0564aSStefano Zampini - x - the approximate solution 6030b0564aSStefano Zampini 6130b0564aSStefano Zampini Output Parameter: 6230b0564aSStefano Zampini . r - location to store the residual 6330b0564aSStefano Zampini 6430b0564aSStefano Zampini Level: developer 6530b0564aSStefano Zampini 66*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetMatResidual()`, `PCMGResidualDefault()` 6730b0564aSStefano Zampini @*/ 689371c9d4SSatish Balay PetscErrorCode PCMGMatResidualDefault(Mat mat, Mat b, Mat x, Mat r) { 6930b0564aSStefano Zampini PetscFunctionBegin; 709566063dSJacob Faibussowitsch PetscCall(MatMatMult(mat, x, MAT_REUSE_MATRIX, PETSC_DEFAULT, &r)); 719566063dSJacob Faibussowitsch PetscCall(MatAYPX(r, -1.0, b, UNKNOWN_NONZERO_PATTERN)); 7230b0564aSStefano Zampini PetscFunctionReturn(0); 7330b0564aSStefano Zampini } 7430b0564aSStefano Zampini 7530b0564aSStefano Zampini /*@C 7630b0564aSStefano Zampini PCMGMatResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system 7730b0564aSStefano Zampini 78*f1580f4eSBarry Smith Collective on mat 7930b0564aSStefano Zampini 8030b0564aSStefano Zampini Input Parameters: 8130b0564aSStefano Zampini + mat - the matrix 8230b0564aSStefano Zampini . b - the right-hand-side 8330b0564aSStefano Zampini - x - the approximate solution 8430b0564aSStefano Zampini 8530b0564aSStefano Zampini Output Parameter: 8630b0564aSStefano Zampini . r - location to store the residual 8730b0564aSStefano Zampini 8830b0564aSStefano Zampini Level: developer 8930b0564aSStefano Zampini 90*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetMatResidualTranspose()` 9130b0564aSStefano Zampini @*/ 929371c9d4SSatish Balay PetscErrorCode PCMGMatResidualTransposeDefault(Mat mat, Mat b, Mat x, Mat r) { 9330b0564aSStefano Zampini PetscFunctionBegin; 949566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(mat, x, MAT_REUSE_MATRIX, PETSC_DEFAULT, &r)); 959566063dSJacob Faibussowitsch PetscCall(MatAYPX(r, -1.0, b, UNKNOWN_NONZERO_PATTERN)); 9630b0564aSStefano Zampini PetscFunctionReturn(0); 9730b0564aSStefano Zampini } 98f39d8e23SSatish Balay /*@ 9997177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 1004b9ad928SBarry Smith 1014b9ad928SBarry Smith Not Collective 1024b9ad928SBarry Smith 1034b9ad928SBarry Smith Input Parameter: 1044b9ad928SBarry Smith . pc - the multigrid context 1054b9ad928SBarry Smith 1064b9ad928SBarry Smith Output Parameter: 1074b9ad928SBarry Smith . ksp - the coarse grid solver context 1084b9ad928SBarry Smith 1094b9ad928SBarry Smith Level: advanced 1104b9ad928SBarry Smith 111*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, `PCMGGetSmoother()` 1124b9ad928SBarry Smith @*/ 1139371c9d4SSatish Balay PetscErrorCode PCMGGetCoarseSolve(PC pc, KSP *ksp) { 114f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 115f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1164b9ad928SBarry Smith 1174b9ad928SBarry Smith PetscFunctionBegin; 118c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 119f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 1204b9ad928SBarry Smith PetscFunctionReturn(0); 1214b9ad928SBarry Smith } 1224b9ad928SBarry Smith 1234b9ad928SBarry Smith /*@C 124*f1580f4eSBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual on the lth level. 1254b9ad928SBarry Smith 126*f1580f4eSBarry Smith Logically Collective on pc 1274b9ad928SBarry Smith 1284b9ad928SBarry Smith Input Parameters: 1294b9ad928SBarry Smith + pc - the multigrid context 1304b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 131157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 132d0e4de75SBarry Smith previous one were provided then a default is used 1334b9ad928SBarry Smith - mat - matrix associated with residual 1344b9ad928SBarry Smith 1354b9ad928SBarry Smith Level: advanced 1364b9ad928SBarry Smith 137*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGResidualDefault()` 1384b9ad928SBarry Smith @*/ 1399371c9d4SSatish Balay PetscErrorCode PCMGSetResidual(PC pc, PetscInt l, PetscErrorCode (*residual)(Mat, Vec, Vec, Vec), Mat mat) { 140f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 141f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1424b9ad928SBarry Smith 1434b9ad928SBarry Smith PetscFunctionBegin; 144c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14528b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 1462fa5cd67SKarl Rupp if (residual) mglevels[l]->residual = residual; 14754b2cd4bSJed Brown if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault; 14830b0564aSStefano Zampini mglevels[l]->matresidual = PCMGMatResidualDefault; 1499566063dSJacob Faibussowitsch if (mat) PetscCall(PetscObjectReference((PetscObject)mat)); 1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->A)); 151f3fbd535SBarry Smith mglevels[l]->A = mat; 1524b9ad928SBarry Smith PetscFunctionReturn(0); 1534b9ad928SBarry Smith } 1544b9ad928SBarry Smith 155fcb023d4SJed Brown /*@C 156fcb023d4SJed Brown PCMGSetResidualTranspose - Sets the function to be used to calculate the residual of the transposed linear system 157fcb023d4SJed Brown on the lth level. 158fcb023d4SJed Brown 159*f1580f4eSBarry Smith Logically Collective on pc 160fcb023d4SJed Brown 161fcb023d4SJed Brown Input Parameters: 162fcb023d4SJed Brown + pc - the multigrid context 163fcb023d4SJed Brown . l - the level (0 is coarsest) to supply 164fcb023d4SJed Brown . residualt - function used to form transpose of residual, if none is provided the previously provide one is used, if no 165fcb023d4SJed Brown previous one were provided then a default is used 166fcb023d4SJed Brown - mat - matrix associated with residual 167fcb023d4SJed Brown 168fcb023d4SJed Brown Level: advanced 169fcb023d4SJed Brown 170*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGResidualTransposeDefault()` 171fcb023d4SJed Brown @*/ 1729371c9d4SSatish Balay PetscErrorCode PCMGSetResidualTranspose(PC pc, PetscInt l, PetscErrorCode (*residualt)(Mat, Vec, Vec, Vec), Mat mat) { 173fcb023d4SJed Brown PC_MG *mg = (PC_MG *)pc->data; 174fcb023d4SJed Brown PC_MG_Levels **mglevels = mg->levels; 175fcb023d4SJed Brown 176fcb023d4SJed Brown PetscFunctionBegin; 177fcb023d4SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 17828b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 179fcb023d4SJed Brown if (residualt) mglevels[l]->residualtranspose = residualt; 180fcb023d4SJed Brown if (!mglevels[l]->residualtranspose) mglevels[l]->residualtranspose = PCMGResidualTransposeDefault; 18130b0564aSStefano Zampini mglevels[l]->matresidualtranspose = PCMGMatResidualTransposeDefault; 1829566063dSJacob Faibussowitsch if (mat) PetscCall(PetscObjectReference((PetscObject)mat)); 1839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->A)); 184fcb023d4SJed Brown mglevels[l]->A = mat; 185fcb023d4SJed Brown PetscFunctionReturn(0); 186fcb023d4SJed Brown } 187fcb023d4SJed Brown 1884b9ad928SBarry Smith /*@ 189aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 190bf5b2e24SBarry Smith interpolation from l-1 to the lth level 1914b9ad928SBarry Smith 192*f1580f4eSBarry Smith Logically Collective on pc 1934b9ad928SBarry Smith 1944b9ad928SBarry Smith Input Parameters: 1954b9ad928SBarry Smith + pc - the multigrid context 1964b9ad928SBarry Smith . mat - the interpolation operator 197bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 1984b9ad928SBarry Smith 1994b9ad928SBarry Smith Level: advanced 2004b9ad928SBarry Smith 2014b9ad928SBarry Smith Notes: 2024b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 2034b9ad928SBarry Smith for the same level. 2044b9ad928SBarry Smith 2054b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 2064b9ad928SBarry Smith out from the matrix size which one it is. 2074b9ad928SBarry Smith 208*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRestriction()` 2094b9ad928SBarry Smith @*/ 2109371c9d4SSatish Balay PetscErrorCode PCMGSetInterpolation(PC pc, PetscInt l, Mat mat) { 211f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 212f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2134b9ad928SBarry Smith 2144b9ad928SBarry Smith PetscFunctionBegin; 215c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 21628b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 21728b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set interpolation routine for coarsest level"); 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 2199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->interpolate)); 2202fa5cd67SKarl Rupp 221f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 2224b9ad928SBarry Smith PetscFunctionReturn(0); 2234b9ad928SBarry Smith } 2244b9ad928SBarry Smith 2258a2c336bSFande Kong /*@ 22695750439SFande Kong PCMGSetOperators - Sets operator and preconditioning matrix for lth level 2278a2c336bSFande Kong 228*f1580f4eSBarry Smith Logically Collective on pc 2298a2c336bSFande Kong 2308a2c336bSFande Kong Input Parameters: 2318a2c336bSFande Kong + pc - the multigrid context 2328a2c336bSFande Kong . Amat - the operator 2338a2c336bSFande Kong . pmat - the preconditioning operator 23495750439SFande Kong - l - the level (0 is the coarsest) to supply 2358a2c336bSFande Kong 2368a2c336bSFande Kong Level: advanced 2378a2c336bSFande Kong 2388a2c336bSFande Kong .keywords: multigrid, set, interpolate, level 2398a2c336bSFande Kong 240*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetGalerkin()`, `PCMGSetRestriction()`, `PCMGSetInterpolation()` 2418a2c336bSFande Kong @*/ 2429371c9d4SSatish Balay PetscErrorCode PCMGSetOperators(PC pc, PetscInt l, Mat Amat, Mat Pmat) { 243360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data; 244360ee056SFande Kong PC_MG_Levels **mglevels = mg->levels; 245360ee056SFande Kong 246360ee056SFande Kong PetscFunctionBegin; 247360ee056SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2488a2c336bSFande Kong PetscValidHeaderSpecific(Amat, MAT_CLASSID, 3); 2498a2c336bSFande Kong PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 4); 25028b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 2519566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[l]->smoothd, Amat, Pmat)); 252360ee056SFande Kong PetscFunctionReturn(0); 253360ee056SFande Kong } 254360ee056SFande Kong 255c88c5224SJed Brown /*@ 256c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 257c88c5224SJed Brown interpolation from l-1 to the lth level 258c88c5224SJed Brown 259*f1580f4eSBarry Smith Logically Collective on pc 260c88c5224SJed Brown 261c88c5224SJed Brown Input Parameters: 262c88c5224SJed Brown + pc - the multigrid context 263c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 264c88c5224SJed Brown 265c88c5224SJed Brown Output Parameter: 2663ad4599aSBarry Smith . mat - the interpolation matrix, can be NULL 267c88c5224SJed Brown 268c88c5224SJed Brown Level: advanced 269c88c5224SJed Brown 270*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()` 271c88c5224SJed Brown @*/ 2729371c9d4SSatish Balay PetscErrorCode PCMGGetInterpolation(PC pc, PetscInt l, Mat *mat) { 273c88c5224SJed Brown PC_MG *mg = (PC_MG *)pc->data; 274c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 275c88c5224SJed Brown 276c88c5224SJed Brown PetscFunctionBegin; 277c88c5224SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2783ad4599aSBarry Smith if (mat) PetscValidPointer(mat, 3); 27928b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 2802472a847SBarry Smith PetscCheck(l > 0 && l < mg->nlevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Level %" PetscInt_FMT " must be in range {1,...,%" PetscInt_FMT "}", l, mg->nlevels - 1); 28148a46eb9SPierre Jolivet if (!mglevels[l]->interpolate && mglevels[l]->restrct) PetscCall(PCMGSetInterpolation(pc, l, mglevels[l]->restrct)); 2823ad4599aSBarry Smith if (mat) *mat = mglevels[l]->interpolate; 283c88c5224SJed Brown PetscFunctionReturn(0); 284c88c5224SJed Brown } 285c88c5224SJed Brown 2868a2c336bSFande Kong /*@ 287eab52d2bSLawrence Mitchell PCMGSetRestriction - Sets the function to be used to restrict dual vectors 2884b9ad928SBarry Smith from level l to l-1. 2894b9ad928SBarry Smith 290*f1580f4eSBarry Smith Logically Collective on pc 2914b9ad928SBarry Smith 2924b9ad928SBarry Smith Input Parameters: 2934b9ad928SBarry Smith + pc - the multigrid context 294c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 295c88c5224SJed Brown - mat - the restriction matrix 2964b9ad928SBarry Smith 2974b9ad928SBarry Smith Level: advanced 2984b9ad928SBarry Smith 2994b9ad928SBarry Smith Notes: 3004b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 3014b9ad928SBarry Smith for the same level. 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 3044b9ad928SBarry Smith out from the matrix size which one it is. 3054b9ad928SBarry Smith 306*f1580f4eSBarry Smith If you do not set this, the transpose of the `Mat` set with `PCMGSetInterpolation()` 307fccaa45eSBarry Smith is used. 308fccaa45eSBarry Smith 309*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInterpolation()` 3104b9ad928SBarry Smith @*/ 3119371c9d4SSatish Balay PetscErrorCode PCMGSetRestriction(PC pc, PetscInt l, Mat mat) { 312f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 313f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 3144b9ad928SBarry Smith 3154b9ad928SBarry Smith PetscFunctionBegin; 316c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 317c88c5224SJed Brown PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 31828b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 31928b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set restriction routine for coarsest level"); 3209566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 3219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->restrct)); 3222fa5cd67SKarl Rupp 323f3fbd535SBarry Smith mglevels[l]->restrct = mat; 3244b9ad928SBarry Smith PetscFunctionReturn(0); 3254b9ad928SBarry Smith } 3264b9ad928SBarry Smith 327c88c5224SJed Brown /*@ 328eab52d2bSLawrence Mitchell PCMGGetRestriction - Gets the function to be used to restrict dual vectors 329c88c5224SJed Brown from level l to l-1. 330c88c5224SJed Brown 331*f1580f4eSBarry Smith Logically Collective on pc 332c88c5224SJed Brown 333c88c5224SJed Brown Input Parameters: 334c88c5224SJed Brown + pc - the multigrid context 335c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 336c88c5224SJed Brown 337c88c5224SJed Brown Output Parameter: 338c88c5224SJed Brown . mat - the restriction matrix 339c88c5224SJed Brown 340c88c5224SJed Brown Level: advanced 341c88c5224SJed Brown 342*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGGetInjection()` 343c88c5224SJed Brown @*/ 3449371c9d4SSatish Balay PetscErrorCode PCMGGetRestriction(PC pc, PetscInt l, Mat *mat) { 345c88c5224SJed Brown PC_MG *mg = (PC_MG *)pc->data; 346c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 347c88c5224SJed Brown 348c88c5224SJed Brown PetscFunctionBegin; 349c88c5224SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3503ad4599aSBarry Smith if (mat) PetscValidPointer(mat, 3); 35128b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 3522472a847SBarry Smith PetscCheck(l > 0 && l < mg->nlevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Level %" PetscInt_FMT " must be in range {1,...,%" PetscInt_FMT "}", l, mg->nlevels - 1); 35348a46eb9SPierre Jolivet if (!mglevels[l]->restrct && mglevels[l]->interpolate) PetscCall(PCMGSetRestriction(pc, l, mglevels[l]->interpolate)); 3543ad4599aSBarry Smith if (mat) *mat = mglevels[l]->restrct; 355c88c5224SJed Brown PetscFunctionReturn(0); 356c88c5224SJed Brown } 357c88c5224SJed Brown 35873250ac0SBarry Smith /*@ 35973250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 36073250ac0SBarry Smith 361*f1580f4eSBarry Smith Logically Collective on pc 362c88c5224SJed Brown 363c88c5224SJed Brown Input Parameters: 364c88c5224SJed Brown + pc - the multigrid context 365c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 366c88c5224SJed Brown . rscale - the scaling 367c88c5224SJed Brown 368c88c5224SJed Brown Level: advanced 369c88c5224SJed Brown 370*f1580f4eSBarry Smith Note: 371*f1580f4eSBarry 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. 372*f1580f4eSBarry Smith It is preferable to use `PCMGSetInjection()` to control moving primal vectors. 373c88c5224SJed Brown 374*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGSetInjection()` 375c88c5224SJed Brown @*/ 3769371c9d4SSatish Balay PetscErrorCode PCMGSetRScale(PC pc, PetscInt l, Vec rscale) { 377c88c5224SJed Brown PC_MG *mg = (PC_MG *)pc->data; 378c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 379c88c5224SJed Brown 380c88c5224SJed Brown PetscFunctionBegin; 381c88c5224SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 38228b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 3832472a847SBarry Smith PetscCheck(l > 0 && l < mg->nlevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Level %" PetscInt_FMT " must be in range {1,...,%" PetscInt_FMT "}", l, mg->nlevels - 1); 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rscale)); 3859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[l]->rscale)); 3862fa5cd67SKarl Rupp 387c88c5224SJed Brown mglevels[l]->rscale = rscale; 388c88c5224SJed Brown PetscFunctionReturn(0); 389c88c5224SJed Brown } 390c88c5224SJed Brown 391c88c5224SJed Brown /*@ 392c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 393c88c5224SJed Brown 394*f1580f4eSBarry Smith Collective on pc 39573250ac0SBarry Smith 39673250ac0SBarry Smith Input Parameters: 39773250ac0SBarry Smith + pc - the multigrid context 39873250ac0SBarry Smith . rscale - the scaling 39973250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 40073250ac0SBarry Smith 40173250ac0SBarry Smith Level: advanced 40273250ac0SBarry Smith 403*f1580f4eSBarry Smith Note: 404*f1580f4eSBarry 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. 405*f1580f4eSBarry Smith It is preferable to use `PCMGGetInjection()` to control moving primal vectors. 40673250ac0SBarry Smith 407*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInterpolation()`, `PCMGGetRestriction()`, `PCMGGetInjection()` 40873250ac0SBarry Smith @*/ 4099371c9d4SSatish Balay PetscErrorCode PCMGGetRScale(PC pc, PetscInt l, Vec *rscale) { 41073250ac0SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 41173250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 41273250ac0SBarry Smith 41373250ac0SBarry Smith PetscFunctionBegin; 41473250ac0SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 41528b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 4162472a847SBarry Smith PetscCheck(l > 0 && l < mg->nlevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Level %" PetscInt_FMT " must be in range {1,...,%" PetscInt_FMT "}", l, mg->nlevels - 1); 417c88c5224SJed Brown if (!mglevels[l]->rscale) { 418c88c5224SJed Brown Mat R; 419c88c5224SJed Brown Vec X, Y, coarse, fine; 420c88c5224SJed Brown PetscInt M, N; 4219566063dSJacob Faibussowitsch PetscCall(PCMGGetRestriction(pc, l, &R)); 4229566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(R, &X, &Y)); 4239566063dSJacob Faibussowitsch PetscCall(MatGetSize(R, &M, &N)); 4242fa5cd67SKarl Rupp if (M < N) { 4252fa5cd67SKarl Rupp fine = X; 4262fa5cd67SKarl Rupp coarse = Y; 4272fa5cd67SKarl Rupp } else if (N < M) { 4289371c9d4SSatish Balay fine = Y; 4299371c9d4SSatish Balay coarse = X; 430ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)R), PETSC_ERR_SUP, "Restriction matrix is square, cannot determine which Vec is coarser"); 4319566063dSJacob Faibussowitsch PetscCall(VecSet(fine, 1.)); 4329566063dSJacob Faibussowitsch PetscCall(MatRestrict(R, fine, coarse)); 4339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fine)); 4349566063dSJacob Faibussowitsch PetscCall(VecReciprocal(coarse)); 435c88c5224SJed Brown mglevels[l]->rscale = coarse; 436c88c5224SJed Brown } 437c88c5224SJed Brown *rscale = mglevels[l]->rscale; 43873250ac0SBarry Smith PetscFunctionReturn(0); 43973250ac0SBarry Smith } 44073250ac0SBarry Smith 441f39d8e23SSatish Balay /*@ 442eab52d2bSLawrence Mitchell PCMGSetInjection - Sets the function to be used to inject primal vectors 443eab52d2bSLawrence Mitchell from level l to l-1. 444eab52d2bSLawrence Mitchell 445*f1580f4eSBarry Smith Logically Collective on pc 446eab52d2bSLawrence Mitchell 447eab52d2bSLawrence Mitchell Input Parameters: 448eab52d2bSLawrence Mitchell + pc - the multigrid context 449eab52d2bSLawrence Mitchell . l - the level (0 is coarsest) to supply [Do not supply 0] 450eab52d2bSLawrence Mitchell - mat - the injection matrix 451eab52d2bSLawrence Mitchell 452eab52d2bSLawrence Mitchell Level: advanced 453eab52d2bSLawrence Mitchell 454*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRestriction()` 455eab52d2bSLawrence Mitchell @*/ 4569371c9d4SSatish Balay PetscErrorCode PCMGSetInjection(PC pc, PetscInt l, Mat mat) { 457eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG *)pc->data; 458eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 459eab52d2bSLawrence Mitchell 460eab52d2bSLawrence Mitchell PetscFunctionBegin; 461eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 462eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 46328b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 46428b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set restriction routine for coarsest level"); 4659566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 4669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->inject)); 467eab52d2bSLawrence Mitchell 468eab52d2bSLawrence Mitchell mglevels[l]->inject = mat; 469eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 470eab52d2bSLawrence Mitchell } 471eab52d2bSLawrence Mitchell 472eab52d2bSLawrence Mitchell /*@ 473eab52d2bSLawrence Mitchell PCMGGetInjection - Gets the function to be used to inject primal vectors 474eab52d2bSLawrence Mitchell from level l to l-1. 475eab52d2bSLawrence Mitchell 476*f1580f4eSBarry Smith Logically Collective on pc 477eab52d2bSLawrence Mitchell 478eab52d2bSLawrence Mitchell Input Parameters: 479eab52d2bSLawrence Mitchell + pc - the multigrid context 480eab52d2bSLawrence Mitchell - l - the level (0 is coarsest) to supply [Do not supply 0] 481eab52d2bSLawrence Mitchell 482eab52d2bSLawrence Mitchell Output Parameter: 48399a38656SLawrence Mitchell . mat - the restriction matrix (may be NULL if no injection is available). 484eab52d2bSLawrence Mitchell 485eab52d2bSLawrence Mitchell Level: advanced 486eab52d2bSLawrence Mitchell 487*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInjection()`, `PCMGetGetRestriction()` 488eab52d2bSLawrence Mitchell @*/ 4899371c9d4SSatish Balay PetscErrorCode PCMGGetInjection(PC pc, PetscInt l, Mat *mat) { 490eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG *)pc->data; 491eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 492eab52d2bSLawrence Mitchell 493eab52d2bSLawrence Mitchell PetscFunctionBegin; 494eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 495eab52d2bSLawrence Mitchell if (mat) PetscValidPointer(mat, 3); 49628b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 4972472a847SBarry Smith PetscCheck(l > 0 && l < mg->nlevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Level %" PetscInt_FMT " must be in range {1,...,%" PetscInt_FMT "}", l, mg->nlevels - 1); 498eab52d2bSLawrence Mitchell if (mat) *mat = mglevels[l]->inject; 499eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 500eab52d2bSLawrence Mitchell } 501eab52d2bSLawrence Mitchell 502eab52d2bSLawrence Mitchell /*@ 503*f1580f4eSBarry Smith PCMGGetSmoother - Gets the `KSP` context to be used as smoother for 504*f1580f4eSBarry Smith both pre- and post-smoothing. Call both `PCMGGetSmootherUp()` and 505*f1580f4eSBarry Smith `PCMGGetSmootherDown()` to use different functions for pre- and 5064b9ad928SBarry Smith post-smoothing. 5074b9ad928SBarry Smith 508*f1580f4eSBarry Smith Not Collective, ksp returned is parallel if pc is 5094b9ad928SBarry Smith 5104b9ad928SBarry Smith Input Parameters: 5114b9ad928SBarry Smith + pc - the multigrid context 5124b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5134b9ad928SBarry Smith 51401d2d390SJose E. Roman Output Parameter: 5154b9ad928SBarry Smith . ksp - the smoother 5164b9ad928SBarry Smith 517*f1580f4eSBarry Smith Note: 518*f1580f4eSBarry Smith Once you have called this routine, you can call `KSPSetOperators()` on the resulting ksp to provide the operators for the smoother for this level. 519*f1580f4eSBarry Smith You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call `KSPGetPC`(ksp,&pc) 520*f1580f4eSBarry Smith and modify PC options for the smoother; for example `PCSetType`(pc,`PCSOR`); to use SOR smoothing. 52157420d5bSBarry Smith 5224b9ad928SBarry Smith Level: advanced 5234b9ad928SBarry Smith 524*f1580f4eSBarry Smith .seealso: PCMG`, ``PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, `PCMGGetCoarseSolve()` 5254b9ad928SBarry Smith @*/ 5269371c9d4SSatish Balay PetscErrorCode PCMGGetSmoother(PC pc, PetscInt l, KSP *ksp) { 527f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 528f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5294b9ad928SBarry Smith 5304b9ad928SBarry Smith PetscFunctionBegin; 531c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 532f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5334b9ad928SBarry Smith PetscFunctionReturn(0); 5344b9ad928SBarry Smith } 5354b9ad928SBarry Smith 536f39d8e23SSatish Balay /*@ 53797177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 5384b9ad928SBarry Smith coarse grid correction (post-smoother). 5394b9ad928SBarry Smith 540*f1580f4eSBarry Smith Not Collective, ksp returned is parallel if pc is 5414b9ad928SBarry Smith 5424b9ad928SBarry Smith Input Parameters: 5434b9ad928SBarry Smith + pc - the multigrid context 5444b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5454b9ad928SBarry Smith 54601d2d390SJose E. Roman Output Parameter: 5474b9ad928SBarry Smith . ksp - the smoother 5484b9ad928SBarry Smith 5494b9ad928SBarry Smith Level: advanced 5504b9ad928SBarry Smith 551*f1580f4eSBarry Smith Note: 552*f1580f4eSBarry Smith Calling this will result in a different pre and post smoother so you may need to set options on the pre smoother also 55389cce641SBarry Smith 554*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()` 5554b9ad928SBarry Smith @*/ 5569371c9d4SSatish Balay PetscErrorCode PCMGGetSmootherUp(PC pc, PetscInt l, KSP *ksp) { 557f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 558f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 559f69a0ea3SMatthew Knepley const char *prefix; 5604b9ad928SBarry Smith MPI_Comm comm; 5614b9ad928SBarry Smith 5624b9ad928SBarry Smith PetscFunctionBegin; 563c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 5644b9ad928SBarry Smith /* 5654b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 5664b9ad928SBarry Smith Thus we check if a different one has already been allocated, 5674b9ad928SBarry Smith if not we allocate it. 5684b9ad928SBarry Smith */ 56928b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "There is no such thing as a up smoother on the coarse grid"); 570f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 57119fd82e9SBarry Smith KSPType ksptype; 57219fd82e9SBarry Smith PCType pctype; 573336babb1SJed Brown PC ipc; 574336babb1SJed Brown PetscReal rtol, abstol, dtol; 575336babb1SJed Brown PetscInt maxits; 576336babb1SJed Brown KSPNormType normtype; 5779566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mglevels[l]->smoothd, &comm)); 5789566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(mglevels[l]->smoothd, &prefix)); 5799566063dSJacob Faibussowitsch PetscCall(KSPGetTolerances(mglevels[l]->smoothd, &rtol, &abstol, &dtol, &maxits)); 5809566063dSJacob Faibussowitsch PetscCall(KSPGetType(mglevels[l]->smoothd, &ksptype)); 5819566063dSJacob Faibussowitsch PetscCall(KSPGetNormType(mglevels[l]->smoothd, &normtype)); 5829566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[l]->smoothd, &ipc)); 5839566063dSJacob Faibussowitsch PetscCall(PCGetType(ipc, &pctype)); 584336babb1SJed Brown 5859566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &mglevels[l]->smoothu)); 5869566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(mglevels[l]->smoothu, pc->erroriffailure)); 5879566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu, (PetscObject)pc, mglevels[0]->levels - l)); 5889566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(mglevels[l]->smoothu, prefix)); 5899566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[l]->smoothu, rtol, abstol, dtol, maxits)); 5909566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[l]->smoothu, ksptype)); 5919566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(mglevels[l]->smoothu, normtype)); 5929566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(mglevels[l]->smoothu, KSPConvergedSkip, NULL, NULL)); 5939566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[l]->smoothu, &ipc)); 5949566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, pctype)); 5959566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)mglevels[l]->smoothu)); 5969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level)); 5974b9ad928SBarry Smith } 598f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 5994b9ad928SBarry Smith PetscFunctionReturn(0); 6004b9ad928SBarry Smith } 6014b9ad928SBarry Smith 602f39d8e23SSatish Balay /*@ 603*f1580f4eSBarry Smith PCMGGetSmootherDown - Gets the `KSP` context to be used as smoother before 6044b9ad928SBarry Smith coarse grid correction (pre-smoother). 6054b9ad928SBarry Smith 606*f1580f4eSBarry Smith Not Collective, ksp returned is parallel if pc is 6074b9ad928SBarry Smith 6084b9ad928SBarry Smith Input Parameters: 6094b9ad928SBarry Smith + pc - the multigrid context 6104b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 6114b9ad928SBarry Smith 61201d2d390SJose E. Roman Output Parameter: 6134b9ad928SBarry Smith . ksp - the smoother 6144b9ad928SBarry Smith 6154b9ad928SBarry Smith Level: advanced 6164b9ad928SBarry Smith 617*f1580f4eSBarry Smith Note: 618*f1580f4eSBarry Smith Calling this will result in a different pre and post smoother so you may need to 61989cce641SBarry Smith set options on the post smoother also 62089cce641SBarry Smith 621*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmoother()` 6224b9ad928SBarry Smith @*/ 6239371c9d4SSatish Balay PetscErrorCode PCMGGetSmootherDown(PC pc, PetscInt l, KSP *ksp) { 624f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 625f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6264b9ad928SBarry Smith 6274b9ad928SBarry Smith PetscFunctionBegin; 628c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6294b9ad928SBarry Smith /* make sure smoother up and down are different */ 6301baa6e33SBarry Smith if (l) PetscCall(PCMGGetSmootherUp(pc, l, NULL)); 631f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 6324b9ad928SBarry Smith PetscFunctionReturn(0); 6334b9ad928SBarry Smith } 6344b9ad928SBarry Smith 6354b9ad928SBarry Smith /*@ 636cab31ae5SJed Brown PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 6374b9ad928SBarry Smith 638*f1580f4eSBarry Smith Logically Collective on pc 6394b9ad928SBarry Smith 6404b9ad928SBarry Smith Input Parameters: 6414b9ad928SBarry Smith + pc - the multigrid context 642c1cbb1deSBarry Smith . l - the level (0 is coarsest) 643*f1580f4eSBarry Smith - c - either `PC_MG_CYCLE_V` or `PC_MG_CYCLE_W` 6444b9ad928SBarry Smith 6454b9ad928SBarry Smith Level: advanced 6464b9ad928SBarry Smith 647*f1580f4eSBarry Smith .seealso: `PCMG`, PCMGCycleType`, `PCMGSetCycleType()` 6484b9ad928SBarry Smith @*/ 6499371c9d4SSatish Balay PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc, PetscInt l, PCMGCycleType c) { 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); 65528b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 656c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, l, 2); 657c51679f6SSatish Balay PetscValidLogicalCollectiveEnum(pc, c, 3); 658f3fbd535SBarry Smith mglevels[l]->cycles = c; 6594b9ad928SBarry Smith PetscFunctionReturn(0); 6604b9ad928SBarry Smith } 6614b9ad928SBarry Smith 6624b9ad928SBarry Smith /*@ 663f3b08a26SMatthew G. Knepley PCMGSetRhs - Sets the vector to be used to store the right-hand side on a particular level. 6644b9ad928SBarry Smith 665*f1580f4eSBarry Smith Logically Collective on pc 6664b9ad928SBarry Smith 6674b9ad928SBarry Smith Input Parameters: 6684b9ad928SBarry Smith + pc - the multigrid context 6694b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 670f3b08a26SMatthew G. Knepley - c - the Vec 6714b9ad928SBarry Smith 6724b9ad928SBarry Smith Level: advanced 6734b9ad928SBarry Smith 674*f1580f4eSBarry Smith Note: 675*f1580f4eSBarry Smith 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. 676fccaa45eSBarry Smith 677*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetX()`, `PCMGSetR()` 6784b9ad928SBarry Smith @*/ 6799371c9d4SSatish Balay PetscErrorCode PCMGSetRhs(PC pc, PetscInt l, Vec c) { 680f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 681f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6824b9ad928SBarry Smith 6834b9ad928SBarry Smith PetscFunctionBegin; 684c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 68528b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 68608401ef6SPierre Jolivet PetscCheck(l != mglevels[0]->levels - 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Do not set rhs for finest level"); 6879566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)c)); 6889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[l]->b)); 6892fa5cd67SKarl Rupp 690f3fbd535SBarry Smith mglevels[l]->b = c; 6914b9ad928SBarry Smith PetscFunctionReturn(0); 6924b9ad928SBarry Smith } 6934b9ad928SBarry Smith 6944b9ad928SBarry Smith /*@ 695f3b08a26SMatthew G. Knepley PCMGSetX - Sets the vector to be used to store the solution on a particular level. 6964b9ad928SBarry Smith 697*f1580f4eSBarry Smith Logically Collective on pc 6984b9ad928SBarry Smith 6994b9ad928SBarry Smith Input Parameters: 7004b9ad928SBarry Smith + pc - the multigrid context 701251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 702f3b08a26SMatthew G. Knepley - c - the Vec 7034b9ad928SBarry Smith 7044b9ad928SBarry Smith Level: advanced 7054b9ad928SBarry Smith 706*f1580f4eSBarry Smith Note: 707*f1580f4eSBarry Smith 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. 708fccaa45eSBarry Smith 709*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRhs()`, `PCMGSetR()` 7104b9ad928SBarry Smith @*/ 7119371c9d4SSatish Balay PetscErrorCode PCMGSetX(PC pc, PetscInt l, Vec c) { 712f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 713f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7144b9ad928SBarry Smith 7154b9ad928SBarry Smith PetscFunctionBegin; 716c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 71728b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 71808401ef6SPierre Jolivet PetscCheck(l != mglevels[0]->levels - 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Do not set x for finest level"); 7199566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)c)); 7209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[l]->x)); 7212fa5cd67SKarl Rupp 722f3fbd535SBarry Smith mglevels[l]->x = c; 7234b9ad928SBarry Smith PetscFunctionReturn(0); 7244b9ad928SBarry Smith } 7254b9ad928SBarry Smith 7264b9ad928SBarry Smith /*@ 727f3b08a26SMatthew G. Knepley PCMGSetR - Sets the vector to be used to store the residual on a particular level. 7284b9ad928SBarry Smith 729*f1580f4eSBarry Smith Logically Collective on pc 7304b9ad928SBarry Smith 7314b9ad928SBarry Smith Input Parameters: 7324b9ad928SBarry Smith + pc - the multigrid context 7334b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 734f3b08a26SMatthew G. Knepley - c - the Vec 7354b9ad928SBarry Smith 7364b9ad928SBarry Smith Level: advanced 7374b9ad928SBarry Smith 738*f1580f4eSBarry Smith Note: 739*f1580f4eSBarry Smith 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. 740fccaa45eSBarry Smith 741*f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRhs()`, `PCMGSetX()` 7424b9ad928SBarry Smith @*/ 7439371c9d4SSatish Balay PetscErrorCode PCMGSetR(PC pc, PetscInt l, Vec c) { 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); 74928b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 75028b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Need not set residual vector for coarse grid"); 7519566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)c)); 7529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[l]->r)); 7532fa5cd67SKarl Rupp 754f3fbd535SBarry Smith mglevels[l]->r = c; 7554b9ad928SBarry Smith PetscFunctionReturn(0); 7564b9ad928SBarry Smith } 757