xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
14b9ad928SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>       /*I "petscksp.h" I*/
34b9ad928SBarry Smith 
4f9426fe0SMark Adams /* ---------------------------------------------------------------------------*/
51f6cc5b2SSatish Balay /*@C
654b2cd4bSJed Brown    PCMGResidualDefault - Default routine to calculate the residual.
74b9ad928SBarry Smith 
8d083f849SBarry Smith    Collective on Mat
94b9ad928SBarry Smith 
104b9ad928SBarry Smith    Input Parameters:
114b9ad928SBarry Smith +  mat - the matrix
124b9ad928SBarry Smith .  b   - the right-hand-side
134b9ad928SBarry Smith -  x   - the approximate solution
144b9ad928SBarry Smith 
154b9ad928SBarry Smith    Output Parameter:
164b9ad928SBarry Smith .  r - location to store the residual
174b9ad928SBarry Smith 
18d0e4de75SBarry Smith    Level: developer
194b9ad928SBarry Smith 
20db781477SPatrick Sanan .seealso: `PCMGSetResidual()`
214b9ad928SBarry Smith @*/
2254b2cd4bSJed Brown PetscErrorCode  PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r)
234b9ad928SBarry Smith {
244b9ad928SBarry Smith   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(MatResidual(mat,b,x,r));
264b9ad928SBarry Smith   PetscFunctionReturn(0);
274b9ad928SBarry Smith }
284b9ad928SBarry Smith 
29fcb023d4SJed Brown /*@C
30fcb023d4SJed Brown    PCMGResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system
31fcb023d4SJed Brown 
32fcb023d4SJed Brown    Collective on Mat
33fcb023d4SJed Brown 
34fcb023d4SJed Brown    Input Parameters:
35fcb023d4SJed Brown +  mat - the matrix
36fcb023d4SJed Brown .  b   - the right-hand-side
37fcb023d4SJed Brown -  x   - the approximate solution
38fcb023d4SJed Brown 
39fcb023d4SJed Brown    Output Parameter:
40fcb023d4SJed Brown .  r - location to store the residual
41fcb023d4SJed Brown 
42fcb023d4SJed Brown    Level: developer
43fcb023d4SJed Brown 
44db781477SPatrick Sanan .seealso: `PCMGSetResidualTranspose()`
45fcb023d4SJed Brown @*/
46fcb023d4SJed Brown PetscErrorCode PCMGResidualTransposeDefault(Mat mat,Vec b,Vec x,Vec r)
47fcb023d4SJed Brown {
48fcb023d4SJed Brown   PetscFunctionBegin;
499566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(mat,x,r));
509566063dSJacob Faibussowitsch   PetscCall(VecAYPX(r,-1.0,b));
51fcb023d4SJed Brown   PetscFunctionReturn(0);
52fcb023d4SJed Brown }
53fcb023d4SJed Brown 
5430b0564aSStefano Zampini /*@C
5530b0564aSStefano Zampini    PCMGMatResidualDefault - Default routine to calculate the residual.
5630b0564aSStefano Zampini 
5730b0564aSStefano Zampini    Collective on Mat
5830b0564aSStefano Zampini 
5930b0564aSStefano Zampini    Input Parameters:
6030b0564aSStefano Zampini +  mat - the matrix
6130b0564aSStefano Zampini .  b   - the right-hand-side
6230b0564aSStefano Zampini -  x   - the approximate solution
6330b0564aSStefano Zampini 
6430b0564aSStefano Zampini    Output Parameter:
6530b0564aSStefano Zampini .  r - location to store the residual
6630b0564aSStefano Zampini 
6730b0564aSStefano Zampini    Level: developer
6830b0564aSStefano Zampini 
69db781477SPatrick Sanan .seealso: `PCMGSetMatResidual()`
7030b0564aSStefano Zampini @*/
7130b0564aSStefano Zampini PetscErrorCode  PCMGMatResidualDefault(Mat mat,Mat b,Mat x,Mat r)
7230b0564aSStefano Zampini {
7330b0564aSStefano Zampini   PetscFunctionBegin;
749566063dSJacob Faibussowitsch   PetscCall(MatMatMult(mat,x,MAT_REUSE_MATRIX,PETSC_DEFAULT,&r));
759566063dSJacob Faibussowitsch   PetscCall(MatAYPX(r,-1.0,b,UNKNOWN_NONZERO_PATTERN));
7630b0564aSStefano Zampini   PetscFunctionReturn(0);
7730b0564aSStefano Zampini }
7830b0564aSStefano Zampini 
7930b0564aSStefano Zampini /*@C
8030b0564aSStefano Zampini    PCMGMatResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system
8130b0564aSStefano Zampini 
8230b0564aSStefano Zampini    Collective on Mat
8330b0564aSStefano Zampini 
8430b0564aSStefano Zampini    Input Parameters:
8530b0564aSStefano Zampini +  mat - the matrix
8630b0564aSStefano Zampini .  b   - the right-hand-side
8730b0564aSStefano Zampini -  x   - the approximate solution
8830b0564aSStefano Zampini 
8930b0564aSStefano Zampini    Output Parameter:
9030b0564aSStefano Zampini .  r - location to store the residual
9130b0564aSStefano Zampini 
9230b0564aSStefano Zampini    Level: developer
9330b0564aSStefano Zampini 
94db781477SPatrick Sanan .seealso: `PCMGSetMatResidualTranspose()`
9530b0564aSStefano Zampini @*/
9630b0564aSStefano Zampini PetscErrorCode PCMGMatResidualTransposeDefault(Mat mat,Mat b,Mat x,Mat r)
9730b0564aSStefano Zampini {
9830b0564aSStefano Zampini   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(MatTransposeMatMult(mat,x,MAT_REUSE_MATRIX,PETSC_DEFAULT,&r));
1009566063dSJacob Faibussowitsch   PetscCall(MatAYPX(r,-1.0,b,UNKNOWN_NONZERO_PATTERN));
10130b0564aSStefano Zampini   PetscFunctionReturn(0);
10230b0564aSStefano Zampini }
103f39d8e23SSatish Balay /*@
10497177400SBarry Smith    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
1054b9ad928SBarry Smith 
1064b9ad928SBarry Smith    Not Collective
1074b9ad928SBarry Smith 
1084b9ad928SBarry Smith    Input Parameter:
1094b9ad928SBarry Smith .  pc - the multigrid context
1104b9ad928SBarry Smith 
1114b9ad928SBarry Smith    Output Parameter:
1124b9ad928SBarry Smith .  ksp - the coarse grid solver context
1134b9ad928SBarry Smith 
1144b9ad928SBarry Smith    Level: advanced
1154b9ad928SBarry Smith 
116db781477SPatrick Sanan .seealso: `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, `PCMGGetSmoother()`
1174b9ad928SBarry Smith @*/
1187087cfbeSBarry Smith PetscErrorCode  PCMGGetCoarseSolve(PC pc,KSP *ksp)
1194b9ad928SBarry Smith {
120f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
121f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
1224b9ad928SBarry Smith 
1234b9ad928SBarry Smith   PetscFunctionBegin;
124c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
125f3fbd535SBarry Smith   *ksp =  mglevels[0]->smoothd;
1264b9ad928SBarry Smith   PetscFunctionReturn(0);
1274b9ad928SBarry Smith }
1284b9ad928SBarry Smith 
1294b9ad928SBarry Smith /*@C
13097177400SBarry Smith    PCMGSetResidual - Sets the function to be used to calculate the residual
1314b9ad928SBarry Smith    on the lth level.
1324b9ad928SBarry Smith 
133d083f849SBarry Smith    Logically Collective on PC
1344b9ad928SBarry Smith 
1354b9ad928SBarry Smith    Input Parameters:
1364b9ad928SBarry Smith +  pc       - the multigrid context
1374b9ad928SBarry Smith .  l        - the level (0 is coarsest) to supply
138157726a2SBarry Smith .  residual - function used to form residual, if none is provided the previously provide one is used, if no
139d0e4de75SBarry Smith               previous one were provided then a default is used
1404b9ad928SBarry Smith -  mat      - matrix associated with residual
1414b9ad928SBarry Smith 
1424b9ad928SBarry Smith    Level: advanced
1434b9ad928SBarry Smith 
144db781477SPatrick Sanan .seealso: `PCMGResidualDefault()`
1454b9ad928SBarry Smith @*/
1467087cfbeSBarry Smith PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
1474b9ad928SBarry Smith {
148f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
149f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1504b9ad928SBarry Smith 
1514b9ad928SBarry Smith   PetscFunctionBegin;
152c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15328b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1542fa5cd67SKarl Rupp   if (residual) mglevels[l]->residual = residual;
15554b2cd4bSJed Brown   if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault;
15630b0564aSStefano Zampini   mglevels[l]->matresidual = PCMGMatResidualDefault;
1579566063dSJacob Faibussowitsch   if (mat) PetscCall(PetscObjectReference((PetscObject)mat));
1589566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->A));
159f3fbd535SBarry Smith   mglevels[l]->A = mat;
1604b9ad928SBarry Smith   PetscFunctionReturn(0);
1614b9ad928SBarry Smith }
1624b9ad928SBarry Smith 
163fcb023d4SJed Brown /*@C
164fcb023d4SJed Brown    PCMGSetResidualTranspose - Sets the function to be used to calculate the residual of the transposed linear system
165fcb023d4SJed Brown    on the lth level.
166fcb023d4SJed Brown 
167fcb023d4SJed Brown    Logically Collective on PC
168fcb023d4SJed Brown 
169fcb023d4SJed Brown    Input Parameters:
170fcb023d4SJed Brown +  pc        - the multigrid context
171fcb023d4SJed Brown .  l         - the level (0 is coarsest) to supply
172fcb023d4SJed Brown .  residualt - function used to form transpose of residual, if none is provided the previously provide one is used, if no
173fcb023d4SJed Brown                previous one were provided then a default is used
174fcb023d4SJed Brown -  mat       - matrix associated with residual
175fcb023d4SJed Brown 
176fcb023d4SJed Brown    Level: advanced
177fcb023d4SJed Brown 
178db781477SPatrick Sanan .seealso: `PCMGResidualTransposeDefault()`
179fcb023d4SJed Brown @*/
180fcb023d4SJed Brown PetscErrorCode  PCMGSetResidualTranspose(PC pc,PetscInt l,PetscErrorCode (*residualt)(Mat,Vec,Vec,Vec),Mat mat)
181fcb023d4SJed Brown {
182fcb023d4SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
183fcb023d4SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
184fcb023d4SJed Brown 
185fcb023d4SJed Brown   PetscFunctionBegin;
186fcb023d4SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
18728b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
188fcb023d4SJed Brown   if (residualt) mglevels[l]->residualtranspose = residualt;
189fcb023d4SJed Brown   if (!mglevels[l]->residualtranspose) mglevels[l]->residualtranspose = PCMGResidualTransposeDefault;
19030b0564aSStefano Zampini   mglevels[l]->matresidualtranspose = PCMGMatResidualTransposeDefault;
1919566063dSJacob Faibussowitsch   if (mat) PetscCall(PetscObjectReference((PetscObject)mat));
1929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->A));
193fcb023d4SJed Brown   mglevels[l]->A = mat;
194fcb023d4SJed Brown   PetscFunctionReturn(0);
195fcb023d4SJed Brown }
196fcb023d4SJed Brown 
1974b9ad928SBarry Smith /*@
198aea2a34eSBarry Smith    PCMGSetInterpolation - Sets the function to be used to calculate the
199bf5b2e24SBarry Smith    interpolation from l-1 to the lth level
2004b9ad928SBarry Smith 
201d083f849SBarry Smith    Logically Collective on PC
2024b9ad928SBarry Smith 
2034b9ad928SBarry Smith    Input Parameters:
2044b9ad928SBarry Smith +  pc  - the multigrid context
2054b9ad928SBarry Smith .  mat - the interpolation operator
206bf5b2e24SBarry Smith -  l   - the level (0 is coarsest) to supply [do not supply 0]
2074b9ad928SBarry Smith 
2084b9ad928SBarry Smith    Level: advanced
2094b9ad928SBarry Smith 
2104b9ad928SBarry Smith    Notes:
2114b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
2124b9ad928SBarry Smith     for the same level.
2134b9ad928SBarry Smith 
2144b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
2154b9ad928SBarry Smith     out from the matrix size which one it is.
2164b9ad928SBarry Smith 
217db781477SPatrick Sanan .seealso: `PCMGSetRestriction()`
2184b9ad928SBarry Smith @*/
2197087cfbeSBarry Smith PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
2204b9ad928SBarry Smith {
221f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
222f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2234b9ad928SBarry Smith 
2244b9ad928SBarry Smith   PetscFunctionBegin;
225c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
22628b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
22728b400f6SJacob Faibussowitsch   PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
2289566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
2299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->interpolate));
2302fa5cd67SKarl Rupp 
231f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
2324b9ad928SBarry Smith   PetscFunctionReturn(0);
2334b9ad928SBarry Smith }
2344b9ad928SBarry Smith 
2358a2c336bSFande Kong /*@
23695750439SFande Kong    PCMGSetOperators - Sets operator and preconditioning matrix for lth level
2378a2c336bSFande Kong 
238d083f849SBarry Smith    Logically Collective on PC
2398a2c336bSFande Kong 
2408a2c336bSFande Kong    Input Parameters:
2418a2c336bSFande Kong +  pc  - the multigrid context
2428a2c336bSFande Kong .  Amat - the operator
2438a2c336bSFande Kong .  pmat - the preconditioning operator
24495750439SFande Kong -  l   - the level (0 is the coarsest) to supply
2458a2c336bSFande Kong 
2468a2c336bSFande Kong    Level: advanced
2478a2c336bSFande Kong 
2488a2c336bSFande Kong .keywords:  multigrid, set, interpolate, level
2498a2c336bSFande Kong 
250db781477SPatrick Sanan .seealso: `PCMGSetRestriction()`, `PCMGSetInterpolation()`
2518a2c336bSFande Kong @*/
2528a2c336bSFande Kong PetscErrorCode  PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat)
253360ee056SFande Kong {
254360ee056SFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
255360ee056SFande Kong   PC_MG_Levels   **mglevels = mg->levels;
256360ee056SFande Kong 
257360ee056SFande Kong   PetscFunctionBegin;
258360ee056SFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2598a2c336bSFande Kong   PetscValidHeaderSpecific(Amat,MAT_CLASSID,3);
2608a2c336bSFande Kong   PetscValidHeaderSpecific(Pmat,MAT_CLASSID,4);
26128b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
2629566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(mglevels[l]->smoothd,Amat,Pmat));
263360ee056SFande Kong   PetscFunctionReturn(0);
264360ee056SFande Kong }
265360ee056SFande Kong 
266c88c5224SJed Brown /*@
267c88c5224SJed Brown    PCMGGetInterpolation - Gets the function to be used to calculate the
268c88c5224SJed Brown    interpolation from l-1 to the lth level
269c88c5224SJed Brown 
270c88c5224SJed Brown    Logically Collective on PC
271c88c5224SJed Brown 
272c88c5224SJed Brown    Input Parameters:
273c88c5224SJed Brown +  pc - the multigrid context
274c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
275c88c5224SJed Brown 
276c88c5224SJed Brown    Output Parameter:
2773ad4599aSBarry Smith .  mat - the interpolation matrix, can be NULL
278c88c5224SJed Brown 
279c88c5224SJed Brown    Level: advanced
280c88c5224SJed Brown 
281db781477SPatrick Sanan .seealso: `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`
282c88c5224SJed Brown @*/
283c88c5224SJed Brown PetscErrorCode  PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat)
284c88c5224SJed Brown {
285c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
286c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
287c88c5224SJed Brown 
288c88c5224SJed Brown   PetscFunctionBegin;
289c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2903ad4599aSBarry Smith   if (mat) PetscValidPointer(mat,3);
29128b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
2922472a847SBarry 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);
2932b3cbbdaSStefano Zampini   if (!mglevels[l]->interpolate && mglevels[l]->restrct) {
2949566063dSJacob Faibussowitsch     PetscCall(PCMGSetInterpolation(pc,l,mglevels[l]->restrct));
295c88c5224SJed Brown   }
2963ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->interpolate;
297c88c5224SJed Brown   PetscFunctionReturn(0);
298c88c5224SJed Brown }
299c88c5224SJed Brown 
3008a2c336bSFande Kong /*@
301eab52d2bSLawrence Mitchell    PCMGSetRestriction - Sets the function to be used to restrict dual vectors
3024b9ad928SBarry Smith    from level l to l-1.
3034b9ad928SBarry Smith 
304d083f849SBarry Smith    Logically Collective on PC
3054b9ad928SBarry Smith 
3064b9ad928SBarry Smith    Input Parameters:
3074b9ad928SBarry Smith +  pc - the multigrid context
308c88c5224SJed Brown .  l - the level (0 is coarsest) to supply [Do not supply 0]
309c88c5224SJed Brown -  mat - the restriction matrix
3104b9ad928SBarry Smith 
3114b9ad928SBarry Smith    Level: advanced
3124b9ad928SBarry Smith 
3134b9ad928SBarry Smith    Notes:
3144b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
3154b9ad928SBarry Smith     for the same level.
3164b9ad928SBarry Smith 
3174b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
3184b9ad928SBarry Smith     out from the matrix size which one it is.
3194b9ad928SBarry Smith 
320aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
321fccaa45eSBarry Smith     is used.
322fccaa45eSBarry Smith 
323db781477SPatrick Sanan .seealso: `PCMGSetInterpolation()`
3244b9ad928SBarry Smith @*/
3257087cfbeSBarry Smith PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
3264b9ad928SBarry Smith {
327f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
328f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith   PetscFunctionBegin;
331c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
332c88c5224SJed Brown   PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
33328b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
33428b400f6SJacob Faibussowitsch   PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
3359566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
3369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->restrct));
3372fa5cd67SKarl Rupp 
338f3fbd535SBarry Smith   mglevels[l]->restrct = mat;
3394b9ad928SBarry Smith   PetscFunctionReturn(0);
3404b9ad928SBarry Smith }
3414b9ad928SBarry Smith 
342c88c5224SJed Brown /*@
343eab52d2bSLawrence Mitchell    PCMGGetRestriction - Gets the function to be used to restrict dual vectors
344c88c5224SJed Brown    from level l to l-1.
345c88c5224SJed Brown 
346d083f849SBarry Smith    Logically Collective on PC
347c88c5224SJed Brown 
348c88c5224SJed Brown    Input Parameters:
349c88c5224SJed Brown +  pc - the multigrid context
350c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
351c88c5224SJed Brown 
352c88c5224SJed Brown    Output Parameter:
353c88c5224SJed Brown .  mat - the restriction matrix
354c88c5224SJed Brown 
355c88c5224SJed Brown    Level: advanced
356c88c5224SJed Brown 
357db781477SPatrick Sanan .seealso: `PCMGGetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGGetInjection()`
358c88c5224SJed Brown @*/
359c88c5224SJed Brown PetscErrorCode  PCMGGetRestriction(PC pc,PetscInt l,Mat *mat)
360c88c5224SJed Brown {
361c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
362c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
363c88c5224SJed Brown 
364c88c5224SJed Brown   PetscFunctionBegin;
365c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3663ad4599aSBarry Smith   if (mat) PetscValidPointer(mat,3);
36728b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
3682472a847SBarry 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);
3692b3cbbdaSStefano Zampini   if (!mglevels[l]->restrct && mglevels[l]->interpolate) {
3709566063dSJacob Faibussowitsch     PetscCall(PCMGSetRestriction(pc,l,mglevels[l]->interpolate));
371c88c5224SJed Brown   }
3723ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->restrct;
373c88c5224SJed Brown   PetscFunctionReturn(0);
374c88c5224SJed Brown }
375c88c5224SJed Brown 
37673250ac0SBarry Smith /*@
37773250ac0SBarry Smith    PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.
37873250ac0SBarry Smith 
379d083f849SBarry Smith    Logically Collective on PC
380c88c5224SJed Brown 
381c88c5224SJed Brown    Input Parameters:
382c88c5224SJed Brown +  pc - the multigrid context
383c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
384c88c5224SJed Brown .  rscale - the scaling
385c88c5224SJed Brown 
386c88c5224SJed Brown    Level: advanced
387c88c5224SJed Brown 
388c88c5224SJed Brown    Notes:
389eab52d2bSLawrence Mitchell        When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R.  It is preferable to use PCMGSetInjection() to control moving primal vectors.
390c88c5224SJed Brown 
391db781477SPatrick Sanan .seealso: `PCMGSetInterpolation()`, `PCMGSetRestriction()`, `PCMGGetRScale()`, `PCMGSetInjection()`
392c88c5224SJed Brown @*/
393c88c5224SJed Brown PetscErrorCode  PCMGSetRScale(PC pc,PetscInt l,Vec rscale)
394c88c5224SJed Brown {
395c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
396c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
397c88c5224SJed Brown 
398c88c5224SJed Brown   PetscFunctionBegin;
399c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
40028b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
4012472a847SBarry 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);
4029566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)rscale));
4039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mglevels[l]->rscale));
4042fa5cd67SKarl Rupp 
405c88c5224SJed Brown   mglevels[l]->rscale = rscale;
406c88c5224SJed Brown   PetscFunctionReturn(0);
407c88c5224SJed Brown }
408c88c5224SJed Brown 
409c88c5224SJed Brown /*@
410c88c5224SJed Brown    PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1.
411c88c5224SJed Brown 
412c88c5224SJed Brown    Collective on PC
41373250ac0SBarry Smith 
41473250ac0SBarry Smith    Input Parameters:
41573250ac0SBarry Smith +  pc - the multigrid context
41673250ac0SBarry Smith .  rscale - the scaling
41773250ac0SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
41873250ac0SBarry Smith 
41973250ac0SBarry Smith    Level: advanced
42073250ac0SBarry Smith 
42173250ac0SBarry Smith    Notes:
422eab52d2bSLawrence Mitchell        When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R.  It is preferable to use PCMGGetInjection() to control moving primal vectors.
42373250ac0SBarry Smith 
424db781477SPatrick Sanan .seealso: `PCMGSetInterpolation()`, `PCMGGetRestriction()`, `PCMGGetInjection()`
42573250ac0SBarry Smith @*/
426c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale)
42773250ac0SBarry Smith {
42873250ac0SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
42973250ac0SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
43073250ac0SBarry Smith 
43173250ac0SBarry Smith   PetscFunctionBegin;
43273250ac0SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
43328b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
4342472a847SBarry 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);
435c88c5224SJed Brown   if (!mglevels[l]->rscale) {
436c88c5224SJed Brown     Mat      R;
437c88c5224SJed Brown     Vec      X,Y,coarse,fine;
438c88c5224SJed Brown     PetscInt M,N;
4399566063dSJacob Faibussowitsch     PetscCall(PCMGGetRestriction(pc,l,&R));
4409566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(R,&X,&Y));
4419566063dSJacob Faibussowitsch     PetscCall(MatGetSize(R,&M,&N));
4422fa5cd67SKarl Rupp     if (M < N) {
4432fa5cd67SKarl Rupp       fine = X;
4442fa5cd67SKarl Rupp       coarse = Y;
4452fa5cd67SKarl Rupp     } else if (N < M) {
4462fa5cd67SKarl Rupp       fine = Y; coarse = X;
447ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser");
4489566063dSJacob Faibussowitsch     PetscCall(VecSet(fine,1.));
4499566063dSJacob Faibussowitsch     PetscCall(MatRestrict(R,fine,coarse));
4509566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&fine));
4519566063dSJacob Faibussowitsch     PetscCall(VecReciprocal(coarse));
452c88c5224SJed Brown     mglevels[l]->rscale = coarse;
453c88c5224SJed Brown   }
454c88c5224SJed Brown   *rscale = mglevels[l]->rscale;
45573250ac0SBarry Smith   PetscFunctionReturn(0);
45673250ac0SBarry Smith }
45773250ac0SBarry Smith 
458f39d8e23SSatish Balay /*@
459eab52d2bSLawrence Mitchell    PCMGSetInjection - Sets the function to be used to inject primal vectors
460eab52d2bSLawrence Mitchell    from level l to l-1.
461eab52d2bSLawrence Mitchell 
462d083f849SBarry Smith    Logically Collective on PC
463eab52d2bSLawrence Mitchell 
464eab52d2bSLawrence Mitchell    Input Parameters:
465eab52d2bSLawrence Mitchell +  pc - the multigrid context
466eab52d2bSLawrence Mitchell .  l - the level (0 is coarsest) to supply [Do not supply 0]
467eab52d2bSLawrence Mitchell -  mat - the injection matrix
468eab52d2bSLawrence Mitchell 
469eab52d2bSLawrence Mitchell    Level: advanced
470eab52d2bSLawrence Mitchell 
471db781477SPatrick Sanan .seealso: `PCMGSetRestriction()`
472eab52d2bSLawrence Mitchell @*/
473eab52d2bSLawrence Mitchell PetscErrorCode  PCMGSetInjection(PC pc,PetscInt l,Mat mat)
474eab52d2bSLawrence Mitchell {
475eab52d2bSLawrence Mitchell   PC_MG          *mg        = (PC_MG*)pc->data;
476eab52d2bSLawrence Mitchell   PC_MG_Levels   **mglevels = mg->levels;
477eab52d2bSLawrence Mitchell 
478eab52d2bSLawrence Mitchell   PetscFunctionBegin;
479eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
480eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
48128b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
48228b400f6SJacob Faibussowitsch   PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
4839566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
4849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mglevels[l]->inject));
485eab52d2bSLawrence Mitchell 
486eab52d2bSLawrence Mitchell   mglevels[l]->inject = mat;
487eab52d2bSLawrence Mitchell   PetscFunctionReturn(0);
488eab52d2bSLawrence Mitchell }
489eab52d2bSLawrence Mitchell 
490eab52d2bSLawrence Mitchell /*@
491eab52d2bSLawrence Mitchell    PCMGGetInjection - Gets the function to be used to inject primal vectors
492eab52d2bSLawrence Mitchell    from level l to l-1.
493eab52d2bSLawrence Mitchell 
494d083f849SBarry Smith    Logically Collective on PC
495eab52d2bSLawrence Mitchell 
496eab52d2bSLawrence Mitchell    Input Parameters:
497eab52d2bSLawrence Mitchell +  pc - the multigrid context
498eab52d2bSLawrence Mitchell -  l - the level (0 is coarsest) to supply [Do not supply 0]
499eab52d2bSLawrence Mitchell 
500eab52d2bSLawrence Mitchell    Output Parameter:
50199a38656SLawrence Mitchell .  mat - the restriction matrix (may be NULL if no injection is available).
502eab52d2bSLawrence Mitchell 
503eab52d2bSLawrence Mitchell    Level: advanced
504eab52d2bSLawrence Mitchell 
505db781477SPatrick Sanan .seealso: `PCMGSetInjection()`, `PCMGetGetRestriction()`
506eab52d2bSLawrence Mitchell @*/
507eab52d2bSLawrence Mitchell PetscErrorCode  PCMGGetInjection(PC pc,PetscInt l,Mat *mat)
508eab52d2bSLawrence Mitchell {
509eab52d2bSLawrence Mitchell   PC_MG          *mg        = (PC_MG*)pc->data;
510eab52d2bSLawrence Mitchell   PC_MG_Levels   **mglevels = mg->levels;
511eab52d2bSLawrence Mitchell 
512eab52d2bSLawrence Mitchell   PetscFunctionBegin;
513eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
514eab52d2bSLawrence Mitchell   if (mat) PetscValidPointer(mat,3);
51528b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
5162472a847SBarry 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);
517eab52d2bSLawrence Mitchell   if (mat) *mat = mglevels[l]->inject;
518eab52d2bSLawrence Mitchell   PetscFunctionReturn(0);
519eab52d2bSLawrence Mitchell }
520eab52d2bSLawrence Mitchell 
521eab52d2bSLawrence Mitchell /*@
52297177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
52397177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
52497177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
5254b9ad928SBarry Smith    post-smoothing.
5264b9ad928SBarry Smith 
5274b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
5284b9ad928SBarry Smith 
5294b9ad928SBarry Smith    Input Parameters:
5304b9ad928SBarry Smith +  pc - the multigrid context
5314b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
5324b9ad928SBarry Smith 
53301d2d390SJose E. Roman    Output Parameter:
5344b9ad928SBarry Smith .  ksp - the smoother
5354b9ad928SBarry Smith 
53657420d5bSBarry Smith    Notes:
53757420d5bSBarry Smith    Once you have called this routine, you can call KSPSetOperators(ksp,...) on the resulting ksp to provide the operators for the smoother for this level.
53857420d5bSBarry Smith    You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc)
53957420d5bSBarry Smith    and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing.
54057420d5bSBarry Smith 
5414b9ad928SBarry Smith    Level: advanced
5424b9ad928SBarry Smith 
543db781477SPatrick Sanan .seealso: `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, `PCMGGetCoarseSolve()`
5444b9ad928SBarry Smith @*/
5457087cfbeSBarry Smith PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
5464b9ad928SBarry Smith {
547f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
548f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
5494b9ad928SBarry Smith 
5504b9ad928SBarry Smith   PetscFunctionBegin;
551c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
552f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
5534b9ad928SBarry Smith   PetscFunctionReturn(0);
5544b9ad928SBarry Smith }
5554b9ad928SBarry Smith 
556f39d8e23SSatish Balay /*@
55797177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
5584b9ad928SBarry Smith    coarse grid correction (post-smoother).
5594b9ad928SBarry Smith 
5604b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
5614b9ad928SBarry Smith 
5624b9ad928SBarry Smith    Input Parameters:
5634b9ad928SBarry Smith +  pc - the multigrid context
5644b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
5654b9ad928SBarry Smith 
56601d2d390SJose E. Roman    Output Parameter:
5674b9ad928SBarry Smith .  ksp - the smoother
5684b9ad928SBarry Smith 
5694b9ad928SBarry Smith    Level: advanced
5704b9ad928SBarry Smith 
57195452b02SPatrick Sanan    Notes:
57295452b02SPatrick Sanan     calling this will result in a different pre and post smoother so you may need to
57389cce641SBarry Smith          set options on the pre smoother also
57489cce641SBarry Smith 
575db781477SPatrick Sanan .seealso: `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`
5764b9ad928SBarry Smith @*/
5777087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
5784b9ad928SBarry Smith {
579f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
580f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
581f69a0ea3SMatthew Knepley   const char     *prefix;
5824b9ad928SBarry Smith   MPI_Comm       comm;
5834b9ad928SBarry Smith 
5844b9ad928SBarry Smith   PetscFunctionBegin;
585c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5864b9ad928SBarry Smith   /*
5874b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
5884b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
5894b9ad928SBarry Smith      if not we allocate it.
5904b9ad928SBarry Smith   */
59128b400f6SJacob Faibussowitsch   PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
592f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
59319fd82e9SBarry Smith     KSPType     ksptype;
59419fd82e9SBarry Smith     PCType      pctype;
595336babb1SJed Brown     PC          ipc;
596336babb1SJed Brown     PetscReal   rtol,abstol,dtol;
597336babb1SJed Brown     PetscInt    maxits;
598336babb1SJed Brown     KSPNormType normtype;
5999566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm));
6009566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix));
6019566063dSJacob Faibussowitsch     PetscCall(KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits));
6029566063dSJacob Faibussowitsch     PetscCall(KSPGetType(mglevels[l]->smoothd,&ksptype));
6039566063dSJacob Faibussowitsch     PetscCall(KSPGetNormType(mglevels[l]->smoothd,&normtype));
6049566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(mglevels[l]->smoothd,&ipc));
6059566063dSJacob Faibussowitsch     PetscCall(PCGetType(ipc,&pctype));
606336babb1SJed Brown 
6079566063dSJacob Faibussowitsch     PetscCall(KSPCreate(comm,&mglevels[l]->smoothu));
6089566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure));
6099566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l));
6109566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix));
6119566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits));
6129566063dSJacob Faibussowitsch     PetscCall(KSPSetType(mglevels[l]->smoothu,ksptype));
6139566063dSJacob Faibussowitsch     PetscCall(KSPSetNormType(mglevels[l]->smoothu,normtype));
6149566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL));
6159566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(mglevels[l]->smoothu,&ipc));
6169566063dSJacob Faibussowitsch     PetscCall(PCSetType(ipc,pctype));
6179566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu));
6189566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level));
6194b9ad928SBarry Smith   }
620f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
6214b9ad928SBarry Smith   PetscFunctionReturn(0);
6224b9ad928SBarry Smith }
6234b9ad928SBarry Smith 
624f39d8e23SSatish Balay /*@
62597177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
6264b9ad928SBarry Smith    coarse grid correction (pre-smoother).
6274b9ad928SBarry Smith 
6284b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
6294b9ad928SBarry Smith 
6304b9ad928SBarry Smith    Input Parameters:
6314b9ad928SBarry Smith +  pc - the multigrid context
6324b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
6334b9ad928SBarry Smith 
63401d2d390SJose E. Roman    Output Parameter:
6354b9ad928SBarry Smith .  ksp - the smoother
6364b9ad928SBarry Smith 
6374b9ad928SBarry Smith    Level: advanced
6384b9ad928SBarry Smith 
63995452b02SPatrick Sanan    Notes:
64095452b02SPatrick Sanan     calling this will result in a different pre and post smoother so you may need to
64189cce641SBarry Smith          set options on the post smoother also
64289cce641SBarry Smith 
643db781477SPatrick Sanan .seealso: `PCMGGetSmootherUp()`, `PCMGGetSmoother()`
6444b9ad928SBarry Smith @*/
6457087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
6464b9ad928SBarry Smith {
647f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
648f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
6494b9ad928SBarry Smith 
6504b9ad928SBarry Smith   PetscFunctionBegin;
651c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
6524b9ad928SBarry Smith   /* make sure smoother up and down are different */
653*1baa6e33SBarry Smith   if (l) PetscCall(PCMGGetSmootherUp(pc,l,NULL));
654f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
6554b9ad928SBarry Smith   PetscFunctionReturn(0);
6564b9ad928SBarry Smith }
6574b9ad928SBarry Smith 
6584b9ad928SBarry Smith /*@
659cab31ae5SJed Brown    PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level.
6604b9ad928SBarry Smith 
661ad4df100SBarry Smith    Logically Collective on PC
6624b9ad928SBarry Smith 
6634b9ad928SBarry Smith    Input Parameters:
6644b9ad928SBarry Smith +  pc - the multigrid context
665c1cbb1deSBarry Smith .  l  - the level (0 is coarsest)
666c1cbb1deSBarry Smith -  c  - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
6674b9ad928SBarry Smith 
6684b9ad928SBarry Smith    Level: advanced
6694b9ad928SBarry Smith 
670db781477SPatrick Sanan .seealso: `PCMGSetCycleType()`
6714b9ad928SBarry Smith @*/
672c1cbb1deSBarry Smith PetscErrorCode  PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c)
6734b9ad928SBarry Smith {
674f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
675f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
6764b9ad928SBarry Smith 
6774b9ad928SBarry Smith   PetscFunctionBegin;
678c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
67928b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
680c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,l,2);
681c51679f6SSatish Balay   PetscValidLogicalCollectiveEnum(pc,c,3);
682f3fbd535SBarry Smith   mglevels[l]->cycles = c;
6834b9ad928SBarry Smith   PetscFunctionReturn(0);
6844b9ad928SBarry Smith }
6854b9ad928SBarry Smith 
6864b9ad928SBarry Smith /*@
687f3b08a26SMatthew G. Knepley   PCMGSetRhs - Sets the vector to be used to store the right-hand side on a particular level.
6884b9ad928SBarry Smith 
689d083f849SBarry Smith    Logically Collective on PC
6904b9ad928SBarry Smith 
6914b9ad928SBarry Smith   Input Parameters:
6924b9ad928SBarry Smith + pc - the multigrid context
6934b9ad928SBarry Smith . l  - the level (0 is coarsest) this is to be used for
694f3b08a26SMatthew G. Knepley - c  - the Vec
6954b9ad928SBarry Smith 
6964b9ad928SBarry Smith   Level: advanced
6974b9ad928SBarry Smith 
69895452b02SPatrick Sanan   Notes:
699f3b08a26SMatthew G. Knepley   If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it.
700fccaa45eSBarry Smith 
701f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, right-hand-side, rhs, level
702db781477SPatrick Sanan .seealso: `PCMGSetX()`, `PCMGSetR()`
7034b9ad928SBarry Smith @*/
7047087cfbeSBarry Smith PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
7054b9ad928SBarry Smith {
706f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
707f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
7084b9ad928SBarry Smith 
7094b9ad928SBarry Smith   PetscFunctionBegin;
710c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
71128b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
71208401ef6SPierre Jolivet   PetscCheck(l != mglevels[0]->levels-1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
7139566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)c));
7149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mglevels[l]->b));
7152fa5cd67SKarl Rupp 
716f3fbd535SBarry Smith   mglevels[l]->b = c;
7174b9ad928SBarry Smith   PetscFunctionReturn(0);
7184b9ad928SBarry Smith }
7194b9ad928SBarry Smith 
7204b9ad928SBarry Smith /*@
721f3b08a26SMatthew G. Knepley   PCMGSetX - Sets the vector to be used to store the solution on a particular level.
7224b9ad928SBarry Smith 
723d083f849SBarry Smith   Logically Collective on PC
7244b9ad928SBarry Smith 
7254b9ad928SBarry Smith   Input Parameters:
7264b9ad928SBarry Smith + pc - the multigrid context
727251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level)
728f3b08a26SMatthew G. Knepley - c - the Vec
7294b9ad928SBarry Smith 
7304b9ad928SBarry Smith   Level: advanced
7314b9ad928SBarry Smith 
73295452b02SPatrick Sanan   Notes:
733f3b08a26SMatthew G. Knepley   If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it.
734fccaa45eSBarry Smith 
735f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, solution, level
736db781477SPatrick Sanan .seealso: `PCMGSetRhs()`, `PCMGSetR()`
7374b9ad928SBarry Smith @*/
7387087cfbeSBarry Smith PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
7394b9ad928SBarry Smith {
740f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
741f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
7424b9ad928SBarry Smith 
7434b9ad928SBarry Smith   PetscFunctionBegin;
744c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
74528b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
74608401ef6SPierre Jolivet   PetscCheck(l != mglevels[0]->levels-1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
7479566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)c));
7489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mglevels[l]->x));
7492fa5cd67SKarl Rupp 
750f3fbd535SBarry Smith   mglevels[l]->x = c;
7514b9ad928SBarry Smith   PetscFunctionReturn(0);
7524b9ad928SBarry Smith }
7534b9ad928SBarry Smith 
7544b9ad928SBarry Smith /*@
755f3b08a26SMatthew G. Knepley   PCMGSetR - Sets the vector to be used to store the residual on a particular level.
7564b9ad928SBarry Smith 
757d083f849SBarry Smith   Logically Collective on PC
7584b9ad928SBarry Smith 
7594b9ad928SBarry Smith   Input Parameters:
7604b9ad928SBarry Smith + pc - the multigrid context
7614b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for
762f3b08a26SMatthew G. Knepley - c - the Vec
7634b9ad928SBarry Smith 
7644b9ad928SBarry Smith   Level: advanced
7654b9ad928SBarry Smith 
76695452b02SPatrick Sanan   Notes:
767f3b08a26SMatthew G. Knepley   If this is not provided PETSc will automatically generate one. You do not need to keep a reference to this vector if you do not need it. PCDestroy() will properly free it.
768fccaa45eSBarry Smith 
769f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, residual, level
770db781477SPatrick Sanan .seealso: `PCMGSetRhs()`, `PCMGSetX()`
7714b9ad928SBarry Smith @*/
7727087cfbeSBarry Smith PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
7734b9ad928SBarry Smith {
774f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
775f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
7764b9ad928SBarry Smith 
7774b9ad928SBarry Smith   PetscFunctionBegin;
778c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
77928b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
78028b400f6SJacob Faibussowitsch   PetscCheck(l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
7819566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)c));
7829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mglevels[l]->r));
7832fa5cd67SKarl Rupp 
784f3fbd535SBarry Smith   mglevels[l]->r = c;
7854b9ad928SBarry Smith   PetscFunctionReturn(0);
7864b9ad928SBarry Smith }
787