xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision c5eb91543d2ee8daf880d93389b892228ddada03)
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;
62*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
63f3fbd535SBarry Smith   *ksp =  mglevels[0]->smoothd;
644b9ad928SBarry Smith   PetscFunctionReturn(0);
654b9ad928SBarry Smith }
664b9ad928SBarry Smith 
674b9ad928SBarry Smith #undef __FUNCT__
68d6337fedSHong Zhang #define __FUNCT__ "PCMGSetResidual"
694b9ad928SBarry Smith /*@C
7097177400SBarry Smith    PCMGSetResidual - Sets the function to be used to calculate the residual
714b9ad928SBarry Smith    on the lth level.
724b9ad928SBarry Smith 
73ad4df100SBarry Smith    Logically Collective on PC and Mat
744b9ad928SBarry Smith 
754b9ad928SBarry Smith    Input Parameters:
764b9ad928SBarry Smith +  pc       - the multigrid context
774b9ad928SBarry Smith .  l        - the level (0 is coarsest) to supply
7897177400SBarry Smith .  residual - function used to form residual (usually PCMGDefaultResidual)
794b9ad928SBarry Smith -  mat      - matrix associated with residual
804b9ad928SBarry Smith 
814b9ad928SBarry Smith    Level: advanced
824b9ad928SBarry Smith 
834b9ad928SBarry Smith .keywords:  MG, set, multigrid, residual, level
844b9ad928SBarry Smith 
8597177400SBarry Smith .seealso: PCMGDefaultResidual()
864b9ad928SBarry Smith @*/
8797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
884b9ad928SBarry Smith {
89f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
90f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
914b9ad928SBarry Smith 
924b9ad928SBarry Smith   PetscFunctionBegin;
93*c5eb9154SBarry 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;
96f3fbd535SBarry Smith   mglevels[l]->A        = mat;
974b9ad928SBarry Smith   PetscFunctionReturn(0);
984b9ad928SBarry Smith }
994b9ad928SBarry Smith 
1004b9ad928SBarry Smith #undef __FUNCT__
101d6337fedSHong Zhang #define __FUNCT__ "PCMGSetInterpolation"
1024b9ad928SBarry Smith /*@
103aea2a34eSBarry Smith    PCMGSetInterpolation - Sets the function to be used to calculate the
104bf5b2e24SBarry Smith    interpolation from l-1 to the lth level
1054b9ad928SBarry Smith 
106ad4df100SBarry Smith    Logically Collective on PC and Mat
1074b9ad928SBarry Smith 
1084b9ad928SBarry Smith    Input Parameters:
1094b9ad928SBarry Smith +  pc  - the multigrid context
1104b9ad928SBarry Smith .  mat - the interpolation operator
111bf5b2e24SBarry Smith -  l   - the level (0 is coarsest) to supply [do not supply 0]
1124b9ad928SBarry Smith 
1134b9ad928SBarry Smith    Level: advanced
1144b9ad928SBarry Smith 
1154b9ad928SBarry Smith    Notes:
1164b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
1174b9ad928SBarry Smith     for the same level.
1184b9ad928SBarry Smith 
1194b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1204b9ad928SBarry Smith     out from the matrix size which one it is.
1214b9ad928SBarry Smith 
1224b9ad928SBarry Smith .keywords:  multigrid, set, interpolate, level
1234b9ad928SBarry Smith 
12497177400SBarry Smith .seealso: PCMGSetRestriction()
1254b9ad928SBarry Smith @*/
126aea2a34eSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
1274b9ad928SBarry Smith {
128f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
129f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
130fccaa45eSBarry Smith   PetscErrorCode ierr;
1314b9ad928SBarry Smith 
1324b9ad928SBarry Smith   PetscFunctionBegin;
133*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
134e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
135e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
136c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
137f3fbd535SBarry Smith   if (mglevels[l]->interpolate) {ierr = MatDestroy(mglevels[l]->interpolate);CHKERRQ(ierr);}
138f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
1394b9ad928SBarry Smith   PetscFunctionReturn(0);
1404b9ad928SBarry Smith }
1414b9ad928SBarry Smith 
1424b9ad928SBarry Smith #undef __FUNCT__
143d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRestriction"
1444b9ad928SBarry Smith /*@
14597177400SBarry Smith    PCMGSetRestriction - Sets the function to be used to restrict vector
1464b9ad928SBarry Smith    from level l to l-1.
1474b9ad928SBarry Smith 
148ad4df100SBarry Smith    Logically Collective on PC and Mat
1494b9ad928SBarry Smith 
1504b9ad928SBarry Smith    Input Parameters:
1514b9ad928SBarry Smith +  pc - the multigrid context
1524b9ad928SBarry Smith .  mat - the restriction matrix
153bf5b2e24SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
1544b9ad928SBarry Smith 
1554b9ad928SBarry Smith    Level: advanced
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith    Notes:
1584b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
1594b9ad928SBarry Smith     for the same level.
1604b9ad928SBarry Smith 
1614b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1624b9ad928SBarry Smith     out from the matrix size which one it is.
1634b9ad928SBarry Smith 
164aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
165fccaa45eSBarry Smith     is used.
166fccaa45eSBarry Smith 
1674b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
1684b9ad928SBarry Smith 
169aea2a34eSBarry Smith .seealso: PCMGSetInterpolation()
1704b9ad928SBarry Smith @*/
17197177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
1724b9ad928SBarry Smith {
173fccaa45eSBarry Smith   PetscErrorCode ierr;
174f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
175f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1764b9ad928SBarry Smith 
1774b9ad928SBarry Smith   PetscFunctionBegin;
178*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
179e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
180e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
181c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
182f3fbd535SBarry Smith   if (mglevels[l]->restrct) {ierr = MatDestroy(mglevels[l]->restrct);CHKERRQ(ierr);}
183f3fbd535SBarry Smith   mglevels[l]->restrct  = mat;
1844b9ad928SBarry Smith   PetscFunctionReturn(0);
1854b9ad928SBarry Smith }
1864b9ad928SBarry Smith 
1874b9ad928SBarry Smith #undef __FUNCT__
188d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmoother"
189f39d8e23SSatish Balay /*@
19097177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
19197177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
19297177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
1934b9ad928SBarry Smith    post-smoothing.
1944b9ad928SBarry Smith 
1954b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
1964b9ad928SBarry Smith 
1974b9ad928SBarry Smith    Input Parameters:
1984b9ad928SBarry Smith +  pc - the multigrid context
1994b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
2004b9ad928SBarry Smith 
2014b9ad928SBarry Smith    Ouput Parameters:
2024b9ad928SBarry Smith .  ksp - the smoother
2034b9ad928SBarry Smith 
2044b9ad928SBarry Smith    Level: advanced
2054b9ad928SBarry Smith 
2064b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
2074b9ad928SBarry Smith 
20897177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2094b9ad928SBarry Smith @*/
21097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
2114b9ad928SBarry Smith {
212f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
213f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2144b9ad928SBarry Smith 
2154b9ad928SBarry Smith   PetscFunctionBegin;
216*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
217f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
2184b9ad928SBarry Smith   PetscFunctionReturn(0);
2194b9ad928SBarry Smith }
2204b9ad928SBarry Smith 
2214b9ad928SBarry Smith #undef __FUNCT__
222d6337fedSHong Zhang #define __FUNCT__ "PCMGGetSmootherUp"
223f39d8e23SSatish Balay /*@
22497177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
2254b9ad928SBarry Smith    coarse grid correction (post-smoother).
2264b9ad928SBarry Smith 
2274b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2284b9ad928SBarry Smith 
2294b9ad928SBarry Smith    Input Parameters:
2304b9ad928SBarry Smith +  pc - the multigrid context
2314b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith    Ouput Parameters:
2344b9ad928SBarry Smith .  ksp - the smoother
2354b9ad928SBarry Smith 
2364b9ad928SBarry Smith    Level: advanced
2374b9ad928SBarry Smith 
2384b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
2394b9ad928SBarry Smith 
24097177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
2414b9ad928SBarry Smith @*/
24297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
2434b9ad928SBarry Smith {
244f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
245f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
246dfbe8321SBarry Smith   PetscErrorCode ierr;
247f69a0ea3SMatthew Knepley   const char     *prefix;
2484b9ad928SBarry Smith   MPI_Comm       comm;
2494b9ad928SBarry Smith 
2504b9ad928SBarry Smith   PetscFunctionBegin;
251*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2524b9ad928SBarry Smith   /*
2534b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
2544b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
2554b9ad928SBarry Smith      if not we allocate it.
2564b9ad928SBarry Smith   */
257e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
258f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
259f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
260f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
261f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
262f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
263f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
264f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
265f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr);
2664b9ad928SBarry Smith   }
267f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
2684b9ad928SBarry Smith   PetscFunctionReturn(0);
2694b9ad928SBarry Smith }
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith #undef __FUNCT__
2729dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetSmootherDown"
273f39d8e23SSatish Balay /*@
27497177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
2754b9ad928SBarry Smith    coarse grid correction (pre-smoother).
2764b9ad928SBarry Smith 
2774b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2784b9ad928SBarry Smith 
2794b9ad928SBarry Smith    Input Parameters:
2804b9ad928SBarry Smith +  pc - the multigrid context
2814b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2824b9ad928SBarry Smith 
2834b9ad928SBarry Smith    Ouput Parameters:
2844b9ad928SBarry Smith .  ksp - the smoother
2854b9ad928SBarry Smith 
2864b9ad928SBarry Smith    Level: advanced
2874b9ad928SBarry Smith 
2884b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
2894b9ad928SBarry Smith 
29097177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
2914b9ad928SBarry Smith @*/
29297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
2934b9ad928SBarry Smith {
294dfbe8321SBarry Smith   PetscErrorCode ierr;
295f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
296f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2974b9ad928SBarry Smith 
2984b9ad928SBarry Smith   PetscFunctionBegin;
299*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3004b9ad928SBarry Smith   /* make sure smoother up and down are different */
301*c5eb9154SBarry Smith   if (l) {
30297177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
303d8148a5aSMatthew Knepley   }
304f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
3054b9ad928SBarry Smith   PetscFunctionReturn(0);
3064b9ad928SBarry Smith }
3074b9ad928SBarry Smith 
3084b9ad928SBarry Smith #undef __FUNCT__
3099dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetCyclesOnLevel"
3104b9ad928SBarry Smith /*@
31197177400SBarry Smith    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
3124b9ad928SBarry Smith 
313ad4df100SBarry Smith    Logically Collective on PC
3144b9ad928SBarry Smith 
3154b9ad928SBarry Smith    Input Parameters:
3164b9ad928SBarry Smith +  pc - the multigrid context
3174b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3184b9ad928SBarry Smith -  n  - the number of cycles
3194b9ad928SBarry Smith 
3204b9ad928SBarry Smith    Level: advanced
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
3234b9ad928SBarry Smith 
32497177400SBarry Smith .seealso: PCMGSetCycles()
3254b9ad928SBarry Smith @*/
32697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
3274b9ad928SBarry Smith {
328f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
329f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith   PetscFunctionBegin;
332*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
333e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
334*c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,l,2);
335*c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,c,3);
336f3fbd535SBarry Smith   mglevels[l]->cycles  = c;
3374b9ad928SBarry Smith   PetscFunctionReturn(0);
3384b9ad928SBarry Smith }
3394b9ad928SBarry Smith 
3404b9ad928SBarry Smith #undef __FUNCT__
341d6337fedSHong Zhang #define __FUNCT__ "PCMGSetRhs"
3424b9ad928SBarry Smith /*@
34397177400SBarry Smith    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
344fccaa45eSBarry Smith    on a particular level.
3454b9ad928SBarry Smith 
346ad4df100SBarry Smith    Logically Collective on PC and Vec
3474b9ad928SBarry Smith 
3484b9ad928SBarry Smith    Input Parameters:
3494b9ad928SBarry Smith +  pc - the multigrid context
3504b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3514b9ad928SBarry Smith -  c  - the space
3524b9ad928SBarry Smith 
3534b9ad928SBarry Smith    Level: advanced
3544b9ad928SBarry Smith 
355fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
356fccaa45eSBarry Smith 
357fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
358fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
359fccaa45eSBarry Smith 
3604b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
3614b9ad928SBarry Smith 
36297177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
3634b9ad928SBarry Smith @*/
36497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetRhs(PC pc,PetscInt l,Vec c)
3654b9ad928SBarry Smith {
366fccaa45eSBarry Smith   PetscErrorCode ierr;
367f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
368f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
3694b9ad928SBarry Smith 
3704b9ad928SBarry Smith   PetscFunctionBegin;
371*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
372e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
373e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
374c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
375f3fbd535SBarry Smith   if (mglevels[l]->b) {ierr = VecDestroy(mglevels[l]->b);CHKERRQ(ierr);}
376f3fbd535SBarry Smith   mglevels[l]->b  = c;
3774b9ad928SBarry Smith   PetscFunctionReturn(0);
3784b9ad928SBarry Smith }
3794b9ad928SBarry Smith 
3804b9ad928SBarry Smith #undef __FUNCT__
381d6337fedSHong Zhang #define __FUNCT__ "PCMGSetX"
3824b9ad928SBarry Smith /*@
38397177400SBarry Smith    PCMGSetX - Sets the vector space to be used to store the solution on a
384fccaa45eSBarry Smith    particular level.
3854b9ad928SBarry Smith 
386ad4df100SBarry Smith    Logically Collective on PC and Vec
3874b9ad928SBarry Smith 
3884b9ad928SBarry Smith    Input Parameters:
3894b9ad928SBarry Smith +  pc - the multigrid context
3904b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
3914b9ad928SBarry Smith -  c - the space
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith    Level: advanced
3944b9ad928SBarry Smith 
395fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
396fccaa45eSBarry Smith 
397fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
398fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
399fccaa45eSBarry Smith 
4004b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
4014b9ad928SBarry Smith 
40297177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
4034b9ad928SBarry Smith @*/
40497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetX(PC pc,PetscInt l,Vec c)
4054b9ad928SBarry Smith {
406fccaa45eSBarry Smith   PetscErrorCode ierr;
407f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
408f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4094b9ad928SBarry Smith 
4104b9ad928SBarry Smith   PetscFunctionBegin;
411*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
412e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
413e7e72b3dSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
414c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
415f3fbd535SBarry Smith   if (mglevels[l]->x) {ierr = VecDestroy(mglevels[l]->x);CHKERRQ(ierr);}
416f3fbd535SBarry Smith   mglevels[l]->x  = c;
4174b9ad928SBarry Smith   PetscFunctionReturn(0);
4184b9ad928SBarry Smith }
4194b9ad928SBarry Smith 
4204b9ad928SBarry Smith #undef __FUNCT__
421d6337fedSHong Zhang #define __FUNCT__ "PCMGSetR"
4224b9ad928SBarry Smith /*@
42397177400SBarry Smith    PCMGSetR - Sets the vector space to be used to store the residual on a
424fccaa45eSBarry Smith    particular level.
4254b9ad928SBarry Smith 
426ad4df100SBarry Smith    Logically Collective on PC and Vec
4274b9ad928SBarry Smith 
4284b9ad928SBarry Smith    Input Parameters:
4294b9ad928SBarry Smith +  pc - the multigrid context
4304b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4314b9ad928SBarry Smith -  c - the space
4324b9ad928SBarry Smith 
4334b9ad928SBarry Smith    Level: advanced
4344b9ad928SBarry Smith 
435fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
436fccaa45eSBarry Smith 
437fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
438fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
439fccaa45eSBarry Smith 
4404b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
4414b9ad928SBarry Smith @*/
44297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetR(PC pc,PetscInt l,Vec c)
4434b9ad928SBarry Smith {
444fccaa45eSBarry Smith   PetscErrorCode ierr;
445f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
446f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
4474b9ad928SBarry Smith 
4484b9ad928SBarry Smith   PetscFunctionBegin;
449*c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
450e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
451e7e72b3dSBarry Smith   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
452c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
453f3fbd535SBarry Smith   if (mglevels[l]->r) {ierr = VecDestroy(mglevels[l]->r);CHKERRQ(ierr);}
454f3fbd535SBarry Smith   mglevels[l]->r  = c;
4554b9ad928SBarry Smith   PetscFunctionReturn(0);
4564b9ad928SBarry Smith }
457