xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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 
2097177400SBarry Smith .seealso: PCMGSetResidual()
214b9ad928SBarry Smith @*/
2254b2cd4bSJed Brown PetscErrorCode  PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r)
234b9ad928SBarry Smith {
24dfbe8321SBarry Smith   PetscErrorCode ierr;
254b9ad928SBarry Smith 
264b9ad928SBarry Smith   PetscFunctionBegin;
27f9426fe0SMark Adams   ierr = MatResidual(mat,b,x,r);CHKERRQ(ierr);
284b9ad928SBarry Smith   PetscFunctionReturn(0);
294b9ad928SBarry Smith }
304b9ad928SBarry Smith 
31fcb023d4SJed Brown /*@C
32fcb023d4SJed Brown    PCMGResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system
33fcb023d4SJed Brown 
34fcb023d4SJed Brown    Collective on Mat
35fcb023d4SJed Brown 
36fcb023d4SJed Brown    Input Parameters:
37fcb023d4SJed Brown +  mat - the matrix
38fcb023d4SJed Brown .  b   - the right-hand-side
39fcb023d4SJed Brown -  x   - the approximate solution
40fcb023d4SJed Brown 
41fcb023d4SJed Brown    Output Parameter:
42fcb023d4SJed Brown .  r - location to store the residual
43fcb023d4SJed Brown 
44fcb023d4SJed Brown    Level: developer
45fcb023d4SJed Brown 
46fcb023d4SJed Brown .seealso: PCMGSetResidualTranspose()
47fcb023d4SJed Brown @*/
48fcb023d4SJed Brown PetscErrorCode PCMGResidualTransposeDefault(Mat mat,Vec b,Vec x,Vec r)
49fcb023d4SJed Brown {
50fcb023d4SJed Brown   PetscErrorCode ierr;
51fcb023d4SJed Brown 
52fcb023d4SJed Brown   PetscFunctionBegin;
53fcb023d4SJed Brown   ierr = MatMultTranspose(mat,x,r);CHKERRQ(ierr);
54fcb023d4SJed Brown   ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
55fcb023d4SJed Brown   PetscFunctionReturn(0);
56fcb023d4SJed Brown }
57fcb023d4SJed Brown 
5830b0564aSStefano Zampini /*@C
5930b0564aSStefano Zampini    PCMGMatResidualDefault - Default routine to calculate the residual.
6030b0564aSStefano Zampini 
6130b0564aSStefano Zampini    Collective on Mat
6230b0564aSStefano Zampini 
6330b0564aSStefano Zampini    Input Parameters:
6430b0564aSStefano Zampini +  mat - the matrix
6530b0564aSStefano Zampini .  b   - the right-hand-side
6630b0564aSStefano Zampini -  x   - the approximate solution
6730b0564aSStefano Zampini 
6830b0564aSStefano Zampini    Output Parameter:
6930b0564aSStefano Zampini .  r - location to store the residual
7030b0564aSStefano Zampini 
7130b0564aSStefano Zampini    Level: developer
7230b0564aSStefano Zampini 
7330b0564aSStefano Zampini .seealso: PCMGSetMatResidual()
7430b0564aSStefano Zampini @*/
7530b0564aSStefano Zampini PetscErrorCode  PCMGMatResidualDefault(Mat mat,Mat b,Mat x,Mat r)
7630b0564aSStefano Zampini {
7730b0564aSStefano Zampini   PetscErrorCode ierr;
7830b0564aSStefano Zampini 
7930b0564aSStefano Zampini   PetscFunctionBegin;
8030b0564aSStefano Zampini   ierr = MatMatMult(mat,x,MAT_REUSE_MATRIX,PETSC_DEFAULT,&r);CHKERRQ(ierr);
8130b0564aSStefano Zampini   ierr = MatAYPX(r,-1.0,b,UNKNOWN_NONZERO_PATTERN);CHKERRQ(ierr);
8230b0564aSStefano Zampini   PetscFunctionReturn(0);
8330b0564aSStefano Zampini }
8430b0564aSStefano Zampini 
8530b0564aSStefano Zampini /*@C
8630b0564aSStefano Zampini    PCMGMatResidualTransposeDefault - Default routine to calculate the residual of the transposed linear system
8730b0564aSStefano Zampini 
8830b0564aSStefano Zampini    Collective on Mat
8930b0564aSStefano Zampini 
9030b0564aSStefano Zampini    Input Parameters:
9130b0564aSStefano Zampini +  mat - the matrix
9230b0564aSStefano Zampini .  b   - the right-hand-side
9330b0564aSStefano Zampini -  x   - the approximate solution
9430b0564aSStefano Zampini 
9530b0564aSStefano Zampini    Output Parameter:
9630b0564aSStefano Zampini .  r - location to store the residual
9730b0564aSStefano Zampini 
9830b0564aSStefano Zampini    Level: developer
9930b0564aSStefano Zampini 
10030b0564aSStefano Zampini .seealso: PCMGSetMatResidualTranspose()
10130b0564aSStefano Zampini @*/
10230b0564aSStefano Zampini PetscErrorCode PCMGMatResidualTransposeDefault(Mat mat,Mat b,Mat x,Mat r)
10330b0564aSStefano Zampini {
10430b0564aSStefano Zampini   PetscErrorCode ierr;
10530b0564aSStefano Zampini 
10630b0564aSStefano Zampini   PetscFunctionBegin;
10730b0564aSStefano Zampini   ierr = MatTransposeMatMult(mat,x,MAT_REUSE_MATRIX,PETSC_DEFAULT,&r);CHKERRQ(ierr);
10830b0564aSStefano Zampini   ierr = MatAYPX(r,-1.0,b,UNKNOWN_NONZERO_PATTERN);CHKERRQ(ierr);
10930b0564aSStefano Zampini   PetscFunctionReturn(0);
11030b0564aSStefano Zampini }
111f39d8e23SSatish Balay /*@
11297177400SBarry Smith    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
1134b9ad928SBarry Smith 
1144b9ad928SBarry Smith    Not Collective
1154b9ad928SBarry Smith 
1164b9ad928SBarry Smith    Input Parameter:
1174b9ad928SBarry Smith .  pc - the multigrid context
1184b9ad928SBarry Smith 
1194b9ad928SBarry Smith    Output Parameter:
1204b9ad928SBarry Smith .  ksp - the coarse grid solver context
1214b9ad928SBarry Smith 
1224b9ad928SBarry Smith    Level: advanced
1234b9ad928SBarry Smith 
12428668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetSmoother()
1254b9ad928SBarry Smith @*/
1267087cfbeSBarry Smith PetscErrorCode  PCMGGetCoarseSolve(PC pc,KSP *ksp)
1274b9ad928SBarry Smith {
128f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
129f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
1304b9ad928SBarry Smith 
1314b9ad928SBarry Smith   PetscFunctionBegin;
132c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
133f3fbd535SBarry Smith   *ksp =  mglevels[0]->smoothd;
1344b9ad928SBarry Smith   PetscFunctionReturn(0);
1354b9ad928SBarry Smith }
1364b9ad928SBarry Smith 
1374b9ad928SBarry Smith /*@C
13897177400SBarry Smith    PCMGSetResidual - Sets the function to be used to calculate the residual
1394b9ad928SBarry Smith    on the lth level.
1404b9ad928SBarry Smith 
141d083f849SBarry Smith    Logically Collective on PC
1424b9ad928SBarry Smith 
1434b9ad928SBarry Smith    Input Parameters:
1444b9ad928SBarry Smith +  pc       - the multigrid context
1454b9ad928SBarry Smith .  l        - the level (0 is coarsest) to supply
146157726a2SBarry Smith .  residual - function used to form residual, if none is provided the previously provide one is used, if no
147d0e4de75SBarry Smith               previous one were provided then a default is used
1484b9ad928SBarry Smith -  mat      - matrix associated with residual
1494b9ad928SBarry Smith 
1504b9ad928SBarry Smith    Level: advanced
1514b9ad928SBarry Smith 
15254b2cd4bSJed Brown .seealso: PCMGResidualDefault()
1534b9ad928SBarry Smith @*/
1547087cfbeSBarry Smith PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
1554b9ad928SBarry Smith {
156f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
157f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
158298cc208SBarry Smith   PetscErrorCode ierr;
1594b9ad928SBarry Smith 
1604b9ad928SBarry Smith   PetscFunctionBegin;
161c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
162*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1632fa5cd67SKarl Rupp   if (residual) mglevels[l]->residual = residual;
16454b2cd4bSJed Brown   if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault;
16530b0564aSStefano Zampini   mglevels[l]->matresidual = PCMGMatResidualDefault;
166f3ae41bdSBarry Smith   if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);}
167298cc208SBarry Smith   ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr);
168f3fbd535SBarry Smith   mglevels[l]->A = mat;
1694b9ad928SBarry Smith   PetscFunctionReturn(0);
1704b9ad928SBarry Smith }
1714b9ad928SBarry Smith 
172fcb023d4SJed Brown /*@C
173fcb023d4SJed Brown    PCMGSetResidualTranspose - Sets the function to be used to calculate the residual of the transposed linear system
174fcb023d4SJed Brown    on the lth level.
175fcb023d4SJed Brown 
176fcb023d4SJed Brown    Logically Collective on PC
177fcb023d4SJed Brown 
178fcb023d4SJed Brown    Input Parameters:
179fcb023d4SJed Brown +  pc        - the multigrid context
180fcb023d4SJed Brown .  l         - the level (0 is coarsest) to supply
181fcb023d4SJed Brown .  residualt - function used to form transpose of residual, if none is provided the previously provide one is used, if no
182fcb023d4SJed Brown                previous one were provided then a default is used
183fcb023d4SJed Brown -  mat       - matrix associated with residual
184fcb023d4SJed Brown 
185fcb023d4SJed Brown    Level: advanced
186fcb023d4SJed Brown 
187fcb023d4SJed Brown .seealso: PCMGResidualTransposeDefault()
188fcb023d4SJed Brown @*/
189fcb023d4SJed Brown PetscErrorCode  PCMGSetResidualTranspose(PC pc,PetscInt l,PetscErrorCode (*residualt)(Mat,Vec,Vec,Vec),Mat mat)
190fcb023d4SJed Brown {
191fcb023d4SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
192fcb023d4SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
193fcb023d4SJed Brown   PetscErrorCode ierr;
194fcb023d4SJed Brown 
195fcb023d4SJed Brown   PetscFunctionBegin;
196fcb023d4SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
197*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
198fcb023d4SJed Brown   if (residualt) mglevels[l]->residualtranspose = residualt;
199fcb023d4SJed Brown   if (!mglevels[l]->residualtranspose) mglevels[l]->residualtranspose = PCMGResidualTransposeDefault;
20030b0564aSStefano Zampini   mglevels[l]->matresidualtranspose = PCMGMatResidualTransposeDefault;
201fcb023d4SJed Brown   if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);}
202fcb023d4SJed Brown   ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr);
203fcb023d4SJed Brown   mglevels[l]->A = mat;
204fcb023d4SJed Brown   PetscFunctionReturn(0);
205fcb023d4SJed Brown }
206fcb023d4SJed Brown 
2074b9ad928SBarry Smith /*@
208aea2a34eSBarry Smith    PCMGSetInterpolation - Sets the function to be used to calculate the
209bf5b2e24SBarry Smith    interpolation from l-1 to the lth level
2104b9ad928SBarry Smith 
211d083f849SBarry Smith    Logically Collective on PC
2124b9ad928SBarry Smith 
2134b9ad928SBarry Smith    Input Parameters:
2144b9ad928SBarry Smith +  pc  - the multigrid context
2154b9ad928SBarry Smith .  mat - the interpolation operator
216bf5b2e24SBarry Smith -  l   - the level (0 is coarsest) to supply [do not supply 0]
2174b9ad928SBarry Smith 
2184b9ad928SBarry Smith    Level: advanced
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith    Notes:
2214b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
2224b9ad928SBarry Smith     for the same level.
2234b9ad928SBarry Smith 
2244b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
2254b9ad928SBarry Smith     out from the matrix size which one it is.
2264b9ad928SBarry Smith 
22797177400SBarry Smith .seealso: PCMGSetRestriction()
2284b9ad928SBarry Smith @*/
2297087cfbeSBarry Smith PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
2304b9ad928SBarry Smith {
231f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
232f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
233fccaa45eSBarry Smith   PetscErrorCode ierr;
2344b9ad928SBarry Smith 
2354b9ad928SBarry Smith   PetscFunctionBegin;
236c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
237*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
238*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
239c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
2406bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr);
2412fa5cd67SKarl Rupp 
242f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
2434b9ad928SBarry Smith   PetscFunctionReturn(0);
2444b9ad928SBarry Smith }
2454b9ad928SBarry Smith 
2468a2c336bSFande Kong /*@
24795750439SFande Kong    PCMGSetOperators - Sets operator and preconditioning matrix for lth level
2488a2c336bSFande Kong 
249d083f849SBarry Smith    Logically Collective on PC
2508a2c336bSFande Kong 
2518a2c336bSFande Kong    Input Parameters:
2528a2c336bSFande Kong +  pc  - the multigrid context
2538a2c336bSFande Kong .  Amat - the operator
2548a2c336bSFande Kong .  pmat - the preconditioning operator
25595750439SFande Kong -  l   - the level (0 is the coarsest) to supply
2568a2c336bSFande Kong 
2578a2c336bSFande Kong    Level: advanced
2588a2c336bSFande Kong 
2598a2c336bSFande Kong .keywords:  multigrid, set, interpolate, level
2608a2c336bSFande Kong 
2618a2c336bSFande Kong .seealso: PCMGSetRestriction(), PCMGSetInterpolation()
2628a2c336bSFande Kong @*/
2638a2c336bSFande Kong PetscErrorCode  PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat)
264360ee056SFande Kong {
265360ee056SFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
266360ee056SFande Kong   PC_MG_Levels   **mglevels = mg->levels;
267360ee056SFande Kong   PetscErrorCode ierr;
268360ee056SFande Kong 
269360ee056SFande Kong   PetscFunctionBegin;
270360ee056SFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2718a2c336bSFande Kong   PetscValidHeaderSpecific(Amat,MAT_CLASSID,3);
2728a2c336bSFande Kong   PetscValidHeaderSpecific(Pmat,MAT_CLASSID,4);
273*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
2748a2c336bSFande Kong   ierr = KSPSetOperators(mglevels[l]->smoothd,Amat,Pmat);CHKERRQ(ierr);
275360ee056SFande Kong   PetscFunctionReturn(0);
276360ee056SFande Kong }
277360ee056SFande Kong 
278c88c5224SJed Brown /*@
279c88c5224SJed Brown    PCMGGetInterpolation - Gets the function to be used to calculate the
280c88c5224SJed Brown    interpolation from l-1 to the lth level
281c88c5224SJed Brown 
282c88c5224SJed Brown    Logically Collective on PC
283c88c5224SJed Brown 
284c88c5224SJed Brown    Input Parameters:
285c88c5224SJed Brown +  pc - the multigrid context
286c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
287c88c5224SJed Brown 
288c88c5224SJed Brown    Output Parameter:
2893ad4599aSBarry Smith .  mat - the interpolation matrix, can be NULL
290c88c5224SJed Brown 
291c88c5224SJed Brown    Level: advanced
292c88c5224SJed Brown 
293c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale()
294c88c5224SJed Brown @*/
295c88c5224SJed Brown PetscErrorCode  PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat)
296c88c5224SJed Brown {
297c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
298c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
299c88c5224SJed Brown   PetscErrorCode ierr;
300c88c5224SJed Brown 
301c88c5224SJed Brown   PetscFunctionBegin;
302c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3033ad4599aSBarry Smith   if (mat) PetscValidPointer(mat,3);
304*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
305*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
306c88c5224SJed Brown   if (!mglevels[l]->interpolate) {
307*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!mglevels[l]->restrct,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()");
308c88c5224SJed Brown     ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr);
309c88c5224SJed Brown   }
3103ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->interpolate;
311c88c5224SJed Brown   PetscFunctionReturn(0);
312c88c5224SJed Brown }
313c88c5224SJed Brown 
3148a2c336bSFande Kong /*@
315eab52d2bSLawrence Mitchell    PCMGSetRestriction - Sets the function to be used to restrict dual vectors
3164b9ad928SBarry Smith    from level l to l-1.
3174b9ad928SBarry Smith 
318d083f849SBarry Smith    Logically Collective on PC
3194b9ad928SBarry Smith 
3204b9ad928SBarry Smith    Input Parameters:
3214b9ad928SBarry Smith +  pc - the multigrid context
322c88c5224SJed Brown .  l - the level (0 is coarsest) to supply [Do not supply 0]
323c88c5224SJed Brown -  mat - the restriction matrix
3244b9ad928SBarry Smith 
3254b9ad928SBarry Smith    Level: advanced
3264b9ad928SBarry Smith 
3274b9ad928SBarry Smith    Notes:
3284b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
3294b9ad928SBarry Smith     for the same level.
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
3324b9ad928SBarry Smith     out from the matrix size which one it is.
3334b9ad928SBarry Smith 
334aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
335fccaa45eSBarry Smith     is used.
336fccaa45eSBarry Smith 
337a4355517SPierre Jolivet .seealso: PCMGSetInterpolation()
3384b9ad928SBarry Smith @*/
3397087cfbeSBarry Smith PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
3404b9ad928SBarry Smith {
341fccaa45eSBarry Smith   PetscErrorCode ierr;
342f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
343f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3444b9ad928SBarry Smith 
3454b9ad928SBarry Smith   PetscFunctionBegin;
346c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
347c88c5224SJed Brown   PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
348*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
349*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
350c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
3516bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr);
3522fa5cd67SKarl Rupp 
353f3fbd535SBarry Smith   mglevels[l]->restrct = mat;
3544b9ad928SBarry Smith   PetscFunctionReturn(0);
3554b9ad928SBarry Smith }
3564b9ad928SBarry Smith 
357c88c5224SJed Brown /*@
358eab52d2bSLawrence Mitchell    PCMGGetRestriction - Gets the function to be used to restrict dual vectors
359c88c5224SJed Brown    from level l to l-1.
360c88c5224SJed Brown 
361d083f849SBarry Smith    Logically Collective on PC
362c88c5224SJed Brown 
363c88c5224SJed Brown    Input Parameters:
364c88c5224SJed Brown +  pc - the multigrid context
365c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
366c88c5224SJed Brown 
367c88c5224SJed Brown    Output Parameter:
368c88c5224SJed Brown .  mat - the restriction matrix
369c88c5224SJed Brown 
370c88c5224SJed Brown    Level: advanced
371c88c5224SJed Brown 
372eab52d2bSLawrence Mitchell .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection()
373c88c5224SJed Brown @*/
374c88c5224SJed Brown PetscErrorCode  PCMGGetRestriction(PC pc,PetscInt l,Mat *mat)
375c88c5224SJed Brown {
376c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
377c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
378c88c5224SJed Brown   PetscErrorCode ierr;
379c88c5224SJed Brown 
380c88c5224SJed Brown   PetscFunctionBegin;
381c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3823ad4599aSBarry Smith   if (mat) PetscValidPointer(mat,3);
383*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
384*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
385c88c5224SJed Brown   if (!mglevels[l]->restrct) {
386*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!mglevels[l]->interpolate,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()");
387c88c5224SJed Brown     ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr);
388c88c5224SJed Brown   }
3893ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->restrct;
390c88c5224SJed Brown   PetscFunctionReturn(0);
391c88c5224SJed Brown }
392c88c5224SJed Brown 
39373250ac0SBarry Smith /*@
39473250ac0SBarry Smith    PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.
39573250ac0SBarry Smith 
396d083f849SBarry Smith    Logically Collective on PC
397c88c5224SJed Brown 
398c88c5224SJed Brown    Input Parameters:
399c88c5224SJed Brown +  pc - the multigrid context
400c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
401c88c5224SJed Brown .  rscale - the scaling
402c88c5224SJed Brown 
403c88c5224SJed Brown    Level: advanced
404c88c5224SJed Brown 
405c88c5224SJed Brown    Notes:
406eab52d2bSLawrence 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.
407c88c5224SJed Brown 
408eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection()
409c88c5224SJed Brown @*/
410c88c5224SJed Brown PetscErrorCode  PCMGSetRScale(PC pc,PetscInt l,Vec rscale)
411c88c5224SJed Brown {
412c88c5224SJed Brown   PetscErrorCode ierr;
413c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
414c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
415c88c5224SJed Brown 
416c88c5224SJed Brown   PetscFunctionBegin;
417c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
418*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
419*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
420c88c5224SJed Brown   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
421c88c5224SJed Brown   ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr);
4222fa5cd67SKarl Rupp 
423c88c5224SJed Brown   mglevels[l]->rscale = rscale;
424c88c5224SJed Brown   PetscFunctionReturn(0);
425c88c5224SJed Brown }
426c88c5224SJed Brown 
427c88c5224SJed Brown /*@
428c88c5224SJed Brown    PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1.
429c88c5224SJed Brown 
430c88c5224SJed Brown    Collective on PC
43173250ac0SBarry Smith 
43273250ac0SBarry Smith    Input Parameters:
43373250ac0SBarry Smith +  pc - the multigrid context
43473250ac0SBarry Smith .  rscale - the scaling
43573250ac0SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
43673250ac0SBarry Smith 
43773250ac0SBarry Smith    Level: advanced
43873250ac0SBarry Smith 
43973250ac0SBarry Smith    Notes:
440eab52d2bSLawrence 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.
44173250ac0SBarry Smith 
442eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection()
44373250ac0SBarry Smith @*/
444c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale)
44573250ac0SBarry Smith {
44673250ac0SBarry Smith   PetscErrorCode ierr;
44773250ac0SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
44873250ac0SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
44973250ac0SBarry Smith 
45073250ac0SBarry Smith   PetscFunctionBegin;
45173250ac0SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
452*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
453*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
454c88c5224SJed Brown   if (!mglevels[l]->rscale) {
455c88c5224SJed Brown     Mat      R;
456c88c5224SJed Brown     Vec      X,Y,coarse,fine;
457c88c5224SJed Brown     PetscInt M,N;
458c88c5224SJed Brown     ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr);
4592a7a6963SBarry Smith     ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr);
460c88c5224SJed Brown     ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr);
4612fa5cd67SKarl Rupp     if (M < N) {
4622fa5cd67SKarl Rupp       fine = X;
4632fa5cd67SKarl Rupp       coarse = Y;
4642fa5cd67SKarl Rupp     } else if (N < M) {
4652fa5cd67SKarl Rupp       fine = Y; coarse = X;
466ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser");
467c88c5224SJed Brown     ierr = VecSet(fine,1.);CHKERRQ(ierr);
468c88c5224SJed Brown     ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr);
469c88c5224SJed Brown     ierr = VecDestroy(&fine);CHKERRQ(ierr);
470c88c5224SJed Brown     ierr = VecReciprocal(coarse);CHKERRQ(ierr);
471c88c5224SJed Brown     mglevels[l]->rscale = coarse;
472c88c5224SJed Brown   }
473c88c5224SJed Brown   *rscale = mglevels[l]->rscale;
47473250ac0SBarry Smith   PetscFunctionReturn(0);
47573250ac0SBarry Smith }
47673250ac0SBarry Smith 
477f39d8e23SSatish Balay /*@
478eab52d2bSLawrence Mitchell    PCMGSetInjection - Sets the function to be used to inject primal vectors
479eab52d2bSLawrence Mitchell    from level l to l-1.
480eab52d2bSLawrence Mitchell 
481d083f849SBarry Smith    Logically Collective on PC
482eab52d2bSLawrence Mitchell 
483eab52d2bSLawrence Mitchell    Input Parameters:
484eab52d2bSLawrence Mitchell +  pc - the multigrid context
485eab52d2bSLawrence Mitchell .  l - the level (0 is coarsest) to supply [Do not supply 0]
486eab52d2bSLawrence Mitchell -  mat - the injection matrix
487eab52d2bSLawrence Mitchell 
488eab52d2bSLawrence Mitchell    Level: advanced
489eab52d2bSLawrence Mitchell 
490eab52d2bSLawrence Mitchell .seealso: PCMGSetRestriction()
491eab52d2bSLawrence Mitchell @*/
492eab52d2bSLawrence Mitchell PetscErrorCode  PCMGSetInjection(PC pc,PetscInt l,Mat mat)
493eab52d2bSLawrence Mitchell {
494eab52d2bSLawrence Mitchell   PetscErrorCode ierr;
495eab52d2bSLawrence Mitchell   PC_MG          *mg        = (PC_MG*)pc->data;
496eab52d2bSLawrence Mitchell   PC_MG_Levels   **mglevels = mg->levels;
497eab52d2bSLawrence Mitchell 
498eab52d2bSLawrence Mitchell   PetscFunctionBegin;
499eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
500eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
501*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
502*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
503eab52d2bSLawrence Mitchell   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
504eab52d2bSLawrence Mitchell   ierr = MatDestroy(&mglevels[l]->inject);CHKERRQ(ierr);
505eab52d2bSLawrence Mitchell 
506eab52d2bSLawrence Mitchell   mglevels[l]->inject = mat;
507eab52d2bSLawrence Mitchell   PetscFunctionReturn(0);
508eab52d2bSLawrence Mitchell }
509eab52d2bSLawrence Mitchell 
510eab52d2bSLawrence Mitchell /*@
511eab52d2bSLawrence Mitchell    PCMGGetInjection - Gets the function to be used to inject primal vectors
512eab52d2bSLawrence Mitchell    from level l to l-1.
513eab52d2bSLawrence Mitchell 
514d083f849SBarry Smith    Logically Collective on PC
515eab52d2bSLawrence Mitchell 
516eab52d2bSLawrence Mitchell    Input Parameters:
517eab52d2bSLawrence Mitchell +  pc - the multigrid context
518eab52d2bSLawrence Mitchell -  l - the level (0 is coarsest) to supply [Do not supply 0]
519eab52d2bSLawrence Mitchell 
520eab52d2bSLawrence Mitchell    Output Parameter:
52199a38656SLawrence Mitchell .  mat - the restriction matrix (may be NULL if no injection is available).
522eab52d2bSLawrence Mitchell 
523eab52d2bSLawrence Mitchell    Level: advanced
524eab52d2bSLawrence Mitchell 
525eab52d2bSLawrence Mitchell .seealso: PCMGSetInjection(), PCMGetGetRestriction()
526eab52d2bSLawrence Mitchell @*/
527eab52d2bSLawrence Mitchell PetscErrorCode  PCMGGetInjection(PC pc,PetscInt l,Mat *mat)
528eab52d2bSLawrence Mitchell {
529eab52d2bSLawrence Mitchell   PC_MG          *mg        = (PC_MG*)pc->data;
530eab52d2bSLawrence Mitchell   PC_MG_Levels   **mglevels = mg->levels;
531eab52d2bSLawrence Mitchell 
532eab52d2bSLawrence Mitchell   PetscFunctionBegin;
533eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
534eab52d2bSLawrence Mitchell   if (mat) PetscValidPointer(mat,3);
535*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
536*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(l <= 0 || mg->nlevels <= l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
537eab52d2bSLawrence Mitchell   if (mat) *mat = mglevels[l]->inject;
538eab52d2bSLawrence Mitchell   PetscFunctionReturn(0);
539eab52d2bSLawrence Mitchell }
540eab52d2bSLawrence Mitchell 
541eab52d2bSLawrence Mitchell /*@
54297177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
54397177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
54497177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
5454b9ad928SBarry Smith    post-smoothing.
5464b9ad928SBarry Smith 
5474b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
5484b9ad928SBarry Smith 
5494b9ad928SBarry Smith    Input Parameters:
5504b9ad928SBarry Smith +  pc - the multigrid context
5514b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
5524b9ad928SBarry Smith 
55301d2d390SJose E. Roman    Output Parameter:
5544b9ad928SBarry Smith .  ksp - the smoother
5554b9ad928SBarry Smith 
55657420d5bSBarry Smith    Notes:
55757420d5bSBarry 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.
55857420d5bSBarry Smith    You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc)
55957420d5bSBarry Smith    and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing.
56057420d5bSBarry Smith 
5614b9ad928SBarry Smith    Level: advanced
5624b9ad928SBarry Smith 
56328668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve()
5644b9ad928SBarry Smith @*/
5657087cfbeSBarry Smith PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
5664b9ad928SBarry Smith {
567f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
568f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
5694b9ad928SBarry Smith 
5704b9ad928SBarry Smith   PetscFunctionBegin;
571c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
572f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
5734b9ad928SBarry Smith   PetscFunctionReturn(0);
5744b9ad928SBarry Smith }
5754b9ad928SBarry Smith 
576f39d8e23SSatish Balay /*@
57797177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
5784b9ad928SBarry Smith    coarse grid correction (post-smoother).
5794b9ad928SBarry Smith 
5804b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
5814b9ad928SBarry Smith 
5824b9ad928SBarry Smith    Input Parameters:
5834b9ad928SBarry Smith +  pc - the multigrid context
5844b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
5854b9ad928SBarry Smith 
58601d2d390SJose E. Roman    Output Parameter:
5874b9ad928SBarry Smith .  ksp - the smoother
5884b9ad928SBarry Smith 
5894b9ad928SBarry Smith    Level: advanced
5904b9ad928SBarry Smith 
59195452b02SPatrick Sanan    Notes:
59295452b02SPatrick Sanan     calling this will result in a different pre and post smoother so you may need to
59389cce641SBarry Smith          set options on the pre smoother also
59489cce641SBarry Smith 
59597177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
5964b9ad928SBarry Smith @*/
5977087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
5984b9ad928SBarry Smith {
599f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
600f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
601dfbe8321SBarry Smith   PetscErrorCode ierr;
602f69a0ea3SMatthew Knepley   const char     *prefix;
6034b9ad928SBarry Smith   MPI_Comm       comm;
6044b9ad928SBarry Smith 
6054b9ad928SBarry Smith   PetscFunctionBegin;
606c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
6074b9ad928SBarry Smith   /*
6084b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
6094b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
6104b9ad928SBarry Smith      if not we allocate it.
6114b9ad928SBarry Smith   */
612*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
613f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
61419fd82e9SBarry Smith     KSPType     ksptype;
61519fd82e9SBarry Smith     PCType      pctype;
616336babb1SJed Brown     PC          ipc;
617336babb1SJed Brown     PetscReal   rtol,abstol,dtol;
618336babb1SJed Brown     PetscInt    maxits;
619336babb1SJed Brown     KSPNormType normtype;
620f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
621f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
622336babb1SJed Brown     ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr);
6233bf036e2SBarry Smith     ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr);
624336babb1SJed Brown     ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr);
625336babb1SJed Brown     ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr);
626336babb1SJed Brown     ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr);
627336babb1SJed Brown 
628f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
629422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr);
630f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
631f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
632336babb1SJed Brown     ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr);
633336babb1SJed Brown     ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr);
634336babb1SJed Brown     ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr);
6350059c7bdSJed Brown     ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
636336babb1SJed Brown     ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr);
637336babb1SJed Brown     ierr = PCSetType(ipc,pctype);CHKERRQ(ierr);
6383bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr);
639ab83eea4SMatthew G. Knepley     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr);
6404b9ad928SBarry Smith   }
641f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
6424b9ad928SBarry Smith   PetscFunctionReturn(0);
6434b9ad928SBarry Smith }
6444b9ad928SBarry Smith 
645f39d8e23SSatish Balay /*@
64697177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
6474b9ad928SBarry Smith    coarse grid correction (pre-smoother).
6484b9ad928SBarry Smith 
6494b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
6504b9ad928SBarry Smith 
6514b9ad928SBarry Smith    Input Parameters:
6524b9ad928SBarry Smith +  pc - the multigrid context
6534b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
6544b9ad928SBarry Smith 
65501d2d390SJose E. Roman    Output Parameter:
6564b9ad928SBarry Smith .  ksp - the smoother
6574b9ad928SBarry Smith 
6584b9ad928SBarry Smith    Level: advanced
6594b9ad928SBarry Smith 
66095452b02SPatrick Sanan    Notes:
66195452b02SPatrick Sanan     calling this will result in a different pre and post smoother so you may need to
66289cce641SBarry Smith          set options on the post smoother also
66389cce641SBarry Smith 
66497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
6654b9ad928SBarry Smith @*/
6667087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
6674b9ad928SBarry Smith {
668dfbe8321SBarry Smith   PetscErrorCode ierr;
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);
6744b9ad928SBarry Smith   /* make sure smoother up and down are different */
675c5eb9154SBarry Smith   if (l) {
6760298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr);
677d8148a5aSMatthew Knepley   }
678f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
6794b9ad928SBarry Smith   PetscFunctionReturn(0);
6804b9ad928SBarry Smith }
6814b9ad928SBarry Smith 
6824b9ad928SBarry Smith /*@
683cab31ae5SJed Brown    PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level.
6844b9ad928SBarry Smith 
685ad4df100SBarry Smith    Logically Collective on PC
6864b9ad928SBarry Smith 
6874b9ad928SBarry Smith    Input Parameters:
6884b9ad928SBarry Smith +  pc - the multigrid context
689c1cbb1deSBarry Smith .  l  - the level (0 is coarsest)
690c1cbb1deSBarry Smith -  c  - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
6914b9ad928SBarry Smith 
6924b9ad928SBarry Smith    Level: advanced
6934b9ad928SBarry Smith 
694c1cbb1deSBarry Smith .seealso: PCMGSetCycleType()
6954b9ad928SBarry Smith @*/
696c1cbb1deSBarry Smith PetscErrorCode  PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c)
6974b9ad928SBarry Smith {
698f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
699f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
7004b9ad928SBarry Smith 
7014b9ad928SBarry Smith   PetscFunctionBegin;
702c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
703*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
704c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,l,2);
705c51679f6SSatish Balay   PetscValidLogicalCollectiveEnum(pc,c,3);
706f3fbd535SBarry Smith   mglevels[l]->cycles = c;
7074b9ad928SBarry Smith   PetscFunctionReturn(0);
7084b9ad928SBarry Smith }
7094b9ad928SBarry Smith 
7104b9ad928SBarry Smith /*@
711f3b08a26SMatthew G. Knepley   PCMGSetRhs - Sets the vector to be used to store the right-hand side on a particular level.
7124b9ad928SBarry Smith 
713d083f849SBarry Smith    Logically Collective on PC
7144b9ad928SBarry Smith 
7154b9ad928SBarry Smith   Input Parameters:
7164b9ad928SBarry Smith + pc - the multigrid context
7174b9ad928SBarry Smith . l  - the level (0 is coarsest) this is to be used for
718f3b08a26SMatthew G. Knepley - c  - the Vec
7194b9ad928SBarry Smith 
7204b9ad928SBarry Smith   Level: advanced
7214b9ad928SBarry Smith 
72295452b02SPatrick Sanan   Notes:
723f3b08a26SMatthew 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.
724fccaa45eSBarry Smith 
725f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, right-hand-side, rhs, level
72697177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
7274b9ad928SBarry Smith @*/
7287087cfbeSBarry Smith PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
7294b9ad928SBarry Smith {
730fccaa45eSBarry Smith   PetscErrorCode ierr;
731f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
732f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
7334b9ad928SBarry Smith 
7344b9ad928SBarry Smith   PetscFunctionBegin;
735c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
736*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
737*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(l == mglevels[0]->levels-1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
738c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
7396bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
7402fa5cd67SKarl Rupp 
741f3fbd535SBarry Smith   mglevels[l]->b = c;
7424b9ad928SBarry Smith   PetscFunctionReturn(0);
7434b9ad928SBarry Smith }
7444b9ad928SBarry Smith 
7454b9ad928SBarry Smith /*@
746f3b08a26SMatthew G. Knepley   PCMGSetX - Sets the vector to be used to store the solution on a particular level.
7474b9ad928SBarry Smith 
748d083f849SBarry Smith   Logically Collective on PC
7494b9ad928SBarry Smith 
7504b9ad928SBarry Smith   Input Parameters:
7514b9ad928SBarry Smith + pc - the multigrid context
752251f4c67SDmitry Karpeev . l - the level (0 is coarsest) this is to be used for (do not supply the finest level)
753f3b08a26SMatthew G. Knepley - c - the Vec
7544b9ad928SBarry Smith 
7554b9ad928SBarry Smith   Level: advanced
7564b9ad928SBarry Smith 
75795452b02SPatrick Sanan   Notes:
758f3b08a26SMatthew 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.
759fccaa45eSBarry Smith 
760f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, solution, level
76197177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
7624b9ad928SBarry Smith @*/
7637087cfbeSBarry Smith PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
7644b9ad928SBarry Smith {
765fccaa45eSBarry Smith   PetscErrorCode ierr;
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);
771*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
772*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(l == mglevels[0]->levels-1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
773c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
7746bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
7752fa5cd67SKarl Rupp 
776f3fbd535SBarry Smith   mglevels[l]->x = c;
7774b9ad928SBarry Smith   PetscFunctionReturn(0);
7784b9ad928SBarry Smith }
7794b9ad928SBarry Smith 
7804b9ad928SBarry Smith /*@
781f3b08a26SMatthew G. Knepley   PCMGSetR - Sets the vector to be used to store the residual on a particular level.
7824b9ad928SBarry Smith 
783d083f849SBarry Smith   Logically Collective on PC
7844b9ad928SBarry Smith 
7854b9ad928SBarry Smith   Input Parameters:
7864b9ad928SBarry Smith + pc - the multigrid context
7874b9ad928SBarry Smith . l - the level (0 is coarsest) this is to be used for
788f3b08a26SMatthew G. Knepley - c - the Vec
7894b9ad928SBarry Smith 
7904b9ad928SBarry Smith   Level: advanced
7914b9ad928SBarry Smith 
79295452b02SPatrick Sanan   Notes:
793f3b08a26SMatthew 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.
794fccaa45eSBarry Smith 
795f3b08a26SMatthew G. Knepley .keywords: MG, multigrid, set, residual, level
796f3b08a26SMatthew G. Knepley .seealso: PCMGSetRhs(), PCMGSetX()
7974b9ad928SBarry Smith @*/
7987087cfbeSBarry Smith PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
7994b9ad928SBarry Smith {
800fccaa45eSBarry Smith   PetscErrorCode ierr;
801f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
802f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
8034b9ad928SBarry Smith 
8044b9ad928SBarry Smith   PetscFunctionBegin;
805c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
806*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
807*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!l,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
808c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
8096bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
8102fa5cd67SKarl Rupp 
811f3fbd535SBarry Smith   mglevels[l]->r = c;
8124b9ad928SBarry Smith   PetscFunctionReturn(0);
8134b9ad928SBarry Smith }
814