xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision feff33ee0b5b037fa8f9f294dede656a2f85cc47)
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 dual vectors
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(), PCMGetSetInjection()
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 dual vectors
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(), PCMGGetInjection()
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.  It is preferable to use PCMGSetInjection() to control moving primal vectors.
271 
272 .keywords: MG, set, multigrid, restriction, level
273 
274 .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection()
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.  It is preferable to use PCMGGetInjection() to control moving primal vectors.
307 
308 .keywords: MG, set, multigrid, restriction, level
309 
310 .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection()
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    PCMGSetInjection - Sets the function to be used to inject primal vectors
347    from level l to l-1.
348 
349    Logically Collective on PC and Mat
350 
351    Input Parameters:
352 +  pc - the multigrid context
353 .  l - the level (0 is coarsest) to supply [Do not supply 0]
354 -  mat - the injection matrix
355 
356    Level: advanced
357 
358 .keywords: MG, set, multigrid, restriction, injection, level
359 
360 .seealso: PCMGSetRestriction()
361 @*/
362 PetscErrorCode  PCMGSetInjection(PC pc,PetscInt l,Mat mat)
363 {
364   PetscErrorCode ierr;
365   PC_MG          *mg        = (PC_MG*)pc->data;
366   PC_MG_Levels   **mglevels = mg->levels;
367 
368   PetscFunctionBegin;
369   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
370   PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
371   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
372   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
373   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
374   ierr = MatDestroy(&mglevels[l]->inject);CHKERRQ(ierr);
375 
376   mglevels[l]->inject = mat;
377   PetscFunctionReturn(0);
378 }
379 
380 /*@
381    PCMGGetInjection - Gets the function to be used to inject primal vectors
382    from level l to l-1.
383 
384    Logically Collective on PC and Mat
385 
386    Input Parameters:
387 +  pc - the multigrid context
388 -  l - the level (0 is coarsest) to supply [Do not supply 0]
389 
390    Output Parameter:
391 .  mat - the restriction matrix (may be NULL if no injection is available).
392 
393    Level: advanced
394 
395 .keywords: MG, get, multigrid, restriction, injection, level
396 
397 .seealso: PCMGSetInjection(), PCMGetGetRestriction()
398 @*/
399 PetscErrorCode  PCMGGetInjection(PC pc,PetscInt l,Mat *mat)
400 {
401   PC_MG          *mg        = (PC_MG*)pc->data;
402   PC_MG_Levels   **mglevels = mg->levels;
403 
404   PetscFunctionBegin;
405   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
406   if (mat) PetscValidPointer(mat,3);
407   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
408   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);
409   if (mat) *mat = mglevels[l]->inject;
410   PetscFunctionReturn(0);
411 }
412 
413 /*@
414    PCMGGetSmoother - Gets the KSP context to be used as smoother for
415    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
416    PCMGGetSmootherDown() to use different functions for pre- and
417    post-smoothing.
418 
419    Not Collective, KSP returned is parallel if PC is
420 
421    Input Parameters:
422 +  pc - the multigrid context
423 -  l - the level (0 is coarsest) to supply
424 
425    Ouput Parameters:
426 .  ksp - the smoother
427 
428    Notes:
429    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.
430    You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc)
431    and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing.
432 
433    Level: advanced
434 
435 .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
436 
437 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
438 @*/
439 PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
440 {
441   PC_MG        *mg        = (PC_MG*)pc->data;
442   PC_MG_Levels **mglevels = mg->levels;
443 
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
446   *ksp = mglevels[l]->smoothd;
447   PetscFunctionReturn(0);
448 }
449 
450 /*@
451    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
452    coarse grid correction (post-smoother).
453 
454    Not Collective, KSP returned is parallel if PC is
455 
456    Input Parameters:
457 +  pc - the multigrid context
458 -  l  - the level (0 is coarsest) to supply
459 
460    Ouput Parameters:
461 .  ksp - the smoother
462 
463    Level: advanced
464 
465    Notes: calling this will result in a different pre and post smoother so you may need to
466          set options on the pre smoother also
467 
468 .keywords: MG, multigrid, get, smoother, up, post-smoother, level
469 
470 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
471 @*/
472 PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
473 {
474   PC_MG          *mg        = (PC_MG*)pc->data;
475   PC_MG_Levels   **mglevels = mg->levels;
476   PetscErrorCode ierr;
477   const char     *prefix;
478   MPI_Comm       comm;
479 
480   PetscFunctionBegin;
481   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
482   /*
483      This is called only if user wants a different pre-smoother from post.
484      Thus we check if a different one has already been allocated,
485      if not we allocate it.
486   */
487   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
488   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
489     KSPType     ksptype;
490     PCType      pctype;
491     PC          ipc;
492     PetscReal   rtol,abstol,dtol;
493     PetscInt    maxits;
494     KSPNormType normtype;
495     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
496     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
497     ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr);
498     ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr);
499     ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr);
500     ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr);
501     ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr);
502 
503     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
504     ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr);
505     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
506     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
507     ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr);
508     ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr);
509     ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr);
510     ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
511     ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr);
512     ierr = PCSetType(ipc,pctype);CHKERRQ(ierr);
513     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr);
514     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr);
515   }
516   if (ksp) *ksp = mglevels[l]->smoothu;
517   PetscFunctionReturn(0);
518 }
519 
520 /*@
521    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
522    coarse grid correction (pre-smoother).
523 
524    Not Collective, KSP returned is parallel if PC is
525 
526    Input Parameters:
527 +  pc - the multigrid context
528 -  l  - the level (0 is coarsest) to supply
529 
530    Ouput Parameters:
531 .  ksp - the smoother
532 
533    Level: advanced
534 
535    Notes: calling this will result in a different pre and post smoother so you may need to
536          set options on the post smoother also
537 
538 .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
539 
540 .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
541 @*/
542 PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
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   /* make sure smoother up and down are different */
551   if (l) {
552     ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr);
553   }
554   *ksp = mglevels[l]->smoothd;
555   PetscFunctionReturn(0);
556 }
557 
558 /*@
559    PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level.
560 
561    Logically Collective on PC
562 
563    Input Parameters:
564 +  pc - the multigrid context
565 .  l  - the level (0 is coarsest)
566 -  c  - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
567 
568    Level: advanced
569 
570 .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
571 
572 .seealso: PCMGSetCycleType()
573 @*/
574 PetscErrorCode  PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c)
575 {
576   PC_MG        *mg        = (PC_MG*)pc->data;
577   PC_MG_Levels **mglevels = mg->levels;
578 
579   PetscFunctionBegin;
580   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
581   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
582   PetscValidLogicalCollectiveInt(pc,l,2);
583   PetscValidLogicalCollectiveEnum(pc,c,3);
584   mglevels[l]->cycles = c;
585   PetscFunctionReturn(0);
586 }
587 
588 /*@
589    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
590    on a particular level.
591 
592    Logically Collective on PC and Vec
593 
594    Input Parameters:
595 +  pc - the multigrid context
596 .  l  - the level (0 is coarsest) this is to be used for
597 -  c  - the space
598 
599    Level: advanced
600 
601    Notes: If this is not provided PETSc will automatically generate one.
602 
603           You do not need to keep a reference to this vector if you do
604           not need it PCDestroy() will properly free it.
605 
606 .keywords: MG, multigrid, set, right-hand-side, rhs, level
607 
608 .seealso: PCMGSetX(), PCMGSetR()
609 @*/
610 PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
611 {
612   PetscErrorCode ierr;
613   PC_MG          *mg        = (PC_MG*)pc->data;
614   PC_MG_Levels   **mglevels = mg->levels;
615 
616   PetscFunctionBegin;
617   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
618   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
619   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
620   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
621   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
622 
623   mglevels[l]->b = c;
624   PetscFunctionReturn(0);
625 }
626 
627 /*@
628    PCMGSetX - Sets the vector space to be used to store the solution on a
629    particular level.
630 
631    Logically Collective on PC and Vec
632 
633    Input Parameters:
634 +  pc - the multigrid context
635 .  l - the level (0 is coarsest) this is to be used for (do not supply the finest level)
636 -  c - the space
637 
638    Level: advanced
639 
640    Notes: If this is not provided PETSc will automatically generate one.
641 
642           You do not need to keep a reference to this vector if you do
643           not need it PCDestroy() will properly free it.
644 
645 .keywords: MG, multigrid, set, solution, level
646 
647 .seealso: PCMGSetRhs(), PCMGSetR()
648 @*/
649 PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
650 {
651   PetscErrorCode ierr;
652   PC_MG          *mg        = (PC_MG*)pc->data;
653   PC_MG_Levels   **mglevels = mg->levels;
654 
655   PetscFunctionBegin;
656   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
657   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
658   if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
659   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
660   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
661 
662   mglevels[l]->x = c;
663   PetscFunctionReturn(0);
664 }
665 
666 /*@
667    PCMGSetR - Sets the vector space to be used to store the residual on a
668    particular level.
669 
670    Logically Collective on PC and Vec
671 
672    Input Parameters:
673 +  pc - the multigrid context
674 .  l - the level (0 is coarsest) this is to be used for
675 -  c - the space
676 
677    Level: advanced
678 
679    Notes: If this is not provided PETSc will automatically generate one.
680 
681           You do not need to keep a reference to this vector if you do
682           not need it PCDestroy() will properly free it.
683 
684 .keywords: MG, multigrid, set, residual, level
685 @*/
686 PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
687 {
688   PetscErrorCode ierr;
689   PC_MG          *mg        = (PC_MG*)pc->data;
690   PC_MG_Levels   **mglevels = mg->levels;
691 
692   PetscFunctionBegin;
693   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
694   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
695   if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
696   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
697   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
698 
699   mglevels[l]->r = c;
700   PetscFunctionReturn(0);
701 }
702