xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision f3f935be36dbbe192ee18859ece7a89379f99cd7)
14b9ad928SBarry Smith 
27c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"       /*I "petscksp.h" I*/
3*f3f935beSJed 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;
904b9ad928SBarry Smith 
914b9ad928SBarry Smith   PetscFunctionBegin;
92c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
93e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
94f3fbd535SBarry Smith   mglevels[l]->residual = residual;
95f3fbd535SBarry Smith   mglevels[l]->A        = mat;
964b9ad928SBarry Smith   PetscFunctionReturn(0);
974b9ad928SBarry Smith }
984b9ad928SBarry Smith 
994b9ad928SBarry Smith #undef __FUNCT__
100d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation"
1014b9ad928SBarry Smith /*@
102aea2a34eSBarry Smith    PCMGSetInterpolation - Sets the function to be used to calculate the
103bf5b2e24SBarry Smith    interpolation from l-1 to the lth level
1044b9ad928SBarry Smith 
105ad4df100SBarry Smith    Logically Collective on PC and Mat
1064b9ad928SBarry Smith 
1074b9ad928SBarry Smith    Input Parameters:
1084b9ad928SBarry Smith +  pc  - the multigrid context
1094b9ad928SBarry Smith .  mat - the interpolation operator
110bf5b2e24SBarry Smith -  l   - the level (0 is coarsest) to supply [do not supply 0]
1114b9ad928SBarry Smith 
1124b9ad928SBarry Smith    Level: advanced
1134b9ad928SBarry Smith 
1144b9ad928SBarry Smith    Notes:
1154b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
1164b9ad928SBarry Smith     for the same level.
1174b9ad928SBarry Smith 
1184b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1194b9ad928SBarry Smith     out from the matrix size which one it is.
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith .keywords:  multigrid, set, interpolate, level
1224b9ad928SBarry Smith 
12397177400SBarry Smith .seealso: PCMGSetRestriction()
1244b9ad928SBarry Smith @*/
1257087cfbeSBarry Smith PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
1264b9ad928SBarry Smith {
127f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
128f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
129fccaa45eSBarry Smith   PetscErrorCode ierr;
1304b9ad928SBarry Smith 
1314b9ad928SBarry Smith   PetscFunctionBegin;
132c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
133e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
134e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
135c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
136f3fbd535SBarry Smith   if (mglevels[l]->interpolate) {ierr = MatDestroy(mglevels[l]->interpolate);CHKERRQ(ierr);}
137f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
1384b9ad928SBarry Smith   PetscFunctionReturn(0);
1394b9ad928SBarry Smith }
1404b9ad928SBarry Smith 
1414b9ad928SBarry Smith #undef __FUNCT__
142d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction"
1434b9ad928SBarry Smith /*@
14497177400SBarry Smith    PCMGSetRestriction - Sets the function to be used to restrict vector
1454b9ad928SBarry Smith    from level l to l-1.
1464b9ad928SBarry Smith 
147ad4df100SBarry Smith    Logically Collective on PC and Mat
1484b9ad928SBarry Smith 
1494b9ad928SBarry Smith    Input Parameters:
1504b9ad928SBarry Smith +  pc - the multigrid context
1514b9ad928SBarry Smith .  mat - the restriction matrix
152bf5b2e24SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
1534b9ad928SBarry Smith 
1544b9ad928SBarry Smith    Level: advanced
1554b9ad928SBarry Smith 
1564b9ad928SBarry Smith    Notes:
1574b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
1584b9ad928SBarry Smith     for the same level.
1594b9ad928SBarry Smith 
1604b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1614b9ad928SBarry Smith     out from the matrix size which one it is.
1624b9ad928SBarry Smith 
163aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
164fccaa45eSBarry Smith     is used.
165fccaa45eSBarry Smith 
1664b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
1674b9ad928SBarry Smith 
168aea2a34eSBarry Smith .seealso: PCMGSetInterpolation()
1694b9ad928SBarry Smith @*/
1707087cfbeSBarry Smith PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
1714b9ad928SBarry Smith {
172fccaa45eSBarry Smith   PetscErrorCode ierr;
173f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
174f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1754b9ad928SBarry Smith 
1764b9ad928SBarry Smith   PetscFunctionBegin;
177c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
178e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
179e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
180c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
181f3fbd535SBarry Smith   if (mglevels[l]->restrct) {ierr = MatDestroy(mglevels[l]->restrct);CHKERRQ(ierr);}
182f3fbd535SBarry Smith   mglevels[l]->restrct  = mat;
1834b9ad928SBarry Smith   PetscFunctionReturn(0);
1844b9ad928SBarry Smith }
1854b9ad928SBarry Smith 
1864b9ad928SBarry Smith #undef __FUNCT__
187d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother"
188f39d8e23SSatish Balay /*@
18997177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
19097177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
19197177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
1924b9ad928SBarry Smith    post-smoothing.
1934b9ad928SBarry Smith 
1944b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
1954b9ad928SBarry Smith 
1964b9ad928SBarry Smith    Input Parameters:
1974b9ad928SBarry Smith +  pc - the multigrid context
1984b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
1994b9ad928SBarry Smith 
2004b9ad928SBarry Smith    Ouput Parameters:
2014b9ad928SBarry Smith .  ksp - the smoother
2024b9ad928SBarry Smith 
2034b9ad928SBarry Smith    Level: advanced
2044b9ad928SBarry Smith 
2054b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
2064b9ad928SBarry Smith 
20797177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2084b9ad928SBarry Smith @*/
2097087cfbeSBarry Smith PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
2104b9ad928SBarry Smith {
211f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
212f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2134b9ad928SBarry Smith 
2144b9ad928SBarry Smith   PetscFunctionBegin;
215c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
216f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
2174b9ad928SBarry Smith   PetscFunctionReturn(0);
2184b9ad928SBarry Smith }
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith #undef __FUNCT__
221d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp"
222f39d8e23SSatish Balay /*@
22397177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
2244b9ad928SBarry Smith    coarse grid correction (post-smoother).
2254b9ad928SBarry Smith 
2264b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2274b9ad928SBarry Smith 
2284b9ad928SBarry Smith    Input Parameters:
2294b9ad928SBarry Smith +  pc - the multigrid context
2304b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2314b9ad928SBarry Smith 
2324b9ad928SBarry Smith    Ouput Parameters:
2334b9ad928SBarry Smith .  ksp - the smoother
2344b9ad928SBarry Smith 
2354b9ad928SBarry Smith    Level: advanced
2364b9ad928SBarry Smith 
2374b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
2384b9ad928SBarry Smith 
23997177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2404b9ad928SBarry Smith @*/
2417087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
2424b9ad928SBarry Smith {
243f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
244f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
245dfbe8321SBarry Smith   PetscErrorCode ierr;
246f69a0ea3SMatthew Knepley   const char     *prefix;
2474b9ad928SBarry Smith   MPI_Comm       comm;
2484b9ad928SBarry Smith 
2494b9ad928SBarry Smith   PetscFunctionBegin;
250c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2514b9ad928SBarry Smith   /*
2524b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
2534b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
2544b9ad928SBarry Smith      if not we allocate it.
2554b9ad928SBarry Smith   */
256e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
257f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
258f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
259f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
260f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
261f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
262f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
263f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
264f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr);
2654b9ad928SBarry Smith   }
266f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
2674b9ad928SBarry Smith   PetscFunctionReturn(0);
2684b9ad928SBarry Smith }
2694b9ad928SBarry Smith 
2704b9ad928SBarry Smith #undef __FUNCT__
2719dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown"
272f39d8e23SSatish Balay /*@
27397177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
2744b9ad928SBarry Smith    coarse grid correction (pre-smoother).
2754b9ad928SBarry Smith 
2764b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2774b9ad928SBarry Smith 
2784b9ad928SBarry Smith    Input Parameters:
2794b9ad928SBarry Smith +  pc - the multigrid context
2804b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2814b9ad928SBarry Smith 
2824b9ad928SBarry Smith    Ouput Parameters:
2834b9ad928SBarry Smith .  ksp - the smoother
2844b9ad928SBarry Smith 
2854b9ad928SBarry Smith    Level: advanced
2864b9ad928SBarry Smith 
2874b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
2884b9ad928SBarry Smith 
28997177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
2904b9ad928SBarry Smith @*/
2917087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
2924b9ad928SBarry Smith {
293dfbe8321SBarry Smith   PetscErrorCode ierr;
294f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
295f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2964b9ad928SBarry Smith 
2974b9ad928SBarry Smith   PetscFunctionBegin;
298c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2994b9ad928SBarry Smith   /* make sure smoother up and down are different */
300c5eb9154SBarry Smith   if (l) {
30197177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
302d8148a5aSMatthew Knepley   }
303f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
3044b9ad928SBarry Smith   PetscFunctionReturn(0);
3054b9ad928SBarry Smith }
3064b9ad928SBarry Smith 
3074b9ad928SBarry Smith #undef __FUNCT__
3089dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel"
3094b9ad928SBarry Smith /*@
31097177400SBarry Smith    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
3114b9ad928SBarry Smith 
312ad4df100SBarry Smith    Logically Collective on PC
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith    Input Parameters:
3154b9ad928SBarry Smith +  pc - the multigrid context
3164b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3174b9ad928SBarry Smith -  n  - the number of cycles
3184b9ad928SBarry Smith 
3194b9ad928SBarry Smith    Level: advanced
3204b9ad928SBarry Smith 
3214b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
3224b9ad928SBarry Smith 
32397177400SBarry Smith .seealso: PCMGSetCycles()
3244b9ad928SBarry Smith @*/
3257087cfbeSBarry Smith PetscErrorCode  PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
3264b9ad928SBarry Smith {
327f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
328f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith   PetscFunctionBegin;
331c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
332e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
333c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,l,2);
334c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,c,3);
335f3fbd535SBarry Smith   mglevels[l]->cycles  = c;
3364b9ad928SBarry Smith   PetscFunctionReturn(0);
3374b9ad928SBarry Smith }
3384b9ad928SBarry Smith 
3394b9ad928SBarry Smith #undef __FUNCT__
340d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs"
3414b9ad928SBarry Smith /*@
34297177400SBarry Smith    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
343fccaa45eSBarry Smith    on a particular level.
3444b9ad928SBarry Smith 
345ad4df100SBarry Smith    Logically Collective on PC and Vec
3464b9ad928SBarry Smith 
3474b9ad928SBarry Smith    Input Parameters:
3484b9ad928SBarry Smith +  pc - the multigrid context
3494b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3504b9ad928SBarry Smith -  c  - the space
3514b9ad928SBarry Smith 
3524b9ad928SBarry Smith    Level: advanced
3534b9ad928SBarry Smith 
354fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
355fccaa45eSBarry Smith 
356fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
357fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
358fccaa45eSBarry Smith 
3594b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
3604b9ad928SBarry Smith 
36197177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
3624b9ad928SBarry Smith @*/
3637087cfbeSBarry Smith PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
3644b9ad928SBarry Smith {
365fccaa45eSBarry Smith   PetscErrorCode ierr;
366f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
367f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith   PetscFunctionBegin;
370c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
371e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
372e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
373c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
374f3fbd535SBarry Smith   if (mglevels[l]->b) {ierr = VecDestroy(mglevels[l]->b);CHKERRQ(ierr);}
375f3fbd535SBarry Smith   mglevels[l]->b  = c;
3764b9ad928SBarry Smith   PetscFunctionReturn(0);
3774b9ad928SBarry Smith }
3784b9ad928SBarry Smith 
3794b9ad928SBarry Smith #undef __FUNCT__
380d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX"
3814b9ad928SBarry Smith /*@
38297177400SBarry Smith    PCMGSetX - Sets the vector space to be used to store the solution on a
383fccaa45eSBarry Smith    particular level.
3844b9ad928SBarry Smith 
385ad4df100SBarry Smith    Logically Collective on PC and Vec
3864b9ad928SBarry Smith 
3874b9ad928SBarry Smith    Input Parameters:
3884b9ad928SBarry Smith +  pc - the multigrid context
3894b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
3904b9ad928SBarry Smith -  c - the space
3914b9ad928SBarry Smith 
3924b9ad928SBarry Smith    Level: advanced
3934b9ad928SBarry Smith 
394fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
395fccaa45eSBarry Smith 
396fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
397fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
398fccaa45eSBarry Smith 
3994b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
4004b9ad928SBarry Smith 
40197177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
4024b9ad928SBarry Smith @*/
4037087cfbeSBarry Smith PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
4044b9ad928SBarry Smith {
405fccaa45eSBarry Smith   PetscErrorCode ierr;
406f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
407f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4084b9ad928SBarry Smith 
4094b9ad928SBarry Smith   PetscFunctionBegin;
410c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
411e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
412e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
413c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
414f3fbd535SBarry Smith   if (mglevels[l]->x) {ierr = VecDestroy(mglevels[l]->x);CHKERRQ(ierr);}
415f3fbd535SBarry Smith   mglevels[l]->x  = c;
4164b9ad928SBarry Smith   PetscFunctionReturn(0);
4174b9ad928SBarry Smith }
4184b9ad928SBarry Smith 
4194b9ad928SBarry Smith #undef __FUNCT__
420d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR"
4214b9ad928SBarry Smith /*@
42297177400SBarry Smith    PCMGSetR - Sets the vector space to be used to store the residual on a
423fccaa45eSBarry Smith    particular level.
4244b9ad928SBarry Smith 
425ad4df100SBarry Smith    Logically Collective on PC and Vec
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith    Input Parameters:
4284b9ad928SBarry Smith +  pc - the multigrid context
4294b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4304b9ad928SBarry Smith -  c - the space
4314b9ad928SBarry Smith 
4324b9ad928SBarry Smith    Level: advanced
4334b9ad928SBarry Smith 
434fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
435fccaa45eSBarry Smith 
436fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
437fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
438fccaa45eSBarry Smith 
4394b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
4404b9ad928SBarry Smith @*/
4417087cfbeSBarry Smith PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
4424b9ad928SBarry Smith {
443fccaa45eSBarry Smith   PetscErrorCode ierr;
444f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
445f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4464b9ad928SBarry Smith 
4474b9ad928SBarry Smith   PetscFunctionBegin;
448c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
449e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
450e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
451c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
452f3fbd535SBarry Smith   if (mglevels[l]->r) {ierr = VecDestroy(mglevels[l]->r);CHKERRQ(ierr);}
453f3fbd535SBarry Smith   mglevels[l]->r  = c;
4544b9ad928SBarry Smith   PetscFunctionReturn(0);
4554b9ad928SBarry Smith }
456