xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision feefa0e191a340680bb02e1467a36facdcb0b150)
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 
7c3339decSBarry Smith   Collective
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));
253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
264b9ad928SBarry Smith }
274b9ad928SBarry Smith 
28fcb023d4SJed Brown /*@C
29fcb023d4SJed Brown   PCMGResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system
30fcb023d4SJed Brown 
31c3339decSBarry Smith   Collective
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));
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51fcb023d4SJed Brown }
52fcb023d4SJed Brown 
5330b0564aSStefano Zampini /*@C
5430b0564aSStefano Zampini   PCMGMatResidualDefault - Default routine to calculate the residual.
5530b0564aSStefano Zampini 
56c3339decSBarry Smith   Collective
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));
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7630b0564aSStefano Zampini }
7730b0564aSStefano Zampini 
7830b0564aSStefano Zampini /*@C
7930b0564aSStefano Zampini   PCMGMatResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system
8030b0564aSStefano Zampini 
81c3339decSBarry Smith   Collective
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));
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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;
1253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 
131c3339decSBarry Smith   Logically Collective
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;
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 
165c3339decSBarry Smith   Logically Collective
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;
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 
199c3339decSBarry Smith   Logically Collective
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;
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2314b9ad928SBarry Smith }
2324b9ad928SBarry Smith 
2338a2c336bSFande Kong /*@
23495750439SFande Kong   PCMGSetOperators - Sets operator and preconditioning matrix for lth level
2358a2c336bSFande Kong 
236c3339decSBarry Smith   Logically Collective
2378a2c336bSFande Kong 
2388a2c336bSFande Kong   Input Parameters:
2398a2c336bSFande Kong + pc   - the multigrid context
2408a2c336bSFande Kong . Amat - the operator
241*feefa0e1SJacob Faibussowitsch . Pmat - the preconditioning operator
24295750439SFande Kong - l    - the level (0 is the coarsest) to supply
2438a2c336bSFande Kong 
2448a2c336bSFande Kong   Level: advanced
2458a2c336bSFande Kong 
246f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetGalerkin()`, `PCMGSetRestriction()`, `PCMGSetInterpolation()`
2478a2c336bSFande Kong @*/
248d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetOperators(PC pc, PetscInt l, Mat Amat, Mat Pmat)
249d71ae5a4SJacob Faibussowitsch {
250360ee056SFande Kong   PC_MG         *mg       = (PC_MG *)pc->data;
251360ee056SFande Kong   PC_MG_Levels **mglevels = mg->levels;
252360ee056SFande Kong 
253360ee056SFande Kong   PetscFunctionBegin;
254360ee056SFande Kong   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2558a2c336bSFande Kong   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 3);
2568a2c336bSFande Kong   PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 4);
25728b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
2589566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(mglevels[l]->smoothd, Amat, Pmat));
2593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
260360ee056SFande Kong }
261360ee056SFande Kong 
262c88c5224SJed Brown /*@
263c88c5224SJed Brown   PCMGGetInterpolation - Gets the function to be used to calculate the
264c88c5224SJed Brown   interpolation from l-1 to the lth level
265c88c5224SJed Brown 
266c3339decSBarry Smith   Logically Collective
267c88c5224SJed Brown 
268c88c5224SJed Brown   Input Parameters:
269c88c5224SJed Brown + pc - the multigrid context
270c88c5224SJed Brown - l  - the level (0 is coarsest) to supply [Do not supply 0]
271c88c5224SJed Brown 
272c88c5224SJed Brown   Output Parameter:
2732fe279fdSBarry Smith . mat - the interpolation matrix, can be `NULL`
274c88c5224SJed Brown 
275c88c5224SJed Brown   Level: advanced
276c88c5224SJed Brown 
277f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`
278c88c5224SJed Brown @*/
279d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetInterpolation(PC pc, PetscInt l, Mat *mat)
280d71ae5a4SJacob Faibussowitsch {
281c88c5224SJed Brown   PC_MG         *mg       = (PC_MG *)pc->data;
282c88c5224SJed Brown   PC_MG_Levels **mglevels = mg->levels;
283c88c5224SJed Brown 
284c88c5224SJed Brown   PetscFunctionBegin;
285c88c5224SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2863ad4599aSBarry Smith   if (mat) PetscValidPointer(mat, 3);
28728b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
2882472a847SBarry 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);
28948a46eb9SPierre Jolivet   if (!mglevels[l]->interpolate && mglevels[l]->restrct) PetscCall(PCMGSetInterpolation(pc, l, mglevels[l]->restrct));
2903ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->interpolate;
2913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
292c88c5224SJed Brown }
293c88c5224SJed Brown 
2948a2c336bSFande Kong /*@
295eab52d2bSLawrence Mitchell   PCMGSetRestriction - Sets the function to be used to restrict dual vectors
2964b9ad928SBarry Smith   from level l to l-1.
2974b9ad928SBarry Smith 
298c3339decSBarry Smith   Logically Collective
2994b9ad928SBarry Smith 
3004b9ad928SBarry Smith   Input Parameters:
3014b9ad928SBarry Smith + pc  - the multigrid context
302c88c5224SJed Brown . l   - the level (0 is coarsest) to supply [Do not supply 0]
303c88c5224SJed Brown - mat - the restriction matrix
3044b9ad928SBarry Smith 
3054b9ad928SBarry Smith   Level: advanced
3064b9ad928SBarry Smith 
3074b9ad928SBarry Smith   Notes:
3084b9ad928SBarry Smith   Usually this is the same matrix used also to set the interpolation
3094b9ad928SBarry Smith   for the same level.
3104b9ad928SBarry Smith 
3114b9ad928SBarry Smith   One can pass in the interpolation matrix or its transpose; PETSc figures
3124b9ad928SBarry Smith   out from the matrix size which one it is.
3134b9ad928SBarry Smith 
314f1580f4eSBarry Smith   If you do not set this, the transpose of the `Mat` set with `PCMGSetInterpolation()`
315fccaa45eSBarry Smith   is used.
316fccaa45eSBarry Smith 
317f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInterpolation()`
3184b9ad928SBarry Smith @*/
319d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetRestriction(PC pc, PetscInt l, Mat mat)
320d71ae5a4SJacob Faibussowitsch {
321f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
322f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith   PetscFunctionBegin;
325c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
326c88c5224SJed Brown   PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
32728b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
32828b400f6SJacob Faibussowitsch   PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set restriction routine for coarsest level");
3299566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
3309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->restrct));
3312fa5cd67SKarl Rupp 
332f3fbd535SBarry Smith   mglevels[l]->restrct = mat;
3333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3344b9ad928SBarry Smith }
3354b9ad928SBarry Smith 
336c88c5224SJed Brown /*@
337eab52d2bSLawrence Mitchell   PCMGGetRestriction - Gets the function to be used to restrict dual vectors
338c88c5224SJed Brown   from level l to l-1.
339c88c5224SJed Brown 
340c3339decSBarry Smith   Logically Collective
341c88c5224SJed Brown 
342c88c5224SJed Brown   Input Parameters:
343c88c5224SJed Brown + pc - the multigrid context
344c88c5224SJed Brown - l  - the level (0 is coarsest) to supply [Do not supply 0]
345c88c5224SJed Brown 
346c88c5224SJed Brown   Output Parameter:
347c88c5224SJed Brown . mat - the restriction matrix
348c88c5224SJed Brown 
349c88c5224SJed Brown   Level: advanced
350c88c5224SJed Brown 
351f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGGetInjection()`
352c88c5224SJed Brown @*/
353d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetRestriction(PC pc, PetscInt l, Mat *mat)
354d71ae5a4SJacob Faibussowitsch {
355c88c5224SJed Brown   PC_MG         *mg       = (PC_MG *)pc->data;
356c88c5224SJed Brown   PC_MG_Levels **mglevels = mg->levels;
357c88c5224SJed Brown 
358c88c5224SJed Brown   PetscFunctionBegin;
359c88c5224SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
3603ad4599aSBarry Smith   if (mat) PetscValidPointer(mat, 3);
36128b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
3622472a847SBarry 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);
36348a46eb9SPierre Jolivet   if (!mglevels[l]->restrct && mglevels[l]->interpolate) PetscCall(PCMGSetRestriction(pc, l, mglevels[l]->interpolate));
3643ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->restrct;
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
366c88c5224SJed Brown }
367c88c5224SJed Brown 
36873250ac0SBarry Smith /*@
36973250ac0SBarry Smith   PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.
37073250ac0SBarry Smith 
371c3339decSBarry Smith   Logically Collective
372c88c5224SJed Brown 
373c88c5224SJed Brown   Input Parameters:
374c88c5224SJed Brown + pc     - the multigrid context
375*feefa0e1SJacob Faibussowitsch . l      - the level (0 is coarsest) to supply [Do not supply 0]
376*feefa0e1SJacob Faibussowitsch - rscale - the scaling
377c88c5224SJed Brown 
378c88c5224SJed Brown   Level: advanced
379c88c5224SJed Brown 
380f1580f4eSBarry Smith   Note:
381f1580f4eSBarry 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.
382f1580f4eSBarry Smith   It is preferable to use `PCMGSetInjection()` to control moving primal vectors.
383c88c5224SJed Brown 
384f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGSetInjection()`
385c88c5224SJed Brown @*/
386d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetRScale(PC pc, PetscInt l, Vec rscale)
387d71ae5a4SJacob Faibussowitsch {
388c88c5224SJed Brown   PC_MG         *mg       = (PC_MG *)pc->data;
389c88c5224SJed Brown   PC_MG_Levels **mglevels = mg->levels;
390c88c5224SJed Brown 
391c88c5224SJed Brown   PetscFunctionBegin;
392c88c5224SJed Brown   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
39328b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
3942472a847SBarry 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);
3959566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)rscale));
3969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mglevels[l]->rscale));
3972fa5cd67SKarl Rupp 
398c88c5224SJed Brown   mglevels[l]->rscale = rscale;
3993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
400c88c5224SJed Brown }
401c88c5224SJed Brown 
402c88c5224SJed Brown /*@
403c88c5224SJed Brown   PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1.
404c88c5224SJed Brown 
405c3339decSBarry Smith   Collective
40673250ac0SBarry Smith 
40773250ac0SBarry Smith   Input Parameters:
40873250ac0SBarry Smith + pc     - the multigrid context
40973250ac0SBarry Smith . rscale - the scaling
41073250ac0SBarry Smith - l      - the level (0 is coarsest) to supply [Do not supply 0]
41173250ac0SBarry Smith 
41273250ac0SBarry Smith   Level: advanced
41373250ac0SBarry Smith 
414f1580f4eSBarry Smith   Note:
415f1580f4eSBarry 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.
416f1580f4eSBarry Smith   It is preferable to use `PCMGGetInjection()` to control moving primal vectors.
41773250ac0SBarry Smith 
418f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInterpolation()`, `PCMGGetRestriction()`, `PCMGGetInjection()`
41973250ac0SBarry Smith @*/
420d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetRScale(PC pc, PetscInt l, Vec *rscale)
421d71ae5a4SJacob Faibussowitsch {
42273250ac0SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
42373250ac0SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
42473250ac0SBarry Smith 
42573250ac0SBarry Smith   PetscFunctionBegin;
42673250ac0SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
42728b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
4282472a847SBarry 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);
429c88c5224SJed Brown   if (!mglevels[l]->rscale) {
430c88c5224SJed Brown     Mat      R;
431c88c5224SJed Brown     Vec      X, Y, coarse, fine;
432c88c5224SJed Brown     PetscInt M, N;
4330fdf79fbSJacob Faibussowitsch 
4349566063dSJacob Faibussowitsch     PetscCall(PCMGGetRestriction(pc, l, &R));
4359566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(R, &X, &Y));
4369566063dSJacob Faibussowitsch     PetscCall(MatGetSize(R, &M, &N));
4370fdf79fbSJacob Faibussowitsch     PetscCheck(N != M, PetscObjectComm((PetscObject)R), PETSC_ERR_SUP, "Restriction matrix is square, cannot determine which Vec is coarser");
4382fa5cd67SKarl Rupp     if (M < N) {
4392fa5cd67SKarl Rupp       fine   = X;
4402fa5cd67SKarl Rupp       coarse = Y;
4410fdf79fbSJacob Faibussowitsch     } else {
4429371c9d4SSatish Balay       fine   = Y;
4439371c9d4SSatish Balay       coarse = X;
4440fdf79fbSJacob Faibussowitsch     }
4459566063dSJacob Faibussowitsch     PetscCall(VecSet(fine, 1.));
4469566063dSJacob Faibussowitsch     PetscCall(MatRestrict(R, fine, coarse));
4479566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&fine));
4489566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(coarse));
449c88c5224SJed Brown     mglevels[l]->rscale = coarse;
450c88c5224SJed Brown   }
451c88c5224SJed Brown   *rscale = mglevels[l]->rscale;
4523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45373250ac0SBarry Smith }
45473250ac0SBarry Smith 
455f39d8e23SSatish Balay /*@
456eab52d2bSLawrence Mitchell   PCMGSetInjection - Sets the function to be used to inject primal vectors
457eab52d2bSLawrence Mitchell   from level l to l-1.
458eab52d2bSLawrence Mitchell 
459c3339decSBarry Smith   Logically Collective
460eab52d2bSLawrence Mitchell 
461eab52d2bSLawrence Mitchell   Input Parameters:
462eab52d2bSLawrence Mitchell + pc  - the multigrid context
463eab52d2bSLawrence Mitchell . l   - the level (0 is coarsest) to supply [Do not supply 0]
464eab52d2bSLawrence Mitchell - mat - the injection matrix
465eab52d2bSLawrence Mitchell 
466eab52d2bSLawrence Mitchell   Level: advanced
467eab52d2bSLawrence Mitchell 
468f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRestriction()`
469eab52d2bSLawrence Mitchell @*/
470d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetInjection(PC pc, PetscInt l, Mat mat)
471d71ae5a4SJacob Faibussowitsch {
472eab52d2bSLawrence Mitchell   PC_MG         *mg       = (PC_MG *)pc->data;
473eab52d2bSLawrence Mitchell   PC_MG_Levels **mglevels = mg->levels;
474eab52d2bSLawrence Mitchell 
475eab52d2bSLawrence Mitchell   PetscFunctionBegin;
476eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
477eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
47828b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
47928b400f6SJacob Faibussowitsch   PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Do not set restriction routine for coarsest level");
4809566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
4819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->inject));
482eab52d2bSLawrence Mitchell 
483eab52d2bSLawrence Mitchell   mglevels[l]->inject = mat;
4843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
485eab52d2bSLawrence Mitchell }
486eab52d2bSLawrence Mitchell 
487eab52d2bSLawrence Mitchell /*@
488eab52d2bSLawrence Mitchell   PCMGGetInjection - Gets the function to be used to inject primal vectors
489eab52d2bSLawrence Mitchell   from level l to l-1.
490eab52d2bSLawrence Mitchell 
491c3339decSBarry Smith   Logically Collective
492eab52d2bSLawrence Mitchell 
493eab52d2bSLawrence Mitchell   Input Parameters:
494eab52d2bSLawrence Mitchell + pc - the multigrid context
495eab52d2bSLawrence Mitchell - l  - the level (0 is coarsest) to supply [Do not supply 0]
496eab52d2bSLawrence Mitchell 
497eab52d2bSLawrence Mitchell   Output Parameter:
49899a38656SLawrence Mitchell . mat - the restriction matrix (may be NULL if no injection is available).
499eab52d2bSLawrence Mitchell 
500eab52d2bSLawrence Mitchell   Level: advanced
501eab52d2bSLawrence Mitchell 
502f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetInjection()`, `PCMGetGetRestriction()`
503eab52d2bSLawrence Mitchell @*/
504d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetInjection(PC pc, PetscInt l, Mat *mat)
505d71ae5a4SJacob Faibussowitsch {
506eab52d2bSLawrence Mitchell   PC_MG         *mg       = (PC_MG *)pc->data;
507eab52d2bSLawrence Mitchell   PC_MG_Levels **mglevels = mg->levels;
508eab52d2bSLawrence Mitchell 
509eab52d2bSLawrence Mitchell   PetscFunctionBegin;
510eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
511eab52d2bSLawrence Mitchell   if (mat) PetscValidPointer(mat, 3);
51228b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
5132472a847SBarry 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);
514eab52d2bSLawrence Mitchell   if (mat) *mat = mglevels[l]->inject;
5153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
516eab52d2bSLawrence Mitchell }
517eab52d2bSLawrence Mitchell 
518eab52d2bSLawrence Mitchell /*@
519f1580f4eSBarry Smith   PCMGGetSmoother - Gets the `KSP` context to be used as smoother for
520f1580f4eSBarry Smith   both pre- and post-smoothing.  Call both `PCMGGetSmootherUp()` and
521f1580f4eSBarry Smith   `PCMGGetSmootherDown()` to use different functions for pre- and
5224b9ad928SBarry Smith   post-smoothing.
5234b9ad928SBarry Smith 
524f1580f4eSBarry Smith   Not Collective, ksp returned is parallel if pc is
5254b9ad928SBarry Smith 
5264b9ad928SBarry Smith   Input Parameters:
5274b9ad928SBarry Smith + pc - the multigrid context
5284b9ad928SBarry Smith - l  - the level (0 is coarsest) to supply
5294b9ad928SBarry Smith 
53001d2d390SJose E. Roman   Output Parameter:
5314b9ad928SBarry Smith . ksp - the smoother
5324b9ad928SBarry Smith 
533f1580f4eSBarry Smith   Note:
534f1580f4eSBarry 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.
535f1580f4eSBarry Smith   You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call `KSPGetPC`(ksp,&pc)
536f1580f4eSBarry Smith   and modify PC options for the smoother; for example `PCSetType`(pc,`PCSOR`); to use SOR smoothing.
53757420d5bSBarry Smith 
5384b9ad928SBarry Smith   Level: advanced
5394b9ad928SBarry Smith 
540*feefa0e1SJacob Faibussowitsch .seealso: `PCMG`, ``PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, `PCMGGetCoarseSolve()`
5414b9ad928SBarry Smith @*/
542d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetSmoother(PC pc, PetscInt l, KSP *ksp)
543d71ae5a4SJacob Faibussowitsch {
544f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
545f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
5464b9ad928SBarry Smith 
5474b9ad928SBarry Smith   PetscFunctionBegin;
548c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
549f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
5503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5514b9ad928SBarry Smith }
5524b9ad928SBarry Smith 
553f39d8e23SSatish Balay /*@
55497177400SBarry Smith   PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
5554b9ad928SBarry Smith   coarse grid correction (post-smoother).
5564b9ad928SBarry Smith 
557f1580f4eSBarry Smith   Not Collective, ksp returned is parallel if pc is
5584b9ad928SBarry Smith 
5594b9ad928SBarry Smith   Input Parameters:
5604b9ad928SBarry Smith + pc - the multigrid context
5614b9ad928SBarry Smith - l  - the level (0 is coarsest) to supply
5624b9ad928SBarry Smith 
56301d2d390SJose E. Roman   Output Parameter:
5644b9ad928SBarry Smith . ksp - the smoother
5654b9ad928SBarry Smith 
5664b9ad928SBarry Smith   Level: advanced
5674b9ad928SBarry Smith 
568f1580f4eSBarry Smith   Note:
569f1580f4eSBarry Smith   Calling this will result in a different pre and post smoother so you may need to set options on the pre smoother also
57089cce641SBarry Smith 
571f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`
5724b9ad928SBarry Smith @*/
573d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetSmootherUp(PC pc, PetscInt l, KSP *ksp)
574d71ae5a4SJacob Faibussowitsch {
575f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
576f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
577f69a0ea3SMatthew Knepley   const char    *prefix;
5784b9ad928SBarry Smith   MPI_Comm       comm;
5794b9ad928SBarry Smith 
5804b9ad928SBarry Smith   PetscFunctionBegin;
581c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5824b9ad928SBarry Smith   /*
5834b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
5844b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
5854b9ad928SBarry Smith      if not we allocate it.
5864b9ad928SBarry Smith   */
58728b400f6SJacob Faibussowitsch   PetscCheck(l, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "There is no such thing as a up smoother on the coarse grid");
588f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
58919fd82e9SBarry Smith     KSPType     ksptype;
59019fd82e9SBarry Smith     PCType      pctype;
591336babb1SJed Brown     PC          ipc;
592336babb1SJed Brown     PetscReal   rtol, abstol, dtol;
593336babb1SJed Brown     PetscInt    maxits;
594336babb1SJed Brown     KSPNormType normtype;
5959566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)mglevels[l]->smoothd, &comm));
5969566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(mglevels[l]->smoothd, &prefix));
5979566063dSJacob Faibussowitsch     PetscCall(KSPGetTolerances(mglevels[l]->smoothd, &rtol, &abstol, &dtol, &maxits));
5989566063dSJacob Faibussowitsch     PetscCall(KSPGetType(mglevels[l]->smoothd, &ksptype));
5999566063dSJacob Faibussowitsch     PetscCall(KSPGetNormType(mglevels[l]->smoothd, &normtype));
6009566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(mglevels[l]->smoothd, &ipc));
6019566063dSJacob Faibussowitsch     PetscCall(PCGetType(ipc, &pctype));
602336babb1SJed Brown 
6039566063dSJacob Faibussowitsch     PetscCall(KSPCreate(comm, &mglevels[l]->smoothu));
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 
638f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetSmootherUp()`, `PCMGGetSmoother()`
6394b9ad928SBarry Smith @*/
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 
665*feefa0e1SJacob Faibussowitsch .seealso: `PCMG`, `PCMGCycleType`, `PCMGSetCycleType()`
6664b9ad928SBarry Smith @*/
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 
696f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetX()`, `PCMGSetR()`
6974b9ad928SBarry Smith @*/
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 
729f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRhs()`, `PCMGSetR()`
7304b9ad928SBarry Smith @*/
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 
762f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetRhs()`, `PCMGSetX()`
7634b9ad928SBarry Smith @*/
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