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