xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 157726a207616ee5b74eb502738d680c1e285726)
14b9ad928SBarry Smith 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h>       /*I "petscksp.h" I*/
3f3f935beSJed Brown                           /*I "petscpcmg.h"   I*/
44b9ad928SBarry Smith 
54b9ad928SBarry Smith #undef __FUNCT__
697177400SBarry Smith #define __FUNCT__ "PCMGDefaultResidual"
71f6cc5b2SSatish Balay /*@C
897177400SBarry Smith    PCMGDefaultResidual - 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 
204b9ad928SBarry Smith    Level: advanced
214b9ad928SBarry Smith 
224b9ad928SBarry Smith .keywords: MG, default, multigrid, residual
234b9ad928SBarry Smith 
2497177400SBarry Smith .seealso: PCMGSetResidual()
254b9ad928SBarry Smith @*/
267087cfbeSBarry Smith PetscErrorCode  PCMGDefaultResidual(Mat mat,Vec b,Vec x,Vec r)
274b9ad928SBarry Smith {
28dfbe8321SBarry Smith   PetscErrorCode ierr;
294b9ad928SBarry Smith 
304b9ad928SBarry Smith   PetscFunctionBegin;
314b9ad928SBarry Smith   ierr = MatMult(mat,x,r);CHKERRQ(ierr);
32efb30889SBarry Smith   ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
334b9ad928SBarry Smith   PetscFunctionReturn(0);
344b9ad928SBarry Smith }
354b9ad928SBarry Smith 
364b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/
374b9ad928SBarry Smith 
384b9ad928SBarry Smith #undef __FUNCT__
39d6337fedSHong Zhang #define __FUNCT__ "PCMGGetCoarseSolve"
40f39d8e23SSatish Balay /*@
4197177400SBarry Smith    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
424b9ad928SBarry Smith 
434b9ad928SBarry Smith    Not Collective
444b9ad928SBarry Smith 
454b9ad928SBarry Smith    Input Parameter:
464b9ad928SBarry Smith .  pc - the multigrid context
474b9ad928SBarry Smith 
484b9ad928SBarry Smith    Output Parameter:
494b9ad928SBarry Smith .  ksp - the coarse grid solver context
504b9ad928SBarry Smith 
514b9ad928SBarry Smith    Level: advanced
524b9ad928SBarry Smith 
534b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid
544b9ad928SBarry Smith @*/
557087cfbeSBarry Smith PetscErrorCode  PCMGGetCoarseSolve(PC pc,KSP *ksp)
564b9ad928SBarry Smith {
57f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
58f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
594b9ad928SBarry Smith 
604b9ad928SBarry Smith   PetscFunctionBegin;
61c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
62f3fbd535SBarry Smith   *ksp =  mglevels[0]->smoothd;
634b9ad928SBarry Smith   PetscFunctionReturn(0);
644b9ad928SBarry Smith }
654b9ad928SBarry Smith 
664b9ad928SBarry Smith #undef __FUNCT__
67d6337fedSHong Zhang #define __FUNCT__ "PCMGSetResidual"
684b9ad928SBarry Smith /*@C
6997177400SBarry Smith    PCMGSetResidual - Sets the function to be used to calculate the residual
704b9ad928SBarry Smith    on the lth level.
714b9ad928SBarry Smith 
72ad4df100SBarry Smith    Logically Collective on PC and Mat
734b9ad928SBarry Smith 
744b9ad928SBarry Smith    Input Parameters:
754b9ad928SBarry Smith +  pc       - the multigrid context
764b9ad928SBarry Smith .  l        - the level (0 is coarsest) to supply
77*157726a2SBarry Smith .  residual - function used to form residual, if none is provided the previously provide one is used, if no
78*157726a2SBarry Smith               previous one were provided then PCMGDefaultResidual() is used
794b9ad928SBarry Smith -  mat      - matrix associated with residual
804b9ad928SBarry Smith 
814b9ad928SBarry Smith    Level: advanced
824b9ad928SBarry Smith 
834b9ad928SBarry Smith .keywords:  MG, set, multigrid, residual, level
844b9ad928SBarry Smith 
8597177400SBarry Smith .seealso: PCMGDefaultResidual()
864b9ad928SBarry Smith @*/
877087cfbeSBarry Smith PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
884b9ad928SBarry Smith {
89f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
90f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
91298cc208SBarry Smith   PetscErrorCode ierr;
924b9ad928SBarry Smith 
934b9ad928SBarry Smith   PetscFunctionBegin;
94c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
95e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
96*157726a2SBarry Smith   if (residual) {
97f3fbd535SBarry Smith     mglevels[l]->residual = residual;
98*157726a2SBarry Smith   } if (!mglevels[l]->residual) {
99*157726a2SBarry Smith     mglevels[l]->residual = PCMGDefaultResidual;
100*157726a2SBarry Smith   }
101f3ae41bdSBarry Smith   if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);}
102298cc208SBarry Smith   ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr);
103f3fbd535SBarry Smith   mglevels[l]->A        = mat;
1044b9ad928SBarry Smith   PetscFunctionReturn(0);
1054b9ad928SBarry Smith }
1064b9ad928SBarry Smith 
1074b9ad928SBarry Smith #undef __FUNCT__
108d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation"
1094b9ad928SBarry Smith /*@
110aea2a34eSBarry Smith    PCMGSetInterpolation - Sets the function to be used to calculate the
111bf5b2e24SBarry Smith    interpolation from l-1 to the lth level
1124b9ad928SBarry Smith 
113ad4df100SBarry Smith    Logically Collective on PC and Mat
1144b9ad928SBarry Smith 
1154b9ad928SBarry Smith    Input Parameters:
1164b9ad928SBarry Smith +  pc  - the multigrid context
1174b9ad928SBarry Smith .  mat - the interpolation operator
118bf5b2e24SBarry Smith -  l   - the level (0 is coarsest) to supply [do not supply 0]
1194b9ad928SBarry Smith 
1204b9ad928SBarry Smith    Level: advanced
1214b9ad928SBarry Smith 
1224b9ad928SBarry Smith    Notes:
1234b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
1244b9ad928SBarry Smith     for the same level.
1254b9ad928SBarry Smith 
1264b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1274b9ad928SBarry Smith     out from the matrix size which one it is.
1284b9ad928SBarry Smith 
1294b9ad928SBarry Smith .keywords:  multigrid, set, interpolate, level
1304b9ad928SBarry Smith 
13197177400SBarry Smith .seealso: PCMGSetRestriction()
1324b9ad928SBarry Smith @*/
1337087cfbeSBarry Smith PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
1344b9ad928SBarry Smith {
135f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
136f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
137fccaa45eSBarry Smith   PetscErrorCode ierr;
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith   PetscFunctionBegin;
140c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
141e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
142e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
143c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1446bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr);
145f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
1464b9ad928SBarry Smith   PetscFunctionReturn(0);
1474b9ad928SBarry Smith }
1484b9ad928SBarry Smith 
1494b9ad928SBarry Smith #undef __FUNCT__
150d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction"
1514b9ad928SBarry Smith /*@
15297177400SBarry Smith    PCMGSetRestriction - Sets the function to be used to restrict vector
1534b9ad928SBarry Smith    from level l to l-1.
1544b9ad928SBarry Smith 
155ad4df100SBarry Smith    Logically Collective on PC and Mat
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith    Input Parameters:
1584b9ad928SBarry Smith +  pc - the multigrid context
1594b9ad928SBarry Smith .  mat - the restriction matrix
160bf5b2e24SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
1614b9ad928SBarry Smith 
1624b9ad928SBarry Smith    Level: advanced
1634b9ad928SBarry Smith 
1644b9ad928SBarry Smith    Notes:
1654b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
1664b9ad928SBarry Smith     for the same level.
1674b9ad928SBarry Smith 
1684b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1694b9ad928SBarry Smith     out from the matrix size which one it is.
1704b9ad928SBarry Smith 
171aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
172fccaa45eSBarry Smith     is used.
173fccaa45eSBarry Smith 
1744b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
1754b9ad928SBarry Smith 
176aea2a34eSBarry Smith .seealso: PCMGSetInterpolation()
1774b9ad928SBarry Smith @*/
1787087cfbeSBarry Smith PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
1794b9ad928SBarry Smith {
180fccaa45eSBarry Smith   PetscErrorCode ierr;
181f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
182f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith   PetscFunctionBegin;
185c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
186e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
187e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
188c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1896bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr);
190f3fbd535SBarry Smith   mglevels[l]->restrct  = mat;
1914b9ad928SBarry Smith   PetscFunctionReturn(0);
1924b9ad928SBarry Smith }
1934b9ad928SBarry Smith 
1944b9ad928SBarry Smith #undef __FUNCT__
19573250ac0SBarry Smith #define __FUNCT__ "PCMGSetRScale"
19673250ac0SBarry Smith /*@
19773250ac0SBarry Smith    PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.
19873250ac0SBarry Smith 
19973250ac0SBarry Smith    Logically Collective on PC and Mat
20073250ac0SBarry Smith 
20173250ac0SBarry Smith    Input Parameters:
20273250ac0SBarry Smith +  pc - the multigrid context
20373250ac0SBarry Smith .  rscale - the scaling
20473250ac0SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
20573250ac0SBarry Smith 
20673250ac0SBarry Smith    Level: advanced
20773250ac0SBarry Smith 
20873250ac0SBarry Smith    Notes:
20973250ac0SBarry 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.
21073250ac0SBarry Smith 
21173250ac0SBarry Smith .keywords: MG, set, multigrid, restriction, level
21273250ac0SBarry Smith 
21373250ac0SBarry Smith .seealso: PCMGSetInterpolation(), PCMGSetRestriction()
21473250ac0SBarry Smith @*/
21573250ac0SBarry Smith PetscErrorCode  PCMGSetRScale(PC pc,PetscInt l,Vec rscale)
21673250ac0SBarry Smith {
21773250ac0SBarry Smith   PetscErrorCode ierr;
21873250ac0SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
21973250ac0SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
22073250ac0SBarry Smith 
22173250ac0SBarry Smith   PetscFunctionBegin;
22273250ac0SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
22373250ac0SBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
22473250ac0SBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
22573250ac0SBarry Smith   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
22673250ac0SBarry Smith   ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr);
22773250ac0SBarry Smith   mglevels[l]->rscale  = rscale;
22873250ac0SBarry Smith   PetscFunctionReturn(0);
22973250ac0SBarry Smith }
23073250ac0SBarry Smith 
23173250ac0SBarry Smith #undef __FUNCT__
232d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother"
233f39d8e23SSatish Balay /*@
23497177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
23597177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
23697177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
2374b9ad928SBarry Smith    post-smoothing.
2384b9ad928SBarry Smith 
2394b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2404b9ad928SBarry Smith 
2414b9ad928SBarry Smith    Input Parameters:
2424b9ad928SBarry Smith +  pc - the multigrid context
2434b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
2444b9ad928SBarry Smith 
2454b9ad928SBarry Smith    Ouput Parameters:
2464b9ad928SBarry Smith .  ksp - the smoother
2474b9ad928SBarry Smith 
2484b9ad928SBarry Smith    Level: advanced
2494b9ad928SBarry Smith 
2504b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
2514b9ad928SBarry Smith 
25297177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2534b9ad928SBarry Smith @*/
2547087cfbeSBarry Smith PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
2554b9ad928SBarry Smith {
256f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
257f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2584b9ad928SBarry Smith 
2594b9ad928SBarry Smith   PetscFunctionBegin;
260c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
261f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
2624b9ad928SBarry Smith   PetscFunctionReturn(0);
2634b9ad928SBarry Smith }
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith #undef __FUNCT__
266d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp"
267f39d8e23SSatish Balay /*@
26897177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
2694b9ad928SBarry Smith    coarse grid correction (post-smoother).
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2724b9ad928SBarry Smith 
2734b9ad928SBarry Smith    Input Parameters:
2744b9ad928SBarry Smith +  pc - the multigrid context
2754b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2764b9ad928SBarry Smith 
2774b9ad928SBarry Smith    Ouput Parameters:
2784b9ad928SBarry Smith .  ksp - the smoother
2794b9ad928SBarry Smith 
2804b9ad928SBarry Smith    Level: advanced
2814b9ad928SBarry Smith 
2824b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
2834b9ad928SBarry Smith 
28497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2854b9ad928SBarry Smith @*/
2867087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
2874b9ad928SBarry Smith {
288f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
289f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
290dfbe8321SBarry Smith   PetscErrorCode ierr;
291f69a0ea3SMatthew Knepley   const char     *prefix;
2924b9ad928SBarry Smith   MPI_Comm       comm;
2934b9ad928SBarry Smith 
2944b9ad928SBarry Smith   PetscFunctionBegin;
295c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2964b9ad928SBarry Smith   /*
2974b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
2984b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
2994b9ad928SBarry Smith      if not we allocate it.
3004b9ad928SBarry Smith   */
301e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
302f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
303f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
304f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
305f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
306f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
307f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
308f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
309f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr);
3104b9ad928SBarry Smith   }
311f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
3124b9ad928SBarry Smith   PetscFunctionReturn(0);
3134b9ad928SBarry Smith }
3144b9ad928SBarry Smith 
3154b9ad928SBarry Smith #undef __FUNCT__
3169dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown"
317f39d8e23SSatish Balay /*@
31897177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
3194b9ad928SBarry Smith    coarse grid correction (pre-smoother).
3204b9ad928SBarry Smith 
3214b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
3224b9ad928SBarry Smith 
3234b9ad928SBarry Smith    Input Parameters:
3244b9ad928SBarry Smith +  pc - the multigrid context
3254b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
3264b9ad928SBarry Smith 
3274b9ad928SBarry Smith    Ouput Parameters:
3284b9ad928SBarry Smith .  ksp - the smoother
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith    Level: advanced
3314b9ad928SBarry Smith 
3324b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
3334b9ad928SBarry Smith 
33497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
3354b9ad928SBarry Smith @*/
3367087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
3374b9ad928SBarry Smith {
338dfbe8321SBarry Smith   PetscErrorCode ierr;
339f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
340f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3414b9ad928SBarry Smith 
3424b9ad928SBarry Smith   PetscFunctionBegin;
343c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3444b9ad928SBarry Smith   /* make sure smoother up and down are different */
345c5eb9154SBarry Smith   if (l) {
34697177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
347d8148a5aSMatthew Knepley   }
348f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
3494b9ad928SBarry Smith   PetscFunctionReturn(0);
3504b9ad928SBarry Smith }
3514b9ad928SBarry Smith 
3524b9ad928SBarry Smith #undef __FUNCT__
3539dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel"
3544b9ad928SBarry Smith /*@
35597177400SBarry Smith    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
3564b9ad928SBarry Smith 
357ad4df100SBarry Smith    Logically Collective on PC
3584b9ad928SBarry Smith 
3594b9ad928SBarry Smith    Input Parameters:
3604b9ad928SBarry Smith +  pc - the multigrid context
3614b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3624b9ad928SBarry Smith -  n  - the number of cycles
3634b9ad928SBarry Smith 
3644b9ad928SBarry Smith    Level: advanced
3654b9ad928SBarry Smith 
3664b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
3674b9ad928SBarry Smith 
36897177400SBarry Smith .seealso: PCMGSetCycles()
3694b9ad928SBarry Smith @*/
3707087cfbeSBarry Smith PetscErrorCode  PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
3714b9ad928SBarry Smith {
372f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
373f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3744b9ad928SBarry Smith 
3754b9ad928SBarry Smith   PetscFunctionBegin;
376c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
377e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
378c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,l,2);
379c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,c,3);
380f3fbd535SBarry Smith   mglevels[l]->cycles  = c;
3814b9ad928SBarry Smith   PetscFunctionReturn(0);
3824b9ad928SBarry Smith }
3834b9ad928SBarry Smith 
3844b9ad928SBarry Smith #undef __FUNCT__
385d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs"
3864b9ad928SBarry Smith /*@
38797177400SBarry Smith    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
388fccaa45eSBarry Smith    on a particular level.
3894b9ad928SBarry Smith 
390ad4df100SBarry Smith    Logically Collective on PC and Vec
3914b9ad928SBarry Smith 
3924b9ad928SBarry Smith    Input Parameters:
3934b9ad928SBarry Smith +  pc - the multigrid context
3944b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3954b9ad928SBarry Smith -  c  - the space
3964b9ad928SBarry Smith 
3974b9ad928SBarry Smith    Level: advanced
3984b9ad928SBarry Smith 
399fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
400fccaa45eSBarry Smith 
401fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
402fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
403fccaa45eSBarry Smith 
4044b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
4054b9ad928SBarry Smith 
40697177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
4074b9ad928SBarry Smith @*/
4087087cfbeSBarry Smith PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
4094b9ad928SBarry Smith {
410fccaa45eSBarry Smith   PetscErrorCode ierr;
411f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
412f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith   PetscFunctionBegin;
415c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
416e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
417e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
418c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
4196bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
420f3fbd535SBarry Smith   mglevels[l]->b  = c;
4214b9ad928SBarry Smith   PetscFunctionReturn(0);
4224b9ad928SBarry Smith }
4234b9ad928SBarry Smith 
4244b9ad928SBarry Smith #undef __FUNCT__
425d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX"
4264b9ad928SBarry Smith /*@
42797177400SBarry Smith    PCMGSetX - Sets the vector space to be used to store the solution on a
428fccaa45eSBarry Smith    particular level.
4294b9ad928SBarry Smith 
430ad4df100SBarry Smith    Logically Collective on PC and Vec
4314b9ad928SBarry Smith 
4324b9ad928SBarry Smith    Input Parameters:
4334b9ad928SBarry Smith +  pc - the multigrid context
4344b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4354b9ad928SBarry Smith -  c - the space
4364b9ad928SBarry Smith 
4374b9ad928SBarry Smith    Level: advanced
4384b9ad928SBarry Smith 
439fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
440fccaa45eSBarry Smith 
441fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
442fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
443fccaa45eSBarry Smith 
4444b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
4454b9ad928SBarry Smith 
44697177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
4474b9ad928SBarry Smith @*/
4487087cfbeSBarry Smith PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
4494b9ad928SBarry Smith {
450fccaa45eSBarry Smith   PetscErrorCode ierr;
451f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
452f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4534b9ad928SBarry Smith 
4544b9ad928SBarry Smith   PetscFunctionBegin;
455c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
456e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
457e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
458c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
4596bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
460f3fbd535SBarry Smith   mglevels[l]->x  = c;
4614b9ad928SBarry Smith   PetscFunctionReturn(0);
4624b9ad928SBarry Smith }
4634b9ad928SBarry Smith 
4644b9ad928SBarry Smith #undef __FUNCT__
465d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR"
4664b9ad928SBarry Smith /*@
46797177400SBarry Smith    PCMGSetR - Sets the vector space to be used to store the residual on a
468fccaa45eSBarry Smith    particular level.
4694b9ad928SBarry Smith 
470ad4df100SBarry Smith    Logically Collective on PC and Vec
4714b9ad928SBarry Smith 
4724b9ad928SBarry Smith    Input Parameters:
4734b9ad928SBarry Smith +  pc - the multigrid context
4744b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4754b9ad928SBarry Smith -  c - the space
4764b9ad928SBarry Smith 
4774b9ad928SBarry Smith    Level: advanced
4784b9ad928SBarry Smith 
479fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
480fccaa45eSBarry Smith 
481fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
482fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
483fccaa45eSBarry Smith 
4844b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
4854b9ad928SBarry Smith @*/
4867087cfbeSBarry Smith PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
4874b9ad928SBarry Smith {
488fccaa45eSBarry Smith   PetscErrorCode ierr;
489f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
490f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4914b9ad928SBarry Smith 
4924b9ad928SBarry Smith   PetscFunctionBegin;
493c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
494e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
495e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
496c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
4976bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
498f3fbd535SBarry Smith   mglevels[l]->r  = c;
4994b9ad928SBarry Smith   PetscFunctionReturn(0);
5004b9ad928SBarry Smith }
501