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