xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 298cc2083debba3835f8dd411d8115b662ff203b)
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
7797177400SBarry Smith .  residual - function used to form residual (usually PCMGDefaultResidual)
784b9ad928SBarry Smith -  mat      - matrix associated with residual
794b9ad928SBarry Smith 
804b9ad928SBarry Smith    Level: advanced
814b9ad928SBarry Smith 
824b9ad928SBarry Smith .keywords:  MG, set, multigrid, residual, level
834b9ad928SBarry Smith 
8497177400SBarry Smith .seealso: PCMGDefaultResidual()
854b9ad928SBarry Smith @*/
867087cfbeSBarry Smith PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
874b9ad928SBarry Smith {
88f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
89f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
90*298cc208SBarry Smith   PetscErrorCode ierr;
914b9ad928SBarry Smith 
924b9ad928SBarry Smith   PetscFunctionBegin;
93c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
94e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
95f3fbd535SBarry Smith   mglevels[l]->residual = residual;
96*298cc208SBarry Smith   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
97*298cc208SBarry Smith   ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr);
98f3fbd535SBarry Smith   mglevels[l]->A        = mat;
994b9ad928SBarry Smith   PetscFunctionReturn(0);
1004b9ad928SBarry Smith }
1014b9ad928SBarry Smith 
1024b9ad928SBarry Smith #undef __FUNCT__
103d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation"
1044b9ad928SBarry Smith /*@
105aea2a34eSBarry Smith    PCMGSetInterpolation - Sets the function to be used to calculate the
106bf5b2e24SBarry Smith    interpolation from l-1 to the lth level
1074b9ad928SBarry Smith 
108ad4df100SBarry Smith    Logically Collective on PC and Mat
1094b9ad928SBarry Smith 
1104b9ad928SBarry Smith    Input Parameters:
1114b9ad928SBarry Smith +  pc  - the multigrid context
1124b9ad928SBarry Smith .  mat - the interpolation operator
113bf5b2e24SBarry Smith -  l   - the level (0 is coarsest) to supply [do not supply 0]
1144b9ad928SBarry Smith 
1154b9ad928SBarry Smith    Level: advanced
1164b9ad928SBarry Smith 
1174b9ad928SBarry Smith    Notes:
1184b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
1194b9ad928SBarry Smith     for the same level.
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1224b9ad928SBarry Smith     out from the matrix size which one it is.
1234b9ad928SBarry Smith 
1244b9ad928SBarry Smith .keywords:  multigrid, set, interpolate, level
1254b9ad928SBarry Smith 
12697177400SBarry Smith .seealso: PCMGSetRestriction()
1274b9ad928SBarry Smith @*/
1287087cfbeSBarry Smith PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
1294b9ad928SBarry Smith {
130f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
131f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
132fccaa45eSBarry Smith   PetscErrorCode ierr;
1334b9ad928SBarry Smith 
1344b9ad928SBarry Smith   PetscFunctionBegin;
135c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
136e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
137e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
138c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1396bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr);
140f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
1414b9ad928SBarry Smith   PetscFunctionReturn(0);
1424b9ad928SBarry Smith }
1434b9ad928SBarry Smith 
1444b9ad928SBarry Smith #undef __FUNCT__
145d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction"
1464b9ad928SBarry Smith /*@
14797177400SBarry Smith    PCMGSetRestriction - Sets the function to be used to restrict vector
1484b9ad928SBarry Smith    from level l to l-1.
1494b9ad928SBarry Smith 
150ad4df100SBarry Smith    Logically Collective on PC and Mat
1514b9ad928SBarry Smith 
1524b9ad928SBarry Smith    Input Parameters:
1534b9ad928SBarry Smith +  pc - the multigrid context
1544b9ad928SBarry Smith .  mat - the restriction matrix
155bf5b2e24SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith    Level: advanced
1584b9ad928SBarry Smith 
1594b9ad928SBarry Smith    Notes:
1604b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
1614b9ad928SBarry Smith     for the same level.
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1644b9ad928SBarry Smith     out from the matrix size which one it is.
1654b9ad928SBarry Smith 
166aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
167fccaa45eSBarry Smith     is used.
168fccaa45eSBarry Smith 
1694b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
1704b9ad928SBarry Smith 
171aea2a34eSBarry Smith .seealso: PCMGSetInterpolation()
1724b9ad928SBarry Smith @*/
1737087cfbeSBarry Smith PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
1744b9ad928SBarry Smith {
175fccaa45eSBarry Smith   PetscErrorCode ierr;
176f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
177f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1784b9ad928SBarry Smith 
1794b9ad928SBarry Smith   PetscFunctionBegin;
180c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
181e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
182e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
183c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1846bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr);
185f3fbd535SBarry Smith   mglevels[l]->restrct  = mat;
1864b9ad928SBarry Smith   PetscFunctionReturn(0);
1874b9ad928SBarry Smith }
1884b9ad928SBarry Smith 
1894b9ad928SBarry Smith #undef __FUNCT__
190d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother"
191f39d8e23SSatish Balay /*@
19297177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
19397177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
19497177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
1954b9ad928SBarry Smith    post-smoothing.
1964b9ad928SBarry Smith 
1974b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
1984b9ad928SBarry Smith 
1994b9ad928SBarry Smith    Input Parameters:
2004b9ad928SBarry Smith +  pc - the multigrid context
2014b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
2024b9ad928SBarry Smith 
2034b9ad928SBarry Smith    Ouput Parameters:
2044b9ad928SBarry Smith .  ksp - the smoother
2054b9ad928SBarry Smith 
2064b9ad928SBarry Smith    Level: advanced
2074b9ad928SBarry Smith 
2084b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
2094b9ad928SBarry Smith 
21097177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2114b9ad928SBarry Smith @*/
2127087cfbeSBarry Smith PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
2134b9ad928SBarry Smith {
214f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
215f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2164b9ad928SBarry Smith 
2174b9ad928SBarry Smith   PetscFunctionBegin;
218c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
219f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
2204b9ad928SBarry Smith   PetscFunctionReturn(0);
2214b9ad928SBarry Smith }
2224b9ad928SBarry Smith 
2234b9ad928SBarry Smith #undef __FUNCT__
224d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp"
225f39d8e23SSatish Balay /*@
22697177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
2274b9ad928SBarry Smith    coarse grid correction (post-smoother).
2284b9ad928SBarry Smith 
2294b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2304b9ad928SBarry Smith 
2314b9ad928SBarry Smith    Input Parameters:
2324b9ad928SBarry Smith +  pc - the multigrid context
2334b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2344b9ad928SBarry Smith 
2354b9ad928SBarry Smith    Ouput Parameters:
2364b9ad928SBarry Smith .  ksp - the smoother
2374b9ad928SBarry Smith 
2384b9ad928SBarry Smith    Level: advanced
2394b9ad928SBarry Smith 
2404b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
2414b9ad928SBarry Smith 
24297177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2434b9ad928SBarry Smith @*/
2447087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
2454b9ad928SBarry Smith {
246f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
247f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
248dfbe8321SBarry Smith   PetscErrorCode ierr;
249f69a0ea3SMatthew Knepley   const char     *prefix;
2504b9ad928SBarry Smith   MPI_Comm       comm;
2514b9ad928SBarry Smith 
2524b9ad928SBarry Smith   PetscFunctionBegin;
253c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2544b9ad928SBarry Smith   /*
2554b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
2564b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
2574b9ad928SBarry Smith      if not we allocate it.
2584b9ad928SBarry Smith   */
259e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
260f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
261f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
262f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
263f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
264f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
265f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
266f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
267f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr);
2684b9ad928SBarry Smith   }
269f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
2704b9ad928SBarry Smith   PetscFunctionReturn(0);
2714b9ad928SBarry Smith }
2724b9ad928SBarry Smith 
2734b9ad928SBarry Smith #undef __FUNCT__
2749dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown"
275f39d8e23SSatish Balay /*@
27697177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
2774b9ad928SBarry Smith    coarse grid correction (pre-smoother).
2784b9ad928SBarry Smith 
2794b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2804b9ad928SBarry Smith 
2814b9ad928SBarry Smith    Input Parameters:
2824b9ad928SBarry Smith +  pc - the multigrid context
2834b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2844b9ad928SBarry Smith 
2854b9ad928SBarry Smith    Ouput Parameters:
2864b9ad928SBarry Smith .  ksp - the smoother
2874b9ad928SBarry Smith 
2884b9ad928SBarry Smith    Level: advanced
2894b9ad928SBarry Smith 
2904b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
2914b9ad928SBarry Smith 
29297177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
2934b9ad928SBarry Smith @*/
2947087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
2954b9ad928SBarry Smith {
296dfbe8321SBarry Smith   PetscErrorCode ierr;
297f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
298f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2994b9ad928SBarry Smith 
3004b9ad928SBarry Smith   PetscFunctionBegin;
301c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3024b9ad928SBarry Smith   /* make sure smoother up and down are different */
303c5eb9154SBarry Smith   if (l) {
30497177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
305d8148a5aSMatthew Knepley   }
306f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
3074b9ad928SBarry Smith   PetscFunctionReturn(0);
3084b9ad928SBarry Smith }
3094b9ad928SBarry Smith 
3104b9ad928SBarry Smith #undef __FUNCT__
3119dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel"
3124b9ad928SBarry Smith /*@
31397177400SBarry Smith    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
3144b9ad928SBarry Smith 
315ad4df100SBarry Smith    Logically Collective on PC
3164b9ad928SBarry Smith 
3174b9ad928SBarry Smith    Input Parameters:
3184b9ad928SBarry Smith +  pc - the multigrid context
3194b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3204b9ad928SBarry Smith -  n  - the number of cycles
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith    Level: advanced
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
3254b9ad928SBarry Smith 
32697177400SBarry Smith .seealso: PCMGSetCycles()
3274b9ad928SBarry Smith @*/
3287087cfbeSBarry Smith PetscErrorCode  PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
3294b9ad928SBarry Smith {
330f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
331f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3324b9ad928SBarry Smith 
3334b9ad928SBarry Smith   PetscFunctionBegin;
334c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
335e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
336c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,l,2);
337c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,c,3);
338f3fbd535SBarry Smith   mglevels[l]->cycles  = c;
3394b9ad928SBarry Smith   PetscFunctionReturn(0);
3404b9ad928SBarry Smith }
3414b9ad928SBarry Smith 
3424b9ad928SBarry Smith #undef __FUNCT__
343d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs"
3444b9ad928SBarry Smith /*@
34597177400SBarry Smith    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
346fccaa45eSBarry Smith    on a particular level.
3474b9ad928SBarry Smith 
348ad4df100SBarry Smith    Logically Collective on PC and Vec
3494b9ad928SBarry Smith 
3504b9ad928SBarry Smith    Input Parameters:
3514b9ad928SBarry Smith +  pc - the multigrid context
3524b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3534b9ad928SBarry Smith -  c  - the space
3544b9ad928SBarry Smith 
3554b9ad928SBarry Smith    Level: advanced
3564b9ad928SBarry Smith 
357fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
358fccaa45eSBarry Smith 
359fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
360fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
361fccaa45eSBarry Smith 
3624b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
3634b9ad928SBarry Smith 
36497177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
3654b9ad928SBarry Smith @*/
3667087cfbeSBarry Smith PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
3674b9ad928SBarry Smith {
368fccaa45eSBarry Smith   PetscErrorCode ierr;
369f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
370f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3714b9ad928SBarry Smith 
3724b9ad928SBarry Smith   PetscFunctionBegin;
373c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
374e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
375e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
376c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
3776bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
378f3fbd535SBarry Smith   mglevels[l]->b  = c;
3794b9ad928SBarry Smith   PetscFunctionReturn(0);
3804b9ad928SBarry Smith }
3814b9ad928SBarry Smith 
3824b9ad928SBarry Smith #undef __FUNCT__
383d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX"
3844b9ad928SBarry Smith /*@
38597177400SBarry Smith    PCMGSetX - Sets the vector space to be used to store the solution on a
386fccaa45eSBarry Smith    particular level.
3874b9ad928SBarry Smith 
388ad4df100SBarry Smith    Logically Collective on PC and Vec
3894b9ad928SBarry Smith 
3904b9ad928SBarry Smith    Input Parameters:
3914b9ad928SBarry Smith +  pc - the multigrid context
3924b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
3934b9ad928SBarry Smith -  c - the space
3944b9ad928SBarry Smith 
3954b9ad928SBarry Smith    Level: advanced
3964b9ad928SBarry Smith 
397fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
398fccaa45eSBarry Smith 
399fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
400fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
401fccaa45eSBarry Smith 
4024b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
4034b9ad928SBarry Smith 
40497177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
4054b9ad928SBarry Smith @*/
4067087cfbeSBarry Smith PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
4074b9ad928SBarry Smith {
408fccaa45eSBarry Smith   PetscErrorCode ierr;
409f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
410f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4114b9ad928SBarry Smith 
4124b9ad928SBarry Smith   PetscFunctionBegin;
413c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
414e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
415e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
416c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
4176bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
418f3fbd535SBarry Smith   mglevels[l]->x  = c;
4194b9ad928SBarry Smith   PetscFunctionReturn(0);
4204b9ad928SBarry Smith }
4214b9ad928SBarry Smith 
4224b9ad928SBarry Smith #undef __FUNCT__
423d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR"
4244b9ad928SBarry Smith /*@
42597177400SBarry Smith    PCMGSetR - Sets the vector space to be used to store the residual on a
426fccaa45eSBarry Smith    particular level.
4274b9ad928SBarry Smith 
428ad4df100SBarry Smith    Logically Collective on PC and Vec
4294b9ad928SBarry Smith 
4304b9ad928SBarry Smith    Input Parameters:
4314b9ad928SBarry Smith +  pc - the multigrid context
4324b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4334b9ad928SBarry Smith -  c - the space
4344b9ad928SBarry Smith 
4354b9ad928SBarry Smith    Level: advanced
4364b9ad928SBarry Smith 
437fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
438fccaa45eSBarry Smith 
439fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
440fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
441fccaa45eSBarry Smith 
4424b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
4434b9ad928SBarry Smith @*/
4447087cfbeSBarry Smith PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
4454b9ad928SBarry Smith {
446fccaa45eSBarry Smith   PetscErrorCode ierr;
447f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
448f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4494b9ad928SBarry Smith 
4504b9ad928SBarry Smith   PetscFunctionBegin;
451c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
452e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
453e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
454c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
4556bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
456f3fbd535SBarry Smith   mglevels[l]->r  = c;
4574b9ad928SBarry Smith   PetscFunctionReturn(0);
4584b9ad928SBarry Smith }
459