xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 89cce641326674dee6eed3a03bd1210b0b453d6c)
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
77157726a2SBarry Smith .  residual - function used to form residual, if none is provided the previously provide one is used, if no
78157726a2SBarry 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");
96157726a2SBarry Smith   if (residual) {
97f3fbd535SBarry Smith     mglevels[l]->residual = residual;
98157726a2SBarry Smith   } if (!mglevels[l]->residual) {
99157726a2SBarry Smith     mglevels[l]->residual = PCMGDefaultResidual;
100157726a2SBarry 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 
282*89cce641SBarry Smith    Notes: calling this will result in a different pre and post smoother so you may need to
283*89cce641SBarry Smith          set options on the pre smoother also
284*89cce641SBarry Smith 
2854b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
2864b9ad928SBarry Smith 
28797177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2884b9ad928SBarry Smith @*/
2897087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
2904b9ad928SBarry Smith {
291f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
292f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
293dfbe8321SBarry Smith   PetscErrorCode ierr;
294f69a0ea3SMatthew Knepley   const char     *prefix;
2954b9ad928SBarry Smith   MPI_Comm       comm;
2964b9ad928SBarry Smith 
2974b9ad928SBarry Smith   PetscFunctionBegin;
298c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2994b9ad928SBarry Smith   /*
3004b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
3014b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
3024b9ad928SBarry Smith      if not we allocate it.
3034b9ad928SBarry Smith   */
304e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
305f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
306f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
307f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
308f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
309f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
310f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
311f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
312f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr);
3134b9ad928SBarry Smith   }
314f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
3154b9ad928SBarry Smith   PetscFunctionReturn(0);
3164b9ad928SBarry Smith }
3174b9ad928SBarry Smith 
3184b9ad928SBarry Smith #undef __FUNCT__
3199dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown"
320f39d8e23SSatish Balay /*@
32197177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
3224b9ad928SBarry Smith    coarse grid correction (pre-smoother).
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
3254b9ad928SBarry Smith 
3264b9ad928SBarry Smith    Input Parameters:
3274b9ad928SBarry Smith +  pc - the multigrid context
3284b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith    Ouput Parameters:
3314b9ad928SBarry Smith .  ksp - the smoother
3324b9ad928SBarry Smith 
3334b9ad928SBarry Smith    Level: advanced
3344b9ad928SBarry Smith 
335*89cce641SBarry Smith    Notes: calling this will result in a different pre and post smoother so you may need to
336*89cce641SBarry Smith          set options on the post smoother also
337*89cce641SBarry Smith 
3384b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
3394b9ad928SBarry Smith 
34097177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
3414b9ad928SBarry Smith @*/
3427087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
3434b9ad928SBarry Smith {
344dfbe8321SBarry Smith   PetscErrorCode ierr;
345f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
346f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3474b9ad928SBarry Smith 
3484b9ad928SBarry Smith   PetscFunctionBegin;
349c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3504b9ad928SBarry Smith   /* make sure smoother up and down are different */
351c5eb9154SBarry Smith   if (l) {
35297177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
353d8148a5aSMatthew Knepley   }
354f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
3554b9ad928SBarry Smith   PetscFunctionReturn(0);
3564b9ad928SBarry Smith }
3574b9ad928SBarry Smith 
3584b9ad928SBarry Smith #undef __FUNCT__
3599dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel"
3604b9ad928SBarry Smith /*@
36197177400SBarry Smith    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
3624b9ad928SBarry Smith 
363ad4df100SBarry Smith    Logically Collective on PC
3644b9ad928SBarry Smith 
3654b9ad928SBarry Smith    Input Parameters:
3664b9ad928SBarry Smith +  pc - the multigrid context
3674b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3684b9ad928SBarry Smith -  n  - the number of cycles
3694b9ad928SBarry Smith 
3704b9ad928SBarry Smith    Level: advanced
3714b9ad928SBarry Smith 
3724b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
3734b9ad928SBarry Smith 
37497177400SBarry Smith .seealso: PCMGSetCycles()
3754b9ad928SBarry Smith @*/
3767087cfbeSBarry Smith PetscErrorCode  PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
3774b9ad928SBarry Smith {
378f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
379f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3804b9ad928SBarry Smith 
3814b9ad928SBarry Smith   PetscFunctionBegin;
382c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
383e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
384c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,l,2);
385c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,c,3);
386f3fbd535SBarry Smith   mglevels[l]->cycles  = c;
3874b9ad928SBarry Smith   PetscFunctionReturn(0);
3884b9ad928SBarry Smith }
3894b9ad928SBarry Smith 
3904b9ad928SBarry Smith #undef __FUNCT__
391d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs"
3924b9ad928SBarry Smith /*@
39397177400SBarry Smith    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
394fccaa45eSBarry Smith    on a particular level.
3954b9ad928SBarry Smith 
396ad4df100SBarry Smith    Logically Collective on PC and Vec
3974b9ad928SBarry Smith 
3984b9ad928SBarry Smith    Input Parameters:
3994b9ad928SBarry Smith +  pc - the multigrid context
4004b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
4014b9ad928SBarry Smith -  c  - the space
4024b9ad928SBarry Smith 
4034b9ad928SBarry Smith    Level: advanced
4044b9ad928SBarry Smith 
405fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
406fccaa45eSBarry Smith 
407fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
408fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
409fccaa45eSBarry Smith 
4104b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
4114b9ad928SBarry Smith 
41297177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
4134b9ad928SBarry Smith @*/
4147087cfbeSBarry Smith PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
4154b9ad928SBarry Smith {
416fccaa45eSBarry Smith   PetscErrorCode ierr;
417f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
418f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4194b9ad928SBarry Smith 
4204b9ad928SBarry Smith   PetscFunctionBegin;
421c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
422e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
423e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
424c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
4256bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
426f3fbd535SBarry Smith   mglevels[l]->b  = c;
4274b9ad928SBarry Smith   PetscFunctionReturn(0);
4284b9ad928SBarry Smith }
4294b9ad928SBarry Smith 
4304b9ad928SBarry Smith #undef __FUNCT__
431d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX"
4324b9ad928SBarry Smith /*@
43397177400SBarry Smith    PCMGSetX - Sets the vector space to be used to store the solution on a
434fccaa45eSBarry Smith    particular level.
4354b9ad928SBarry Smith 
436ad4df100SBarry Smith    Logically Collective on PC and Vec
4374b9ad928SBarry Smith 
4384b9ad928SBarry Smith    Input Parameters:
4394b9ad928SBarry Smith +  pc - the multigrid context
4404b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4414b9ad928SBarry Smith -  c - the space
4424b9ad928SBarry Smith 
4434b9ad928SBarry Smith    Level: advanced
4444b9ad928SBarry Smith 
445fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
446fccaa45eSBarry Smith 
447fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
448fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
449fccaa45eSBarry Smith 
4504b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
4514b9ad928SBarry Smith 
45297177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
4534b9ad928SBarry Smith @*/
4547087cfbeSBarry Smith PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
4554b9ad928SBarry Smith {
456fccaa45eSBarry Smith   PetscErrorCode ierr;
457f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
458f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4594b9ad928SBarry Smith 
4604b9ad928SBarry Smith   PetscFunctionBegin;
461c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
462e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
463e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
464c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
4656bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
466f3fbd535SBarry Smith   mglevels[l]->x  = c;
4674b9ad928SBarry Smith   PetscFunctionReturn(0);
4684b9ad928SBarry Smith }
4694b9ad928SBarry Smith 
4704b9ad928SBarry Smith #undef __FUNCT__
471d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR"
4724b9ad928SBarry Smith /*@
47397177400SBarry Smith    PCMGSetR - Sets the vector space to be used to store the residual on a
474fccaa45eSBarry Smith    particular level.
4754b9ad928SBarry Smith 
476ad4df100SBarry Smith    Logically Collective on PC and Vec
4774b9ad928SBarry Smith 
4784b9ad928SBarry Smith    Input Parameters:
4794b9ad928SBarry Smith +  pc - the multigrid context
4804b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4814b9ad928SBarry Smith -  c - the space
4824b9ad928SBarry Smith 
4834b9ad928SBarry Smith    Level: advanced
4844b9ad928SBarry Smith 
485fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
486fccaa45eSBarry Smith 
487fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
488fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
489fccaa45eSBarry Smith 
4904b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
4914b9ad928SBarry Smith @*/
4927087cfbeSBarry Smith PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
4934b9ad928SBarry Smith {
494fccaa45eSBarry Smith   PetscErrorCode ierr;
495f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
496f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4974b9ad928SBarry Smith 
4984b9ad928SBarry Smith   PetscFunctionBegin;
499c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
500e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
501e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
502c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
5036bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
504f3fbd535SBarry Smith   mglevels[l]->r  = c;
5054b9ad928SBarry Smith   PetscFunctionReturn(0);
5064b9ad928SBarry Smith }
507