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