xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 226f8a8a5081bc6ad7227cd631662400f0d6e2a0)
1af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/
24b9ad928SBarry Smith 
3cc4c1da9SBarry Smith /*@
454b2cd4bSJed Brown   PCMGResidualDefault - Default routine to calculate the residual.
54b9ad928SBarry Smith 
6c3339decSBarry Smith   Collective
74b9ad928SBarry Smith 
84b9ad928SBarry Smith   Input Parameters:
94b9ad928SBarry Smith + mat - the matrix
10dd8e379bSPierre Jolivet . b   - the right-hand side
114b9ad928SBarry Smith - x   - the approximate solution
124b9ad928SBarry Smith 
134b9ad928SBarry Smith   Output Parameter:
144b9ad928SBarry Smith . r - location to store the residual
154b9ad928SBarry Smith 
16d0e4de75SBarry Smith   Level: developer
174b9ad928SBarry Smith 
18562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetResidual()`, `PCMGSetMatResidual()`
194b9ad928SBarry Smith @*/
PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r)20d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGResidualDefault(Mat mat, Vec b, Vec x, Vec r)
21d71ae5a4SJacob Faibussowitsch {
224b9ad928SBarry Smith   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(MatResidual(mat, b, x, r));
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
254b9ad928SBarry Smith }
264b9ad928SBarry Smith 
27cc4c1da9SBarry Smith /*@
28fcb023d4SJed Brown   PCMGResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system
29fcb023d4SJed Brown 
30c3339decSBarry Smith   Collective
31fcb023d4SJed Brown 
32fcb023d4SJed Brown   Input Parameters:
33fcb023d4SJed Brown + mat - the matrix
34dd8e379bSPierre Jolivet . 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 
42562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetResidualTranspose()`, `PCMGMatResidualTransposeDefault()`
43fcb023d4SJed Brown @*/
PCMGResidualTransposeDefault(Mat mat,Vec b,Vec x,Vec r)44d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGResidualTransposeDefault(Mat mat, Vec b, Vec x, Vec r)
45d71ae5a4SJacob Faibussowitsch {
46fcb023d4SJed Brown   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(mat, x, r));
489566063dSJacob Faibussowitsch   PetscCall(VecAYPX(r, -1.0, b));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50fcb023d4SJed Brown }
51fcb023d4SJed Brown 
52cc4c1da9SBarry Smith /*@
5330b0564aSStefano Zampini   PCMGMatResidualDefault - Default routine to calculate the residual.
5430b0564aSStefano Zampini 
55c3339decSBarry Smith   Collective
5630b0564aSStefano Zampini 
5730b0564aSStefano Zampini   Input Parameters:
5830b0564aSStefano Zampini + mat - the matrix
59dd8e379bSPierre Jolivet . b   - the right-hand side
6030b0564aSStefano Zampini - x   - the approximate solution
6130b0564aSStefano Zampini 
6230b0564aSStefano Zampini   Output Parameter:
6330b0564aSStefano Zampini . r - location to store the residual
6430b0564aSStefano Zampini 
6530b0564aSStefano Zampini   Level: developer
6630b0564aSStefano Zampini 
67562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetMatResidual()`, `PCMGResidualDefault()`
6830b0564aSStefano Zampini @*/
PCMGMatResidualDefault(Mat mat,Mat b,Mat x,Mat r)69d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMatResidualDefault(Mat mat, Mat b, Mat x, Mat r)
70d71ae5a4SJacob Faibussowitsch {
7130b0564aSStefano Zampini   PetscFunctionBegin;
72fb842aefSJose E. Roman   PetscCall(MatMatMult(mat, x, MAT_REUSE_MATRIX, PETSC_CURRENT, &r));
739566063dSJacob Faibussowitsch   PetscCall(MatAYPX(r, -1.0, b, UNKNOWN_NONZERO_PATTERN));
743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7530b0564aSStefano Zampini }
7630b0564aSStefano Zampini 
77cc4c1da9SBarry Smith /*@
7830b0564aSStefano Zampini   PCMGMatResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system
7930b0564aSStefano Zampini 
80c3339decSBarry Smith   Collective
8130b0564aSStefano Zampini 
8230b0564aSStefano Zampini   Input Parameters:
8330b0564aSStefano Zampini + mat - the matrix
84dd8e379bSPierre Jolivet . b   - the right-hand side
8530b0564aSStefano Zampini - x   - the approximate solution
8630b0564aSStefano Zampini 
8730b0564aSStefano Zampini   Output Parameter:
8830b0564aSStefano Zampini . r - location to store the residual
8930b0564aSStefano Zampini 
9030b0564aSStefano Zampini   Level: developer
9130b0564aSStefano Zampini 
92562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetMatResidualTranspose()`
9330b0564aSStefano Zampini @*/
PCMGMatResidualTransposeDefault(Mat mat,Mat b,Mat x,Mat r)94d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMatResidualTransposeDefault(Mat mat, Mat b, Mat x, Mat r)
95d71ae5a4SJacob Faibussowitsch {
9630b0564aSStefano Zampini   PetscFunctionBegin;
97fb842aefSJose E. Roman   PetscCall(MatTransposeMatMult(mat, x, MAT_REUSE_MATRIX, PETSC_CURRENT, &r));
989566063dSJacob Faibussowitsch   PetscCall(MatAYPX(r, -1.0, b, UNKNOWN_NONZERO_PATTERN));
993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10030b0564aSStefano Zampini }
101f39d8e23SSatish Balay /*@
10297177400SBarry Smith   PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
1034b9ad928SBarry Smith 
1044b9ad928SBarry Smith   Not Collective
1054b9ad928SBarry Smith 
1064b9ad928SBarry Smith   Input Parameter:
1074b9ad928SBarry Smith . pc - the multigrid context
1084b9ad928SBarry Smith 
1094b9ad928SBarry Smith   Output Parameter:
1104b9ad928SBarry Smith . ksp - the coarse grid solver context
1114b9ad928SBarry Smith 
1124b9ad928SBarry Smith   Level: advanced
1134b9ad928SBarry Smith 
114562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, `PCMGGetSmoother()`
1154b9ad928SBarry Smith @*/
PCMGGetCoarseSolve(PC pc,KSP * ksp)116d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetCoarseSolve(PC pc, KSP *ksp)
117d71ae5a4SJacob Faibussowitsch {
118f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
119f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith   PetscFunctionBegin;
122c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
123f3fbd535SBarry Smith   *ksp = mglevels[0]->smoothd;
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1254b9ad928SBarry Smith }
1264b9ad928SBarry Smith 
1274b9ad928SBarry Smith /*@C
128f1580f4eSBarry Smith   PCMGSetResidual - Sets the function to be used to calculate the residual on the lth level.
1294b9ad928SBarry Smith 
130c3339decSBarry Smith   Logically Collective
1314b9ad928SBarry Smith 
1324b9ad928SBarry Smith   Input Parameters:
1334b9ad928SBarry Smith + pc       - the multigrid context
1344b9ad928SBarry Smith . l        - the level (0 is coarsest) to supply
135157726a2SBarry Smith . residual - function used to form residual, if none is provided the previously provide one is used, if no
136d0e4de75SBarry Smith               previous one were provided then a default is used
1374b9ad928SBarry Smith - mat      - matrix associated with residual
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith   Level: advanced
1404b9ad928SBarry Smith 
141562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGResidualDefault()`
1424b9ad928SBarry Smith @*/
PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (* residual)(Mat,Vec,Vec,Vec),Mat mat)143d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetResidual(PC pc, PetscInt l, PetscErrorCode (*residual)(Mat, Vec, Vec, Vec), Mat mat)
144d71ae5a4SJacob Faibussowitsch {
145f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
146f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
1474b9ad928SBarry Smith 
1484b9ad928SBarry Smith   PetscFunctionBegin;
149c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15028b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
1512fa5cd67SKarl Rupp   if (residual) mglevels[l]->residual = residual;
15254b2cd4bSJed Brown   if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault;
15330b0564aSStefano Zampini   mglevels[l]->matresidual = PCMGMatResidualDefault;
1549566063dSJacob Faibussowitsch   if (mat) PetscCall(PetscObjectReference((PetscObject)mat));
1559566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->A));
156f3fbd535SBarry Smith   mglevels[l]->A = mat;
1573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1584b9ad928SBarry Smith }
1594b9ad928SBarry Smith 
160fcb023d4SJed Brown /*@C
161fcb023d4SJed Brown   PCMGSetResidualTranspose - Sets the function to be used to calculate the residual of the transposed linear system
162fcb023d4SJed Brown   on the lth level.
163fcb023d4SJed Brown 
164c3339decSBarry Smith   Logically Collective
165fcb023d4SJed Brown 
166fcb023d4SJed Brown   Input Parameters:
167fcb023d4SJed Brown + pc        - the multigrid context
168fcb023d4SJed Brown . l         - the level (0 is coarsest) to supply
169fcb023d4SJed Brown . residualt - function used to form transpose of residual, if none is provided the previously provide one is used, if no
170fcb023d4SJed Brown                previous one were provided then a default is used
171fcb023d4SJed Brown - mat       - matrix associated with residual
172fcb023d4SJed Brown 
173fcb023d4SJed Brown   Level: advanced
174fcb023d4SJed Brown 
175562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGResidualTransposeDefault()`
176fcb023d4SJed Brown @*/
PCMGSetResidualTranspose(PC pc,PetscInt l,PetscErrorCode (* residualt)(Mat,Vec,Vec,Vec),Mat mat)177d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetResidualTranspose(PC pc, PetscInt l, PetscErrorCode (*residualt)(Mat, Vec, Vec, Vec), Mat mat)
178d71ae5a4SJacob Faibussowitsch {
179fcb023d4SJed Brown   PC_MG         *mg       = (PC_MG *)pc->data;
180fcb023d4SJed Brown   PC_MG_Levels **mglevels = mg->levels;
181fcb023d4SJed Brown 
182fcb023d4SJed Brown   PetscFunctionBegin;
183fcb023d4SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
18428b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
185fcb023d4SJed Brown   if (residualt) mglevels[l]->residualtranspose = residualt;
186fcb023d4SJed Brown   if (!mglevels[l]->residualtranspose) mglevels[l]->residualtranspose = PCMGResidualTransposeDefault;
18730b0564aSStefano Zampini   mglevels[l]->matresidualtranspose = PCMGMatResidualTransposeDefault;
1889566063dSJacob Faibussowitsch   if (mat) PetscCall(PetscObjectReference((PetscObject)mat));
1899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->A));
190fcb023d4SJed Brown   mglevels[l]->A = mat;
1913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
192fcb023d4SJed Brown }
193fcb023d4SJed Brown 
1944b9ad928SBarry Smith /*@
195aea2a34eSBarry Smith   PCMGSetInterpolation - Sets the function to be used to calculate the
196bf5b2e24SBarry Smith   interpolation from l-1 to the lth level
1974b9ad928SBarry Smith 
198c3339decSBarry Smith   Logically Collective
1994b9ad928SBarry Smith 
2004b9ad928SBarry Smith   Input Parameters:
2014b9ad928SBarry Smith + pc  - the multigrid context
2024b9ad928SBarry Smith . mat - the interpolation operator
203bf5b2e24SBarry Smith - l   - the level (0 is coarsest) to supply [do not supply 0]
2044b9ad928SBarry Smith 
2054b9ad928SBarry Smith   Level: advanced
2064b9ad928SBarry Smith 
2074b9ad928SBarry Smith   Notes:
2084b9ad928SBarry Smith   Usually this is the same matrix used also to set the restriction
2094b9ad928SBarry Smith   for the same level.
2104b9ad928SBarry Smith 
2114b9ad928SBarry Smith   One can pass in the interpolation matrix or its transpose; PETSc figures
2124b9ad928SBarry Smith   out from the matrix size which one it is.
2134b9ad928SBarry Smith 
214562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetRestriction()`
2154b9ad928SBarry Smith @*/
PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)216d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetInterpolation(PC pc, PetscInt l, Mat mat)
217d71ae5a4SJacob Faibussowitsch {
218f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
219f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
2204b9ad928SBarry Smith 
2214b9ad928SBarry Smith   PetscFunctionBegin;
222c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
22328b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
22428b400f6SJacob Faibussowitsch   PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set interpolation routine for coarsest level");
2259566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
2269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->interpolate));
2272fa5cd67SKarl Rupp 
228f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
2293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2304b9ad928SBarry Smith }
2314b9ad928SBarry Smith 
2328a2c336bSFande Kong /*@
233*7addb90fSBarry Smith   PCMGSetOperators - Sets operator and matrix from which to construct a preconditioner for lth level
2348a2c336bSFande Kong 
235c3339decSBarry Smith   Logically Collective
2368a2c336bSFande Kong 
2378a2c336bSFande Kong   Input Parameters:
2388a2c336bSFande Kong + pc   - the multigrid context
2398a2c336bSFande Kong . Amat - the operator
240feefa0e1SJacob Faibussowitsch . Pmat - the preconditioning operator
24195750439SFande Kong - l    - the level (0 is the coarsest) to supply
2428a2c336bSFande Kong 
2438a2c336bSFande Kong   Level: advanced
2448a2c336bSFande Kong 
245562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetGalerkin()`, `PCMGSetRestriction()`, `PCMGSetInterpolation()`
2468a2c336bSFande Kong @*/
PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat)247d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetOperators(PC pc, PetscInt l, Mat Amat, Mat Pmat)
248d71ae5a4SJacob Faibussowitsch {
249360ee056SFande Kong   PC_MG         *mg       = (PC_MG *)pc->data;
250360ee056SFande Kong   PC_MG_Levels **mglevels = mg->levels;
251360ee056SFande Kong 
252360ee056SFande Kong   PetscFunctionBegin;
253360ee056SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2548a2c336bSFande Kong   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 3);
2558a2c336bSFande Kong   PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 4);
25628b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
2579566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(mglevels[l]->smoothd, Amat, Pmat));
2583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
259360ee056SFande Kong }
260360ee056SFande Kong 
261c88c5224SJed Brown /*@
262c88c5224SJed Brown   PCMGGetInterpolation - Gets the function to be used to calculate the
263c88c5224SJed Brown   interpolation from l-1 to the lth level
264c88c5224SJed Brown 
265c3339decSBarry Smith   Logically Collective
266c88c5224SJed Brown 
267c88c5224SJed Brown   Input Parameters:
268c88c5224SJed Brown + pc - the multigrid context
269c88c5224SJed Brown - l  - the level (0 is coarsest) to supply [Do not supply 0]
270c88c5224SJed Brown 
271c88c5224SJed Brown   Output Parameter:
2722fe279fdSBarry Smith . mat - the interpolation matrix, can be `NULL`
273c88c5224SJed Brown 
274c88c5224SJed Brown   Level: advanced
275c88c5224SJed Brown 
276562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`
277c88c5224SJed Brown @*/
PCMGGetInterpolation(PC pc,PetscInt l,Mat * mat)278d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetInterpolation(PC pc, PetscInt l, Mat *mat)
279d71ae5a4SJacob Faibussowitsch {
280c88c5224SJed Brown   PC_MG         *mg       = (PC_MG *)pc->data;
281c88c5224SJed Brown   PC_MG_Levels **mglevels = mg->levels;
282c88c5224SJed Brown 
283c88c5224SJed Brown   PetscFunctionBegin;
284c88c5224SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2854f572ea9SToby Isaac   if (mat) PetscAssertPointer(mat, 3);
28628b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
2872472a847SBarry 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);
28848a46eb9SPierre Jolivet   if (!mglevels[l]->interpolate && mglevels[l]->restrct) PetscCall(PCMGSetInterpolation(pc, l, mglevels[l]->restrct));
2893ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->interpolate;
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
291c88c5224SJed Brown }
292c88c5224SJed Brown 
2938a2c336bSFande Kong /*@
294eab52d2bSLawrence Mitchell   PCMGSetRestriction - Sets the function to be used to restrict dual vectors
2954b9ad928SBarry Smith   from level l to l-1.
2964b9ad928SBarry Smith 
297c3339decSBarry Smith   Logically Collective
2984b9ad928SBarry Smith 
2994b9ad928SBarry Smith   Input Parameters:
3004b9ad928SBarry Smith + pc  - the multigrid context
301c88c5224SJed Brown . l   - the level (0 is coarsest) to supply [Do not supply 0]
302c88c5224SJed Brown - mat - the restriction matrix
3034b9ad928SBarry Smith 
3044b9ad928SBarry Smith   Level: advanced
3054b9ad928SBarry Smith 
3064b9ad928SBarry Smith   Notes:
3074b9ad928SBarry Smith   Usually this is the same matrix used also to set the interpolation
3084b9ad928SBarry Smith   for the same level.
3094b9ad928SBarry Smith 
3104b9ad928SBarry Smith   One can pass in the interpolation matrix or its transpose; PETSc figures
3114b9ad928SBarry Smith   out from the matrix size which one it is.
3124b9ad928SBarry Smith 
313f1580f4eSBarry Smith   If you do not set this, the transpose of the `Mat` set with `PCMGSetInterpolation()`
314fccaa45eSBarry Smith   is used.
315fccaa45eSBarry Smith 
316562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetInterpolation()`
3174b9ad928SBarry Smith @*/
PCMGSetRestriction(PC pc,PetscInt l,Mat mat)318d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetRestriction(PC pc, PetscInt l, Mat mat)
319d71ae5a4SJacob Faibussowitsch {
320f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
321f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
3224b9ad928SBarry Smith 
3234b9ad928SBarry Smith   PetscFunctionBegin;
324c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
325c88c5224SJed Brown   PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
32628b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
32728b400f6SJacob Faibussowitsch   PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set restriction routine for coarsest level");
3289566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
3299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->restrct));
3302fa5cd67SKarl Rupp 
331f3fbd535SBarry Smith   mglevels[l]->restrct = mat;
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3334b9ad928SBarry Smith }
3344b9ad928SBarry Smith 
335c88c5224SJed Brown /*@
336ffbd2f08SBarry Smith   PCMGGetRestriction - Gets the function to be used to restrict dual (i.e. residual) vectors
337c88c5224SJed Brown   from level l to l-1.
338c88c5224SJed Brown 
339c3339decSBarry Smith   Logically Collective
340c88c5224SJed Brown 
341c88c5224SJed Brown   Input Parameters:
342c88c5224SJed Brown + pc - the multigrid context
343c88c5224SJed Brown - l  - the level (0 is coarsest) to supply [Do not supply 0]
344c88c5224SJed Brown 
345c88c5224SJed Brown   Output Parameter:
346c88c5224SJed Brown . mat - the restriction matrix
347c88c5224SJed Brown 
348c88c5224SJed Brown   Level: advanced
349c88c5224SJed Brown 
350562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGGetInjection()`
351c88c5224SJed Brown @*/
PCMGGetRestriction(PC pc,PetscInt l,Mat * mat)352d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetRestriction(PC pc, PetscInt l, Mat *mat)
353d71ae5a4SJacob Faibussowitsch {
354c88c5224SJed Brown   PC_MG         *mg       = (PC_MG *)pc->data;
355c88c5224SJed Brown   PC_MG_Levels **mglevels = mg->levels;
356c88c5224SJed Brown 
357c88c5224SJed Brown   PetscFunctionBegin;
358c88c5224SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3594f572ea9SToby Isaac   if (mat) PetscAssertPointer(mat, 3);
36028b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
3612472a847SBarry 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);
36248a46eb9SPierre Jolivet   if (!mglevels[l]->restrct && mglevels[l]->interpolate) PetscCall(PCMGSetRestriction(pc, l, mglevels[l]->interpolate));
3633ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->restrct;
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
365c88c5224SJed Brown }
366c88c5224SJed Brown 
36773250ac0SBarry Smith /*@
36873250ac0SBarry Smith   PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.
36973250ac0SBarry Smith 
370c3339decSBarry Smith   Logically Collective
371c88c5224SJed Brown 
372c88c5224SJed Brown   Input Parameters:
373c88c5224SJed Brown + pc     - the multigrid context
374feefa0e1SJacob Faibussowitsch . l      - the level (0 is coarsest) to supply [Do not supply 0]
375feefa0e1SJacob Faibussowitsch - rscale - the scaling
376c88c5224SJed Brown 
377c88c5224SJed Brown   Level: advanced
378c88c5224SJed Brown 
379f1580f4eSBarry Smith   Note:
380f1580f4eSBarry 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.
381f1580f4eSBarry Smith   It is preferable to use `PCMGSetInjection()` to control moving primal vectors.
382c88c5224SJed Brown 
383562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGSetInjection()`
384c88c5224SJed Brown @*/
PCMGSetRScale(PC pc,PetscInt l,Vec rscale)385d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetRScale(PC pc, PetscInt l, Vec rscale)
386d71ae5a4SJacob Faibussowitsch {
387c88c5224SJed Brown   PC_MG         *mg       = (PC_MG *)pc->data;
388c88c5224SJed Brown   PC_MG_Levels **mglevels = mg->levels;
389c88c5224SJed Brown 
390c88c5224SJed Brown   PetscFunctionBegin;
391c88c5224SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
39228b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
3932472a847SBarry 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);
3949566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)rscale));
3959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mglevels[l]->rscale));
3962fa5cd67SKarl Rupp 
397c88c5224SJed Brown   mglevels[l]->rscale = rscale;
3983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
399c88c5224SJed Brown }
400c88c5224SJed Brown 
401c88c5224SJed Brown /*@
402c88c5224SJed Brown   PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1.
403c88c5224SJed Brown 
404c3339decSBarry Smith   Collective
40573250ac0SBarry Smith 
40673250ac0SBarry Smith   Input Parameters:
40773250ac0SBarry Smith + pc     - the multigrid context
40873250ac0SBarry Smith . rscale - the scaling
40973250ac0SBarry Smith - l      - the level (0 is coarsest) to supply [Do not supply 0]
41073250ac0SBarry Smith 
41173250ac0SBarry Smith   Level: advanced
41273250ac0SBarry Smith 
413f1580f4eSBarry Smith   Note:
414f1580f4eSBarry 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.
415f1580f4eSBarry Smith   It is preferable to use `PCMGGetInjection()` to control moving primal vectors.
41673250ac0SBarry Smith 
417562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetInterpolation()`, `PCMGGetRestriction()`, `PCMGGetInjection()`
41873250ac0SBarry Smith @*/
PCMGGetRScale(PC pc,PetscInt l,Vec * rscale)419d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetRScale(PC pc, PetscInt l, Vec *rscale)
420d71ae5a4SJacob Faibussowitsch {
42173250ac0SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
42273250ac0SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
42373250ac0SBarry Smith 
42473250ac0SBarry Smith   PetscFunctionBegin;
42573250ac0SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
42628b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
4272472a847SBarry 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);
428c88c5224SJed Brown   if (!mglevels[l]->rscale) {
429c88c5224SJed Brown     Mat      R;
430c88c5224SJed Brown     Vec      X, Y, coarse, fine;
431c88c5224SJed Brown     PetscInt M, N;
4320fdf79fbSJacob Faibussowitsch 
4339566063dSJacob Faibussowitsch     PetscCall(PCMGGetRestriction(pc, l, &R));
4349566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(R, &X, &Y));
4359566063dSJacob Faibussowitsch     PetscCall(MatGetSize(R, &M, &N));
4360fdf79fbSJacob Faibussowitsch     PetscCheck(N != M, PetscObjectComm((PetscObject)R), PETSC_ERR_SUP, "Restriction matrix is square, cannot determine which Vec is coarser");
4372fa5cd67SKarl Rupp     if (M < N) {
4382fa5cd67SKarl Rupp       fine   = X;
4392fa5cd67SKarl Rupp       coarse = Y;
4400fdf79fbSJacob Faibussowitsch     } else {
4419371c9d4SSatish Balay       fine   = Y;
4429371c9d4SSatish Balay       coarse = X;
4430fdf79fbSJacob Faibussowitsch     }
4449566063dSJacob Faibussowitsch     PetscCall(VecSet(fine, 1.));
4459566063dSJacob Faibussowitsch     PetscCall(MatRestrict(R, fine, coarse));
4469566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&fine));
4479566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(coarse));
448c88c5224SJed Brown     mglevels[l]->rscale = coarse;
449c88c5224SJed Brown   }
450c88c5224SJed Brown   *rscale = mglevels[l]->rscale;
4513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45273250ac0SBarry Smith }
45373250ac0SBarry Smith 
454f39d8e23SSatish Balay /*@
455ffbd2f08SBarry Smith   PCMGSetInjection - Sets the function to be used to inject primal (i.e. solution) vectors
456eab52d2bSLawrence Mitchell   from level l to l-1.
457eab52d2bSLawrence Mitchell 
458c3339decSBarry Smith   Logically Collective
459eab52d2bSLawrence Mitchell 
460eab52d2bSLawrence Mitchell   Input Parameters:
461eab52d2bSLawrence Mitchell + pc  - the multigrid context
462eab52d2bSLawrence Mitchell . l   - the level (0 is coarsest) to supply [Do not supply 0]
463eab52d2bSLawrence Mitchell - mat - the injection matrix
464eab52d2bSLawrence Mitchell 
465eab52d2bSLawrence Mitchell   Level: advanced
466eab52d2bSLawrence Mitchell 
467562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetRestriction()`
468eab52d2bSLawrence Mitchell @*/
PCMGSetInjection(PC pc,PetscInt l,Mat mat)469d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetInjection(PC pc, PetscInt l, Mat mat)
470d71ae5a4SJacob Faibussowitsch {
471eab52d2bSLawrence Mitchell   PC_MG         *mg       = (PC_MG *)pc->data;
472eab52d2bSLawrence Mitchell   PC_MG_Levels **mglevels = mg->levels;
473eab52d2bSLawrence Mitchell 
474eab52d2bSLawrence Mitchell   PetscFunctionBegin;
475eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
476eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
47728b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
47828b400f6SJacob Faibussowitsch   PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set restriction routine for coarsest level");
4799566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
4809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->inject));
481eab52d2bSLawrence Mitchell 
482eab52d2bSLawrence Mitchell   mglevels[l]->inject = mat;
4833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
484eab52d2bSLawrence Mitchell }
485eab52d2bSLawrence Mitchell 
486eab52d2bSLawrence Mitchell /*@
487ffbd2f08SBarry Smith   PCMGGetInjection - Gets the function to be used to inject primal vectors (i.e. solutions)
488eab52d2bSLawrence Mitchell   from level l to l-1.
489eab52d2bSLawrence Mitchell 
490c3339decSBarry Smith   Logically Collective
491eab52d2bSLawrence Mitchell 
492eab52d2bSLawrence Mitchell   Input Parameters:
493eab52d2bSLawrence Mitchell + pc - the multigrid context
494eab52d2bSLawrence Mitchell - l  - the level (0 is coarsest) to supply [Do not supply 0]
495eab52d2bSLawrence Mitchell 
496eab52d2bSLawrence Mitchell   Output Parameter:
497dc9a610eSPierre Jolivet . mat - the restriction matrix (may be `NULL` if no injection is available).
498eab52d2bSLawrence Mitchell 
499eab52d2bSLawrence Mitchell   Level: advanced
500eab52d2bSLawrence Mitchell 
501562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetInjection()`, `PCMGetGetRestriction()`
502eab52d2bSLawrence Mitchell @*/
PCMGGetInjection(PC pc,PetscInt l,Mat * mat)503d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetInjection(PC pc, PetscInt l, Mat *mat)
504d71ae5a4SJacob Faibussowitsch {
505eab52d2bSLawrence Mitchell   PC_MG         *mg       = (PC_MG *)pc->data;
506eab52d2bSLawrence Mitchell   PC_MG_Levels **mglevels = mg->levels;
507eab52d2bSLawrence Mitchell 
508eab52d2bSLawrence Mitchell   PetscFunctionBegin;
509eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5104f572ea9SToby Isaac   if (mat) PetscAssertPointer(mat, 3);
51128b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
5122472a847SBarry 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);
513eab52d2bSLawrence Mitchell   if (mat) *mat = mglevels[l]->inject;
5143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
515eab52d2bSLawrence Mitchell }
516eab52d2bSLawrence Mitchell 
517eab52d2bSLawrence Mitchell /*@
518f1580f4eSBarry Smith   PCMGGetSmoother - Gets the `KSP` context to be used as smoother for
519f1580f4eSBarry Smith   both pre- and post-smoothing.  Call both `PCMGGetSmootherUp()` and
520f1580f4eSBarry Smith   `PCMGGetSmootherDown()` to use different functions for pre- and
5214b9ad928SBarry Smith   post-smoothing.
5224b9ad928SBarry Smith 
523f1580f4eSBarry Smith   Not Collective, ksp returned is parallel if pc is
5244b9ad928SBarry Smith 
5254b9ad928SBarry Smith   Input Parameters:
5264b9ad928SBarry Smith + pc - the multigrid context
5274b9ad928SBarry Smith - l  - the level (0 is coarsest) to supply
5284b9ad928SBarry Smith 
52901d2d390SJose E. Roman   Output Parameter:
5304b9ad928SBarry Smith . ksp - the smoother
5314b9ad928SBarry Smith 
532f1580f4eSBarry Smith   Note:
533f1580f4eSBarry 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.
534f1580f4eSBarry Smith   You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call `KSPGetPC`(ksp,&pc)
535f1580f4eSBarry Smith   and modify PC options for the smoother; for example `PCSetType`(pc,`PCSOR`); to use SOR smoothing.
53657420d5bSBarry Smith 
5374b9ad928SBarry Smith   Level: advanced
5384b9ad928SBarry Smith 
539a94f484eSPierre Jolivet .seealso: [](ch_ksp), `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, `PCMGGetCoarseSolve()`
5404b9ad928SBarry Smith @*/
PCMGGetSmoother(PC pc,PetscInt l,KSP * ksp)541d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetSmoother(PC pc, PetscInt l, KSP *ksp)
542d71ae5a4SJacob Faibussowitsch {
543f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
544f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
5454b9ad928SBarry Smith 
5464b9ad928SBarry Smith   PetscFunctionBegin;
547c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
548f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
5493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5504b9ad928SBarry Smith }
5514b9ad928SBarry Smith 
552f39d8e23SSatish Balay /*@
55397177400SBarry Smith   PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
5544b9ad928SBarry Smith   coarse grid correction (post-smoother).
5554b9ad928SBarry Smith 
556f1580f4eSBarry Smith   Not Collective, ksp returned is parallel if pc is
5574b9ad928SBarry Smith 
5584b9ad928SBarry Smith   Input Parameters:
5594b9ad928SBarry Smith + pc - the multigrid context
5604b9ad928SBarry Smith - l  - the level (0 is coarsest) to supply
5614b9ad928SBarry Smith 
56201d2d390SJose E. Roman   Output Parameter:
5634b9ad928SBarry Smith . ksp - the smoother
5644b9ad928SBarry Smith 
5654b9ad928SBarry Smith   Level: advanced
5664b9ad928SBarry Smith 
567f1580f4eSBarry Smith   Note:
568f1580f4eSBarry Smith   Calling this will result in a different pre and post smoother so you may need to set options on the pre smoother also
56989cce641SBarry Smith 
570562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetSmootherDown()`
5714b9ad928SBarry Smith @*/
PCMGGetSmootherUp(PC pc,PetscInt l,KSP * ksp)572d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetSmootherUp(PC pc, PetscInt l, KSP *ksp)
573d71ae5a4SJacob Faibussowitsch {
574f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
575f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
576f69a0ea3SMatthew Knepley   const char    *prefix;
5774b9ad928SBarry Smith   MPI_Comm       comm;
5784b9ad928SBarry Smith 
5794b9ad928SBarry Smith   PetscFunctionBegin;
580c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5814b9ad928SBarry Smith   /*
5824b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
5834b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
5844b9ad928SBarry Smith      if not we allocate it.
5854b9ad928SBarry Smith   */
58628b400f6SJacob Faibussowitsch   PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "There is no such thing as a up smoother on the coarse grid");
587f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
58819fd82e9SBarry Smith     KSPType     ksptype;
58919fd82e9SBarry Smith     PCType      pctype;
590336babb1SJed Brown     PC          ipc;
591336babb1SJed Brown     PetscReal   rtol, abstol, dtol;
592336babb1SJed Brown     PetscInt    maxits;
593336babb1SJed Brown     KSPNormType normtype;
5949566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)mglevels[l]->smoothd, &comm));
5959566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(mglevels[l]->smoothd, &prefix));
5969566063dSJacob Faibussowitsch     PetscCall(KSPGetTolerances(mglevels[l]->smoothd, &rtol, &abstol, &dtol, &maxits));
5979566063dSJacob Faibussowitsch     PetscCall(KSPGetType(mglevels[l]->smoothd, &ksptype));
5989566063dSJacob Faibussowitsch     PetscCall(KSPGetNormType(mglevels[l]->smoothd, &normtype));
5999566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(mglevels[l]->smoothd, &ipc));
6009566063dSJacob Faibussowitsch     PetscCall(PCGetType(ipc, &pctype));
601336babb1SJed Brown 
6029566063dSJacob Faibussowitsch     PetscCall(KSPCreate(comm, &mglevels[l]->smoothu));
6033821be0aSBarry Smith     PetscCall(KSPSetNestLevel(mglevels[l]->smoothu, pc->kspnestlevel));
6049566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(mglevels[l]->smoothu, pc->erroriffailure));
6059566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu, (PetscObject)pc, mglevels[0]->levels - l));
6069566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(mglevels[l]->smoothu, prefix));
6079566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[l]->smoothu, rtol, abstol, dtol, maxits));
6089566063dSJacob Faibussowitsch     PetscCall(KSPSetType(mglevels[l]->smoothu, ksptype));
6099566063dSJacob Faibussowitsch     PetscCall(KSPSetNormType(mglevels[l]->smoothu, normtype));
6109566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(mglevels[l]->smoothu, KSPConvergedSkip, NULL, NULL));
6119566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(mglevels[l]->smoothu, &ipc));
6129566063dSJacob Faibussowitsch     PetscCall(PCSetType(ipc, pctype));
6139566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level));
6144b9ad928SBarry Smith   }
615f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6174b9ad928SBarry Smith }
6184b9ad928SBarry Smith 
619f39d8e23SSatish Balay /*@
620f1580f4eSBarry Smith   PCMGGetSmootherDown - Gets the `KSP` context to be used as smoother before
6214b9ad928SBarry Smith   coarse grid correction (pre-smoother).
6224b9ad928SBarry Smith 
623f1580f4eSBarry Smith   Not Collective, ksp returned is parallel if pc is
6244b9ad928SBarry Smith 
6254b9ad928SBarry Smith   Input Parameters:
6264b9ad928SBarry Smith + pc - the multigrid context
6274b9ad928SBarry Smith - l  - the level (0 is coarsest) to supply
6284b9ad928SBarry Smith 
62901d2d390SJose E. Roman   Output Parameter:
6304b9ad928SBarry Smith . ksp - the smoother
6314b9ad928SBarry Smith 
6324b9ad928SBarry Smith   Level: advanced
6334b9ad928SBarry Smith 
634f1580f4eSBarry Smith   Note:
635f1580f4eSBarry Smith   Calling this will result in a different pre and post smoother so you may need to
63689cce641SBarry Smith   set options on the post smoother also
63789cce641SBarry Smith 
638562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmoother()`
6394b9ad928SBarry Smith @*/
PCMGGetSmootherDown(PC pc,PetscInt l,KSP * ksp)640d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetSmootherDown(PC pc, PetscInt l, KSP *ksp)
641d71ae5a4SJacob Faibussowitsch {
642f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
643f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
6444b9ad928SBarry Smith 
6454b9ad928SBarry Smith   PetscFunctionBegin;
646c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
6474b9ad928SBarry Smith   /* make sure smoother up and down are different */
6481baa6e33SBarry Smith   if (l) PetscCall(PCMGGetSmootherUp(pc, l, NULL));
649f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
6503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6514b9ad928SBarry Smith }
6524b9ad928SBarry Smith 
6534b9ad928SBarry Smith /*@
654cab31ae5SJed Brown   PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level.
6554b9ad928SBarry Smith 
656c3339decSBarry Smith   Logically Collective
6574b9ad928SBarry Smith 
6584b9ad928SBarry Smith   Input Parameters:
6594b9ad928SBarry Smith + pc - the multigrid context
660c1cbb1deSBarry Smith . l  - the level (0 is coarsest)
661f1580f4eSBarry Smith - c  - either `PC_MG_CYCLE_V` or `PC_MG_CYCLE_W`
6624b9ad928SBarry Smith 
6634b9ad928SBarry Smith   Level: advanced
6644b9ad928SBarry Smith 
665562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGCycleType`, `PCMGSetCycleType()`
6664b9ad928SBarry Smith @*/
PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c)667d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc, PetscInt l, PCMGCycleType c)
668d71ae5a4SJacob Faibussowitsch {
669f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
670f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
6714b9ad928SBarry Smith 
6724b9ad928SBarry Smith   PetscFunctionBegin;
673c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
67428b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
675c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc, l, 2);
676c51679f6SSatish Balay   PetscValidLogicalCollectiveEnum(pc, c, 3);
677f3fbd535SBarry Smith   mglevels[l]->cycles = c;
6783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6794b9ad928SBarry Smith }
6804b9ad928SBarry Smith 
6814b9ad928SBarry Smith /*@
682f3b08a26SMatthew G. Knepley   PCMGSetRhs - Sets the vector to be used to store the right-hand side on a particular level.
6834b9ad928SBarry Smith 
684c3339decSBarry Smith   Logically Collective
6854b9ad928SBarry Smith 
6864b9ad928SBarry Smith   Input Parameters:
6874b9ad928SBarry Smith + pc - the multigrid context
6884b9ad928SBarry Smith . l  - the level (0 is coarsest) this is to be used for
689f3b08a26SMatthew G. Knepley - c  - the Vec
6904b9ad928SBarry Smith 
6914b9ad928SBarry Smith   Level: advanced
6924b9ad928SBarry Smith 
693f1580f4eSBarry Smith   Note:
694f1580f4eSBarry 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.
695fccaa45eSBarry Smith 
696562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetX()`, `PCMGSetR()`
6974b9ad928SBarry Smith @*/
PCMGSetRhs(PC pc,PetscInt l,Vec c)698d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetRhs(PC pc, PetscInt l, Vec c)
699d71ae5a4SJacob Faibussowitsch {
700f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
701f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
7024b9ad928SBarry Smith 
7034b9ad928SBarry Smith   PetscFunctionBegin;
704c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
70528b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
70608401ef6SPierre Jolivet   PetscCheck(l != mglevels[0]->levels - 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Do not set rhs for finest level");
7079566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)c));
7089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mglevels[l]->b));
7092fa5cd67SKarl Rupp 
710f3fbd535SBarry Smith   mglevels[l]->b = c;
7113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7124b9ad928SBarry Smith }
7134b9ad928SBarry Smith 
7144b9ad928SBarry Smith /*@
715f3b08a26SMatthew G. Knepley   PCMGSetX - Sets the vector to be used to store the solution on a particular level.
7164b9ad928SBarry Smith 
717c3339decSBarry Smith   Logically Collective
7184b9ad928SBarry Smith 
7194b9ad928SBarry Smith   Input Parameters:
7204b9ad928SBarry Smith + pc - the multigrid context
721251f4c67SDmitry Karpeev . l  - the level (0 is coarsest) this is to be used for (do not supply the finest level)
722f3b08a26SMatthew G. Knepley - c  - the Vec
7234b9ad928SBarry Smith 
7244b9ad928SBarry Smith   Level: advanced
7254b9ad928SBarry Smith 
726f1580f4eSBarry Smith   Note:
727f1580f4eSBarry 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.
728fccaa45eSBarry Smith 
729562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetRhs()`, `PCMGSetR()`
7304b9ad928SBarry Smith @*/
PCMGSetX(PC pc,PetscInt l,Vec c)731d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetX(PC pc, PetscInt l, Vec c)
732d71ae5a4SJacob Faibussowitsch {
733f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
734f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
7354b9ad928SBarry Smith 
7364b9ad928SBarry Smith   PetscFunctionBegin;
737c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
73828b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
73908401ef6SPierre Jolivet   PetscCheck(l != mglevels[0]->levels - 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Do not set x for finest level");
7409566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)c));
7419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mglevels[l]->x));
7422fa5cd67SKarl Rupp 
743f3fbd535SBarry Smith   mglevels[l]->x = c;
7443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7454b9ad928SBarry Smith }
7464b9ad928SBarry Smith 
7474b9ad928SBarry Smith /*@
748f3b08a26SMatthew G. Knepley   PCMGSetR - Sets the vector to be used to store the residual on a particular level.
7494b9ad928SBarry Smith 
750c3339decSBarry Smith   Logically Collective
7514b9ad928SBarry Smith 
7524b9ad928SBarry Smith   Input Parameters:
7534b9ad928SBarry Smith + pc - the multigrid context
7544b9ad928SBarry Smith . l  - the level (0 is coarsest) this is to be used for
755f3b08a26SMatthew G. Knepley - c  - the Vec
7564b9ad928SBarry Smith 
7574b9ad928SBarry Smith   Level: advanced
7584b9ad928SBarry Smith 
759f1580f4eSBarry Smith   Note:
760f1580f4eSBarry 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.
761fccaa45eSBarry Smith 
762562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetRhs()`, `PCMGSetX()`
7634b9ad928SBarry Smith @*/
PCMGSetR(PC pc,PetscInt l,Vec c)764d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetR(PC pc, PetscInt l, Vec c)
765d71ae5a4SJacob Faibussowitsch {
766f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
767f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
7684b9ad928SBarry Smith 
7694b9ad928SBarry Smith   PetscFunctionBegin;
770c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
77128b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
77228b400f6SJacob Faibussowitsch   PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Need not set residual vector for coarse grid");
7739566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)c));
7749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mglevels[l]->r));
7752fa5cd67SKarl Rupp 
776f3fbd535SBarry Smith   mglevels[l]->r = c;
7773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7784b9ad928SBarry Smith }
779