xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision ad4df1005a67cf9575a5c4294f11fd96229b2bb4)
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 {
58f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
59f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
604b9ad928SBarry Smith 
614b9ad928SBarry Smith   PetscFunctionBegin;
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 
72*ad4df100SBarry 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 @*/
8697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT 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;
92e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
93f3fbd535SBarry Smith   mglevels[l]->residual = residual;
94f3fbd535SBarry Smith   mglevels[l]->A        = mat;
954b9ad928SBarry Smith   PetscFunctionReturn(0);
964b9ad928SBarry Smith }
974b9ad928SBarry Smith 
984b9ad928SBarry Smith #undef __FUNCT__
99d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation"
1004b9ad928SBarry Smith /*@
101aea2a34eSBarry Smith    PCMGSetInterpolation - Sets the function to be used to calculate the
102bf5b2e24SBarry Smith    interpolation from l-1 to the lth level
1034b9ad928SBarry Smith 
104*ad4df100SBarry Smith    Logically Collective on PC and Mat
1054b9ad928SBarry Smith 
1064b9ad928SBarry Smith    Input Parameters:
1074b9ad928SBarry Smith +  pc  - the multigrid context
1084b9ad928SBarry Smith .  mat - the interpolation operator
109bf5b2e24SBarry Smith -  l   - the level (0 is coarsest) to supply [do not supply 0]
1104b9ad928SBarry Smith 
1114b9ad928SBarry Smith    Level: advanced
1124b9ad928SBarry Smith 
1134b9ad928SBarry Smith    Notes:
1144b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
1154b9ad928SBarry Smith     for the same level.
1164b9ad928SBarry Smith 
1174b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1184b9ad928SBarry Smith     out from the matrix size which one it is.
1194b9ad928SBarry Smith 
1204b9ad928SBarry Smith .keywords:  multigrid, set, interpolate, level
1214b9ad928SBarry Smith 
12297177400SBarry Smith .seealso: PCMGSetRestriction()
1234b9ad928SBarry Smith @*/
124aea2a34eSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
1254b9ad928SBarry Smith {
126f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
127f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
128fccaa45eSBarry Smith   PetscErrorCode ierr;
1294b9ad928SBarry Smith 
1304b9ad928SBarry Smith   PetscFunctionBegin;
131e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
132e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
133c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
134f3fbd535SBarry Smith   if (mglevels[l]->interpolate) {ierr = MatDestroy(mglevels[l]->interpolate);CHKERRQ(ierr);}
135f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
1364b9ad928SBarry Smith   PetscFunctionReturn(0);
1374b9ad928SBarry Smith }
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith #undef __FUNCT__
140d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction"
1414b9ad928SBarry Smith /*@
14297177400SBarry Smith    PCMGSetRestriction - Sets the function to be used to restrict vector
1434b9ad928SBarry Smith    from level l to l-1.
1444b9ad928SBarry Smith 
145*ad4df100SBarry Smith    Logically Collective on PC and Mat
1464b9ad928SBarry Smith 
1474b9ad928SBarry Smith    Input Parameters:
1484b9ad928SBarry Smith +  pc - the multigrid context
1494b9ad928SBarry Smith .  mat - the restriction matrix
150bf5b2e24SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
1514b9ad928SBarry Smith 
1524b9ad928SBarry Smith    Level: advanced
1534b9ad928SBarry Smith 
1544b9ad928SBarry Smith    Notes:
1554b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
1564b9ad928SBarry Smith     for the same level.
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1594b9ad928SBarry Smith     out from the matrix size which one it is.
1604b9ad928SBarry Smith 
161aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
162fccaa45eSBarry Smith     is used.
163fccaa45eSBarry Smith 
1644b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
1654b9ad928SBarry Smith 
166aea2a34eSBarry Smith .seealso: PCMGSetInterpolation()
1674b9ad928SBarry Smith @*/
16897177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
1694b9ad928SBarry Smith {
170fccaa45eSBarry Smith   PetscErrorCode ierr;
171f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
172f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1734b9ad928SBarry Smith 
1744b9ad928SBarry Smith   PetscFunctionBegin;
175e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
176e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
177c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
178f3fbd535SBarry Smith   if (mglevels[l]->restrct) {ierr = MatDestroy(mglevels[l]->restrct);CHKERRQ(ierr);}
179f3fbd535SBarry Smith   mglevels[l]->restrct  = mat;
1804b9ad928SBarry Smith   PetscFunctionReturn(0);
1814b9ad928SBarry Smith }
1824b9ad928SBarry Smith 
1834b9ad928SBarry Smith #undef __FUNCT__
184d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother"
185f39d8e23SSatish Balay /*@
18697177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
18797177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
18897177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
1894b9ad928SBarry Smith    post-smoothing.
1904b9ad928SBarry Smith 
1914b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
1924b9ad928SBarry Smith 
1934b9ad928SBarry Smith    Input Parameters:
1944b9ad928SBarry Smith +  pc - the multigrid context
1954b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
1964b9ad928SBarry Smith 
1974b9ad928SBarry Smith    Ouput Parameters:
1984b9ad928SBarry Smith .  ksp - the smoother
1994b9ad928SBarry Smith 
2004b9ad928SBarry Smith    Level: advanced
2014b9ad928SBarry Smith 
2024b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
2034b9ad928SBarry Smith 
20497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2054b9ad928SBarry Smith @*/
20697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
2074b9ad928SBarry Smith {
208f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
209f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2104b9ad928SBarry Smith 
2114b9ad928SBarry Smith   PetscFunctionBegin;
212f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
2134b9ad928SBarry Smith   PetscFunctionReturn(0);
2144b9ad928SBarry Smith }
2154b9ad928SBarry Smith 
2164b9ad928SBarry Smith #undef __FUNCT__
217d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp"
218f39d8e23SSatish Balay /*@
21997177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
2204b9ad928SBarry Smith    coarse grid correction (post-smoother).
2214b9ad928SBarry Smith 
2224b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2234b9ad928SBarry Smith 
2244b9ad928SBarry Smith    Input Parameters:
2254b9ad928SBarry Smith +  pc - the multigrid context
2264b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2274b9ad928SBarry Smith 
2284b9ad928SBarry Smith    Ouput Parameters:
2294b9ad928SBarry Smith .  ksp - the smoother
2304b9ad928SBarry Smith 
2314b9ad928SBarry Smith    Level: advanced
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
2344b9ad928SBarry Smith 
23597177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2364b9ad928SBarry Smith @*/
23797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
2384b9ad928SBarry Smith {
239f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
240f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
241dfbe8321SBarry Smith   PetscErrorCode ierr;
242f69a0ea3SMatthew Knepley   const char     *prefix;
2434b9ad928SBarry Smith   MPI_Comm       comm;
2444b9ad928SBarry Smith 
2454b9ad928SBarry Smith   PetscFunctionBegin;
2464b9ad928SBarry Smith   /*
2474b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
2484b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
2494b9ad928SBarry Smith      if not we allocate it.
2504b9ad928SBarry Smith   */
251e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
252f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
253f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
254f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
255f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
256f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
257f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
258f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
259f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr);
2604b9ad928SBarry Smith   }
261f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
2624b9ad928SBarry Smith   PetscFunctionReturn(0);
2634b9ad928SBarry Smith }
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith #undef __FUNCT__
2669dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown"
267f39d8e23SSatish Balay /*@
26897177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
2694b9ad928SBarry Smith    coarse grid correction (pre-smoother).
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2724b9ad928SBarry Smith 
2734b9ad928SBarry Smith    Input Parameters:
2744b9ad928SBarry Smith +  pc - the multigrid context
2754b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2764b9ad928SBarry Smith 
2774b9ad928SBarry Smith    Ouput Parameters:
2784b9ad928SBarry Smith .  ksp - the smoother
2794b9ad928SBarry Smith 
2804b9ad928SBarry Smith    Level: advanced
2814b9ad928SBarry Smith 
2824b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
2834b9ad928SBarry Smith 
28497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
2854b9ad928SBarry Smith @*/
28697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
2874b9ad928SBarry Smith {
288dfbe8321SBarry Smith   PetscErrorCode ierr;
289f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
290f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2914b9ad928SBarry Smith 
2924b9ad928SBarry Smith   PetscFunctionBegin;
2934b9ad928SBarry Smith   /* make sure smoother up and down are different */
294d8148a5aSMatthew Knepley   if (l != 0) {
29597177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
296d8148a5aSMatthew Knepley   }
297f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
2984b9ad928SBarry Smith   PetscFunctionReturn(0);
2994b9ad928SBarry Smith }
3004b9ad928SBarry Smith 
3014b9ad928SBarry Smith #undef __FUNCT__
3029dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel"
3034b9ad928SBarry Smith /*@
30497177400SBarry Smith    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
3054b9ad928SBarry Smith 
306*ad4df100SBarry Smith    Logically Collective on PC
3074b9ad928SBarry Smith 
3084b9ad928SBarry Smith    Input Parameters:
3094b9ad928SBarry Smith +  pc - the multigrid context
3104b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3114b9ad928SBarry Smith -  n  - the number of cycles
3124b9ad928SBarry Smith 
3134b9ad928SBarry Smith    Level: advanced
3144b9ad928SBarry Smith 
3154b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
3164b9ad928SBarry Smith 
31797177400SBarry Smith .seealso: PCMGSetCycles()
3184b9ad928SBarry Smith @*/
31997177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
3204b9ad928SBarry Smith {
321f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
322f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith   PetscFunctionBegin;
325e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
326f3fbd535SBarry Smith   mglevels[l]->cycles  = c;
3274b9ad928SBarry Smith   PetscFunctionReturn(0);
3284b9ad928SBarry Smith }
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith #undef __FUNCT__
331d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs"
3324b9ad928SBarry Smith /*@
33397177400SBarry Smith    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
334fccaa45eSBarry Smith    on a particular level.
3354b9ad928SBarry Smith 
336*ad4df100SBarry Smith    Logically Collective on PC and Vec
3374b9ad928SBarry Smith 
3384b9ad928SBarry Smith    Input Parameters:
3394b9ad928SBarry Smith +  pc - the multigrid context
3404b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3414b9ad928SBarry Smith -  c  - the space
3424b9ad928SBarry Smith 
3434b9ad928SBarry Smith    Level: advanced
3444b9ad928SBarry Smith 
345fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
346fccaa45eSBarry Smith 
347fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
348fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
349fccaa45eSBarry Smith 
3504b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
3514b9ad928SBarry Smith 
35297177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
3534b9ad928SBarry Smith @*/
35497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRhs(PC pc,PetscInt l,Vec c)
3554b9ad928SBarry Smith {
356fccaa45eSBarry Smith   PetscErrorCode ierr;
357f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
358f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3594b9ad928SBarry Smith 
3604b9ad928SBarry Smith   PetscFunctionBegin;
361e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
362e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
363c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
364f3fbd535SBarry Smith   if (mglevels[l]->b) {ierr = VecDestroy(mglevels[l]->b);CHKERRQ(ierr);}
365f3fbd535SBarry Smith   mglevels[l]->b  = c;
3664b9ad928SBarry Smith   PetscFunctionReturn(0);
3674b9ad928SBarry Smith }
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith #undef __FUNCT__
370d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX"
3714b9ad928SBarry Smith /*@
37297177400SBarry Smith    PCMGSetX - Sets the vector space to be used to store the solution on a
373fccaa45eSBarry Smith    particular level.
3744b9ad928SBarry Smith 
375*ad4df100SBarry Smith    Logically Collective on PC and Vec
3764b9ad928SBarry Smith 
3774b9ad928SBarry Smith    Input Parameters:
3784b9ad928SBarry Smith +  pc - the multigrid context
3794b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
3804b9ad928SBarry Smith -  c - the space
3814b9ad928SBarry Smith 
3824b9ad928SBarry Smith    Level: advanced
3834b9ad928SBarry Smith 
384fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
385fccaa45eSBarry Smith 
386fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
387fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
388fccaa45eSBarry Smith 
3894b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
3904b9ad928SBarry Smith 
39197177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
3924b9ad928SBarry Smith @*/
39397177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetX(PC pc,PetscInt l,Vec c)
3944b9ad928SBarry Smith {
395fccaa45eSBarry Smith   PetscErrorCode ierr;
396f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
397f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3984b9ad928SBarry Smith 
3994b9ad928SBarry Smith   PetscFunctionBegin;
400e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
401e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
402c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
403f3fbd535SBarry Smith   if (mglevels[l]->x) {ierr = VecDestroy(mglevels[l]->x);CHKERRQ(ierr);}
404f3fbd535SBarry Smith   mglevels[l]->x  = c;
4054b9ad928SBarry Smith   PetscFunctionReturn(0);
4064b9ad928SBarry Smith }
4074b9ad928SBarry Smith 
4084b9ad928SBarry Smith #undef __FUNCT__
409d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR"
4104b9ad928SBarry Smith /*@
41197177400SBarry Smith    PCMGSetR - Sets the vector space to be used to store the residual on a
412fccaa45eSBarry Smith    particular level.
4134b9ad928SBarry Smith 
414*ad4df100SBarry Smith    Logically Collective on PC and Vec
4154b9ad928SBarry Smith 
4164b9ad928SBarry Smith    Input Parameters:
4174b9ad928SBarry Smith +  pc - the multigrid context
4184b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4194b9ad928SBarry Smith -  c - the space
4204b9ad928SBarry Smith 
4214b9ad928SBarry Smith    Level: advanced
4224b9ad928SBarry Smith 
423fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
424fccaa45eSBarry Smith 
425fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
426fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
427fccaa45eSBarry Smith 
4284b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
4294b9ad928SBarry Smith @*/
43097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetR(PC pc,PetscInt l,Vec c)
4314b9ad928SBarry Smith {
432fccaa45eSBarry Smith   PetscErrorCode ierr;
433f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
434f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4354b9ad928SBarry Smith 
4364b9ad928SBarry Smith   PetscFunctionBegin;
437e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
438e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
439c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
440f3fbd535SBarry Smith   if (mglevels[l]->r) {ierr = VecDestroy(mglevels[l]->r);CHKERRQ(ierr);}
441f3fbd535SBarry Smith   mglevels[l]->r  = c;
4424b9ad928SBarry Smith   PetscFunctionReturn(0);
4434b9ad928SBarry Smith }
444