xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 28668aa566ec83483fc98fc3129fc8da840155ee)
14b9ad928SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>       /*I "petscksp.h" I*/
34b9ad928SBarry Smith 
4f9426fe0SMark Adams /* ---------------------------------------------------------------------------*/
51f6cc5b2SSatish Balay /*@C
654b2cd4bSJed Brown    PCMGResidualDefault - Default routine to calculate the residual.
74b9ad928SBarry Smith 
84b9ad928SBarry Smith    Collective on Mat and Vec
94b9ad928SBarry Smith 
104b9ad928SBarry Smith    Input Parameters:
114b9ad928SBarry Smith +  mat - the matrix
124b9ad928SBarry Smith .  b   - the right-hand-side
134b9ad928SBarry Smith -  x   - the approximate solution
144b9ad928SBarry Smith 
154b9ad928SBarry Smith    Output Parameter:
164b9ad928SBarry Smith .  r - location to store the residual
174b9ad928SBarry Smith 
18d0e4de75SBarry Smith    Level: developer
194b9ad928SBarry Smith 
204b9ad928SBarry Smith .keywords: MG, default, multigrid, residual
214b9ad928SBarry Smith 
2297177400SBarry Smith .seealso: PCMGSetResidual()
234b9ad928SBarry Smith @*/
2454b2cd4bSJed Brown PetscErrorCode  PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r)
254b9ad928SBarry Smith {
26dfbe8321SBarry Smith   PetscErrorCode ierr;
274b9ad928SBarry Smith 
284b9ad928SBarry Smith   PetscFunctionBegin;
29f9426fe0SMark Adams   ierr = MatResidual(mat,b,x,r);CHKERRQ(ierr);
304b9ad928SBarry Smith   PetscFunctionReturn(0);
314b9ad928SBarry Smith }
324b9ad928SBarry Smith 
33f39d8e23SSatish Balay /*@
3497177400SBarry Smith    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
354b9ad928SBarry Smith 
364b9ad928SBarry Smith    Not Collective
374b9ad928SBarry Smith 
384b9ad928SBarry Smith    Input Parameter:
394b9ad928SBarry Smith .  pc - the multigrid context
404b9ad928SBarry Smith 
414b9ad928SBarry Smith    Output Parameter:
424b9ad928SBarry Smith .  ksp - the coarse grid solver context
434b9ad928SBarry Smith 
444b9ad928SBarry Smith    Level: advanced
454b9ad928SBarry Smith 
464b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid
47*28668aa5SBarry Smith 
48*28668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetSmoother()
494b9ad928SBarry Smith @*/
507087cfbeSBarry Smith PetscErrorCode  PCMGGetCoarseSolve(PC pc,KSP *ksp)
514b9ad928SBarry Smith {
52f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
53f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
544b9ad928SBarry Smith 
554b9ad928SBarry Smith   PetscFunctionBegin;
56c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
57f3fbd535SBarry Smith   *ksp =  mglevels[0]->smoothd;
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
614b9ad928SBarry Smith /*@C
6297177400SBarry Smith    PCMGSetResidual - Sets the function to be used to calculate the residual
634b9ad928SBarry Smith    on the lth level.
644b9ad928SBarry Smith 
65ad4df100SBarry Smith    Logically Collective on PC and Mat
664b9ad928SBarry Smith 
674b9ad928SBarry Smith    Input Parameters:
684b9ad928SBarry Smith +  pc       - the multigrid context
694b9ad928SBarry Smith .  l        - the level (0 is coarsest) to supply
70157726a2SBarry Smith .  residual - function used to form residual, if none is provided the previously provide one is used, if no
71d0e4de75SBarry Smith               previous one were provided then a default is used
724b9ad928SBarry Smith -  mat      - matrix associated with residual
734b9ad928SBarry Smith 
744b9ad928SBarry Smith    Level: advanced
754b9ad928SBarry Smith 
764b9ad928SBarry Smith .keywords:  MG, set, multigrid, residual, level
774b9ad928SBarry Smith 
7854b2cd4bSJed Brown .seealso: PCMGResidualDefault()
794b9ad928SBarry Smith @*/
807087cfbeSBarry Smith PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
814b9ad928SBarry Smith {
82f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
83f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
84298cc208SBarry Smith   PetscErrorCode ierr;
854b9ad928SBarry Smith 
864b9ad928SBarry Smith   PetscFunctionBegin;
87c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
88ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
892fa5cd67SKarl Rupp   if (residual) mglevels[l]->residual = residual;
9054b2cd4bSJed Brown   if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault;
91f3ae41bdSBarry Smith   if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);}
92298cc208SBarry Smith   ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr);
93f3fbd535SBarry Smith   mglevels[l]->A = mat;
944b9ad928SBarry Smith   PetscFunctionReturn(0);
954b9ad928SBarry Smith }
964b9ad928SBarry Smith 
974b9ad928SBarry Smith /*@
98aea2a34eSBarry Smith    PCMGSetInterpolation - Sets the function to be used to calculate the
99bf5b2e24SBarry Smith    interpolation from l-1 to the lth level
1004b9ad928SBarry Smith 
101ad4df100SBarry Smith    Logically Collective on PC and Mat
1024b9ad928SBarry Smith 
1034b9ad928SBarry Smith    Input Parameters:
1044b9ad928SBarry Smith +  pc  - the multigrid context
1054b9ad928SBarry Smith .  mat - the interpolation operator
106bf5b2e24SBarry Smith -  l   - the level (0 is coarsest) to supply [do not supply 0]
1074b9ad928SBarry Smith 
1084b9ad928SBarry Smith    Level: advanced
1094b9ad928SBarry Smith 
1104b9ad928SBarry Smith    Notes:
1114b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
1124b9ad928SBarry Smith     for the same level.
1134b9ad928SBarry Smith 
1144b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1154b9ad928SBarry Smith     out from the matrix size which one it is.
1164b9ad928SBarry Smith 
1174b9ad928SBarry Smith .keywords:  multigrid, set, interpolate, level
1184b9ad928SBarry Smith 
11997177400SBarry Smith .seealso: PCMGSetRestriction()
1204b9ad928SBarry Smith @*/
1217087cfbeSBarry Smith PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
1224b9ad928SBarry Smith {
123f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
124f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
125fccaa45eSBarry Smith   PetscErrorCode ierr;
1264b9ad928SBarry Smith 
1274b9ad928SBarry Smith   PetscFunctionBegin;
128c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
129ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
130ce94432eSBarry Smith   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
131c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1326bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr);
1332fa5cd67SKarl Rupp 
134f3fbd535SBarry Smith   mglevels[l]->interpolate = mat;
1354b9ad928SBarry Smith   PetscFunctionReturn(0);
1364b9ad928SBarry Smith }
1374b9ad928SBarry Smith 
138c88c5224SJed Brown /*@
139c88c5224SJed Brown    PCMGGetInterpolation - Gets the function to be used to calculate the
140c88c5224SJed Brown    interpolation from l-1 to the lth level
141c88c5224SJed Brown 
142c88c5224SJed Brown    Logically Collective on PC
143c88c5224SJed Brown 
144c88c5224SJed Brown    Input Parameters:
145c88c5224SJed Brown +  pc - the multigrid context
146c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
147c88c5224SJed Brown 
148c88c5224SJed Brown    Output Parameter:
1493ad4599aSBarry Smith .  mat - the interpolation matrix, can be NULL
150c88c5224SJed Brown 
151c88c5224SJed Brown    Level: advanced
152c88c5224SJed Brown 
153c88c5224SJed Brown .keywords: MG, get, multigrid, interpolation, level
154c88c5224SJed Brown 
155c88c5224SJed Brown .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale()
156c88c5224SJed Brown @*/
157c88c5224SJed Brown PetscErrorCode  PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat)
158c88c5224SJed Brown {
159c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
160c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
161c88c5224SJed Brown   PetscErrorCode ierr;
162c88c5224SJed Brown 
163c88c5224SJed Brown   PetscFunctionBegin;
164c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1653ad4599aSBarry Smith   if (mat) PetscValidPointer(mat,3);
166ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
167ce94432eSBarry Smith   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
168c88c5224SJed Brown   if (!mglevels[l]->interpolate) {
1695aa31b60SBarry Smith     if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()");
170c88c5224SJed Brown     ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr);
171c88c5224SJed Brown   }
1723ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->interpolate;
173c88c5224SJed Brown   PetscFunctionReturn(0);
174c88c5224SJed Brown }
175c88c5224SJed Brown 
1764b9ad928SBarry Smith /*@
177eab52d2bSLawrence Mitchell    PCMGSetRestriction - Sets the function to be used to restrict dual vectors
1784b9ad928SBarry Smith    from level l to l-1.
1794b9ad928SBarry Smith 
180ad4df100SBarry Smith    Logically Collective on PC and Mat
1814b9ad928SBarry Smith 
1824b9ad928SBarry Smith    Input Parameters:
1834b9ad928SBarry Smith +  pc - the multigrid context
184c88c5224SJed Brown .  l - the level (0 is coarsest) to supply [Do not supply 0]
185c88c5224SJed Brown -  mat - the restriction matrix
1864b9ad928SBarry Smith 
1874b9ad928SBarry Smith    Level: advanced
1884b9ad928SBarry Smith 
1894b9ad928SBarry Smith    Notes:
1904b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
1914b9ad928SBarry Smith     for the same level.
1924b9ad928SBarry Smith 
1934b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
1944b9ad928SBarry Smith     out from the matrix size which one it is.
1954b9ad928SBarry Smith 
196aea2a34eSBarry Smith          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
197fccaa45eSBarry Smith     is used.
198fccaa45eSBarry Smith 
1994b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
2004b9ad928SBarry Smith 
201eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGetSetInjection()
2024b9ad928SBarry Smith @*/
2037087cfbeSBarry Smith PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
2044b9ad928SBarry Smith {
205fccaa45eSBarry Smith   PetscErrorCode ierr;
206f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
207f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
2084b9ad928SBarry Smith 
2094b9ad928SBarry Smith   PetscFunctionBegin;
210c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
211c88c5224SJed Brown   PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
212ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
213ce94432eSBarry Smith   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
214c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
2156bf464f9SBarry Smith   ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr);
2162fa5cd67SKarl Rupp 
217f3fbd535SBarry Smith   mglevels[l]->restrct = mat;
2184b9ad928SBarry Smith   PetscFunctionReturn(0);
2194b9ad928SBarry Smith }
2204b9ad928SBarry Smith 
221c88c5224SJed Brown /*@
222eab52d2bSLawrence Mitchell    PCMGGetRestriction - Gets the function to be used to restrict dual vectors
223c88c5224SJed Brown    from level l to l-1.
224c88c5224SJed Brown 
225c88c5224SJed Brown    Logically Collective on PC and Mat
226c88c5224SJed Brown 
227c88c5224SJed Brown    Input Parameters:
228c88c5224SJed Brown +  pc - the multigrid context
229c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
230c88c5224SJed Brown 
231c88c5224SJed Brown    Output Parameter:
232c88c5224SJed Brown .  mat - the restriction matrix
233c88c5224SJed Brown 
234c88c5224SJed Brown    Level: advanced
235c88c5224SJed Brown 
236c88c5224SJed Brown .keywords: MG, get, multigrid, restriction, level
237c88c5224SJed Brown 
238eab52d2bSLawrence Mitchell .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection()
239c88c5224SJed Brown @*/
240c88c5224SJed Brown PetscErrorCode  PCMGGetRestriction(PC pc,PetscInt l,Mat *mat)
241c88c5224SJed Brown {
242c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
243c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
244c88c5224SJed Brown   PetscErrorCode ierr;
245c88c5224SJed Brown 
246c88c5224SJed Brown   PetscFunctionBegin;
247c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2483ad4599aSBarry Smith   if (mat) PetscValidPointer(mat,3);
249ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
250ce94432eSBarry Smith   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
251c88c5224SJed Brown   if (!mglevels[l]->restrct) {
252ce94432eSBarry Smith     if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()");
253c88c5224SJed Brown     ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr);
254c88c5224SJed Brown   }
2553ad4599aSBarry Smith   if (mat) *mat = mglevels[l]->restrct;
256c88c5224SJed Brown   PetscFunctionReturn(0);
257c88c5224SJed Brown }
258c88c5224SJed Brown 
25973250ac0SBarry Smith /*@
26073250ac0SBarry Smith    PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.
26173250ac0SBarry Smith 
262c88c5224SJed Brown    Logically Collective on PC and Vec
263c88c5224SJed Brown 
264c88c5224SJed Brown    Input Parameters:
265c88c5224SJed Brown +  pc - the multigrid context
266c88c5224SJed Brown -  l - the level (0 is coarsest) to supply [Do not supply 0]
267c88c5224SJed Brown .  rscale - the scaling
268c88c5224SJed Brown 
269c88c5224SJed Brown    Level: advanced
270c88c5224SJed Brown 
271c88c5224SJed Brown    Notes:
272eab52d2bSLawrence Mitchell        When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R.  It is preferable to use PCMGSetInjection() to control moving primal vectors.
273c88c5224SJed Brown 
274c88c5224SJed Brown .keywords: MG, set, multigrid, restriction, level
275c88c5224SJed Brown 
276eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection()
277c88c5224SJed Brown @*/
278c88c5224SJed Brown PetscErrorCode  PCMGSetRScale(PC pc,PetscInt l,Vec rscale)
279c88c5224SJed Brown {
280c88c5224SJed Brown   PetscErrorCode ierr;
281c88c5224SJed Brown   PC_MG          *mg        = (PC_MG*)pc->data;
282c88c5224SJed Brown   PC_MG_Levels   **mglevels = mg->levels;
283c88c5224SJed Brown 
284c88c5224SJed Brown   PetscFunctionBegin;
285c88c5224SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
286ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
287ce94432eSBarry Smith   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
288c88c5224SJed Brown   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
289c88c5224SJed Brown   ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr);
2902fa5cd67SKarl Rupp 
291c88c5224SJed Brown   mglevels[l]->rscale = rscale;
292c88c5224SJed Brown   PetscFunctionReturn(0);
293c88c5224SJed Brown }
294c88c5224SJed Brown 
295c88c5224SJed Brown /*@
296c88c5224SJed Brown    PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1.
297c88c5224SJed Brown 
298c88c5224SJed Brown    Collective on PC
29973250ac0SBarry Smith 
30073250ac0SBarry Smith    Input Parameters:
30173250ac0SBarry Smith +  pc - the multigrid context
30273250ac0SBarry Smith .  rscale - the scaling
30373250ac0SBarry Smith -  l - the level (0 is coarsest) to supply [Do not supply 0]
30473250ac0SBarry Smith 
30573250ac0SBarry Smith    Level: advanced
30673250ac0SBarry Smith 
30773250ac0SBarry Smith    Notes:
308eab52d2bSLawrence Mitchell        When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R.  It is preferable to use PCMGGetInjection() to control moving primal vectors.
30973250ac0SBarry Smith 
31073250ac0SBarry Smith .keywords: MG, set, multigrid, restriction, level
31173250ac0SBarry Smith 
312eab52d2bSLawrence Mitchell .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection()
31373250ac0SBarry Smith @*/
314c88c5224SJed Brown PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale)
31573250ac0SBarry Smith {
31673250ac0SBarry Smith   PetscErrorCode ierr;
31773250ac0SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
31873250ac0SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
31973250ac0SBarry Smith 
32073250ac0SBarry Smith   PetscFunctionBegin;
32173250ac0SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
322ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
323ce94432eSBarry Smith   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
324c88c5224SJed Brown   if (!mglevels[l]->rscale) {
325c88c5224SJed Brown     Mat      R;
326c88c5224SJed Brown     Vec      X,Y,coarse,fine;
327c88c5224SJed Brown     PetscInt M,N;
328c88c5224SJed Brown     ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr);
3292a7a6963SBarry Smith     ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr);
330c88c5224SJed Brown     ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr);
3312fa5cd67SKarl Rupp     if (M < N) {
3322fa5cd67SKarl Rupp       fine = X;
3332fa5cd67SKarl Rupp       coarse = Y;
3342fa5cd67SKarl Rupp     } else if (N < M) {
3352fa5cd67SKarl Rupp       fine = Y; coarse = X;
336ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser");
337c88c5224SJed Brown     ierr = VecSet(fine,1.);CHKERRQ(ierr);
338c88c5224SJed Brown     ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr);
339c88c5224SJed Brown     ierr = VecDestroy(&fine);CHKERRQ(ierr);
340c88c5224SJed Brown     ierr = VecReciprocal(coarse);CHKERRQ(ierr);
341c88c5224SJed Brown     mglevels[l]->rscale = coarse;
342c88c5224SJed Brown   }
343c88c5224SJed Brown   *rscale = mglevels[l]->rscale;
34473250ac0SBarry Smith   PetscFunctionReturn(0);
34573250ac0SBarry Smith }
34673250ac0SBarry Smith 
347f39d8e23SSatish Balay /*@
348eab52d2bSLawrence Mitchell    PCMGSetInjection - Sets the function to be used to inject primal vectors
349eab52d2bSLawrence Mitchell    from level l to l-1.
350eab52d2bSLawrence Mitchell 
351eab52d2bSLawrence Mitchell    Logically Collective on PC and Mat
352eab52d2bSLawrence Mitchell 
353eab52d2bSLawrence Mitchell    Input Parameters:
354eab52d2bSLawrence Mitchell +  pc - the multigrid context
355eab52d2bSLawrence Mitchell .  l - the level (0 is coarsest) to supply [Do not supply 0]
356eab52d2bSLawrence Mitchell -  mat - the injection matrix
357eab52d2bSLawrence Mitchell 
358eab52d2bSLawrence Mitchell    Level: advanced
359eab52d2bSLawrence Mitchell 
360eab52d2bSLawrence Mitchell .keywords: MG, set, multigrid, restriction, injection, level
361eab52d2bSLawrence Mitchell 
362eab52d2bSLawrence Mitchell .seealso: PCMGSetRestriction()
363eab52d2bSLawrence Mitchell @*/
364eab52d2bSLawrence Mitchell PetscErrorCode  PCMGSetInjection(PC pc,PetscInt l,Mat mat)
365eab52d2bSLawrence Mitchell {
366eab52d2bSLawrence Mitchell   PetscErrorCode ierr;
367eab52d2bSLawrence Mitchell   PC_MG          *mg        = (PC_MG*)pc->data;
368eab52d2bSLawrence Mitchell   PC_MG_Levels   **mglevels = mg->levels;
369eab52d2bSLawrence Mitchell 
370eab52d2bSLawrence Mitchell   PetscFunctionBegin;
371eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
372eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
373eab52d2bSLawrence Mitchell   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
374eab52d2bSLawrence Mitchell   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
375eab52d2bSLawrence Mitchell   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
376eab52d2bSLawrence Mitchell   ierr = MatDestroy(&mglevels[l]->inject);CHKERRQ(ierr);
377eab52d2bSLawrence Mitchell 
378eab52d2bSLawrence Mitchell   mglevels[l]->inject = mat;
379eab52d2bSLawrence Mitchell   PetscFunctionReturn(0);
380eab52d2bSLawrence Mitchell }
381eab52d2bSLawrence Mitchell 
382eab52d2bSLawrence Mitchell /*@
383eab52d2bSLawrence Mitchell    PCMGGetInjection - Gets the function to be used to inject primal vectors
384eab52d2bSLawrence Mitchell    from level l to l-1.
385eab52d2bSLawrence Mitchell 
386eab52d2bSLawrence Mitchell    Logically Collective on PC and Mat
387eab52d2bSLawrence Mitchell 
388eab52d2bSLawrence Mitchell    Input Parameters:
389eab52d2bSLawrence Mitchell +  pc - the multigrid context
390eab52d2bSLawrence Mitchell -  l - the level (0 is coarsest) to supply [Do not supply 0]
391eab52d2bSLawrence Mitchell 
392eab52d2bSLawrence Mitchell    Output Parameter:
39399a38656SLawrence Mitchell .  mat - the restriction matrix (may be NULL if no injection is available).
394eab52d2bSLawrence Mitchell 
395eab52d2bSLawrence Mitchell    Level: advanced
396eab52d2bSLawrence Mitchell 
397eab52d2bSLawrence Mitchell .keywords: MG, get, multigrid, restriction, injection, level
398eab52d2bSLawrence Mitchell 
399eab52d2bSLawrence Mitchell .seealso: PCMGSetInjection(), PCMGetGetRestriction()
400eab52d2bSLawrence Mitchell @*/
401eab52d2bSLawrence Mitchell PetscErrorCode  PCMGGetInjection(PC pc,PetscInt l,Mat *mat)
402eab52d2bSLawrence Mitchell {
403eab52d2bSLawrence Mitchell   PC_MG          *mg        = (PC_MG*)pc->data;
404eab52d2bSLawrence Mitchell   PC_MG_Levels   **mglevels = mg->levels;
405eab52d2bSLawrence Mitchell 
406eab52d2bSLawrence Mitchell   PetscFunctionBegin;
407eab52d2bSLawrence Mitchell   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
408eab52d2bSLawrence Mitchell   if (mat) PetscValidPointer(mat,3);
409eab52d2bSLawrence Mitchell   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
410eab52d2bSLawrence Mitchell   if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1);
411eab52d2bSLawrence Mitchell   if (mat) *mat = mglevels[l]->inject;
412eab52d2bSLawrence Mitchell   PetscFunctionReturn(0);
413eab52d2bSLawrence Mitchell }
414eab52d2bSLawrence Mitchell 
415eab52d2bSLawrence Mitchell /*@
41697177400SBarry Smith    PCMGGetSmoother - Gets the KSP context to be used as smoother for
41797177400SBarry Smith    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
41897177400SBarry Smith    PCMGGetSmootherDown() to use different functions for pre- and
4194b9ad928SBarry Smith    post-smoothing.
4204b9ad928SBarry Smith 
4214b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
4224b9ad928SBarry Smith 
4234b9ad928SBarry Smith    Input Parameters:
4244b9ad928SBarry Smith +  pc - the multigrid context
4254b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith    Ouput Parameters:
4284b9ad928SBarry Smith .  ksp - the smoother
4294b9ad928SBarry Smith 
43057420d5bSBarry Smith    Notes:
43157420d5bSBarry Smith    Once you have called this routine, you can call KSPSetOperators(ksp,...) on the resulting ksp to provide the operators for the smoother for this level.
43257420d5bSBarry Smith    You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc)
43357420d5bSBarry Smith    and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing.
43457420d5bSBarry Smith 
4354b9ad928SBarry Smith    Level: advanced
4364b9ad928SBarry Smith 
4374b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
4384b9ad928SBarry Smith 
439*28668aa5SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve()
4404b9ad928SBarry Smith @*/
4417087cfbeSBarry Smith PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
4424b9ad928SBarry Smith {
443f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
444f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
4454b9ad928SBarry Smith 
4464b9ad928SBarry Smith   PetscFunctionBegin;
447c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
448f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
4494b9ad928SBarry Smith   PetscFunctionReturn(0);
4504b9ad928SBarry Smith }
4514b9ad928SBarry Smith 
452f39d8e23SSatish Balay /*@
45397177400SBarry Smith    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
4544b9ad928SBarry Smith    coarse grid correction (post-smoother).
4554b9ad928SBarry Smith 
4564b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
4574b9ad928SBarry Smith 
4584b9ad928SBarry Smith    Input Parameters:
4594b9ad928SBarry Smith +  pc - the multigrid context
4604b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
4614b9ad928SBarry Smith 
4624b9ad928SBarry Smith    Ouput Parameters:
4634b9ad928SBarry Smith .  ksp - the smoother
4644b9ad928SBarry Smith 
4654b9ad928SBarry Smith    Level: advanced
4664b9ad928SBarry Smith 
46795452b02SPatrick Sanan    Notes:
46895452b02SPatrick Sanan     calling this will result in a different pre and post smoother so you may need to
46989cce641SBarry Smith          set options on the pre smoother also
47089cce641SBarry Smith 
4714b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
4724b9ad928SBarry Smith 
47397177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
4744b9ad928SBarry Smith @*/
4757087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
4764b9ad928SBarry Smith {
477f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
478f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
479dfbe8321SBarry Smith   PetscErrorCode ierr;
480f69a0ea3SMatthew Knepley   const char     *prefix;
4814b9ad928SBarry Smith   MPI_Comm       comm;
4824b9ad928SBarry Smith 
4834b9ad928SBarry Smith   PetscFunctionBegin;
484c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4854b9ad928SBarry Smith   /*
4864b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
4874b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
4884b9ad928SBarry Smith      if not we allocate it.
4894b9ad928SBarry Smith   */
490ce94432eSBarry Smith   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
491f3fbd535SBarry Smith   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
49219fd82e9SBarry Smith     KSPType     ksptype;
49319fd82e9SBarry Smith     PCType      pctype;
494336babb1SJed Brown     PC          ipc;
495336babb1SJed Brown     PetscReal   rtol,abstol,dtol;
496336babb1SJed Brown     PetscInt    maxits;
497336babb1SJed Brown     KSPNormType normtype;
498f3fbd535SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
499f3fbd535SBarry Smith     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
500336babb1SJed Brown     ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr);
5013bf036e2SBarry Smith     ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr);
502336babb1SJed Brown     ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr);
503336babb1SJed Brown     ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr);
504336babb1SJed Brown     ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr);
505336babb1SJed Brown 
506f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
507422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr);
508f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
509f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
510336babb1SJed Brown     ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr);
511336babb1SJed Brown     ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr);
512336babb1SJed Brown     ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr);
5130059c7bdSJed Brown     ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
514336babb1SJed Brown     ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr);
515336babb1SJed Brown     ierr = PCSetType(ipc,pctype);CHKERRQ(ierr);
5163bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr);
517ab83eea4SMatthew G. Knepley     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr);
5184b9ad928SBarry Smith   }
519f3fbd535SBarry Smith   if (ksp) *ksp = mglevels[l]->smoothu;
5204b9ad928SBarry Smith   PetscFunctionReturn(0);
5214b9ad928SBarry Smith }
5224b9ad928SBarry Smith 
523f39d8e23SSatish Balay /*@
52497177400SBarry Smith    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
5254b9ad928SBarry Smith    coarse grid correction (pre-smoother).
5264b9ad928SBarry Smith 
5274b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
5284b9ad928SBarry Smith 
5294b9ad928SBarry Smith    Input Parameters:
5304b9ad928SBarry Smith +  pc - the multigrid context
5314b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
5324b9ad928SBarry Smith 
5334b9ad928SBarry Smith    Ouput Parameters:
5344b9ad928SBarry Smith .  ksp - the smoother
5354b9ad928SBarry Smith 
5364b9ad928SBarry Smith    Level: advanced
5374b9ad928SBarry Smith 
53895452b02SPatrick Sanan    Notes:
53995452b02SPatrick Sanan     calling this will result in a different pre and post smoother so you may need to
54089cce641SBarry Smith          set options on the post smoother also
54189cce641SBarry Smith 
5424b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
5434b9ad928SBarry Smith 
54497177400SBarry Smith .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
5454b9ad928SBarry Smith @*/
5467087cfbeSBarry Smith PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
5474b9ad928SBarry Smith {
548dfbe8321SBarry Smith   PetscErrorCode ierr;
549f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
550f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
5514b9ad928SBarry Smith 
5524b9ad928SBarry Smith   PetscFunctionBegin;
553c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
5544b9ad928SBarry Smith   /* make sure smoother up and down are different */
555c5eb9154SBarry Smith   if (l) {
5560298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr);
557d8148a5aSMatthew Knepley   }
558f3fbd535SBarry Smith   *ksp = mglevels[l]->smoothd;
5594b9ad928SBarry Smith   PetscFunctionReturn(0);
5604b9ad928SBarry Smith }
5614b9ad928SBarry Smith 
5624b9ad928SBarry Smith /*@
563cab31ae5SJed Brown    PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level.
5644b9ad928SBarry Smith 
565ad4df100SBarry Smith    Logically Collective on PC
5664b9ad928SBarry Smith 
5674b9ad928SBarry Smith    Input Parameters:
5684b9ad928SBarry Smith +  pc - the multigrid context
569c1cbb1deSBarry Smith .  l  - the level (0 is coarsest)
570c1cbb1deSBarry Smith -  c  - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
5714b9ad928SBarry Smith 
5724b9ad928SBarry Smith    Level: advanced
5734b9ad928SBarry Smith 
5744b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
5754b9ad928SBarry Smith 
576c1cbb1deSBarry Smith .seealso: PCMGSetCycleType()
5774b9ad928SBarry Smith @*/
578c1cbb1deSBarry Smith PetscErrorCode  PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c)
5794b9ad928SBarry Smith {
580f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
581f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
5824b9ad928SBarry Smith 
5834b9ad928SBarry Smith   PetscFunctionBegin;
584c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
585ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
586c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,l,2);
587c51679f6SSatish Balay   PetscValidLogicalCollectiveEnum(pc,c,3);
588f3fbd535SBarry Smith   mglevels[l]->cycles = c;
5894b9ad928SBarry Smith   PetscFunctionReturn(0);
5904b9ad928SBarry Smith }
5914b9ad928SBarry Smith 
5924b9ad928SBarry Smith /*@
59397177400SBarry Smith    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
594fccaa45eSBarry Smith    on a particular level.
5954b9ad928SBarry Smith 
596ad4df100SBarry Smith    Logically Collective on PC and Vec
5974b9ad928SBarry Smith 
5984b9ad928SBarry Smith    Input Parameters:
5994b9ad928SBarry Smith +  pc - the multigrid context
6004b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
6014b9ad928SBarry Smith -  c  - the space
6024b9ad928SBarry Smith 
6034b9ad928SBarry Smith    Level: advanced
6044b9ad928SBarry Smith 
60595452b02SPatrick Sanan    Notes:
60695452b02SPatrick Sanan     If this is not provided PETSc will automatically generate one.
607fccaa45eSBarry Smith 
608fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
609fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
610fccaa45eSBarry Smith 
6114b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
6124b9ad928SBarry Smith 
61397177400SBarry Smith .seealso: PCMGSetX(), PCMGSetR()
6144b9ad928SBarry Smith @*/
6157087cfbeSBarry Smith PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
6164b9ad928SBarry Smith {
617fccaa45eSBarry Smith   PetscErrorCode ierr;
618f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
619f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
6204b9ad928SBarry Smith 
6214b9ad928SBarry Smith   PetscFunctionBegin;
622c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
623ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
624ce94432eSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
625c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
6266bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
6272fa5cd67SKarl Rupp 
628f3fbd535SBarry Smith   mglevels[l]->b = c;
6294b9ad928SBarry Smith   PetscFunctionReturn(0);
6304b9ad928SBarry Smith }
6314b9ad928SBarry Smith 
6324b9ad928SBarry Smith /*@
63397177400SBarry Smith    PCMGSetX - Sets the vector space to be used to store the solution on a
634fccaa45eSBarry Smith    particular level.
6354b9ad928SBarry Smith 
636ad4df100SBarry Smith    Logically Collective on PC and Vec
6374b9ad928SBarry Smith 
6384b9ad928SBarry Smith    Input Parameters:
6394b9ad928SBarry Smith +  pc - the multigrid context
640251f4c67SDmitry Karpeev .  l - the level (0 is coarsest) this is to be used for (do not supply the finest level)
6414b9ad928SBarry Smith -  c - the space
6424b9ad928SBarry Smith 
6434b9ad928SBarry Smith    Level: advanced
6444b9ad928SBarry Smith 
64595452b02SPatrick Sanan    Notes:
64695452b02SPatrick Sanan     If this is not provided PETSc will automatically generate one.
647fccaa45eSBarry Smith 
648fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
649fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
650fccaa45eSBarry Smith 
6514b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
6524b9ad928SBarry Smith 
65397177400SBarry Smith .seealso: PCMGSetRhs(), PCMGSetR()
6544b9ad928SBarry Smith @*/
6557087cfbeSBarry Smith PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
6564b9ad928SBarry Smith {
657fccaa45eSBarry Smith   PetscErrorCode ierr;
658f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
659f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
6604b9ad928SBarry Smith 
6614b9ad928SBarry Smith   PetscFunctionBegin;
662c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
663ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
664ce94432eSBarry Smith   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
665c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
6666bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
6672fa5cd67SKarl Rupp 
668f3fbd535SBarry Smith   mglevels[l]->x = c;
6694b9ad928SBarry Smith   PetscFunctionReturn(0);
6704b9ad928SBarry Smith }
6714b9ad928SBarry Smith 
6724b9ad928SBarry Smith /*@
67397177400SBarry Smith    PCMGSetR - Sets the vector space to be used to store the residual on a
674fccaa45eSBarry Smith    particular level.
6754b9ad928SBarry Smith 
676ad4df100SBarry Smith    Logically Collective on PC and Vec
6774b9ad928SBarry Smith 
6784b9ad928SBarry Smith    Input Parameters:
6794b9ad928SBarry Smith +  pc - the multigrid context
6804b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
6814b9ad928SBarry Smith -  c - the space
6824b9ad928SBarry Smith 
6834b9ad928SBarry Smith    Level: advanced
6844b9ad928SBarry Smith 
68595452b02SPatrick Sanan    Notes:
68695452b02SPatrick Sanan     If this is not provided PETSc will automatically generate one.
687fccaa45eSBarry Smith 
688fccaa45eSBarry Smith           You do not need to keep a reference to this vector if you do
689fccaa45eSBarry Smith           not need it PCDestroy() will properly free it.
690fccaa45eSBarry Smith 
6914b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
6924b9ad928SBarry Smith @*/
6937087cfbeSBarry Smith PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
6944b9ad928SBarry Smith {
695fccaa45eSBarry Smith   PetscErrorCode ierr;
696f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
697f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
6984b9ad928SBarry Smith 
6994b9ad928SBarry Smith   PetscFunctionBegin;
700c5eb9154SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
701ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
702ce94432eSBarry Smith   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
703c3122656SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
7046bf464f9SBarry Smith   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
7052fa5cd67SKarl Rupp 
706f3fbd535SBarry Smith   mglevels[l]->r = c;
7074b9ad928SBarry Smith   PetscFunctionReturn(0);
7084b9ad928SBarry Smith }
709