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