xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision b41ce5d507ea9a58bfa83cf403107a702e77a67d)
1 
2 #include <petsc/private/pcmgimpl.h>       /*I "petscksp.h" I*/
3 
4 /* ---------------------------------------------------------------------------*/
5 /*@C
6    PCMGResidualDefault - Default routine to calculate the residual.
7 
8    Collective on Mat and Vec
9 
10    Input Parameters:
11 +  mat - the matrix
12 .  b   - the right-hand-side
13 -  x   - the approximate solution
14 
15    Output Parameter:
16 .  r - location to store the residual
17 
18    Level: developer
19 
20 .keywords: MG, default, multigrid, residual
21 
22 .seealso: PCMGSetResidual()
23 @*/
24 PetscErrorCode  PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r)
25 {
26   PetscErrorCode ierr;
27 
28   PetscFunctionBegin;
29   ierr = MatResidual(mat,b,x,r);CHKERRQ(ierr);
30   PetscFunctionReturn(0);
31 }
32 
33 /*@
34    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
35 
36    Not Collective
37 
38    Input Parameter:
39 .  pc - the multigrid context
40 
41    Output Parameter:
42 .  ksp - the coarse grid solver context
43 
44    Level: advanced
45 
46 .keywords: MG, multigrid, get, coarse grid
47 @*/
48 PetscErrorCode  PCMGGetCoarseSolve(PC pc,KSP *ksp)
49 {
50   PC_MG        *mg        = (PC_MG*)pc->data;
51   PC_MG_Levels **mglevels = mg->levels;
52 
53   PetscFunctionBegin;
54   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
55   *ksp =  mglevels[0]->smoothd;
56   PetscFunctionReturn(0);
57 }
58 
59 /*@C
60    PCMGSetResidual - Sets the function to be used to calculate the residual
61    on the lth level.
62 
63    Logically Collective on PC and Mat
64 
65    Input Parameters:
66 +  pc       - the multigrid context
67 .  l        - the level (0 is coarsest) to supply
68 .  residual - function used to form residual, if none is provided the previously provide one is used, if no
69               previous one were provided then a default is used
70 -  mat      - matrix associated with residual
71 
72    Level: advanced
73 
74 .keywords:  MG, set, multigrid, residual, level
75 
76 .seealso: PCMGResidualDefault()
77 @*/
78 PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
79 {
80   PC_MG          *mg        = (PC_MG*)pc->data;
81   PC_MG_Levels   **mglevels = mg->levels;
82   PetscErrorCode ierr;
83 
84   PetscFunctionBegin;
85   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
86   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
87   if (residual) mglevels[l]->residual = residual;
88   if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault;
89   if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);}
90   ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr);
91   mglevels[l]->A = mat;
92   PetscFunctionReturn(0);
93 }
94 
95 /*@
96    PCMGSetInterpolation - Sets the function to be used to calculate the
97    interpolation from l-1 to the lth level
98 
99    Logically Collective on PC and Mat
100 
101    Input Parameters:
102 +  pc  - the multigrid context
103 .  mat - the interpolation operator
104 -  l   - the level (0 is coarsest) to supply [do not supply 0]
105 
106    Level: advanced
107 
108    Notes:
109           Usually this is the same matrix used also to set the restriction
110     for the same level.
111 
112           One can pass in the interpolation matrix or its transpose; PETSc figures
113     out from the matrix size which one it is.
114 
115 .keywords:  multigrid, set, interpolate, level
116 
117 .seealso: PCMGSetRestriction()
118 @*/
119 PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
120 {
121   PC_MG          *mg        = (PC_MG*)pc->data;
122   PC_MG_Levels   **mglevels = mg->levels;
123   PetscErrorCode ierr;
124 
125   PetscFunctionBegin;
126   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
127   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
128   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
129   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
130   ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr);
131 
132   mglevels[l]->interpolate = mat;
133   PetscFunctionReturn(0);
134 }
135 
136 /*@
137    PCMGGetInterpolation - Gets the function to be used to calculate the
138    interpolation from l-1 to the lth level
139 
140    Logically Collective on PC
141 
142    Input Parameters:
143 +  pc - the multigrid context
144 -  l - the level (0 is coarsest) to supply [Do not supply 0]
145 
146    Output Parameter:
147 .  mat - the interpolation matrix, can be NULL
148 
149    Level: advanced
150 
151 .keywords: MG, get, multigrid, interpolation, level
152 
153 .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale()
154 @*/
155 PetscErrorCode  PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat)
156 {
157   PC_MG          *mg        = (PC_MG*)pc->data;
158   PC_MG_Levels   **mglevels = mg->levels;
159   PetscErrorCode ierr;
160 
161   PetscFunctionBegin;
162   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
163   if (mat) PetscValidPointer(mat,3);
164   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
165   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);
166   if (!mglevels[l]->interpolate) {
167     if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()");
168     ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr);
169   }
170   if (mat) *mat = mglevels[l]->interpolate;
171   PetscFunctionReturn(0);
172 }
173 
174 /*@
175    PCMGSetRestriction - Sets the function to be used to restrict vector
176    from level l to l-1.
177 
178    Logically Collective on PC and Mat
179 
180    Input Parameters:
181 +  pc - the multigrid context
182 .  l - the level (0 is coarsest) to supply [Do not supply 0]
183 -  mat - the restriction matrix
184 
185    Level: advanced
186 
187    Notes:
188           Usually this is the same matrix used also to set the interpolation
189     for the same level.
190 
191           One can pass in the interpolation matrix or its transpose; PETSc figures
192     out from the matrix size which one it is.
193 
194          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
195     is used.
196 
197 .keywords: MG, set, multigrid, restriction, level
198 
199 .seealso: PCMGSetInterpolation()
200 @*/
201 PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
202 {
203   PetscErrorCode ierr;
204   PC_MG          *mg        = (PC_MG*)pc->data;
205   PC_MG_Levels   **mglevels = mg->levels;
206 
207   PetscFunctionBegin;
208   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
209   PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
210   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
211   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
212   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
213   ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr);
214 
215   mglevels[l]->restrct = mat;
216   PetscFunctionReturn(0);
217 }
218 
219 /*@
220    PCMGGetRestriction - Gets the function to be used to restrict vector
221    from level l to l-1.
222 
223    Logically Collective on PC and Mat
224 
225    Input Parameters:
226 +  pc - the multigrid context
227 -  l - the level (0 is coarsest) to supply [Do not supply 0]
228 
229    Output Parameter:
230 .  mat - the restriction matrix
231 
232    Level: advanced
233 
234 .keywords: MG, get, multigrid, restriction, level
235 
236 .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale()
237 @*/
238 PetscErrorCode  PCMGGetRestriction(PC pc,PetscInt l,Mat *mat)
239 {
240   PC_MG          *mg        = (PC_MG*)pc->data;
241   PC_MG_Levels   **mglevels = mg->levels;
242   PetscErrorCode ierr;
243 
244   PetscFunctionBegin;
245   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
246   if (mat) PetscValidPointer(mat,3);
247   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
248   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);
249   if (!mglevels[l]->restrct) {
250     if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()");
251     ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr);
252   }
253   if (mat) *mat = mglevels[l]->restrct;
254   PetscFunctionReturn(0);
255 }
256 
257 /*@
258    PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.
259 
260    Logically Collective on PC and Vec
261 
262    Input Parameters:
263 +  pc - the multigrid context
264 -  l - the level (0 is coarsest) to supply [Do not supply 0]
265 .  rscale - the scaling
266 
267    Level: advanced
268 
269    Notes:
270        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.
271 
272 .keywords: MG, set, multigrid, restriction, level
273 
274 .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale()
275 @*/
276 PetscErrorCode  PCMGSetRScale(PC pc,PetscInt l,Vec rscale)
277 {
278   PetscErrorCode ierr;
279   PC_MG          *mg        = (PC_MG*)pc->data;
280   PC_MG_Levels   **mglevels = mg->levels;
281 
282   PetscFunctionBegin;
283   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
284   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
285   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);
286   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
287   ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr);
288 
289   mglevels[l]->rscale = rscale;
290   PetscFunctionReturn(0);
291 }
292 
293 /*@
294    PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1.
295 
296    Collective on PC
297 
298    Input Parameters:
299 +  pc - the multigrid context
300 .  rscale - the scaling
301 -  l - the level (0 is coarsest) to supply [Do not supply 0]
302 
303    Level: advanced
304 
305    Notes:
306        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.
307 
308 .keywords: MG, set, multigrid, restriction, level
309 
310 .seealso: PCMGSetInterpolation(), PCMGGetRestriction()
311 @*/
312 PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale)
313 {
314   PetscErrorCode ierr;
315   PC_MG          *mg        = (PC_MG*)pc->data;
316   PC_MG_Levels   **mglevels = mg->levels;
317 
318   PetscFunctionBegin;
319   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
320   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
321   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);
322   if (!mglevels[l]->rscale) {
323     Mat      R;
324     Vec      X,Y,coarse,fine;
325     PetscInt M,N;
326     ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr);
327     ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr);
328     ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr);
329     if (M < N) {
330       fine = X;
331       coarse = Y;
332     } else if (N < M) {
333       fine = Y; coarse = X;
334     } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser");
335     ierr = VecSet(fine,1.);CHKERRQ(ierr);
336     ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr);
337     ierr = VecDestroy(&fine);CHKERRQ(ierr);
338     ierr = VecReciprocal(coarse);CHKERRQ(ierr);
339     mglevels[l]->rscale = coarse;
340   }
341   *rscale = mglevels[l]->rscale;
342   PetscFunctionReturn(0);
343 }
344 
345 /*@
346    PCMGGetSmoother - Gets the KSP context to be used as smoother for
347    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
348    PCMGGetSmootherDown() to use different functions for pre- and
349    post-smoothing.
350 
351    Not Collective, KSP returned is parallel if PC is
352 
353    Input Parameters:
354 +  pc - the multigrid context
355 -  l - the level (0 is coarsest) to supply
356 
357    Ouput Parameters:
358 .  ksp - the smoother
359 
360    Notes:
361    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.
362    You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc)
363    and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing.
364 
365    Level: advanced
366 
367 .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
368 
369 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
370 @*/
371 PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
372 {
373   PC_MG        *mg        = (PC_MG*)pc->data;
374   PC_MG_Levels **mglevels = mg->levels;
375 
376   PetscFunctionBegin;
377   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
378   *ksp = mglevels[l]->smoothd;
379   PetscFunctionReturn(0);
380 }
381 
382 /*@
383    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
384    coarse grid correction (post-smoother).
385 
386    Not Collective, KSP returned is parallel if PC is
387 
388    Input Parameters:
389 +  pc - the multigrid context
390 -  l  - the level (0 is coarsest) to supply
391 
392    Ouput Parameters:
393 .  ksp - the smoother
394 
395    Level: advanced
396 
397    Notes: calling this will result in a different pre and post smoother so you may need to
398          set options on the pre smoother also
399 
400 .keywords: MG, multigrid, get, smoother, up, post-smoother, level
401 
402 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
403 @*/
404 PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
405 {
406   PC_MG          *mg        = (PC_MG*)pc->data;
407   PC_MG_Levels   **mglevels = mg->levels;
408   PetscErrorCode ierr;
409   const char     *prefix;
410   MPI_Comm       comm;
411 
412   PetscFunctionBegin;
413   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
414   /*
415      This is called only if user wants a different pre-smoother from post.
416      Thus we check if a different one has already been allocated,
417      if not we allocate it.
418   */
419   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
420   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
421     KSPType     ksptype;
422     PCType      pctype;
423     PC          ipc;
424     PetscReal   rtol,abstol,dtol;
425     PetscInt    maxits;
426     KSPNormType normtype;
427     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
428     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
429     ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr);
430     ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr);
431     ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr);
432     ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr);
433     ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr);
434 
435     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
436     ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr);
437     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
438     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
439     ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr);
440     ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr);
441     ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr);
442     ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
443     ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr);
444     ierr = PCSetType(ipc,pctype);CHKERRQ(ierr);
445     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr);
446     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr);
447   }
448   if (ksp) *ksp = mglevels[l]->smoothu;
449   PetscFunctionReturn(0);
450 }
451 
452 /*@
453    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
454    coarse grid correction (pre-smoother).
455 
456    Not Collective, KSP returned is parallel if PC is
457 
458    Input Parameters:
459 +  pc - the multigrid context
460 -  l  - the level (0 is coarsest) to supply
461 
462    Ouput Parameters:
463 .  ksp - the smoother
464 
465    Level: advanced
466 
467    Notes: calling this will result in a different pre and post smoother so you may need to
468          set options on the post smoother also
469 
470 .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
471 
472 .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
473 @*/
474 PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
475 {
476   PetscErrorCode ierr;
477   PC_MG          *mg        = (PC_MG*)pc->data;
478   PC_MG_Levels   **mglevels = mg->levels;
479 
480   PetscFunctionBegin;
481   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
482   /* make sure smoother up and down are different */
483   if (l) {
484     ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr);
485   }
486   *ksp = mglevels[l]->smoothd;
487   PetscFunctionReturn(0);
488 }
489 
490 /*@
491    PCMGSetCycleTypeOnLevel - Sets the number of cycles to run on this level.
492 
493    Logically Collective on PC
494 
495    Input Parameters:
496 +  pc - the multigrid context
497 .  l  - the level (0 is coarsest)
498 -  c  - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
499 
500    Level: advanced
501 
502 .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
503 
504 .seealso: PCMGSetCycleType()
505 @*/
506 PetscErrorCode  PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c)
507 {
508   PC_MG        *mg        = (PC_MG*)pc->data;
509   PC_MG_Levels **mglevels = mg->levels;
510 
511   PetscFunctionBegin;
512   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
513   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
514   PetscValidLogicalCollectiveInt(pc,l,2);
515   PetscValidLogicalCollectiveEnum(pc,c,3);
516   mglevels[l]->cycles = c;
517   PetscFunctionReturn(0);
518 }
519 
520 /*@
521    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
522    on a particular level.
523 
524    Logically Collective on PC and Vec
525 
526    Input Parameters:
527 +  pc - the multigrid context
528 .  l  - the level (0 is coarsest) this is to be used for
529 -  c  - the space
530 
531    Level: advanced
532 
533    Notes: If this is not provided PETSc will automatically generate one.
534 
535           You do not need to keep a reference to this vector if you do
536           not need it PCDestroy() will properly free it.
537 
538 .keywords: MG, multigrid, set, right-hand-side, rhs, level
539 
540 .seealso: PCMGSetX(), PCMGSetR()
541 @*/
542 PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
543 {
544   PetscErrorCode ierr;
545   PC_MG          *mg        = (PC_MG*)pc->data;
546   PC_MG_Levels   **mglevels = mg->levels;
547 
548   PetscFunctionBegin;
549   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
550   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
551   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
552   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
553   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
554 
555   mglevels[l]->b = c;
556   PetscFunctionReturn(0);
557 }
558 
559 /*@
560    PCMGSetX - Sets the vector space to be used to store the solution on a
561    particular level.
562 
563    Logically Collective on PC and Vec
564 
565    Input Parameters:
566 +  pc - the multigrid context
567 .  l - the level (0 is coarsest) this is to be used for (do not supply the finest level)
568 -  c - the space
569 
570    Level: advanced
571 
572    Notes: If this is not provided PETSc will automatically generate one.
573 
574           You do not need to keep a reference to this vector if you do
575           not need it PCDestroy() will properly free it.
576 
577 .keywords: MG, multigrid, set, solution, level
578 
579 .seealso: PCMGSetRhs(), PCMGSetR()
580 @*/
581 PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
582 {
583   PetscErrorCode ierr;
584   PC_MG          *mg        = (PC_MG*)pc->data;
585   PC_MG_Levels   **mglevels = mg->levels;
586 
587   PetscFunctionBegin;
588   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
589   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
590   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
591   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
592   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
593 
594   mglevels[l]->x = c;
595   PetscFunctionReturn(0);
596 }
597 
598 /*@
599    PCMGSetR - Sets the vector space to be used to store the residual on a
600    particular level.
601 
602    Logically Collective on PC and Vec
603 
604    Input Parameters:
605 +  pc - the multigrid context
606 .  l - the level (0 is coarsest) this is to be used for
607 -  c - the space
608 
609    Level: advanced
610 
611    Notes: If this is not provided PETSc will automatically generate one.
612 
613           You do not need to keep a reference to this vector if you do
614           not need it PCDestroy() will properly free it.
615 
616 .keywords: MG, multigrid, set, residual, level
617 @*/
618 PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
619 {
620   PetscErrorCode ierr;
621   PC_MG          *mg        = (PC_MG*)pc->data;
622   PC_MG_Levels   **mglevels = mg->levels;
623 
624   PetscFunctionBegin;
625   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
626   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
627   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
628   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
629   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
630 
631   mglevels[l]->r = c;
632   PetscFunctionReturn(0);
633 }
634