xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 4b9ad92859ccb93b5e851e53cb8c4c04ac10e155)
1*4b9ad928SBarry Smith /*$Id: mgfunc.c,v 1.41 2001/08/07 03:03:36 balay Exp $*/
2*4b9ad928SBarry Smith 
3*4b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h"       /*I "petscksp.h" I*/
4*4b9ad928SBarry Smith                           /*I "petscmg.h"   I*/
5*4b9ad928SBarry Smith 
6*4b9ad928SBarry Smith #undef __FUNCT__
7*4b9ad928SBarry Smith #define __FUNCT__ "MGDefaultResidual"
8*4b9ad928SBarry Smith /*@C
9*4b9ad928SBarry Smith    MGDefaultResidual - Default routine to calculate the residual.
10*4b9ad928SBarry Smith 
11*4b9ad928SBarry Smith    Collective on Mat and Vec
12*4b9ad928SBarry Smith 
13*4b9ad928SBarry Smith    Input Parameters:
14*4b9ad928SBarry Smith +  mat - the matrix
15*4b9ad928SBarry Smith .  b   - the right-hand-side
16*4b9ad928SBarry Smith -  x   - the approximate solution
17*4b9ad928SBarry Smith 
18*4b9ad928SBarry Smith    Output Parameter:
19*4b9ad928SBarry Smith .  r - location to store the residual
20*4b9ad928SBarry Smith 
21*4b9ad928SBarry Smith    Level: advanced
22*4b9ad928SBarry Smith 
23*4b9ad928SBarry Smith .keywords: MG, default, multigrid, residual
24*4b9ad928SBarry Smith 
25*4b9ad928SBarry Smith .seealso: MGSetResidual()
26*4b9ad928SBarry Smith @*/
27*4b9ad928SBarry Smith int MGDefaultResidual(Mat mat,Vec b,Vec x,Vec r)
28*4b9ad928SBarry Smith {
29*4b9ad928SBarry Smith   int    ierr;
30*4b9ad928SBarry Smith   PetscScalar mone = -1.0;
31*4b9ad928SBarry Smith 
32*4b9ad928SBarry Smith   PetscFunctionBegin;
33*4b9ad928SBarry Smith   ierr = MatMult(mat,x,r);CHKERRQ(ierr);
34*4b9ad928SBarry Smith   ierr = VecAYPX(&mone,b,r);CHKERRQ(ierr);
35*4b9ad928SBarry Smith   PetscFunctionReturn(0);
36*4b9ad928SBarry Smith }
37*4b9ad928SBarry Smith 
38*4b9ad928SBarry Smith /* ---------------------------------------------------------------------------*/
39*4b9ad928SBarry Smith 
40*4b9ad928SBarry Smith #undef __FUNCT__
41*4b9ad928SBarry Smith #define __FUNCT__ "MGGetCoarseSolve"
42*4b9ad928SBarry Smith /*@C
43*4b9ad928SBarry Smith    MGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
44*4b9ad928SBarry Smith 
45*4b9ad928SBarry Smith    Not Collective
46*4b9ad928SBarry Smith 
47*4b9ad928SBarry Smith    Input Parameter:
48*4b9ad928SBarry Smith .  pc - the multigrid context
49*4b9ad928SBarry Smith 
50*4b9ad928SBarry Smith    Output Parameter:
51*4b9ad928SBarry Smith .  ksp - the coarse grid solver context
52*4b9ad928SBarry Smith 
53*4b9ad928SBarry Smith    Level: advanced
54*4b9ad928SBarry Smith 
55*4b9ad928SBarry Smith .keywords: MG, multigrid, get, coarse grid
56*4b9ad928SBarry Smith @*/
57*4b9ad928SBarry Smith int MGGetCoarseSolve(PC pc,KSP *ksp)
58*4b9ad928SBarry Smith {
59*4b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
60*4b9ad928SBarry Smith 
61*4b9ad928SBarry Smith   PetscFunctionBegin;
62*4b9ad928SBarry Smith   *ksp =  mg[0]->smoothd;
63*4b9ad928SBarry Smith   PetscFunctionReturn(0);
64*4b9ad928SBarry Smith }
65*4b9ad928SBarry Smith 
66*4b9ad928SBarry Smith #undef __FUNCT__
67*4b9ad928SBarry Smith #define __FUNCT__ "MGSetResidual"
68*4b9ad928SBarry Smith /*@C
69*4b9ad928SBarry Smith    MGSetResidual - Sets the function to be used to calculate the residual
70*4b9ad928SBarry Smith    on the lth level.
71*4b9ad928SBarry Smith 
72*4b9ad928SBarry Smith    Collective on PC and Mat
73*4b9ad928SBarry Smith 
74*4b9ad928SBarry Smith    Input Parameters:
75*4b9ad928SBarry Smith +  pc       - the multigrid context
76*4b9ad928SBarry Smith .  l        - the level (0 is coarsest) to supply
77*4b9ad928SBarry Smith .  residual - function used to form residual (usually MGDefaultResidual)
78*4b9ad928SBarry Smith -  mat      - matrix associated with residual
79*4b9ad928SBarry Smith 
80*4b9ad928SBarry Smith    Level: advanced
81*4b9ad928SBarry Smith 
82*4b9ad928SBarry Smith .keywords:  MG, set, multigrid, residual, level
83*4b9ad928SBarry Smith 
84*4b9ad928SBarry Smith .seealso: MGDefaultResidual()
85*4b9ad928SBarry Smith @*/
86*4b9ad928SBarry Smith int MGSetResidual(PC pc,int l,int (*residual)(Mat,Vec,Vec,Vec),Mat mat)
87*4b9ad928SBarry Smith {
88*4b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
89*4b9ad928SBarry Smith 
90*4b9ad928SBarry Smith   PetscFunctionBegin;
91*4b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
92*4b9ad928SBarry Smith 
93*4b9ad928SBarry Smith   mg[l]->residual = residual;
94*4b9ad928SBarry Smith   mg[l]->A        = mat;
95*4b9ad928SBarry Smith   PetscFunctionReturn(0);
96*4b9ad928SBarry Smith }
97*4b9ad928SBarry Smith 
98*4b9ad928SBarry Smith #undef __FUNCT__
99*4b9ad928SBarry Smith #define __FUNCT__ "MGSetInterpolate"
100*4b9ad928SBarry Smith /*@
101*4b9ad928SBarry Smith    MGSetInterpolate - Sets the function to be used to calculate the
102*4b9ad928SBarry Smith    interpolation on the lth level.
103*4b9ad928SBarry Smith 
104*4b9ad928SBarry Smith    Collective on PC and Mat
105*4b9ad928SBarry Smith 
106*4b9ad928SBarry Smith    Input Parameters:
107*4b9ad928SBarry Smith +  pc  - the multigrid context
108*4b9ad928SBarry Smith .  mat - the interpolation operator
109*4b9ad928SBarry Smith -  l   - the level (0 is coarsest) to supply
110*4b9ad928SBarry Smith 
111*4b9ad928SBarry Smith    Level: advanced
112*4b9ad928SBarry Smith 
113*4b9ad928SBarry Smith    Notes:
114*4b9ad928SBarry Smith           Usually this is the same matrix used also to set the restriction
115*4b9ad928SBarry Smith     for the same level.
116*4b9ad928SBarry Smith 
117*4b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
118*4b9ad928SBarry Smith     out from the matrix size which one it is.
119*4b9ad928SBarry Smith 
120*4b9ad928SBarry Smith .keywords:  multigrid, set, interpolate, level
121*4b9ad928SBarry Smith 
122*4b9ad928SBarry Smith .seealso: MGSetRestriction()
123*4b9ad928SBarry Smith @*/
124*4b9ad928SBarry Smith int MGSetInterpolate(PC pc,int l,Mat mat)
125*4b9ad928SBarry Smith {
126*4b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
127*4b9ad928SBarry Smith 
128*4b9ad928SBarry Smith   PetscFunctionBegin;
129*4b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
130*4b9ad928SBarry Smith   mg[l]->interpolate = mat;
131*4b9ad928SBarry Smith   PetscFunctionReturn(0);
132*4b9ad928SBarry Smith }
133*4b9ad928SBarry Smith 
134*4b9ad928SBarry Smith #undef __FUNCT__
135*4b9ad928SBarry Smith #define __FUNCT__ "MGSetRestriction"
136*4b9ad928SBarry Smith /*@
137*4b9ad928SBarry Smith    MGSetRestriction - Sets the function to be used to restrict vector
138*4b9ad928SBarry Smith    from level l to l-1.
139*4b9ad928SBarry Smith 
140*4b9ad928SBarry Smith    Collective on PC and Mat
141*4b9ad928SBarry Smith 
142*4b9ad928SBarry Smith    Input Parameters:
143*4b9ad928SBarry Smith +  pc - the multigrid context
144*4b9ad928SBarry Smith .  mat - the restriction matrix
145*4b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
146*4b9ad928SBarry Smith 
147*4b9ad928SBarry Smith    Level: advanced
148*4b9ad928SBarry Smith 
149*4b9ad928SBarry Smith    Notes:
150*4b9ad928SBarry Smith           Usually this is the same matrix used also to set the interpolation
151*4b9ad928SBarry Smith     for the same level.
152*4b9ad928SBarry Smith 
153*4b9ad928SBarry Smith           One can pass in the interpolation matrix or its transpose; PETSc figures
154*4b9ad928SBarry Smith     out from the matrix size which one it is.
155*4b9ad928SBarry Smith 
156*4b9ad928SBarry Smith .keywords: MG, set, multigrid, restriction, level
157*4b9ad928SBarry Smith 
158*4b9ad928SBarry Smith .seealso: MGSetInterpolate()
159*4b9ad928SBarry Smith @*/
160*4b9ad928SBarry Smith int MGSetRestriction(PC pc,int l,Mat mat)
161*4b9ad928SBarry Smith {
162*4b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
163*4b9ad928SBarry Smith 
164*4b9ad928SBarry Smith   PetscFunctionBegin;
165*4b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
166*4b9ad928SBarry Smith   mg[l]->restrct  = mat;
167*4b9ad928SBarry Smith   PetscFunctionReturn(0);
168*4b9ad928SBarry Smith }
169*4b9ad928SBarry Smith 
170*4b9ad928SBarry Smith #undef __FUNCT__
171*4b9ad928SBarry Smith #define __FUNCT__ "MGGetSmoother"
172*4b9ad928SBarry Smith /*@C
173*4b9ad928SBarry Smith    MGGetSmoother - Gets the KSP context to be used as smoother for
174*4b9ad928SBarry Smith    both pre- and post-smoothing.  Call both MGGetSmootherUp() and
175*4b9ad928SBarry Smith    MGGetSmootherDown() to use different functions for pre- and
176*4b9ad928SBarry Smith    post-smoothing.
177*4b9ad928SBarry Smith 
178*4b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
179*4b9ad928SBarry Smith 
180*4b9ad928SBarry Smith    Input Parameters:
181*4b9ad928SBarry Smith +  pc - the multigrid context
182*4b9ad928SBarry Smith -  l - the level (0 is coarsest) to supply
183*4b9ad928SBarry Smith 
184*4b9ad928SBarry Smith    Ouput Parameters:
185*4b9ad928SBarry Smith .  ksp - the smoother
186*4b9ad928SBarry Smith 
187*4b9ad928SBarry Smith    Level: advanced
188*4b9ad928SBarry Smith 
189*4b9ad928SBarry Smith .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
190*4b9ad928SBarry Smith 
191*4b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmootherDown()
192*4b9ad928SBarry Smith @*/
193*4b9ad928SBarry Smith int MGGetSmoother(PC pc,int l,KSP *ksp)
194*4b9ad928SBarry Smith {
195*4b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
196*4b9ad928SBarry Smith 
197*4b9ad928SBarry Smith   PetscFunctionBegin;
198*4b9ad928SBarry Smith   *ksp = mg[l]->smoothd;
199*4b9ad928SBarry Smith   PetscFunctionReturn(0);
200*4b9ad928SBarry Smith }
201*4b9ad928SBarry Smith 
202*4b9ad928SBarry Smith #undef __FUNCT__
203*4b9ad928SBarry Smith #define __FUNCT__ "MGGetSmootherUp"
204*4b9ad928SBarry Smith /*@C
205*4b9ad928SBarry Smith    MGGetSmootherUp - Gets the KSP context to be used as smoother after
206*4b9ad928SBarry Smith    coarse grid correction (post-smoother).
207*4b9ad928SBarry Smith 
208*4b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
209*4b9ad928SBarry Smith 
210*4b9ad928SBarry Smith    Input Parameters:
211*4b9ad928SBarry Smith +  pc - the multigrid context
212*4b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
213*4b9ad928SBarry Smith 
214*4b9ad928SBarry Smith    Ouput Parameters:
215*4b9ad928SBarry Smith .  ksp - the smoother
216*4b9ad928SBarry Smith 
217*4b9ad928SBarry Smith    Level: advanced
218*4b9ad928SBarry Smith 
219*4b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, up, post-smoother, level
220*4b9ad928SBarry Smith 
221*4b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmootherDown()
222*4b9ad928SBarry Smith @*/
223*4b9ad928SBarry Smith int MGGetSmootherUp(PC pc,int l,KSP *ksp)
224*4b9ad928SBarry Smith {
225*4b9ad928SBarry Smith   MG       *mg = (MG*)pc->data;
226*4b9ad928SBarry Smith   int      ierr;
227*4b9ad928SBarry Smith   char     *prefix;
228*4b9ad928SBarry Smith   MPI_Comm comm;
229*4b9ad928SBarry Smith 
230*4b9ad928SBarry Smith   PetscFunctionBegin;
231*4b9ad928SBarry Smith   /*
232*4b9ad928SBarry Smith      This is called only if user wants a different pre-smoother from post.
233*4b9ad928SBarry Smith      Thus we check if a different one has already been allocated,
234*4b9ad928SBarry Smith      if not we allocate it.
235*4b9ad928SBarry Smith   */
236*4b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
237*4b9ad928SBarry Smith 
238*4b9ad928SBarry Smith   if (mg[l]->smoothu == mg[l]->smoothd) {
239*4b9ad928SBarry Smith     ierr = PetscObjectGetComm((PetscObject)mg[l]->smoothd,&comm);CHKERRQ(ierr);
240*4b9ad928SBarry Smith     ierr = KSPCreate(comm,&mg[l]->smoothu);CHKERRQ(ierr);
241*4b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[l]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
242*4b9ad928SBarry Smith     ierr = KSPSetOptionsPrefix(mg[l]->smoothu,prefix);CHKERRQ(ierr);
243*4b9ad928SBarry Smith     ierr = KSPAppendOptionsPrefix(mg[l]->smoothd,"mg_levels_");CHKERRQ(ierr);
244*4b9ad928SBarry Smith     PetscLogObjectParent(pc,mg[l]->smoothu);
245*4b9ad928SBarry Smith   }
246*4b9ad928SBarry Smith   if (ksp) *ksp = mg[l]->smoothu;
247*4b9ad928SBarry Smith   PetscFunctionReturn(0);
248*4b9ad928SBarry Smith }
249*4b9ad928SBarry Smith 
250*4b9ad928SBarry Smith #undef __FUNCT__
251*4b9ad928SBarry Smith #define __FUNCT__ "MGGetSmootherDown"
252*4b9ad928SBarry Smith /*@C
253*4b9ad928SBarry Smith    MGGetSmootherDown - Gets the KSP context to be used as smoother before
254*4b9ad928SBarry Smith    coarse grid correction (pre-smoother).
255*4b9ad928SBarry Smith 
256*4b9ad928SBarry Smith    Not Collective, KSP returned is parallel if PC is
257*4b9ad928SBarry Smith 
258*4b9ad928SBarry Smith    Input Parameters:
259*4b9ad928SBarry Smith +  pc - the multigrid context
260*4b9ad928SBarry Smith -  l  - the level (0 is coarsest) to supply
261*4b9ad928SBarry Smith 
262*4b9ad928SBarry Smith    Ouput Parameters:
263*4b9ad928SBarry Smith .  ksp - the smoother
264*4b9ad928SBarry Smith 
265*4b9ad928SBarry Smith    Level: advanced
266*4b9ad928SBarry Smith 
267*4b9ad928SBarry Smith .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
268*4b9ad928SBarry Smith 
269*4b9ad928SBarry Smith .seealso: MGGetSmootherUp(), MGGetSmoother()
270*4b9ad928SBarry Smith @*/
271*4b9ad928SBarry Smith int MGGetSmootherDown(PC pc,int l,KSP *ksp)
272*4b9ad928SBarry Smith {
273*4b9ad928SBarry Smith   int ierr;
274*4b9ad928SBarry Smith   MG  *mg = (MG*)pc->data;
275*4b9ad928SBarry Smith 
276*4b9ad928SBarry Smith   PetscFunctionBegin;
277*4b9ad928SBarry Smith   /* make sure smoother up and down are different */
278*4b9ad928SBarry Smith   ierr = MGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
279*4b9ad928SBarry Smith   *ksp = mg[l]->smoothd;
280*4b9ad928SBarry Smith   PetscFunctionReturn(0);
281*4b9ad928SBarry Smith }
282*4b9ad928SBarry Smith 
283*4b9ad928SBarry Smith #undef __FUNCT__
284*4b9ad928SBarry Smith #define __FUNCT__ "MGSetCyclesOnLevel"
285*4b9ad928SBarry Smith /*@
286*4b9ad928SBarry Smith    MGSetCyclesOnLevel - Sets the number of cycles to run on this level.
287*4b9ad928SBarry Smith 
288*4b9ad928SBarry Smith    Collective on PC
289*4b9ad928SBarry Smith 
290*4b9ad928SBarry Smith    Input Parameters:
291*4b9ad928SBarry Smith +  pc - the multigrid context
292*4b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
293*4b9ad928SBarry Smith -  n  - the number of cycles
294*4b9ad928SBarry Smith 
295*4b9ad928SBarry Smith    Level: advanced
296*4b9ad928SBarry Smith 
297*4b9ad928SBarry Smith .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
298*4b9ad928SBarry Smith 
299*4b9ad928SBarry Smith .seealso: MGSetCycles()
300*4b9ad928SBarry Smith @*/
301*4b9ad928SBarry Smith int MGSetCyclesOnLevel(PC pc,int l,int c)
302*4b9ad928SBarry Smith {
303*4b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
304*4b9ad928SBarry Smith 
305*4b9ad928SBarry Smith   PetscFunctionBegin;
306*4b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
307*4b9ad928SBarry Smith   mg[l]->cycles  = c;
308*4b9ad928SBarry Smith   PetscFunctionReturn(0);
309*4b9ad928SBarry Smith }
310*4b9ad928SBarry Smith 
311*4b9ad928SBarry Smith #undef __FUNCT__
312*4b9ad928SBarry Smith #define __FUNCT__ "MGSetRhs"
313*4b9ad928SBarry Smith /*@
314*4b9ad928SBarry Smith    MGSetRhs - Sets the vector space to be used to store the right-hand side
315*4b9ad928SBarry Smith    on a particular level.  The user should free this space at the conclusion
316*4b9ad928SBarry Smith    of multigrid use.
317*4b9ad928SBarry Smith 
318*4b9ad928SBarry Smith    Collective on PC and Vec
319*4b9ad928SBarry Smith 
320*4b9ad928SBarry Smith    Input Parameters:
321*4b9ad928SBarry Smith +  pc - the multigrid context
322*4b9ad928SBarry Smith .  l  - the level (0 is coarsest) this is to be used for
323*4b9ad928SBarry Smith -  c  - the space
324*4b9ad928SBarry Smith 
325*4b9ad928SBarry Smith    Level: advanced
326*4b9ad928SBarry Smith 
327*4b9ad928SBarry Smith .keywords: MG, multigrid, set, right-hand-side, rhs, level
328*4b9ad928SBarry Smith 
329*4b9ad928SBarry Smith .seealso: MGSetX(), MGSetR()
330*4b9ad928SBarry Smith @*/
331*4b9ad928SBarry Smith int MGSetRhs(PC pc,int l,Vec c)
332*4b9ad928SBarry Smith {
333*4b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
334*4b9ad928SBarry Smith 
335*4b9ad928SBarry Smith   PetscFunctionBegin;
336*4b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
337*4b9ad928SBarry Smith   mg[l]->b  = c;
338*4b9ad928SBarry Smith   PetscFunctionReturn(0);
339*4b9ad928SBarry Smith }
340*4b9ad928SBarry Smith 
341*4b9ad928SBarry Smith #undef __FUNCT__
342*4b9ad928SBarry Smith #define __FUNCT__ "MGSetX"
343*4b9ad928SBarry Smith /*@
344*4b9ad928SBarry Smith    MGSetX - Sets the vector space to be used to store the solution on a
345*4b9ad928SBarry Smith    particular level.  The user should free this space at the conclusion
346*4b9ad928SBarry Smith    of multigrid use.
347*4b9ad928SBarry Smith 
348*4b9ad928SBarry Smith    Collective on PC and Vec
349*4b9ad928SBarry Smith 
350*4b9ad928SBarry Smith    Input Parameters:
351*4b9ad928SBarry Smith +  pc - the multigrid context
352*4b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
353*4b9ad928SBarry Smith -  c - the space
354*4b9ad928SBarry Smith 
355*4b9ad928SBarry Smith    Level: advanced
356*4b9ad928SBarry Smith 
357*4b9ad928SBarry Smith .keywords: MG, multigrid, set, solution, level
358*4b9ad928SBarry Smith 
359*4b9ad928SBarry Smith .seealso: MGSetRhs(), MGSetR()
360*4b9ad928SBarry Smith @*/
361*4b9ad928SBarry Smith int MGSetX(PC pc,int l,Vec c)
362*4b9ad928SBarry Smith {
363*4b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
364*4b9ad928SBarry Smith 
365*4b9ad928SBarry Smith   PetscFunctionBegin;
366*4b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
367*4b9ad928SBarry Smith   mg[l]->x  = c;
368*4b9ad928SBarry Smith   PetscFunctionReturn(0);
369*4b9ad928SBarry Smith }
370*4b9ad928SBarry Smith 
371*4b9ad928SBarry Smith #undef __FUNCT__
372*4b9ad928SBarry Smith #define __FUNCT__ "MGSetR"
373*4b9ad928SBarry Smith /*@
374*4b9ad928SBarry Smith    MGSetR - Sets the vector space to be used to store the residual on a
375*4b9ad928SBarry Smith    particular level.  The user should free this space at the conclusion of
376*4b9ad928SBarry Smith    multigrid use.
377*4b9ad928SBarry Smith 
378*4b9ad928SBarry Smith    Collective on PC and Vec
379*4b9ad928SBarry Smith 
380*4b9ad928SBarry Smith    Input Parameters:
381*4b9ad928SBarry Smith +  pc - the multigrid context
382*4b9ad928SBarry Smith .  l - the level (0 is coarsest) this is to be used for
383*4b9ad928SBarry Smith -  c - the space
384*4b9ad928SBarry Smith 
385*4b9ad928SBarry Smith    Level: advanced
386*4b9ad928SBarry Smith 
387*4b9ad928SBarry Smith .keywords: MG, multigrid, set, residual, level
388*4b9ad928SBarry Smith @*/
389*4b9ad928SBarry Smith int MGSetR(PC pc,int l,Vec c)
390*4b9ad928SBarry Smith {
391*4b9ad928SBarry Smith   MG *mg = (MG*)pc->data;
392*4b9ad928SBarry Smith 
393*4b9ad928SBarry Smith   PetscFunctionBegin;
394*4b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
395*4b9ad928SBarry Smith   mg[l]->r  = c;
396*4b9ad928SBarry Smith   PetscFunctionReturn(0);
397*4b9ad928SBarry Smith }
398*4b9ad928SBarry Smith 
399*4b9ad928SBarry Smith 
400*4b9ad928SBarry Smith 
401*4b9ad928SBarry Smith 
402*4b9ad928SBarry Smith 
403*4b9ad928SBarry Smith 
404*4b9ad928SBarry Smith 
405*4b9ad928SBarry Smith 
406