xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 3ad4599a42e9c85bdaf3228acb15d805000068aa)
14b9ad928SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>       /*I "petscksp.h" I*/
34b9ad928SBarry Smith 
4f9426fe0SMark Adams /* ---------------------------------------------------------------------------*/
54b9ad928SBarry Smith #undef __FUNCT__
654b2cd4bSJed Brown #define __FUNCT__ "PCMGResidualDefault"
71f6cc5b2SSatish Balay /*@C
854b2cd4bSJed Brown    PCMGResidualDefault - Default routine to calculate the residual.
94b9ad928SBarry Smith 
104b9ad928SBarry Smith    Collective on Mat and Vec
114b9ad928SBarry Smith 
124b9ad928SBarry Smith    Input Parameters:
134b9ad928SBarry Smith +  mat - the matrix
144b9ad928SBarry Smith .  b   - the right-hand-side
154b9ad928SBarry Smith -  x   - the approximate solution
164b9ad928SBarry Smith 
174b9ad928SBarry Smith    Output Parameter:
184b9ad928SBarry Smith .  r - location to store the residual
194b9ad928SBarry Smith 
20d0e4de75SBarry Smith    Level: developer
214b9ad928SBarry Smith 
224b9ad928SBarry Smith .keywords: MG, default, multigrid, residual
234b9ad928SBarry Smith 
2497177400SBarry Smith .seealso: PCMGSetResidual()
254b9ad928SBarry Smith @*/
2654b2cd4bSJed Brown PetscErrorCode  PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r)
274b9ad928SBarry Smith {
28dfbe8321SBarry Smith   PetscErrorCode ierr;
294b9ad928SBarry Smith 
304b9ad928SBarry Smith   PetscFunctionBegin;
31f9426fe0SMark Adams   ierr = MatResidual(mat,b,x,r);CHKERRQ(ierr);
324b9ad928SBarry Smith   PetscFunctionReturn(0);
334b9ad928SBarry Smith }
344b9ad928SBarry Smith 
354b9ad928SBarry Smith #undef __FUNCT__
36d6337fedSHong Zhang #define __FUNCT__ "PCMGGetCoarseSolve"
37f39d8e23SSatish Balay /*@
3897177400SBarry Smith    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
394b9ad928SBarry Smith 
404b9ad928SBarry Smith    Not Collective
414b9ad928SBarry Smith 
424b9ad928SBarry Smith    Input Parameter:
434b9ad928SBarry Smith .  pc - the multigrid context
444b9ad928SBarry Smith 
454b9ad928SBarry Smith    Output Parameter:
464b9ad928SBarry Smith .  ksp - the coarse grid solver context
474b9ad928SBarry Smith 
484b9ad928SBarry Smith    Level: advanced
494b9ad928SBarry Smith 
504b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid
514b9ad928SBarry Smith @*/
527087cfbeSBarry Smith PetscErrorCode  PCMGGetCoarseSolve(PC pc,KSP *ksp)
534b9ad928SBarry Smith {
54f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
55f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
564b9ad928SBarry Smith 
574b9ad928SBarry Smith   PetscFunctionBegin;
58c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
59f3fbd535SBarry Smith   *ksp =  mglevels[0]->smoothd;
604b9ad928SBarry Smith   PetscFunctionReturn(0);
614b9ad928SBarry Smith }
624b9ad928SBarry Smith 
634b9ad928SBarry Smith #undef __FUNCT__
64d6337fedSHong Zhang #define __FUNCT__ "PCMGSetResidual"
654b9ad928SBarry Smith /*@C
6697177400SBarry Smith    PCMGSetResidual - Sets the function to be used to calculate the residual
674b9ad928SBarry Smith    on the lth level.
684b9ad928SBarry Smith 
69ad4df100SBarry Smith    Logically Collective on PC and Mat
704b9ad928SBarry Smith 
714b9ad928SBarry Smith    Input Parameters:
724b9ad928SBarry Smith +  pc       - the multigrid context
734b9ad928SBarry Smith .  l        - the level (0 is coarsest) to supply
74157726a2SBarry Smith .  residual - function used to form residual, if none is provided the previously provide one is used, if no
75d0e4de75SBarry Smith               previous one were provided then a default is used
764b9ad928SBarry Smith -  mat      - matrix associated with residual
774b9ad928SBarry Smith 
784b9ad928SBarry Smith    Level: advanced
794b9ad928SBarry Smith 
804b9ad928SBarry Smith .keywords:  MG, set, multigrid, residual, level
814b9ad928SBarry Smith 
8254b2cd4bSJed Brown .seealso: PCMGResidualDefault()
834b9ad928SBarry Smith @*/
847087cfbeSBarry Smith PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
854b9ad928SBarry Smith {
86f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
87f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
88298cc208SBarry Smith   PetscErrorCode ierr;
894b9ad928SBarry Smith 
904b9ad928SBarry Smith   PetscFunctionBegin;
91c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
92ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
932fa5cd67SKarl Rupp   if (residual) mglevels[l]->residual = residual;
9454b2cd4bSJed Brown   if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault;
95f3ae41bdSBarry Smith   if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);}
96298cc208SBarry Smith   ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr);
97f3fbd535SBarry Smith   mglevels[l]->A = mat;
984b9ad928SBarry Smith   PetscFunctionReturn(0);
994b9ad928SBarry Smith }
1004b9ad928SBarry Smith 
1014b9ad928SBarry Smith #undef __FUNCT__
102d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation"
1034b9ad928SBarry Smith /*@
104aea2a34eSBarry Smith    PCMGSetInterpolation - Sets the function to be used to calculate the
105bf5b2e24SBarry Smith    interpolation from l-1 to the lth level
1064b9ad928SBarry Smith 
107ad4df100SBarry Smith    Logically Collective on PC and Mat
1084b9ad928SBarry Smith 
1094b9ad928SBarry Smith    Input Parameters:
1104b9ad928SBarry Smith +  pc  - the multigrid context
1114b9ad928SBarry Smith .  mat - the interpolation operator
112bf5b2e24SBarry Smith -  l   - the level (0 is coarsest) to supply [do not supply 0]
1134b9ad928SBarry Smith 
1144b9ad928SBarry Smith    Level: advanced
1154b9ad928SBarry Smith 
1164b9ad928SBarry Smith    Notes:
1174b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
1184b9ad928SBarry Smith     for the same level.
1194b9ad928SBarry Smith 
1204b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1214b9ad928SBarry Smith     out from the matrix size which one it is.
1224b9ad928SBarry Smith 
1234b9ad928SBarry Smith .keywords:  multigrid, set, interpolate, level
1244b9ad928SBarry Smith 
12597177400SBarry Smith .seealso: PCMGSetRestriction()
1264b9ad928SBarry Smith @*/
1277087cfbeSBarry Smith PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
1284b9ad928SBarry Smith {
129f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
130f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
131fccaa45eSBarry Smith   PetscErrorCode ierr;
1324b9ad928SBarry Smith 
1334b9ad928SBarry Smith   PetscFunctionBegin;
134c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
135ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
136ce94432eSBarry Smith   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
137c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1386bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr);
1392fa5cd67SKarl Rupp 
140f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
1414b9ad928SBarry Smith   PetscFunctionReturn(0);
1424b9ad928SBarry Smith }
1434b9ad928SBarry Smith 
1444b9ad928SBarry Smith #undef __FUNCT__
145c88c5224SJed Brown #define __FUNCT__ "PCMGGetInterpolation"
146c88c5224SJed Brown /*@
147c88c5224SJed Brown    PCMGGetInterpolation - Gets the function to be used to calculate the
148c88c5224SJed Brown    interpolation from l-1 to the lth level
149c88c5224SJed Brown 
150c88c5224SJed Brown    Logically Collective on PC
151c88c5224SJed Brown 
152c88c5224SJed Brown    Input Parameters:
153c88c5224SJed Brown +  pc - the multigrid context
154c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
155c88c5224SJed Brown 
156c88c5224SJed Brown    Output Parameter:
157*3ad4599aSBarry Smith .  mat - the interpolation matrix, can be NULL
158c88c5224SJed Brown 
159c88c5224SJed Brown    Level: advanced
160c88c5224SJed Brown 
161c88c5224SJed Brown .keywords: MG, get, multigrid, interpolation, level
162c88c5224SJed Brown 
163c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale()
164c88c5224SJed Brown @*/
165c88c5224SJed Brown PetscErrorCode  PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat)
166c88c5224SJed Brown {
167c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
168c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
169c88c5224SJed Brown   PetscErrorCode ierr;
170c88c5224SJed Brown 
171c88c5224SJed Brown   PetscFunctionBegin;
172c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
173*3ad4599aSBarry Smith   if (mat) PetscValidPointer(mat,3);
174ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
175ce94432eSBarry Smith   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
176c88c5224SJed Brown   if (!mglevels[l]->interpolate) {
1775aa31b60SBarry Smith     if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()");
178c88c5224SJed Brown     ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr);
179c88c5224SJed Brown   }
180*3ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->interpolate;
181c88c5224SJed Brown   PetscFunctionReturn(0);
182c88c5224SJed Brown }
183c88c5224SJed Brown 
184c88c5224SJed Brown #undef __FUNCT__
185d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction"
1864b9ad928SBarry Smith /*@
18797177400SBarry Smith    PCMGSetRestriction - Sets the function to be used to restrict vector
1884b9ad928SBarry Smith    from level l to l-1.
1894b9ad928SBarry Smith 
190ad4df100SBarry Smith    Logically Collective on PC and Mat
1914b9ad928SBarry Smith 
1924b9ad928SBarry Smith    Input Parameters:
1934b9ad928SBarry Smith +  pc - the multigrid context
194c88c5224SJed Brown .  l - the level (0 is coarsest) to supply [Do not supply 0]
195c88c5224SJed Brown -  mat - the restriction matrix
1964b9ad928SBarry Smith 
1974b9ad928SBarry Smith    Level: advanced
1984b9ad928SBarry Smith 
1994b9ad928SBarry Smith    Notes:
2004b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
2014b9ad928SBarry Smith     for the same level.
2024b9ad928SBarry Smith 
2034b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
2044b9ad928SBarry Smith     out from the matrix size which one it is.
2054b9ad928SBarry Smith 
206aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
207fccaa45eSBarry Smith     is used.
208fccaa45eSBarry Smith 
2094b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
2104b9ad928SBarry Smith 
211aea2a34eSBarry Smith .seealso: PCMGSetInterpolation()
2124b9ad928SBarry Smith @*/
2137087cfbeSBarry Smith PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
2144b9ad928SBarry Smith {
215fccaa45eSBarry Smith   PetscErrorCode ierr;
216f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
217f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2184b9ad928SBarry Smith 
2194b9ad928SBarry Smith   PetscFunctionBegin;
220c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
221c88c5224SJed Brown   PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
222ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
223ce94432eSBarry Smith   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
224c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
2256bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr);
2262fa5cd67SKarl Rupp 
227f3fbd535SBarry Smith   mglevels[l]->restrct = mat;
2284b9ad928SBarry Smith   PetscFunctionReturn(0);
2294b9ad928SBarry Smith }
2304b9ad928SBarry Smith 
2314b9ad928SBarry Smith #undef __FUNCT__
232c88c5224SJed Brown #define __FUNCT__ "PCMGGetRestriction"
233c88c5224SJed Brown /*@
234c88c5224SJed Brown    PCMGGetRestriction - Gets the function to be used to restrict vector
235c88c5224SJed Brown    from level l to l-1.
236c88c5224SJed Brown 
237c88c5224SJed Brown    Logically Collective on PC and Mat
238c88c5224SJed Brown 
239c88c5224SJed Brown    Input Parameters:
240c88c5224SJed Brown +  pc - the multigrid context
241c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
242c88c5224SJed Brown 
243c88c5224SJed Brown    Output Parameter:
244c88c5224SJed Brown .  mat - the restriction matrix
245c88c5224SJed Brown 
246c88c5224SJed Brown    Level: advanced
247c88c5224SJed Brown 
248c88c5224SJed Brown .keywords: MG, get, multigrid, restriction, level
249c88c5224SJed Brown 
250c88c5224SJed Brown .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale()
251c88c5224SJed Brown @*/
252c88c5224SJed Brown PetscErrorCode  PCMGGetRestriction(PC pc,PetscInt l,Mat *mat)
253c88c5224SJed Brown {
254c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
255c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
256c88c5224SJed Brown   PetscErrorCode ierr;
257c88c5224SJed Brown 
258c88c5224SJed Brown   PetscFunctionBegin;
259c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
260*3ad4599aSBarry Smith   if (mat) PetscValidPointer(mat,3);
261ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
262ce94432eSBarry Smith   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
263c88c5224SJed Brown   if (!mglevels[l]->restrct) {
264ce94432eSBarry Smith     if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()");
265c88c5224SJed Brown     ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr);
266c88c5224SJed Brown   }
267*3ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->restrct;
268c88c5224SJed Brown   PetscFunctionReturn(0);
269c88c5224SJed Brown }
270c88c5224SJed Brown 
271c88c5224SJed Brown #undef __FUNCT__
27273250ac0SBarry Smith #define __FUNCT__ "PCMGSetRScale"
27373250ac0SBarry Smith /*@
27473250ac0SBarry Smith    PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.
27573250ac0SBarry Smith 
276c88c5224SJed Brown    Logically Collective on PC and Vec
277c88c5224SJed Brown 
278c88c5224SJed Brown    Input Parameters:
279c88c5224SJed Brown +  pc - the multigrid context
280c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
281c88c5224SJed Brown .  rscale - the scaling
282c88c5224SJed Brown 
283c88c5224SJed Brown    Level: advanced
284c88c5224SJed Brown 
285c88c5224SJed Brown    Notes:
286c88c5224SJed Brown        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.
287c88c5224SJed Brown 
288c88c5224SJed Brown .keywords: MG, set, multigrid, restriction, level
289c88c5224SJed Brown 
290c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale()
291c88c5224SJed Brown @*/
292c88c5224SJed Brown PetscErrorCode  PCMGSetRScale(PC pc,PetscInt l,Vec rscale)
293c88c5224SJed Brown {
294c88c5224SJed Brown   PetscErrorCode ierr;
295c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
296c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
297c88c5224SJed Brown 
298c88c5224SJed Brown   PetscFunctionBegin;
299c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
300ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
301ce94432eSBarry Smith   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
302c88c5224SJed Brown   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
303c88c5224SJed Brown   ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr);
3042fa5cd67SKarl Rupp 
305c88c5224SJed Brown   mglevels[l]->rscale = rscale;
306c88c5224SJed Brown   PetscFunctionReturn(0);
307c88c5224SJed Brown }
308c88c5224SJed Brown 
309c88c5224SJed Brown #undef __FUNCT__
310c88c5224SJed Brown #define __FUNCT__ "PCMGGetRScale"
311c88c5224SJed Brown /*@
312c88c5224SJed Brown    PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1.
313c88c5224SJed Brown 
314c88c5224SJed Brown    Collective on PC
31573250ac0SBarry Smith 
31673250ac0SBarry Smith    Input Parameters:
31773250ac0SBarry Smith +  pc - the multigrid context
31873250ac0SBarry Smith .  rscale - the scaling
31973250ac0SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
32073250ac0SBarry Smith 
32173250ac0SBarry Smith    Level: advanced
32273250ac0SBarry Smith 
32373250ac0SBarry Smith    Notes:
32473250ac0SBarry Smith        When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R.
32573250ac0SBarry Smith 
32673250ac0SBarry Smith .keywords: MG, set, multigrid, restriction, level
32773250ac0SBarry Smith 
328c88c5224SJed Brown .seealso: PCMGSetInterpolation(), PCMGGetRestriction()
32973250ac0SBarry Smith @*/
330c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale)
33173250ac0SBarry Smith {
33273250ac0SBarry Smith   PetscErrorCode ierr;
33373250ac0SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
33473250ac0SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
33573250ac0SBarry Smith 
33673250ac0SBarry Smith   PetscFunctionBegin;
33773250ac0SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
338ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
339ce94432eSBarry Smith   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
340c88c5224SJed Brown   if (!mglevels[l]->rscale) {
341c88c5224SJed Brown     Mat      R;
342c88c5224SJed Brown     Vec      X,Y,coarse,fine;
343c88c5224SJed Brown     PetscInt M,N;
344c88c5224SJed Brown     ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr);
3452a7a6963SBarry Smith     ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr);
346c88c5224SJed Brown     ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr);
3472fa5cd67SKarl Rupp     if (M < N) {
3482fa5cd67SKarl Rupp       fine = X;
3492fa5cd67SKarl Rupp       coarse = Y;
3502fa5cd67SKarl Rupp     } else if (N < M) {
3512fa5cd67SKarl Rupp       fine = Y; coarse = X;
352ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser");
353c88c5224SJed Brown     ierr = VecSet(fine,1.);CHKERRQ(ierr);
354c88c5224SJed Brown     ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr);
355c88c5224SJed Brown     ierr = VecDestroy(&fine);CHKERRQ(ierr);
356c88c5224SJed Brown     ierr = VecReciprocal(coarse);CHKERRQ(ierr);
357c88c5224SJed Brown     mglevels[l]->rscale = coarse;
358c88c5224SJed Brown   }
359c88c5224SJed Brown   *rscale = mglevels[l]->rscale;
36073250ac0SBarry Smith   PetscFunctionReturn(0);
36173250ac0SBarry Smith }
36273250ac0SBarry Smith 
36373250ac0SBarry Smith #undef __FUNCT__
364d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother"
365f39d8e23SSatish Balay /*@
36697177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
36797177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
36897177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
3694b9ad928SBarry Smith    post-smoothing.
3704b9ad928SBarry Smith 
3714b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
3724b9ad928SBarry Smith 
3734b9ad928SBarry Smith    Input Parameters:
3744b9ad928SBarry Smith +  pc - the multigrid context
3754b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
3764b9ad928SBarry Smith 
3774b9ad928SBarry Smith    Ouput Parameters:
3784b9ad928SBarry Smith .  ksp - the smoother
3794b9ad928SBarry Smith 
38057420d5bSBarry Smith    Notes:
38157420d5bSBarry 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.
38257420d5bSBarry Smith    You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc)
38357420d5bSBarry Smith    and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing.
38457420d5bSBarry Smith 
3854b9ad928SBarry Smith    Level: advanced
3864b9ad928SBarry Smith 
3874b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
3884b9ad928SBarry Smith 
38997177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
3904b9ad928SBarry Smith @*/
3917087cfbeSBarry Smith PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
3924b9ad928SBarry Smith {
393f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
394f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
3954b9ad928SBarry Smith 
3964b9ad928SBarry Smith   PetscFunctionBegin;
397c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
398f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
3994b9ad928SBarry Smith   PetscFunctionReturn(0);
4004b9ad928SBarry Smith }
4014b9ad928SBarry Smith 
4024b9ad928SBarry Smith #undef __FUNCT__
403d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp"
404f39d8e23SSatish Balay /*@
40597177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
4064b9ad928SBarry Smith    coarse grid correction (post-smoother).
4074b9ad928SBarry Smith 
4084b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
4094b9ad928SBarry Smith 
4104b9ad928SBarry Smith    Input Parameters:
4114b9ad928SBarry Smith +  pc - the multigrid context
4124b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith    Ouput Parameters:
4154b9ad928SBarry Smith .  ksp - the smoother
4164b9ad928SBarry Smith 
4174b9ad928SBarry Smith    Level: advanced
4184b9ad928SBarry Smith 
41989cce641SBarry Smith    Notes: calling this will result in a different pre and post smoother so you may need to
42089cce641SBarry Smith          set options on the pre smoother also
42189cce641SBarry Smith 
4224b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
4234b9ad928SBarry Smith 
42497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
4254b9ad928SBarry Smith @*/
4267087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
4274b9ad928SBarry Smith {
428f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
429f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
430dfbe8321SBarry Smith   PetscErrorCode ierr;
431f69a0ea3SMatthew Knepley   const char     *prefix;
4324b9ad928SBarry Smith   MPI_Comm       comm;
4334b9ad928SBarry Smith 
4344b9ad928SBarry Smith   PetscFunctionBegin;
435c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4364b9ad928SBarry Smith   /*
4374b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
4384b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
4394b9ad928SBarry Smith      if not we allocate it.
4404b9ad928SBarry Smith   */
441ce94432eSBarry Smith   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
442f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
44319fd82e9SBarry Smith     KSPType     ksptype;
44419fd82e9SBarry Smith     PCType      pctype;
445336babb1SJed Brown     PC          ipc;
446336babb1SJed Brown     PetscReal   rtol,abstol,dtol;
447336babb1SJed Brown     PetscInt    maxits;
448336babb1SJed Brown     KSPNormType normtype;
449f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
450f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
451336babb1SJed Brown     ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr);
4523bf036e2SBarry Smith     ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr);
453336babb1SJed Brown     ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr);
454336babb1SJed Brown     ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr);
455336babb1SJed Brown     ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr);
456336babb1SJed Brown 
457f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
458422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr);
459f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
460f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
461336babb1SJed Brown     ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr);
462336babb1SJed Brown     ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr);
463336babb1SJed Brown     ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr);
4640059c7bdSJed Brown     ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
465336babb1SJed Brown     ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr);
466336babb1SJed Brown     ierr = PCSetType(ipc,pctype);CHKERRQ(ierr);
4673bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr);
468ab83eea4SMatthew G. Knepley     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr);
4694b9ad928SBarry Smith   }
470f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
4714b9ad928SBarry Smith   PetscFunctionReturn(0);
4724b9ad928SBarry Smith }
4734b9ad928SBarry Smith 
4744b9ad928SBarry Smith #undef __FUNCT__
4759dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown"
476f39d8e23SSatish Balay /*@
47797177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
4784b9ad928SBarry Smith    coarse grid correction (pre-smoother).
4794b9ad928SBarry Smith 
4804b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
4814b9ad928SBarry Smith 
4824b9ad928SBarry Smith    Input Parameters:
4834b9ad928SBarry Smith +  pc - the multigrid context
4844b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
4854b9ad928SBarry Smith 
4864b9ad928SBarry Smith    Ouput Parameters:
4874b9ad928SBarry Smith .  ksp - the smoother
4884b9ad928SBarry Smith 
4894b9ad928SBarry Smith    Level: advanced
4904b9ad928SBarry Smith 
49189cce641SBarry Smith    Notes: calling this will result in a different pre and post smoother so you may need to
49289cce641SBarry Smith          set options on the post smoother also
49389cce641SBarry Smith 
4944b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
4954b9ad928SBarry Smith 
49697177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
4974b9ad928SBarry Smith @*/
4987087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
4994b9ad928SBarry Smith {
500dfbe8321SBarry Smith   PetscErrorCode ierr;
501f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
502f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
5034b9ad928SBarry Smith 
5044b9ad928SBarry Smith   PetscFunctionBegin;
505c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5064b9ad928SBarry Smith   /* make sure smoother up and down are different */
507c5eb9154SBarry Smith   if (l) {
5080298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr);
509d8148a5aSMatthew Knepley   }
510f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
5114b9ad928SBarry Smith   PetscFunctionReturn(0);
5124b9ad928SBarry Smith }
5134b9ad928SBarry Smith 
5144b9ad928SBarry Smith #undef __FUNCT__
5159dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel"
5164b9ad928SBarry Smith /*@
51797177400SBarry Smith    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
5184b9ad928SBarry Smith 
519ad4df100SBarry Smith    Logically Collective on PC
5204b9ad928SBarry Smith 
5214b9ad928SBarry Smith    Input Parameters:
5224b9ad928SBarry Smith +  pc - the multigrid context
5234b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
5244b9ad928SBarry Smith -  n  - the number of cycles
5254b9ad928SBarry Smith 
5264b9ad928SBarry Smith    Level: advanced
5274b9ad928SBarry Smith 
5284b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
5294b9ad928SBarry Smith 
53097177400SBarry Smith .seealso: PCMGSetCycles()
5314b9ad928SBarry Smith @*/
5327087cfbeSBarry Smith PetscErrorCode  PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
5334b9ad928SBarry Smith {
534f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
535f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
5364b9ad928SBarry Smith 
5374b9ad928SBarry Smith   PetscFunctionBegin;
538c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
539ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
540c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,l,2);
541c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,c,3);
542f3fbd535SBarry Smith   mglevels[l]->cycles = c;
5434b9ad928SBarry Smith   PetscFunctionReturn(0);
5444b9ad928SBarry Smith }
5454b9ad928SBarry Smith 
5464b9ad928SBarry Smith #undef __FUNCT__
547d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs"
5484b9ad928SBarry Smith /*@
54997177400SBarry Smith    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
550fccaa45eSBarry Smith    on a particular level.
5514b9ad928SBarry Smith 
552ad4df100SBarry Smith    Logically Collective on PC and Vec
5534b9ad928SBarry Smith 
5544b9ad928SBarry Smith    Input Parameters:
5554b9ad928SBarry Smith +  pc - the multigrid context
5564b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
5574b9ad928SBarry Smith -  c  - the space
5584b9ad928SBarry Smith 
5594b9ad928SBarry Smith    Level: advanced
5604b9ad928SBarry Smith 
561fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
562fccaa45eSBarry Smith 
563fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
564fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
565fccaa45eSBarry Smith 
5664b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
5674b9ad928SBarry Smith 
56897177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
5694b9ad928SBarry Smith @*/
5707087cfbeSBarry Smith PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
5714b9ad928SBarry Smith {
572fccaa45eSBarry Smith   PetscErrorCode ierr;
573f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
574f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
5754b9ad928SBarry Smith 
5764b9ad928SBarry Smith   PetscFunctionBegin;
577c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
578ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
579ce94432eSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
580c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
5816bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
5822fa5cd67SKarl Rupp 
583f3fbd535SBarry Smith   mglevels[l]->b = c;
5844b9ad928SBarry Smith   PetscFunctionReturn(0);
5854b9ad928SBarry Smith }
5864b9ad928SBarry Smith 
5874b9ad928SBarry Smith #undef __FUNCT__
588d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX"
5894b9ad928SBarry Smith /*@
59097177400SBarry Smith    PCMGSetX - Sets the vector space to be used to store the solution on a
591fccaa45eSBarry Smith    particular level.
5924b9ad928SBarry Smith 
593ad4df100SBarry Smith    Logically Collective on PC and Vec
5944b9ad928SBarry Smith 
5954b9ad928SBarry Smith    Input Parameters:
5964b9ad928SBarry Smith +  pc - the multigrid context
597251f4c67SDmitry Karpeev .  l - the level (0 is coarsest) this is to be used for (do not supply the finest level)
5984b9ad928SBarry Smith -  c - the space
5994b9ad928SBarry Smith 
6004b9ad928SBarry Smith    Level: advanced
6014b9ad928SBarry Smith 
602fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
603fccaa45eSBarry Smith 
604fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
605fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
606fccaa45eSBarry Smith 
6074b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
6084b9ad928SBarry Smith 
60997177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
6104b9ad928SBarry Smith @*/
6117087cfbeSBarry Smith PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
6124b9ad928SBarry Smith {
613fccaa45eSBarry Smith   PetscErrorCode ierr;
614f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
615f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
6164b9ad928SBarry Smith 
6174b9ad928SBarry Smith   PetscFunctionBegin;
618c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
619ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
620ce94432eSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
621c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
6226bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
6232fa5cd67SKarl Rupp 
624f3fbd535SBarry Smith   mglevels[l]->x = c;
6254b9ad928SBarry Smith   PetscFunctionReturn(0);
6264b9ad928SBarry Smith }
6274b9ad928SBarry Smith 
6284b9ad928SBarry Smith #undef __FUNCT__
629d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR"
6304b9ad928SBarry Smith /*@
63197177400SBarry Smith    PCMGSetR - Sets the vector space to be used to store the residual on a
632fccaa45eSBarry Smith    particular level.
6334b9ad928SBarry Smith 
634ad4df100SBarry Smith    Logically Collective on PC and Vec
6354b9ad928SBarry Smith 
6364b9ad928SBarry Smith    Input Parameters:
6374b9ad928SBarry Smith +  pc - the multigrid context
6384b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
6394b9ad928SBarry Smith -  c - the space
6404b9ad928SBarry Smith 
6414b9ad928SBarry Smith    Level: advanced
6424b9ad928SBarry Smith 
643fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
644fccaa45eSBarry Smith 
645fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
646fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
647fccaa45eSBarry Smith 
6484b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
6494b9ad928SBarry Smith @*/
6507087cfbeSBarry Smith PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
6514b9ad928SBarry Smith {
652fccaa45eSBarry Smith   PetscErrorCode ierr;
653f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
654f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
6554b9ad928SBarry Smith 
6564b9ad928SBarry Smith   PetscFunctionBegin;
657c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
658ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
659ce94432eSBarry Smith   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
660c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
6616bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
6622fa5cd67SKarl Rupp 
663f3fbd535SBarry Smith   mglevels[l]->r = c;
6644b9ad928SBarry Smith   PetscFunctionReturn(0);
6654b9ad928SBarry Smith }
666