xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision fccaa45e966e97669b5dcb5b296f1f67be099458)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
24b9ad928SBarry Smith 
34b9ad928SBarry 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__
74b9ad928SBarry Smith #define __FUNCT__ "MGDefaultResidual"
84b9ad928SBarry Smith /*@C
94b9ad928SBarry Smith    MGDefaultResidual - 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 
254b9ad928SBarry Smith .seealso: MGSetResidual()
264b9ad928SBarry Smith @*/
27dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGDefaultResidual(Mat mat,Vec b,Vec x,Vec r)
284b9ad928SBarry Smith {
29dfbe8321SBarry Smith   PetscErrorCode ierr;
304b9ad928SBarry Smith   PetscScalar    mone = -1.0;
314b9ad928SBarry Smith 
324b9ad928SBarry Smith   PetscFunctionBegin;
334b9ad928SBarry Smith   ierr = MatMult(mat,x,r);CHKERRQ(ierr);
344b9ad928SBarry Smith   ierr = VecAYPX(&mone,b,r);CHKERRQ(ierr);
354b9ad928SBarry Smith   PetscFunctionReturn(0);
364b9ad928SBarry Smith }
374b9ad928SBarry Smith 
384b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/
394b9ad928SBarry Smith 
404b9ad928SBarry Smith #undef __FUNCT__
414b9ad928SBarry Smith #define __FUNCT__ "MGGetCoarseSolve"
424b9ad928SBarry Smith /*@C
434b9ad928SBarry Smith    MGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
444b9ad928SBarry Smith 
454b9ad928SBarry Smith    Not Collective
464b9ad928SBarry Smith 
474b9ad928SBarry Smith    Input Parameter:
484b9ad928SBarry Smith .  pc - the multigrid context
494b9ad928SBarry Smith 
504b9ad928SBarry Smith    Output Parameter:
514b9ad928SBarry Smith .  ksp - the coarse grid solver context
524b9ad928SBarry Smith 
534b9ad928SBarry Smith    Level: advanced
544b9ad928SBarry Smith 
554b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid
564b9ad928SBarry Smith @*/
57dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGGetCoarseSolve(PC pc,KSP *ksp)
584b9ad928SBarry Smith {
594b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
604b9ad928SBarry Smith 
614b9ad928SBarry Smith   PetscFunctionBegin;
624b9ad928SBarry Smith   *ksp =  mg[0]->smoothd;
634b9ad928SBarry Smith   PetscFunctionReturn(0);
644b9ad928SBarry Smith }
654b9ad928SBarry Smith 
664b9ad928SBarry Smith #undef __FUNCT__
674b9ad928SBarry Smith #define __FUNCT__ "MGSetResidual"
684b9ad928SBarry Smith /*@C
694b9ad928SBarry Smith    MGSetResidual - 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
774b9ad928SBarry Smith .  residual - function used to form residual (usually MGDefaultResidual)
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 
844b9ad928SBarry Smith .seealso: MGDefaultResidual()
854b9ad928SBarry Smith @*/
86dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
874b9ad928SBarry Smith {
884b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
894b9ad928SBarry Smith 
904b9ad928SBarry Smith   PetscFunctionBegin;
914b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
924b9ad928SBarry Smith 
934b9ad928SBarry Smith   mg[l]->residual = residual;
944b9ad928SBarry Smith   mg[l]->A        = mat;
954b9ad928SBarry Smith   PetscFunctionReturn(0);
964b9ad928SBarry Smith }
974b9ad928SBarry Smith 
984b9ad928SBarry Smith #undef __FUNCT__
994b9ad928SBarry Smith #define __FUNCT__ "MGSetInterpolate"
1004b9ad928SBarry Smith /*@
1014b9ad928SBarry Smith    MGSetInterpolate - Sets the function to be used to calculate the
1024b9ad928SBarry Smith    interpolation on the lth level.
1034b9ad928SBarry Smith 
1044b9ad928SBarry Smith    Collective on PC and Mat
1054b9ad928SBarry Smith 
1064b9ad928SBarry Smith    Input Parameters:
1074b9ad928SBarry Smith +  pc  - the multigrid context
1084b9ad928SBarry Smith .  mat - the interpolation operator
1094b9ad928SBarry Smith -  l   - the level (0 is coarsest) to supply
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 
120*fccaa45eSBarry Smith          If you do not set this, the transpose of the Mat set with MGSetRestriction()
121*fccaa45eSBarry Smith     is used.
122*fccaa45eSBarry Smith 
1234b9ad928SBarry Smith .keywords:  multigrid, set, interpolate, level
1244b9ad928SBarry Smith 
1254b9ad928SBarry Smith .seealso: MGSetRestriction()
1264b9ad928SBarry Smith @*/
127dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetInterpolate(PC pc,PetscInt l,Mat mat)
1284b9ad928SBarry Smith {
1294b9ad928SBarry Smith   MG             *mg = (MG*)pc->data;
130*fccaa45eSBarry Smith   PetscErrorCode ierr;
1314b9ad928SBarry Smith 
1324b9ad928SBarry Smith   PetscFunctionBegin;
1334b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1344b9ad928SBarry Smith   mg[l]->interpolate = mat;
135*fccaa45eSBarry Smith   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1364b9ad928SBarry Smith   PetscFunctionReturn(0);
1374b9ad928SBarry Smith }
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith #undef __FUNCT__
1404b9ad928SBarry Smith #define __FUNCT__ "MGSetRestriction"
1414b9ad928SBarry Smith /*@
1424b9ad928SBarry Smith    MGSetRestriction - Sets the function to be used to restrict vector
1434b9ad928SBarry Smith    from level l to l-1.
1444b9ad928SBarry Smith 
1454b9ad928SBarry Smith    Collective on PC and Mat
1464b9ad928SBarry Smith 
1474b9ad928SBarry Smith    Input Parameters:
1484b9ad928SBarry Smith +  pc - the multigrid context
1494b9ad928SBarry Smith .  mat - the restriction matrix
1504b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
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 
161*fccaa45eSBarry Smith          If you do not set this, the transpose of the Mat set with MGSetInterpolate()
162*fccaa45eSBarry Smith     is used.
163*fccaa45eSBarry Smith 
1644b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
1654b9ad928SBarry Smith 
1664b9ad928SBarry Smith .seealso: MGSetInterpolate()
1674b9ad928SBarry Smith @*/
168dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetRestriction(PC pc,PetscInt l,Mat mat)
1694b9ad928SBarry Smith {
170*fccaa45eSBarry Smith   PetscErrorCode ierr;
1714b9ad928SBarry Smith   MG             *mg = (MG*)pc->data;
1724b9ad928SBarry Smith 
1734b9ad928SBarry Smith   PetscFunctionBegin;
1744b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1754b9ad928SBarry Smith   mg[l]->restrct  = mat;
176*fccaa45eSBarry Smith   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1774b9ad928SBarry Smith   PetscFunctionReturn(0);
1784b9ad928SBarry Smith }
1794b9ad928SBarry Smith 
1804b9ad928SBarry Smith #undef __FUNCT__
1814b9ad928SBarry Smith #define __FUNCT__ "MGGetSmoother"
1824b9ad928SBarry Smith /*@C
1834b9ad928SBarry Smith    MGGetSmoother - Gets the KSP context to be used as smoother for
1844b9ad928SBarry Smith    both pre- and post-smoothing.  Call both MGGetSmootherUp() and
1854b9ad928SBarry Smith    MGGetSmootherDown() to use different functions for pre- and
1864b9ad928SBarry Smith    post-smoothing.
1874b9ad928SBarry Smith 
1884b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
1894b9ad928SBarry Smith 
1904b9ad928SBarry Smith    Input Parameters:
1914b9ad928SBarry Smith +  pc - the multigrid context
1924b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
1934b9ad928SBarry Smith 
1944b9ad928SBarry Smith    Ouput Parameters:
1954b9ad928SBarry Smith .  ksp - the smoother
1964b9ad928SBarry Smith 
1974b9ad928SBarry Smith    Level: advanced
1984b9ad928SBarry Smith 
1994b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
2004b9ad928SBarry Smith 
2014b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmootherDown()
2024b9ad928SBarry Smith @*/
203dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGGetSmoother(PC pc,PetscInt l,KSP *ksp)
2044b9ad928SBarry Smith {
2054b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
2064b9ad928SBarry Smith 
2074b9ad928SBarry Smith   PetscFunctionBegin;
2084b9ad928SBarry Smith   *ksp = mg[l]->smoothd;
2094b9ad928SBarry Smith   PetscFunctionReturn(0);
2104b9ad928SBarry Smith }
2114b9ad928SBarry Smith 
2124b9ad928SBarry Smith #undef __FUNCT__
2134b9ad928SBarry Smith #define __FUNCT__ "MGGetSmootherUp"
2144b9ad928SBarry Smith /*@C
2154b9ad928SBarry Smith    MGGetSmootherUp - Gets the KSP context to be used as smoother after
2164b9ad928SBarry Smith    coarse grid correction (post-smoother).
2174b9ad928SBarry Smith 
2184b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith    Input Parameters:
2214b9ad928SBarry Smith +  pc - the multigrid context
2224b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2234b9ad928SBarry Smith 
2244b9ad928SBarry Smith    Ouput Parameters:
2254b9ad928SBarry Smith .  ksp - the smoother
2264b9ad928SBarry Smith 
2274b9ad928SBarry Smith    Level: advanced
2284b9ad928SBarry Smith 
2294b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
2304b9ad928SBarry Smith 
2314b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmootherDown()
2324b9ad928SBarry Smith @*/
233dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
2344b9ad928SBarry Smith {
2354b9ad928SBarry Smith   MG             *mg = (MG*)pc->data;
236dfbe8321SBarry Smith   PetscErrorCode ierr;
2374b9ad928SBarry Smith   char           *prefix;
2384b9ad928SBarry Smith   MPI_Comm       comm;
2394b9ad928SBarry Smith 
2404b9ad928SBarry Smith   PetscFunctionBegin;
2414b9ad928SBarry Smith   /*
2424b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
2434b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
2444b9ad928SBarry Smith      if not we allocate it.
2454b9ad928SBarry Smith   */
246b03c7568SBarry Smith   if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
2474b9ad928SBarry Smith   if (mg[l]->smoothu == mg[l]->smoothd) {
2484b9ad928SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mg[l]->smoothd,&comm);CHKERRQ(ierr);
2492e70eadcSBarry Smith     ierr = KSPGetOptionsPrefix(mg[l]->smoothd,&prefix);CHKERRQ(ierr);
2504b9ad928SBarry Smith     ierr = KSPCreate(comm,&mg[l]->smoothu);CHKERRQ(ierr);
251d7d82daaSBarry Smith     ierr = KSPSetTolerances(mg[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
2524b9ad928SBarry Smith     ierr = KSPSetOptionsPrefix(mg[l]->smoothu,prefix);CHKERRQ(ierr);
25352e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,mg[l]->smoothu);CHKERRQ(ierr);
2544b9ad928SBarry Smith   }
2554b9ad928SBarry Smith   if (ksp) *ksp = mg[l]->smoothu;
2564b9ad928SBarry Smith   PetscFunctionReturn(0);
2574b9ad928SBarry Smith }
2584b9ad928SBarry Smith 
2594b9ad928SBarry Smith #undef __FUNCT__
2604b9ad928SBarry Smith #define __FUNCT__ "MGGetSmootherDown"
2614b9ad928SBarry Smith /*@C
2624b9ad928SBarry Smith    MGGetSmootherDown - Gets the KSP context to be used as smoother before
2634b9ad928SBarry Smith    coarse grid correction (pre-smoother).
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
2664b9ad928SBarry Smith 
2674b9ad928SBarry Smith    Input Parameters:
2684b9ad928SBarry Smith +  pc - the multigrid context
2694b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith    Ouput Parameters:
2724b9ad928SBarry Smith .  ksp - the smoother
2734b9ad928SBarry Smith 
2744b9ad928SBarry Smith    Level: advanced
2754b9ad928SBarry Smith 
2764b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
2774b9ad928SBarry Smith 
2784b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmoother()
2794b9ad928SBarry Smith @*/
280dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
2814b9ad928SBarry Smith {
282dfbe8321SBarry Smith   PetscErrorCode ierr;
2834b9ad928SBarry Smith   MG             *mg = (MG*)pc->data;
2844b9ad928SBarry Smith 
2854b9ad928SBarry Smith   PetscFunctionBegin;
2864b9ad928SBarry Smith   /* make sure smoother up and down are different */
2874b9ad928SBarry Smith   ierr = MGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
2884b9ad928SBarry Smith   *ksp = mg[l]->smoothd;
2894b9ad928SBarry Smith   PetscFunctionReturn(0);
2904b9ad928SBarry Smith }
2914b9ad928SBarry Smith 
2924b9ad928SBarry Smith #undef __FUNCT__
2934b9ad928SBarry Smith #define __FUNCT__ "MGSetCyclesOnLevel"
2944b9ad928SBarry Smith /*@
2954b9ad928SBarry Smith    MGSetCyclesOnLevel - Sets the number of cycles to run on this level.
2964b9ad928SBarry Smith 
2974b9ad928SBarry Smith    Collective on PC
2984b9ad928SBarry Smith 
2994b9ad928SBarry Smith    Input Parameters:
3004b9ad928SBarry Smith +  pc - the multigrid context
3014b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3024b9ad928SBarry Smith -  n  - the number of cycles
3034b9ad928SBarry Smith 
3044b9ad928SBarry Smith    Level: advanced
3054b9ad928SBarry Smith 
3064b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
3074b9ad928SBarry Smith 
3084b9ad928SBarry Smith .seealso: MGSetCycles()
3094b9ad928SBarry Smith @*/
310dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
3114b9ad928SBarry Smith {
3124b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith   PetscFunctionBegin;
3154b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
3164b9ad928SBarry Smith   mg[l]->cycles  = c;
3174b9ad928SBarry Smith   PetscFunctionReturn(0);
3184b9ad928SBarry Smith }
3194b9ad928SBarry Smith 
3204b9ad928SBarry Smith #undef __FUNCT__
3214b9ad928SBarry Smith #define __FUNCT__ "MGSetRhs"
3224b9ad928SBarry Smith /*@
3234b9ad928SBarry Smith    MGSetRhs - Sets the vector space to be used to store the right-hand side
324*fccaa45eSBarry Smith    on a particular level.
3254b9ad928SBarry Smith 
3264b9ad928SBarry Smith    Collective on PC and Vec
3274b9ad928SBarry Smith 
3284b9ad928SBarry Smith    Input Parameters:
3294b9ad928SBarry Smith +  pc - the multigrid context
3304b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
3314b9ad928SBarry Smith -  c  - the space
3324b9ad928SBarry Smith 
3334b9ad928SBarry Smith    Level: advanced
3344b9ad928SBarry Smith 
335*fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
336*fccaa45eSBarry Smith 
337*fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
338*fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
339*fccaa45eSBarry Smith 
3404b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
3414b9ad928SBarry Smith 
3424b9ad928SBarry Smith .seealso: MGSetX(), MGSetR()
3434b9ad928SBarry Smith @*/
344dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetRhs(PC pc,PetscInt l,Vec c)
3454b9ad928SBarry Smith {
346*fccaa45eSBarry Smith   PetscErrorCode ierr;
3474b9ad928SBarry Smith   MG             *mg = (MG*)pc->data;
3484b9ad928SBarry Smith 
3494b9ad928SBarry Smith   PetscFunctionBegin;
3504b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
351*fccaa45eSBarry Smith   if (l == mg[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
3524b9ad928SBarry Smith   mg[l]->b  = c;
353*fccaa45eSBarry Smith   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
3544b9ad928SBarry Smith   PetscFunctionReturn(0);
3554b9ad928SBarry Smith }
3564b9ad928SBarry Smith 
3574b9ad928SBarry Smith #undef __FUNCT__
3584b9ad928SBarry Smith #define __FUNCT__ "MGSetX"
3594b9ad928SBarry Smith /*@
3604b9ad928SBarry Smith    MGSetX - Sets the vector space to be used to store the solution on a
361*fccaa45eSBarry Smith    particular level.
3624b9ad928SBarry Smith 
3634b9ad928SBarry Smith    Collective on PC and Vec
3644b9ad928SBarry Smith 
3654b9ad928SBarry Smith    Input Parameters:
3664b9ad928SBarry Smith +  pc - the multigrid context
3674b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
3684b9ad928SBarry Smith -  c - the space
3694b9ad928SBarry Smith 
3704b9ad928SBarry Smith    Level: advanced
3714b9ad928SBarry Smith 
372*fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
373*fccaa45eSBarry Smith 
374*fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
375*fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
376*fccaa45eSBarry Smith 
3774b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
3784b9ad928SBarry Smith 
3794b9ad928SBarry Smith .seealso: MGSetRhs(), MGSetR()
3804b9ad928SBarry Smith @*/
381dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetX(PC pc,PetscInt l,Vec c)
3824b9ad928SBarry Smith {
383*fccaa45eSBarry Smith   PetscErrorCode ierr;
3844b9ad928SBarry Smith   MG             *mg = (MG*)pc->data;
3854b9ad928SBarry Smith 
3864b9ad928SBarry Smith   PetscFunctionBegin;
3874b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
388*fccaa45eSBarry Smith   if (l == mg[0]->levels-1) SETERRQ(PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
3894b9ad928SBarry Smith   mg[l]->x  = c;
390*fccaa45eSBarry Smith   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
3914b9ad928SBarry Smith   PetscFunctionReturn(0);
3924b9ad928SBarry Smith }
3934b9ad928SBarry Smith 
3944b9ad928SBarry Smith #undef __FUNCT__
3954b9ad928SBarry Smith #define __FUNCT__ "MGSetR"
3964b9ad928SBarry Smith /*@
3974b9ad928SBarry Smith    MGSetR - Sets the vector space to be used to store the residual on a
398*fccaa45eSBarry Smith    particular level.
3994b9ad928SBarry Smith 
4004b9ad928SBarry Smith    Collective on PC and Vec
4014b9ad928SBarry Smith 
4024b9ad928SBarry Smith    Input Parameters:
4034b9ad928SBarry Smith +  pc - the multigrid context
4044b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
4054b9ad928SBarry Smith -  c - the space
4064b9ad928SBarry Smith 
4074b9ad928SBarry Smith    Level: advanced
4084b9ad928SBarry Smith 
409*fccaa45eSBarry Smith    Notes: If this is not provided PETSc will automatically generate one.
410*fccaa45eSBarry Smith 
411*fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
412*fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
413*fccaa45eSBarry Smith 
4144b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
4154b9ad928SBarry Smith @*/
416dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT MGSetR(PC pc,PetscInt l,Vec c)
4174b9ad928SBarry Smith {
418*fccaa45eSBarry Smith   PetscErrorCode ierr;
4194b9ad928SBarry Smith   MG             *mg = (MG*)pc->data;
4204b9ad928SBarry Smith 
4214b9ad928SBarry Smith   PetscFunctionBegin;
4224b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
4234b9ad928SBarry Smith   mg[l]->r  = c;
424*fccaa45eSBarry Smith   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
4254b9ad928SBarry Smith   PetscFunctionReturn(0);
4264b9ad928SBarry Smith }
4274b9ad928SBarry Smith 
4284b9ad928SBarry Smith 
4294b9ad928SBarry Smith 
4304b9ad928SBarry Smith 
4314b9ad928SBarry Smith 
4324b9ad928SBarry Smith 
4334b9ad928SBarry Smith 
4344b9ad928SBarry Smith 
435