xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision f3fbd5355936bed48b7956a540fdcc016dd18626)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
24b9ad928SBarry Smith 
37c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"       /*I "petscksp.h" I*/
44b9ad928SBarry Smith                           /*I "petscmg.h"   I*/
54b9ad928SBarry Smith 
64b9ad928SBarry Smith #undef __FUNCT__
797177400SBarry Smith #define __FUNCT__ "PCMGDefaultResidual"
81f6cc5b2SSatish Balay /*@C
997177400SBarry Smith    PCMGDefaultResidual - Default routine to calculate the residual.
104b9ad928SBarry Smith 
114b9ad928SBarry Smith    Collective on Mat and Vec
124b9ad928SBarry Smith 
134b9ad928SBarry Smith    Input Parameters:
144b9ad928SBarry Smith +  mat - the matrix
154b9ad928SBarry Smith .  b   - the right-hand-side
164b9ad928SBarry Smith -  x   - the approximate solution
174b9ad928SBarry Smith 
184b9ad928SBarry Smith    Output Parameter:
194b9ad928SBarry Smith .  r - location to store the residual
204b9ad928SBarry Smith 
214b9ad928SBarry Smith    Level: advanced
224b9ad928SBarry Smith 
234b9ad928SBarry Smith .keywords: MG, default, multigrid, residual
244b9ad928SBarry Smith 
2597177400SBarry Smith .seealso: PCMGSetResidual()
264b9ad928SBarry Smith @*/
2797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGDefaultResidual(Mat mat,Vec b,Vec x,Vec r)
284b9ad928SBarry Smith {
29dfbe8321SBarry Smith   PetscErrorCode ierr;
304b9ad928SBarry Smith 
314b9ad928SBarry Smith   PetscFunctionBegin;
324b9ad928SBarry Smith   ierr = MatMult(mat,x,r);CHKERRQ(ierr);
33efb30889SBarry Smith   ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
344b9ad928SBarry Smith   PetscFunctionReturn(0);
354b9ad928SBarry Smith }
364b9ad928SBarry Smith 
374b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/
384b9ad928SBarry Smith 
394b9ad928SBarry Smith #undef __FUNCT__
40d6337fedSHong Zhang #define __FUNCT__ "PCMGGetCoarseSolve"
41f39d8e23SSatish Balay /*@
4297177400SBarry Smith    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
434b9ad928SBarry Smith 
444b9ad928SBarry Smith    Not Collective
454b9ad928SBarry Smith 
464b9ad928SBarry Smith    Input Parameter:
474b9ad928SBarry Smith .  pc - the multigrid context
484b9ad928SBarry Smith 
494b9ad928SBarry Smith    Output Parameter:
504b9ad928SBarry Smith .  ksp - the coarse grid solver context
514b9ad928SBarry Smith 
524b9ad928SBarry Smith    Level: advanced
534b9ad928SBarry Smith 
544b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid
554b9ad928SBarry Smith @*/
5697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetCoarseSolve(PC pc,KSP *ksp)
574b9ad928SBarry Smith {
58*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
59*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
604b9ad928SBarry Smith 
614b9ad928SBarry Smith   PetscFunctionBegin;
62*f3fbd535SBarry 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 
724b9ad928SBarry Smith    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 @*/
8697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
874b9ad928SBarry Smith {
88*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
89*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
904b9ad928SBarry Smith 
914b9ad928SBarry Smith   PetscFunctionBegin;
92*f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
934b9ad928SBarry Smith 
94*f3fbd535SBarry Smith   mglevels[l]->residual = residual;
95*f3fbd535SBarry 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 
1054b9ad928SBarry Smith    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 @*/
125aea2a34eSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
1264b9ad928SBarry Smith {
127*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
128*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
129fccaa45eSBarry Smith   PetscErrorCode ierr;
1304b9ad928SBarry Smith 
1314b9ad928SBarry Smith   PetscFunctionBegin;
132*f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1336d9142f8SRichard Katz   if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
134c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
135*f3fbd535SBarry Smith   if (mglevels[l]->interpolate) {ierr = MatDestroy(mglevels[l]->interpolate);CHKERRQ(ierr);}
136*f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
1374b9ad928SBarry Smith   PetscFunctionReturn(0);
1384b9ad928SBarry Smith }
1394b9ad928SBarry Smith 
1404b9ad928SBarry Smith #undef __FUNCT__
141d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction"
1424b9ad928SBarry Smith /*@
14397177400SBarry Smith    PCMGSetRestriction - Sets the function to be used to restrict vector
1444b9ad928SBarry Smith    from level l to l-1.
1454b9ad928SBarry Smith 
1464b9ad928SBarry Smith    Collective on PC and Mat
1474b9ad928SBarry Smith 
1484b9ad928SBarry Smith    Input Parameters:
1494b9ad928SBarry Smith +  pc - the multigrid context
1504b9ad928SBarry Smith .  mat - the restriction matrix
151bf5b2e24SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
1524b9ad928SBarry Smith 
1534b9ad928SBarry Smith    Level: advanced
1544b9ad928SBarry Smith 
1554b9ad928SBarry Smith    Notes:
1564b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
1574b9ad928SBarry Smith     for the same level.
1584b9ad928SBarry Smith 
1594b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1604b9ad928SBarry Smith     out from the matrix size which one it is.
1614b9ad928SBarry Smith 
162aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
163fccaa45eSBarry Smith     is used.
164fccaa45eSBarry Smith 
1654b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
1664b9ad928SBarry Smith 
167aea2a34eSBarry Smith .seealso: PCMGSetInterpolation()
1684b9ad928SBarry Smith @*/
16997177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
1704b9ad928SBarry Smith {
171fccaa45eSBarry Smith   PetscErrorCode ierr;
172*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
173*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1744b9ad928SBarry Smith 
1754b9ad928SBarry Smith   PetscFunctionBegin;
176*f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1776d9142f8SRichard Katz   if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
178c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
179*f3fbd535SBarry Smith   if (mglevels[l]->restrct) {ierr = MatDestroy(mglevels[l]->restrct);CHKERRQ(ierr);}
180*f3fbd535SBarry Smith   mglevels[l]->restrct  = mat;
1814b9ad928SBarry Smith   PetscFunctionReturn(0);
1824b9ad928SBarry Smith }
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith #undef __FUNCT__
185d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother"
186f39d8e23SSatish Balay /*@
18797177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
18897177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
18997177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
1904b9ad928SBarry Smith    post-smoothing.
1914b9ad928SBarry Smith 
1924b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
1934b9ad928SBarry Smith 
1944b9ad928SBarry Smith    Input Parameters:
1954b9ad928SBarry Smith +  pc - the multigrid context
1964b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
1974b9ad928SBarry Smith 
1984b9ad928SBarry Smith    Ouput Parameters:
1994b9ad928SBarry Smith .  ksp - the smoother
2004b9ad928SBarry Smith 
2014b9ad928SBarry Smith    Level: advanced
2024b9ad928SBarry Smith 
2034b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
2044b9ad928SBarry Smith 
20597177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2064b9ad928SBarry Smith @*/
20797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
2084b9ad928SBarry Smith {
209*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
210*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2114b9ad928SBarry Smith 
2124b9ad928SBarry Smith   PetscFunctionBegin;
213*f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
2144b9ad928SBarry Smith   PetscFunctionReturn(0);
2154b9ad928SBarry Smith }
2164b9ad928SBarry Smith 
2174b9ad928SBarry Smith #undef __FUNCT__
218d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp"
219f39d8e23SSatish Balay /*@
22097177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
2214b9ad928SBarry Smith    coarse grid correction (post-smoother).
2224b9ad928SBarry Smith 
2234b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2244b9ad928SBarry Smith 
2254b9ad928SBarry Smith    Input Parameters:
2264b9ad928SBarry Smith +  pc - the multigrid context
2274b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2284b9ad928SBarry Smith 
2294b9ad928SBarry Smith    Ouput Parameters:
2304b9ad928SBarry Smith .  ksp - the smoother
2314b9ad928SBarry Smith 
2324b9ad928SBarry Smith    Level: advanced
2334b9ad928SBarry Smith 
2344b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
2354b9ad928SBarry Smith 
23697177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2374b9ad928SBarry Smith @*/
23897177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
2394b9ad928SBarry Smith {
240*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
241*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
242dfbe8321SBarry Smith   PetscErrorCode ierr;
243f69a0ea3SMatthew Knepley   const char     *prefix;
2444b9ad928SBarry Smith   MPI_Comm       comm;
2454b9ad928SBarry Smith 
2464b9ad928SBarry Smith   PetscFunctionBegin;
2474b9ad928SBarry Smith   /*
2484b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
2494b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
2504b9ad928SBarry Smith      if not we allocate it.
2514b9ad928SBarry Smith   */
252b03c7568SBarry Smith   if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
253*f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
254*f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
255*f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
256*f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
257*f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
258*f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
259*f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
260*f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr);
2614b9ad928SBarry Smith   }
262*f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
2634b9ad928SBarry Smith   PetscFunctionReturn(0);
2644b9ad928SBarry Smith }
2654b9ad928SBarry Smith 
2664b9ad928SBarry Smith #undef __FUNCT__
2679dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown"
268f39d8e23SSatish Balay /*@
26997177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
2704b9ad928SBarry Smith    coarse grid correction (pre-smoother).
2714b9ad928SBarry Smith 
2724b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2734b9ad928SBarry Smith 
2744b9ad928SBarry Smith    Input Parameters:
2754b9ad928SBarry Smith +  pc - the multigrid context
2764b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2774b9ad928SBarry Smith 
2784b9ad928SBarry Smith    Ouput Parameters:
2794b9ad928SBarry Smith .  ksp - the smoother
2804b9ad928SBarry Smith 
2814b9ad928SBarry Smith    Level: advanced
2824b9ad928SBarry Smith 
2834b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
2844b9ad928SBarry Smith 
28597177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
2864b9ad928SBarry Smith @*/
28797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
2884b9ad928SBarry Smith {
289dfbe8321SBarry Smith   PetscErrorCode ierr;
290*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
291*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2924b9ad928SBarry Smith 
2934b9ad928SBarry Smith   PetscFunctionBegin;
2944b9ad928SBarry Smith   /* make sure smoother up and down are different */
295d8148a5aSMatthew Knepley   if (l != 0) {
29697177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
297d8148a5aSMatthew Knepley   }
298*f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
2994b9ad928SBarry Smith   PetscFunctionReturn(0);
3004b9ad928SBarry Smith }
3014b9ad928SBarry Smith 
3024b9ad928SBarry Smith #undef __FUNCT__
3039dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel"
3044b9ad928SBarry Smith /*@
30597177400SBarry Smith    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
3064b9ad928SBarry Smith 
3074b9ad928SBarry Smith    Collective on PC
3084b9ad928SBarry Smith 
3094b9ad928SBarry Smith    Input Parameters:
3104b9ad928SBarry Smith +  pc - the multigrid context
3114b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3124b9ad928SBarry Smith -  n  - the number of cycles
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith    Level: advanced
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
3174b9ad928SBarry Smith 
31897177400SBarry Smith .seealso: PCMGSetCycles()
3194b9ad928SBarry Smith @*/
32097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
3214b9ad928SBarry Smith {
322*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
323*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3244b9ad928SBarry Smith 
3254b9ad928SBarry Smith   PetscFunctionBegin;
326*f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
327*f3fbd535SBarry Smith   mglevels[l]->cycles  = c;
3284b9ad928SBarry Smith   PetscFunctionReturn(0);
3294b9ad928SBarry Smith }
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith #undef __FUNCT__
332d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs"
3334b9ad928SBarry Smith /*@
33497177400SBarry Smith    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
335fccaa45eSBarry Smith    on a particular level.
3364b9ad928SBarry Smith 
3374b9ad928SBarry Smith    Collective on PC and Vec
3384b9ad928SBarry Smith 
3394b9ad928SBarry Smith    Input Parameters:
3404b9ad928SBarry Smith +  pc - the multigrid context
3414b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3424b9ad928SBarry Smith -  c  - the space
3434b9ad928SBarry Smith 
3444b9ad928SBarry Smith    Level: advanced
3454b9ad928SBarry Smith 
346fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
347fccaa45eSBarry Smith 
348fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
349fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
350fccaa45eSBarry Smith 
3514b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
3524b9ad928SBarry Smith 
35397177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
3544b9ad928SBarry Smith @*/
35597177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRhs(PC pc,PetscInt l,Vec c)
3564b9ad928SBarry Smith {
357fccaa45eSBarry Smith   PetscErrorCode ierr;
358*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
359*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3604b9ad928SBarry Smith 
3614b9ad928SBarry Smith   PetscFunctionBegin;
362*f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
363*f3fbd535SBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
364c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
365*f3fbd535SBarry Smith   if (mglevels[l]->b) {ierr = VecDestroy(mglevels[l]->b);CHKERRQ(ierr);}
366*f3fbd535SBarry Smith   mglevels[l]->b  = c;
3674b9ad928SBarry Smith   PetscFunctionReturn(0);
3684b9ad928SBarry Smith }
3694b9ad928SBarry Smith 
3704b9ad928SBarry Smith #undef __FUNCT__
371d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX"
3724b9ad928SBarry Smith /*@
37397177400SBarry Smith    PCMGSetX - Sets the vector space to be used to store the solution on a
374fccaa45eSBarry Smith    particular level.
3754b9ad928SBarry Smith 
3764b9ad928SBarry Smith    Collective on PC and Vec
3774b9ad928SBarry Smith 
3784b9ad928SBarry Smith    Input Parameters:
3794b9ad928SBarry Smith +  pc - the multigrid context
3804b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
3814b9ad928SBarry Smith -  c - the space
3824b9ad928SBarry Smith 
3834b9ad928SBarry Smith    Level: advanced
3844b9ad928SBarry Smith 
385fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
386fccaa45eSBarry Smith 
387fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
388fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
389fccaa45eSBarry Smith 
3904b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
3914b9ad928SBarry Smith 
39297177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
3934b9ad928SBarry Smith @*/
39497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetX(PC pc,PetscInt l,Vec c)
3954b9ad928SBarry Smith {
396fccaa45eSBarry Smith   PetscErrorCode ierr;
397*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
398*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3994b9ad928SBarry Smith 
4004b9ad928SBarry Smith   PetscFunctionBegin;
401*f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
402*f3fbd535SBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
403c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
404*f3fbd535SBarry Smith   if (mglevels[l]->x) {ierr = VecDestroy(mglevels[l]->x);CHKERRQ(ierr);}
405*f3fbd535SBarry Smith   mglevels[l]->x  = c;
4064b9ad928SBarry Smith   PetscFunctionReturn(0);
4074b9ad928SBarry Smith }
4084b9ad928SBarry Smith 
4094b9ad928SBarry Smith #undef __FUNCT__
410d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR"
4114b9ad928SBarry Smith /*@
41297177400SBarry Smith    PCMGSetR - Sets the vector space to be used to store the residual on a
413fccaa45eSBarry Smith    particular level.
4144b9ad928SBarry Smith 
4154b9ad928SBarry Smith    Collective on PC and Vec
4164b9ad928SBarry Smith 
4174b9ad928SBarry Smith    Input Parameters:
4184b9ad928SBarry Smith +  pc - the multigrid context
4194b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4204b9ad928SBarry Smith -  c - the space
4214b9ad928SBarry Smith 
4224b9ad928SBarry Smith    Level: advanced
4234b9ad928SBarry Smith 
424fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
425fccaa45eSBarry Smith 
426fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
427fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
428fccaa45eSBarry Smith 
4294b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
4304b9ad928SBarry Smith @*/
43197177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetR(PC pc,PetscInt l,Vec c)
4324b9ad928SBarry Smith {
433fccaa45eSBarry Smith   PetscErrorCode ierr;
434*f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
435*f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4364b9ad928SBarry Smith 
4374b9ad928SBarry Smith   PetscFunctionBegin;
438*f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
439630ba2eeSBarry Smith   if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
440c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
441*f3fbd535SBarry Smith   if (mglevels[l]->r) {ierr = VecDestroy(mglevels[l]->r);CHKERRQ(ierr);}
442*f3fbd535SBarry Smith   mglevels[l]->r  = c;
4434b9ad928SBarry Smith   PetscFunctionReturn(0);
4444b9ad928SBarry Smith }
445