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