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 7f1580f4eSBarry 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 19f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetResidual()`, `PCMGSetMatResidual()` 204b9ad928SBarry Smith @*/ 21d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGResidualDefault(Mat mat, Vec b, Vec x, Vec r) 22d71ae5a4SJacob Faibussowitsch { 234b9ad928SBarry Smith PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(MatResidual(mat, b, x, r)); 254b9ad928SBarry Smith PetscFunctionReturn(0); 264b9ad928SBarry Smith } 274b9ad928SBarry Smith 28fcb023d4SJed Brown /*@C 29fcb023d4SJed Brown PCMGResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system 30fcb023d4SJed Brown 31f1580f4eSBarry Smith Collective on mat 32fcb023d4SJed Brown 33fcb023d4SJed Brown Input Parameters: 34fcb023d4SJed Brown + mat - the matrix 35fcb023d4SJed Brown . b - the right-hand-side 36fcb023d4SJed Brown - x - the approximate solution 37fcb023d4SJed Brown 38fcb023d4SJed Brown Output Parameter: 39fcb023d4SJed Brown . r - location to store the residual 40fcb023d4SJed Brown 41fcb023d4SJed Brown Level: developer 42fcb023d4SJed Brown 43f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetResidualTranspose()`, `PCMGMatResidualTransposeDefault()` 44fcb023d4SJed Brown @*/ 45d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGResidualTransposeDefault(Mat mat, Vec b, Vec x, Vec r) 46d71ae5a4SJacob Faibussowitsch { 47fcb023d4SJed Brown PetscFunctionBegin; 489566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat, x, r)); 499566063dSJacob Faibussowitsch PetscCall(VecAYPX(r, -1.0, b)); 50fcb023d4SJed Brown PetscFunctionReturn(0); 51fcb023d4SJed Brown } 52fcb023d4SJed Brown 5330b0564aSStefano Zampini /*@C 5430b0564aSStefano Zampini PCMGMatResidualDefault - Default routine to calculate the residual. 5530b0564aSStefano Zampini 56f1580f4eSBarry Smith Collective on mat 5730b0564aSStefano Zampini 5830b0564aSStefano Zampini Input Parameters: 5930b0564aSStefano Zampini + mat - the matrix 6030b0564aSStefano Zampini . b - the right-hand-side 6130b0564aSStefano Zampini - x - the approximate solution 6230b0564aSStefano Zampini 6330b0564aSStefano Zampini Output Parameter: 6430b0564aSStefano Zampini . r - location to store the residual 6530b0564aSStefano Zampini 6630b0564aSStefano Zampini Level: developer 6730b0564aSStefano Zampini 68f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetMatResidual()`, `PCMGResidualDefault()` 6930b0564aSStefano Zampini @*/ 70d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMatResidualDefault(Mat mat, Mat b, Mat x, Mat r) 71d71ae5a4SJacob Faibussowitsch { 7230b0564aSStefano Zampini PetscFunctionBegin; 739566063dSJacob Faibussowitsch PetscCall(MatMatMult(mat, x, MAT_REUSE_MATRIX, PETSC_DEFAULT, &r)); 749566063dSJacob Faibussowitsch PetscCall(MatAYPX(r, -1.0, b, UNKNOWN_NONZERO_PATTERN)); 7530b0564aSStefano Zampini PetscFunctionReturn(0); 7630b0564aSStefano Zampini } 7730b0564aSStefano Zampini 7830b0564aSStefano Zampini /*@C 7930b0564aSStefano Zampini PCMGMatResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system 8030b0564aSStefano Zampini 81f1580f4eSBarry Smith Collective on mat 8230b0564aSStefano Zampini 8330b0564aSStefano Zampini Input Parameters: 8430b0564aSStefano Zampini + mat - the matrix 8530b0564aSStefano Zampini . b - the right-hand-side 8630b0564aSStefano Zampini - x - the approximate solution 8730b0564aSStefano Zampini 8830b0564aSStefano Zampini Output Parameter: 8930b0564aSStefano Zampini . r - location to store the residual 9030b0564aSStefano Zampini 9130b0564aSStefano Zampini Level: developer 9230b0564aSStefano Zampini 93f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetMatResidualTranspose()` 9430b0564aSStefano Zampini @*/ 95d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMatResidualTransposeDefault(Mat mat, Mat b, Mat x, Mat r) 96d71ae5a4SJacob Faibussowitsch { 9730b0564aSStefano Zampini PetscFunctionBegin; 989566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(mat, x, MAT_REUSE_MATRIX, PETSC_DEFAULT, &r)); 999566063dSJacob Faibussowitsch PetscCall(MatAYPX(r, -1.0, b, UNKNOWN_NONZERO_PATTERN)); 10030b0564aSStefano Zampini PetscFunctionReturn(0); 10130b0564aSStefano Zampini } 102f39d8e23SSatish Balay /*@ 10397177400SBarry Smith PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 1044b9ad928SBarry Smith 1054b9ad928SBarry Smith Not Collective 1064b9ad928SBarry Smith 1074b9ad928SBarry Smith Input Parameter: 1084b9ad928SBarry Smith . pc - the multigrid context 1094b9ad928SBarry Smith 1104b9ad928SBarry Smith Output Parameter: 1114b9ad928SBarry Smith . ksp - the coarse grid solver context 1124b9ad928SBarry Smith 1134b9ad928SBarry Smith Level: advanced 1144b9ad928SBarry Smith 115f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, `PCMGGetSmoother()` 1164b9ad928SBarry Smith @*/ 117d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetCoarseSolve(PC pc, KSP *ksp) 118d71ae5a4SJacob Faibussowitsch { 119f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 120f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith PetscFunctionBegin; 123c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 124f3fbd535SBarry Smith *ksp = mglevels[0]->smoothd; 1254b9ad928SBarry Smith PetscFunctionReturn(0); 1264b9ad928SBarry Smith } 1274b9ad928SBarry Smith 1284b9ad928SBarry Smith /*@C 129f1580f4eSBarry Smith PCMGSetResidual - Sets the function to be used to calculate the residual on the lth level. 1304b9ad928SBarry Smith 131f1580f4eSBarry Smith Logically Collective on pc 1324b9ad928SBarry Smith 1334b9ad928SBarry Smith Input Parameters: 1344b9ad928SBarry Smith + pc - the multigrid context 1354b9ad928SBarry Smith . l - the level (0 is coarsest) to supply 136157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no 137d0e4de75SBarry Smith previous one were provided then a default is used 1384b9ad928SBarry Smith - mat - matrix associated with residual 1394b9ad928SBarry Smith 1404b9ad928SBarry Smith Level: advanced 1414b9ad928SBarry Smith 142f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGResidualDefault()` 1434b9ad928SBarry Smith @*/ 144d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetResidual(PC pc, PetscInt l, PetscErrorCode (*residual)(Mat, Vec, Vec, Vec), Mat mat) 145d71ae5a4SJacob Faibussowitsch { 146f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 147f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1484b9ad928SBarry Smith 1494b9ad928SBarry Smith PetscFunctionBegin; 150c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15128b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 1522fa5cd67SKarl Rupp if (residual) mglevels[l]->residual = residual; 15354b2cd4bSJed Brown if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault; 15430b0564aSStefano Zampini mglevels[l]->matresidual = PCMGMatResidualDefault; 1559566063dSJacob Faibussowitsch if (mat) PetscCall(PetscObjectReference((PetscObject)mat)); 1569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->A)); 157f3fbd535SBarry Smith mglevels[l]->A = mat; 1584b9ad928SBarry Smith PetscFunctionReturn(0); 1594b9ad928SBarry Smith } 1604b9ad928SBarry Smith 161fcb023d4SJed Brown /*@C 162fcb023d4SJed Brown PCMGSetResidualTranspose - Sets the function to be used to calculate the residual of the transposed linear system 163fcb023d4SJed Brown on the lth level. 164fcb023d4SJed Brown 165f1580f4eSBarry Smith Logically Collective on pc 166fcb023d4SJed Brown 167fcb023d4SJed Brown Input Parameters: 168fcb023d4SJed Brown + pc - the multigrid context 169fcb023d4SJed Brown . l - the level (0 is coarsest) to supply 170fcb023d4SJed Brown . residualt - function used to form transpose of residual, if none is provided the previously provide one is used, if no 171fcb023d4SJed Brown previous one were provided then a default is used 172fcb023d4SJed Brown - mat - matrix associated with residual 173fcb023d4SJed Brown 174fcb023d4SJed Brown Level: advanced 175fcb023d4SJed Brown 176f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGResidualTransposeDefault()` 177fcb023d4SJed Brown @*/ 178d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetResidualTranspose(PC pc, PetscInt l, PetscErrorCode (*residualt)(Mat, Vec, Vec, Vec), Mat mat) 179d71ae5a4SJacob Faibussowitsch { 180fcb023d4SJed Brown PC_MG *mg = (PC_MG *)pc->data; 181fcb023d4SJed Brown PC_MG_Levels **mglevels = mg->levels; 182fcb023d4SJed Brown 183fcb023d4SJed Brown PetscFunctionBegin; 184fcb023d4SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 18528b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 186fcb023d4SJed Brown if (residualt) mglevels[l]->residualtranspose = residualt; 187fcb023d4SJed Brown if (!mglevels[l]->residualtranspose) mglevels[l]->residualtranspose = PCMGResidualTransposeDefault; 18830b0564aSStefano Zampini mglevels[l]->matresidualtranspose = PCMGMatResidualTransposeDefault; 1899566063dSJacob Faibussowitsch if (mat) PetscCall(PetscObjectReference((PetscObject)mat)); 1909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->A)); 191fcb023d4SJed Brown mglevels[l]->A = mat; 192fcb023d4SJed Brown PetscFunctionReturn(0); 193fcb023d4SJed Brown } 194fcb023d4SJed Brown 1954b9ad928SBarry Smith /*@ 196aea2a34eSBarry Smith PCMGSetInterpolation - Sets the function to be used to calculate the 197bf5b2e24SBarry Smith interpolation from l-1 to the lth level 1984b9ad928SBarry Smith 199f1580f4eSBarry Smith Logically Collective on pc 2004b9ad928SBarry Smith 2014b9ad928SBarry Smith Input Parameters: 2024b9ad928SBarry Smith + pc - the multigrid context 2034b9ad928SBarry Smith . mat - the interpolation operator 204bf5b2e24SBarry Smith - l - the level (0 is coarsest) to supply [do not supply 0] 2054b9ad928SBarry Smith 2064b9ad928SBarry Smith Level: advanced 2074b9ad928SBarry Smith 2084b9ad928SBarry Smith Notes: 2094b9ad928SBarry Smith Usually this is the same matrix used also to set the restriction 2104b9ad928SBarry Smith for the same level. 2114b9ad928SBarry Smith 2124b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 2134b9ad928SBarry Smith out from the matrix size which one it is. 2144b9ad928SBarry Smith 215f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRestriction()` 2164b9ad928SBarry Smith @*/ 217d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetInterpolation(PC pc, PetscInt l, Mat mat) 218d71ae5a4SJacob Faibussowitsch { 219f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 220f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2214b9ad928SBarry Smith 2224b9ad928SBarry Smith PetscFunctionBegin; 223c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 22428b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 22528b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set interpolation routine for coarsest level"); 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->interpolate)); 2282fa5cd67SKarl Rupp 229f3fbd535SBarry Smith mglevels[l]->interpolate = mat; 2304b9ad928SBarry Smith PetscFunctionReturn(0); 2314b9ad928SBarry Smith } 2324b9ad928SBarry Smith 2338a2c336bSFande Kong /*@ 23495750439SFande Kong PCMGSetOperators - Sets operator and preconditioning matrix for lth level 2358a2c336bSFande Kong 236f1580f4eSBarry Smith Logically Collective on pc 2378a2c336bSFande Kong 2388a2c336bSFande Kong Input Parameters: 2398a2c336bSFande Kong + pc - the multigrid context 2408a2c336bSFande Kong . Amat - the operator 2418a2c336bSFande Kong . pmat - the preconditioning operator 24295750439SFande Kong - l - the level (0 is the coarsest) to supply 2438a2c336bSFande Kong 2448a2c336bSFande Kong Level: advanced 2458a2c336bSFande Kong 2468a2c336bSFande Kong .keywords: multigrid, set, interpolate, level 2478a2c336bSFande Kong 248f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetGalerkin()`, `PCMGSetRestriction()`, `PCMGSetInterpolation()` 2498a2c336bSFande Kong @*/ 250d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetOperators(PC pc, PetscInt l, Mat Amat, Mat Pmat) 251d71ae5a4SJacob Faibussowitsch { 252360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data; 253360ee056SFande Kong PC_MG_Levels **mglevels = mg->levels; 254360ee056SFande Kong 255360ee056SFande Kong PetscFunctionBegin; 256360ee056SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2578a2c336bSFande Kong PetscValidHeaderSpecific(Amat, MAT_CLASSID, 3); 2588a2c336bSFande Kong PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 4); 25928b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 2609566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(mglevels[l]->smoothd, Amat, Pmat)); 261360ee056SFande Kong PetscFunctionReturn(0); 262360ee056SFande Kong } 263360ee056SFande Kong 264c88c5224SJed Brown /*@ 265c88c5224SJed Brown PCMGGetInterpolation - Gets the function to be used to calculate the 266c88c5224SJed Brown interpolation from l-1 to the lth level 267c88c5224SJed Brown 268f1580f4eSBarry Smith Logically Collective on pc 269c88c5224SJed Brown 270c88c5224SJed Brown Input Parameters: 271c88c5224SJed Brown + pc - the multigrid context 272c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 273c88c5224SJed Brown 274c88c5224SJed Brown Output Parameter: 2753ad4599aSBarry Smith . mat - the interpolation matrix, can be NULL 276c88c5224SJed Brown 277c88c5224SJed Brown Level: advanced 278c88c5224SJed Brown 279f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()` 280c88c5224SJed Brown @*/ 281d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetInterpolation(PC pc, PetscInt l, Mat *mat) 282d71ae5a4SJacob Faibussowitsch { 283c88c5224SJed Brown PC_MG *mg = (PC_MG *)pc->data; 284c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 285c88c5224SJed Brown 286c88c5224SJed Brown PetscFunctionBegin; 287c88c5224SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2883ad4599aSBarry Smith if (mat) PetscValidPointer(mat, 3); 28928b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 2902472a847SBarry 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); 29148a46eb9SPierre Jolivet if (!mglevels[l]->interpolate && mglevels[l]->restrct) PetscCall(PCMGSetInterpolation(pc, l, mglevels[l]->restrct)); 2923ad4599aSBarry Smith if (mat) *mat = mglevels[l]->interpolate; 293c88c5224SJed Brown PetscFunctionReturn(0); 294c88c5224SJed Brown } 295c88c5224SJed Brown 2968a2c336bSFande Kong /*@ 297eab52d2bSLawrence Mitchell PCMGSetRestriction - Sets the function to be used to restrict dual vectors 2984b9ad928SBarry Smith from level l to l-1. 2994b9ad928SBarry Smith 300f1580f4eSBarry Smith Logically Collective on pc 3014b9ad928SBarry Smith 3024b9ad928SBarry Smith Input Parameters: 3034b9ad928SBarry Smith + pc - the multigrid context 304c88c5224SJed Brown . l - the level (0 is coarsest) to supply [Do not supply 0] 305c88c5224SJed Brown - mat - the restriction matrix 3064b9ad928SBarry Smith 3074b9ad928SBarry Smith Level: advanced 3084b9ad928SBarry Smith 3094b9ad928SBarry Smith Notes: 3104b9ad928SBarry Smith Usually this is the same matrix used also to set the interpolation 3114b9ad928SBarry Smith for the same level. 3124b9ad928SBarry Smith 3134b9ad928SBarry Smith One can pass in the interpolation matrix or its transpose; PETSc figures 3144b9ad928SBarry Smith out from the matrix size which one it is. 3154b9ad928SBarry Smith 316f1580f4eSBarry Smith If you do not set this, the transpose of the `Mat` set with `PCMGSetInterpolation()` 317fccaa45eSBarry Smith is used. 318fccaa45eSBarry Smith 319f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInterpolation()` 3204b9ad928SBarry Smith @*/ 321d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetRestriction(PC pc, PetscInt l, Mat mat) 322d71ae5a4SJacob Faibussowitsch { 323f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 324f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 3254b9ad928SBarry Smith 3264b9ad928SBarry Smith PetscFunctionBegin; 327c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 328c88c5224SJed Brown PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 32928b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 33028b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set restriction routine for coarsest level"); 3319566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 3329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->restrct)); 3332fa5cd67SKarl Rupp 334f3fbd535SBarry Smith mglevels[l]->restrct = mat; 3354b9ad928SBarry Smith PetscFunctionReturn(0); 3364b9ad928SBarry Smith } 3374b9ad928SBarry Smith 338c88c5224SJed Brown /*@ 339eab52d2bSLawrence Mitchell PCMGGetRestriction - Gets the function to be used to restrict dual vectors 340c88c5224SJed Brown from level l to l-1. 341c88c5224SJed Brown 342f1580f4eSBarry Smith Logically Collective on pc 343c88c5224SJed Brown 344c88c5224SJed Brown Input Parameters: 345c88c5224SJed Brown + pc - the multigrid context 346c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 347c88c5224SJed Brown 348c88c5224SJed Brown Output Parameter: 349c88c5224SJed Brown . mat - the restriction matrix 350c88c5224SJed Brown 351c88c5224SJed Brown Level: advanced 352c88c5224SJed Brown 353f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGGetInjection()` 354c88c5224SJed Brown @*/ 355d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetRestriction(PC pc, PetscInt l, Mat *mat) 356d71ae5a4SJacob Faibussowitsch { 357c88c5224SJed Brown PC_MG *mg = (PC_MG *)pc->data; 358c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 359c88c5224SJed Brown 360c88c5224SJed Brown PetscFunctionBegin; 361c88c5224SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3623ad4599aSBarry Smith if (mat) PetscValidPointer(mat, 3); 36328b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 3642472a847SBarry 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); 36548a46eb9SPierre Jolivet if (!mglevels[l]->restrct && mglevels[l]->interpolate) PetscCall(PCMGSetRestriction(pc, l, mglevels[l]->interpolate)); 3663ad4599aSBarry Smith if (mat) *mat = mglevels[l]->restrct; 367c88c5224SJed Brown PetscFunctionReturn(0); 368c88c5224SJed Brown } 369c88c5224SJed Brown 37073250ac0SBarry Smith /*@ 37173250ac0SBarry Smith PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 37273250ac0SBarry Smith 373f1580f4eSBarry Smith Logically Collective on pc 374c88c5224SJed Brown 375c88c5224SJed Brown Input Parameters: 376c88c5224SJed Brown + pc - the multigrid context 377c88c5224SJed Brown - l - the level (0 is coarsest) to supply [Do not supply 0] 378c88c5224SJed Brown . rscale - the scaling 379c88c5224SJed Brown 380c88c5224SJed Brown Level: advanced 381c88c5224SJed Brown 382f1580f4eSBarry Smith Note: 383f1580f4eSBarry 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. 384f1580f4eSBarry Smith It is preferable to use `PCMGSetInjection()` to control moving primal vectors. 385c88c5224SJed Brown 386f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGSetInjection()` 387c88c5224SJed Brown @*/ 388d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetRScale(PC pc, PetscInt l, Vec rscale) 389d71ae5a4SJacob Faibussowitsch { 390c88c5224SJed Brown PC_MG *mg = (PC_MG *)pc->data; 391c88c5224SJed Brown PC_MG_Levels **mglevels = mg->levels; 392c88c5224SJed Brown 393c88c5224SJed Brown PetscFunctionBegin; 394c88c5224SJed Brown PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 39528b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 3962472a847SBarry 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); 3979566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rscale)); 3989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[l]->rscale)); 3992fa5cd67SKarl Rupp 400c88c5224SJed Brown mglevels[l]->rscale = rscale; 401c88c5224SJed Brown PetscFunctionReturn(0); 402c88c5224SJed Brown } 403c88c5224SJed Brown 404c88c5224SJed Brown /*@ 405c88c5224SJed Brown PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 406c88c5224SJed Brown 407f1580f4eSBarry Smith Collective on pc 40873250ac0SBarry Smith 40973250ac0SBarry Smith Input Parameters: 41073250ac0SBarry Smith + pc - the multigrid context 41173250ac0SBarry Smith . rscale - the scaling 41273250ac0SBarry Smith - l - the level (0 is coarsest) to supply [Do not supply 0] 41373250ac0SBarry Smith 41473250ac0SBarry Smith Level: advanced 41573250ac0SBarry Smith 416f1580f4eSBarry Smith Note: 417f1580f4eSBarry 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. 418f1580f4eSBarry Smith It is preferable to use `PCMGGetInjection()` to control moving primal vectors. 41973250ac0SBarry Smith 420f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInterpolation()`, `PCMGGetRestriction()`, `PCMGGetInjection()` 42173250ac0SBarry Smith @*/ 422d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetRScale(PC pc, PetscInt l, Vec *rscale) 423d71ae5a4SJacob Faibussowitsch { 42473250ac0SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 42573250ac0SBarry Smith PC_MG_Levels **mglevels = mg->levels; 42673250ac0SBarry Smith 42773250ac0SBarry Smith PetscFunctionBegin; 42873250ac0SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 42928b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 4302472a847SBarry 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); 431c88c5224SJed Brown if (!mglevels[l]->rscale) { 432c88c5224SJed Brown Mat R; 433c88c5224SJed Brown Vec X, Y, coarse, fine; 434c88c5224SJed Brown PetscInt M, N; 435*0fdf79fbSJacob Faibussowitsch 4369566063dSJacob Faibussowitsch PetscCall(PCMGGetRestriction(pc, l, &R)); 4379566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(R, &X, &Y)); 4389566063dSJacob Faibussowitsch PetscCall(MatGetSize(R, &M, &N)); 439*0fdf79fbSJacob Faibussowitsch PetscCheck(N != M, PetscObjectComm((PetscObject)R), PETSC_ERR_SUP, "Restriction matrix is square, cannot determine which Vec is coarser"); 4402fa5cd67SKarl Rupp if (M < N) { 4412fa5cd67SKarl Rupp fine = X; 4422fa5cd67SKarl Rupp coarse = Y; 443*0fdf79fbSJacob Faibussowitsch } else { 4449371c9d4SSatish Balay fine = Y; 4459371c9d4SSatish Balay coarse = X; 446*0fdf79fbSJacob Faibussowitsch } 4479566063dSJacob Faibussowitsch PetscCall(VecSet(fine, 1.)); 4489566063dSJacob Faibussowitsch PetscCall(MatRestrict(R, fine, coarse)); 4499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fine)); 4509566063dSJacob Faibussowitsch PetscCall(VecReciprocal(coarse)); 451c88c5224SJed Brown mglevels[l]->rscale = coarse; 452c88c5224SJed Brown } 453c88c5224SJed Brown *rscale = mglevels[l]->rscale; 45473250ac0SBarry Smith PetscFunctionReturn(0); 45573250ac0SBarry Smith } 45673250ac0SBarry Smith 457f39d8e23SSatish Balay /*@ 458eab52d2bSLawrence Mitchell PCMGSetInjection - Sets the function to be used to inject primal vectors 459eab52d2bSLawrence Mitchell from level l to l-1. 460eab52d2bSLawrence Mitchell 461f1580f4eSBarry Smith Logically Collective on pc 462eab52d2bSLawrence Mitchell 463eab52d2bSLawrence Mitchell Input Parameters: 464eab52d2bSLawrence Mitchell + pc - the multigrid context 465eab52d2bSLawrence Mitchell . l - the level (0 is coarsest) to supply [Do not supply 0] 466eab52d2bSLawrence Mitchell - mat - the injection matrix 467eab52d2bSLawrence Mitchell 468eab52d2bSLawrence Mitchell Level: advanced 469eab52d2bSLawrence Mitchell 470f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRestriction()` 471eab52d2bSLawrence Mitchell @*/ 472d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetInjection(PC pc, PetscInt l, Mat mat) 473d71ae5a4SJacob Faibussowitsch { 474eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG *)pc->data; 475eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 476eab52d2bSLawrence Mitchell 477eab52d2bSLawrence Mitchell PetscFunctionBegin; 478eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 479eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 48028b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 48128b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set restriction routine for coarsest level"); 4829566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 4839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mglevels[l]->inject)); 484eab52d2bSLawrence Mitchell 485eab52d2bSLawrence Mitchell mglevels[l]->inject = mat; 486eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 487eab52d2bSLawrence Mitchell } 488eab52d2bSLawrence Mitchell 489eab52d2bSLawrence Mitchell /*@ 490eab52d2bSLawrence Mitchell PCMGGetInjection - Gets the function to be used to inject primal vectors 491eab52d2bSLawrence Mitchell from level l to l-1. 492eab52d2bSLawrence Mitchell 493f1580f4eSBarry Smith Logically Collective on pc 494eab52d2bSLawrence Mitchell 495eab52d2bSLawrence Mitchell Input Parameters: 496eab52d2bSLawrence Mitchell + pc - the multigrid context 497eab52d2bSLawrence Mitchell - l - the level (0 is coarsest) to supply [Do not supply 0] 498eab52d2bSLawrence Mitchell 499eab52d2bSLawrence Mitchell Output Parameter: 50099a38656SLawrence Mitchell . mat - the restriction matrix (may be NULL if no injection is available). 501eab52d2bSLawrence Mitchell 502eab52d2bSLawrence Mitchell Level: advanced 503eab52d2bSLawrence Mitchell 504f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInjection()`, `PCMGetGetRestriction()` 505eab52d2bSLawrence Mitchell @*/ 506d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetInjection(PC pc, PetscInt l, Mat *mat) 507d71ae5a4SJacob Faibussowitsch { 508eab52d2bSLawrence Mitchell PC_MG *mg = (PC_MG *)pc->data; 509eab52d2bSLawrence Mitchell PC_MG_Levels **mglevels = mg->levels; 510eab52d2bSLawrence Mitchell 511eab52d2bSLawrence Mitchell PetscFunctionBegin; 512eab52d2bSLawrence Mitchell PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 513eab52d2bSLawrence Mitchell if (mat) PetscValidPointer(mat, 3); 51428b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 5152472a847SBarry 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); 516eab52d2bSLawrence Mitchell if (mat) *mat = mglevels[l]->inject; 517eab52d2bSLawrence Mitchell PetscFunctionReturn(0); 518eab52d2bSLawrence Mitchell } 519eab52d2bSLawrence Mitchell 520eab52d2bSLawrence Mitchell /*@ 521f1580f4eSBarry Smith PCMGGetSmoother - Gets the `KSP` context to be used as smoother for 522f1580f4eSBarry Smith both pre- and post-smoothing. Call both `PCMGGetSmootherUp()` and 523f1580f4eSBarry Smith `PCMGGetSmootherDown()` to use different functions for pre- and 5244b9ad928SBarry Smith post-smoothing. 5254b9ad928SBarry Smith 526f1580f4eSBarry Smith Not Collective, ksp returned is parallel if pc is 5274b9ad928SBarry Smith 5284b9ad928SBarry Smith Input Parameters: 5294b9ad928SBarry Smith + pc - the multigrid context 5304b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5314b9ad928SBarry Smith 53201d2d390SJose E. Roman Output Parameter: 5334b9ad928SBarry Smith . ksp - the smoother 5344b9ad928SBarry Smith 535f1580f4eSBarry Smith Note: 536f1580f4eSBarry 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. 537f1580f4eSBarry Smith You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call `KSPGetPC`(ksp,&pc) 538f1580f4eSBarry Smith and modify PC options for the smoother; for example `PCSetType`(pc,`PCSOR`); to use SOR smoothing. 53957420d5bSBarry Smith 5404b9ad928SBarry Smith Level: advanced 5414b9ad928SBarry Smith 542f1580f4eSBarry Smith .seealso: PCMG`, ``PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, `PCMGGetCoarseSolve()` 5434b9ad928SBarry Smith @*/ 544d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetSmoother(PC pc, PetscInt l, KSP *ksp) 545d71ae5a4SJacob Faibussowitsch { 546f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 547f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 5484b9ad928SBarry Smith 5494b9ad928SBarry Smith PetscFunctionBegin; 550c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 551f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 5524b9ad928SBarry Smith PetscFunctionReturn(0); 5534b9ad928SBarry Smith } 5544b9ad928SBarry Smith 555f39d8e23SSatish Balay /*@ 55697177400SBarry Smith PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 5574b9ad928SBarry Smith coarse grid correction (post-smoother). 5584b9ad928SBarry Smith 559f1580f4eSBarry Smith Not Collective, ksp returned is parallel if pc is 5604b9ad928SBarry Smith 5614b9ad928SBarry Smith Input Parameters: 5624b9ad928SBarry Smith + pc - the multigrid context 5634b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 5644b9ad928SBarry Smith 56501d2d390SJose E. Roman Output Parameter: 5664b9ad928SBarry Smith . ksp - the smoother 5674b9ad928SBarry Smith 5684b9ad928SBarry Smith Level: advanced 5694b9ad928SBarry Smith 570f1580f4eSBarry Smith Note: 571f1580f4eSBarry Smith Calling this will result in a different pre and post smoother so you may need to set options on the pre smoother also 57289cce641SBarry Smith 573f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()` 5744b9ad928SBarry Smith @*/ 575d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetSmootherUp(PC pc, PetscInt l, KSP *ksp) 576d71ae5a4SJacob Faibussowitsch { 577f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 578f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 579f69a0ea3SMatthew Knepley const char *prefix; 5804b9ad928SBarry Smith MPI_Comm comm; 5814b9ad928SBarry Smith 5824b9ad928SBarry Smith PetscFunctionBegin; 583c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 5844b9ad928SBarry Smith /* 5854b9ad928SBarry Smith This is called only if user wants a different pre-smoother from post. 5864b9ad928SBarry Smith Thus we check if a different one has already been allocated, 5874b9ad928SBarry Smith if not we allocate it. 5884b9ad928SBarry Smith */ 58928b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "There is no such thing as a up smoother on the coarse grid"); 590f3fbd535SBarry Smith if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 59119fd82e9SBarry Smith KSPType ksptype; 59219fd82e9SBarry Smith PCType pctype; 593336babb1SJed Brown PC ipc; 594336babb1SJed Brown PetscReal rtol, abstol, dtol; 595336babb1SJed Brown PetscInt maxits; 596336babb1SJed Brown KSPNormType normtype; 5979566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mglevels[l]->smoothd, &comm)); 5989566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(mglevels[l]->smoothd, &prefix)); 5999566063dSJacob Faibussowitsch PetscCall(KSPGetTolerances(mglevels[l]->smoothd, &rtol, &abstol, &dtol, &maxits)); 6009566063dSJacob Faibussowitsch PetscCall(KSPGetType(mglevels[l]->smoothd, &ksptype)); 6019566063dSJacob Faibussowitsch PetscCall(KSPGetNormType(mglevels[l]->smoothd, &normtype)); 6029566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[l]->smoothd, &ipc)); 6039566063dSJacob Faibussowitsch PetscCall(PCGetType(ipc, &pctype)); 604336babb1SJed Brown 6059566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &mglevels[l]->smoothu)); 6069566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(mglevels[l]->smoothu, pc->erroriffailure)); 6079566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu, (PetscObject)pc, mglevels[0]->levels - l)); 6089566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(mglevels[l]->smoothu, prefix)); 6099566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(mglevels[l]->smoothu, rtol, abstol, dtol, maxits)); 6109566063dSJacob Faibussowitsch PetscCall(KSPSetType(mglevels[l]->smoothu, ksptype)); 6119566063dSJacob Faibussowitsch PetscCall(KSPSetNormType(mglevels[l]->smoothu, normtype)); 6129566063dSJacob Faibussowitsch PetscCall(KSPSetConvergenceTest(mglevels[l]->smoothu, KSPConvergedSkip, NULL, NULL)); 6139566063dSJacob Faibussowitsch PetscCall(KSPGetPC(mglevels[l]->smoothu, &ipc)); 6149566063dSJacob Faibussowitsch PetscCall(PCSetType(ipc, pctype)); 6159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level)); 6164b9ad928SBarry Smith } 617f3fbd535SBarry Smith if (ksp) *ksp = mglevels[l]->smoothu; 6184b9ad928SBarry Smith PetscFunctionReturn(0); 6194b9ad928SBarry Smith } 6204b9ad928SBarry Smith 621f39d8e23SSatish Balay /*@ 622f1580f4eSBarry Smith PCMGGetSmootherDown - Gets the `KSP` context to be used as smoother before 6234b9ad928SBarry Smith coarse grid correction (pre-smoother). 6244b9ad928SBarry Smith 625f1580f4eSBarry Smith Not Collective, ksp returned is parallel if pc is 6264b9ad928SBarry Smith 6274b9ad928SBarry Smith Input Parameters: 6284b9ad928SBarry Smith + pc - the multigrid context 6294b9ad928SBarry Smith - l - the level (0 is coarsest) to supply 6304b9ad928SBarry Smith 63101d2d390SJose E. Roman Output Parameter: 6324b9ad928SBarry Smith . ksp - the smoother 6334b9ad928SBarry Smith 6344b9ad928SBarry Smith Level: advanced 6354b9ad928SBarry Smith 636f1580f4eSBarry Smith Note: 637f1580f4eSBarry Smith Calling this will result in a different pre and post smoother so you may need to 63889cce641SBarry Smith set options on the post smoother also 63989cce641SBarry Smith 640f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmoother()` 6414b9ad928SBarry Smith @*/ 642d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetSmootherDown(PC pc, PetscInt l, KSP *ksp) 643d71ae5a4SJacob Faibussowitsch { 644f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 645f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6464b9ad928SBarry Smith 6474b9ad928SBarry Smith PetscFunctionBegin; 648c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6494b9ad928SBarry Smith /* make sure smoother up and down are different */ 6501baa6e33SBarry Smith if (l) PetscCall(PCMGGetSmootherUp(pc, l, NULL)); 651f3fbd535SBarry Smith *ksp = mglevels[l]->smoothd; 6524b9ad928SBarry Smith PetscFunctionReturn(0); 6534b9ad928SBarry Smith } 6544b9ad928SBarry Smith 6554b9ad928SBarry Smith /*@ 656cab31ae5SJed Brown PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 6574b9ad928SBarry Smith 658f1580f4eSBarry Smith Logically Collective on pc 6594b9ad928SBarry Smith 6604b9ad928SBarry Smith Input Parameters: 6614b9ad928SBarry Smith + pc - the multigrid context 662c1cbb1deSBarry Smith . l - the level (0 is coarsest) 663f1580f4eSBarry Smith - c - either `PC_MG_CYCLE_V` or `PC_MG_CYCLE_W` 6644b9ad928SBarry Smith 6654b9ad928SBarry Smith Level: advanced 6664b9ad928SBarry Smith 667f1580f4eSBarry Smith .seealso: `PCMG`, PCMGCycleType`, `PCMGSetCycleType()` 6684b9ad928SBarry Smith @*/ 669d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc, PetscInt l, PCMGCycleType c) 670d71ae5a4SJacob Faibussowitsch { 671f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 672f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 6734b9ad928SBarry Smith 6744b9ad928SBarry Smith PetscFunctionBegin; 675c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 67628b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 677c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, l, 2); 678c51679f6SSatish Balay PetscValidLogicalCollectiveEnum(pc, c, 3); 679f3fbd535SBarry Smith mglevels[l]->cycles = c; 6804b9ad928SBarry Smith PetscFunctionReturn(0); 6814b9ad928SBarry Smith } 6824b9ad928SBarry Smith 6834b9ad928SBarry Smith /*@ 684f3b08a26SMatthew G. Knepley PCMGSetRhs - Sets the vector to be used to store the right-hand side on a particular level. 6854b9ad928SBarry Smith 686f1580f4eSBarry Smith Logically Collective on pc 6874b9ad928SBarry Smith 6884b9ad928SBarry Smith Input Parameters: 6894b9ad928SBarry Smith + pc - the multigrid context 6904b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 691f3b08a26SMatthew G. Knepley - c - the Vec 6924b9ad928SBarry Smith 6934b9ad928SBarry Smith Level: advanced 6944b9ad928SBarry Smith 695f1580f4eSBarry Smith Note: 696f1580f4eSBarry 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. 697fccaa45eSBarry Smith 698f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetX()`, `PCMGSetR()` 6994b9ad928SBarry Smith @*/ 700d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetRhs(PC pc, PetscInt l, Vec c) 701d71ae5a4SJacob Faibussowitsch { 702f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 703f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7044b9ad928SBarry Smith 7054b9ad928SBarry Smith PetscFunctionBegin; 706c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 70728b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 70808401ef6SPierre Jolivet PetscCheck(l != mglevels[0]->levels - 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Do not set rhs for finest level"); 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)c)); 7109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[l]->b)); 7112fa5cd67SKarl Rupp 712f3fbd535SBarry Smith mglevels[l]->b = c; 7134b9ad928SBarry Smith PetscFunctionReturn(0); 7144b9ad928SBarry Smith } 7154b9ad928SBarry Smith 7164b9ad928SBarry Smith /*@ 717f3b08a26SMatthew G. Knepley PCMGSetX - Sets the vector to be used to store the solution on a particular level. 7184b9ad928SBarry Smith 719f1580f4eSBarry Smith Logically Collective on pc 7204b9ad928SBarry Smith 7214b9ad928SBarry Smith Input Parameters: 7224b9ad928SBarry Smith + pc - the multigrid context 723251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 724f3b08a26SMatthew G. Knepley - c - the Vec 7254b9ad928SBarry Smith 7264b9ad928SBarry Smith Level: advanced 7274b9ad928SBarry Smith 728f1580f4eSBarry Smith Note: 729f1580f4eSBarry 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. 730fccaa45eSBarry Smith 731f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRhs()`, `PCMGSetR()` 7324b9ad928SBarry Smith @*/ 733d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetX(PC pc, PetscInt l, Vec c) 734d71ae5a4SJacob Faibussowitsch { 735f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 736f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7374b9ad928SBarry Smith 7384b9ad928SBarry Smith PetscFunctionBegin; 739c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 74028b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 74108401ef6SPierre Jolivet PetscCheck(l != mglevels[0]->levels - 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Do not set x for finest level"); 7429566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)c)); 7439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[l]->x)); 7442fa5cd67SKarl Rupp 745f3fbd535SBarry Smith mglevels[l]->x = c; 7464b9ad928SBarry Smith PetscFunctionReturn(0); 7474b9ad928SBarry Smith } 7484b9ad928SBarry Smith 7494b9ad928SBarry Smith /*@ 750f3b08a26SMatthew G. Knepley PCMGSetR - Sets the vector to be used to store the residual on a particular level. 7514b9ad928SBarry Smith 752f1580f4eSBarry Smith Logically Collective on pc 7534b9ad928SBarry Smith 7544b9ad928SBarry Smith Input Parameters: 7554b9ad928SBarry Smith + pc - the multigrid context 7564b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for 757f3b08a26SMatthew G. Knepley - c - the Vec 7584b9ad928SBarry Smith 7594b9ad928SBarry Smith Level: advanced 7604b9ad928SBarry Smith 761f1580f4eSBarry Smith Note: 762f1580f4eSBarry 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. 763fccaa45eSBarry Smith 764f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRhs()`, `PCMGSetX()` 7654b9ad928SBarry Smith @*/ 766d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetR(PC pc, PetscInt l, Vec c) 767d71ae5a4SJacob Faibussowitsch { 768f3fbd535SBarry Smith PC_MG *mg = (PC_MG *)pc->data; 769f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7704b9ad928SBarry Smith 7714b9ad928SBarry Smith PetscFunctionBegin; 772c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 77328b400f6SJacob Faibussowitsch PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 77428b400f6SJacob Faibussowitsch PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Need not set residual vector for coarse grid"); 7759566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)c)); 7769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mglevels[l]->r)); 7772fa5cd67SKarl Rupp 778f3fbd535SBarry Smith mglevels[l]->r = c; 7794b9ad928SBarry Smith PetscFunctionReturn(0); 7804b9ad928SBarry Smith } 781