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