xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision b90f4e8e2d1ae90880ffa0fd3113e4b748be9ce0)
1 
2 #include <../src/ksp/pc/impls/mg/mgimpl.h>       /*I "petscksp.h" I*/
3                           /*I "petscpcmg.h"   I*/
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "PCMGDefaultResidual"
7 /*@C
8    PCMGDefaultResidual - Default routine to calculate the residual.
9 
10    Collective on Mat and Vec
11 
12    Input Parameters:
13 +  mat - the matrix
14 .  b   - the right-hand-side
15 -  x   - the approximate solution
16 
17    Output Parameter:
18 .  r - location to store the residual
19 
20    Level: advanced
21 
22 .keywords: MG, default, multigrid, residual
23 
24 .seealso: PCMGSetResidual()
25 @*/
26 PetscErrorCode  PCMGDefaultResidual(Mat mat,Vec b,Vec x,Vec r)
27 {
28   PetscErrorCode ierr;
29 
30   PetscFunctionBegin;
31   ierr = MatMult(mat,x,r);CHKERRQ(ierr);
32   ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
33   PetscFunctionReturn(0);
34 }
35 
36 /* ---------------------------------------------------------------------------*/
37 
38 #undef __FUNCT__
39 #define __FUNCT__ "PCMGGetCoarseSolve"
40 /*@
41    PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
42 
43    Not Collective
44 
45    Input Parameter:
46 .  pc - the multigrid context
47 
48    Output Parameter:
49 .  ksp - the coarse grid solver context
50 
51    Level: advanced
52 
53 .keywords: MG, multigrid, get, coarse grid
54 @*/
55 PetscErrorCode  PCMGGetCoarseSolve(PC pc,KSP *ksp)
56 {
57   PC_MG          *mg = (PC_MG*)pc->data;
58   PC_MG_Levels   **mglevels = mg->levels;
59 
60   PetscFunctionBegin;
61   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
62   *ksp =  mglevels[0]->smoothd;
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "PCMGSetResidual"
68 /*@C
69    PCMGSetResidual - Sets the function to be used to calculate the residual
70    on the lth level.
71 
72    Logically Collective on PC and Mat
73 
74    Input Parameters:
75 +  pc       - the multigrid context
76 .  l        - the level (0 is coarsest) to supply
77 .  residual - function used to form residual, if none is provided the previously provide one is used, if no
78               previous one were provided then PCMGDefaultResidual() is used
79 -  mat      - matrix associated with residual
80 
81    Level: advanced
82 
83 .keywords:  MG, set, multigrid, residual, level
84 
85 .seealso: PCMGDefaultResidual()
86 @*/
87 PetscErrorCode  PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
88 {
89   PC_MG          *mg = (PC_MG*)pc->data;
90   PC_MG_Levels   **mglevels = mg->levels;
91   PetscErrorCode ierr;
92 
93   PetscFunctionBegin;
94   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
95   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
96   if (residual) {
97     mglevels[l]->residual = residual;
98   } if (!mglevels[l]->residual) {
99     mglevels[l]->residual = PCMGDefaultResidual;
100   }
101   if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);}
102   ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr);
103   mglevels[l]->A        = mat;
104   PetscFunctionReturn(0);
105 }
106 
107 #undef __FUNCT__
108 #define __FUNCT__ "PCMGSetInterpolation"
109 /*@
110    PCMGSetInterpolation - Sets the function to be used to calculate the
111    interpolation from l-1 to the lth level
112 
113    Logically Collective on PC and Mat
114 
115    Input Parameters:
116 +  pc  - the multigrid context
117 .  mat - the interpolation operator
118 -  l   - the level (0 is coarsest) to supply [do not supply 0]
119 
120    Level: advanced
121 
122    Notes:
123           Usually this is the same matrix used also to set the restriction
124     for the same level.
125 
126           One can pass in the interpolation matrix or its transpose; PETSc figures
127     out from the matrix size which one it is.
128 
129 .keywords:  multigrid, set, interpolate, level
130 
131 .seealso: PCMGSetRestriction()
132 @*/
133 PetscErrorCode  PCMGSetInterpolation(PC pc,PetscInt l,Mat mat)
134 {
135   PC_MG          *mg = (PC_MG*)pc->data;
136   PC_MG_Levels   **mglevels = mg->levels;
137   PetscErrorCode ierr;
138 
139   PetscFunctionBegin;
140   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
141   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
142   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level");
143   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
144   ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr);
145   mglevels[l]->interpolate = mat;
146   PetscFunctionReturn(0);
147 }
148 
149 #undef __FUNCT__
150 #define __FUNCT__ "PCMGSetRestriction"
151 /*@
152    PCMGSetRestriction - Sets the function to be used to restrict vector
153    from level l to l-1.
154 
155    Logically Collective on PC and Mat
156 
157    Input Parameters:
158 +  pc - the multigrid context
159 .  mat - the restriction matrix
160 -  l - the level (0 is coarsest) to supply [Do not supply 0]
161 
162    Level: advanced
163 
164    Notes:
165           Usually this is the same matrix used also to set the interpolation
166     for the same level.
167 
168           One can pass in the interpolation matrix or its transpose; PETSc figures
169     out from the matrix size which one it is.
170 
171          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
172     is used.
173 
174 .keywords: MG, set, multigrid, restriction, level
175 
176 .seealso: PCMGSetInterpolation()
177 @*/
178 PetscErrorCode  PCMGSetRestriction(PC pc,PetscInt l,Mat mat)
179 {
180   PetscErrorCode ierr;
181   PC_MG          *mg = (PC_MG*)pc->data;
182   PC_MG_Levels   **mglevels = mg->levels;
183 
184   PetscFunctionBegin;
185   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
186   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
187   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
188   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
189   ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr);
190   mglevels[l]->restrct  = mat;
191   PetscFunctionReturn(0);
192 }
193 
194 #undef __FUNCT__
195 #define __FUNCT__ "PCMGSetRScale"
196 /*@
197    PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1.
198 
199    Logically Collective on PC and Mat
200 
201    Input Parameters:
202 +  pc - the multigrid context
203 .  rscale - the scaling
204 -  l - the level (0 is coarsest) to supply [Do not supply 0]
205 
206    Level: advanced
207 
208    Notes:
209        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.
210 
211 .keywords: MG, set, multigrid, restriction, level
212 
213 .seealso: PCMGSetInterpolation(), PCMGSetRestriction()
214 @*/
215 PetscErrorCode  PCMGSetRScale(PC pc,PetscInt l,Vec rscale)
216 {
217   PetscErrorCode ierr;
218   PC_MG          *mg = (PC_MG*)pc->data;
219   PC_MG_Levels   **mglevels = mg->levels;
220 
221   PetscFunctionBegin;
222   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
223   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
224   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level");
225   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
226   ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr);
227   mglevels[l]->rscale  = rscale;
228   PetscFunctionReturn(0);
229 }
230 
231 #undef __FUNCT__
232 #define __FUNCT__ "PCMGGetSmoother"
233 /*@
234    PCMGGetSmoother - Gets the KSP context to be used as smoother for
235    both pre- and post-smoothing.  Call both PCMGGetSmootherUp() and
236    PCMGGetSmootherDown() to use different functions for pre- and
237    post-smoothing.
238 
239    Not Collective, KSP returned is parallel if PC is
240 
241    Input Parameters:
242 +  pc - the multigrid context
243 -  l - the level (0 is coarsest) to supply
244 
245    Ouput Parameters:
246 .  ksp - the smoother
247 
248    Level: advanced
249 
250 .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
251 
252 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
253 @*/
254 PetscErrorCode  PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp)
255 {
256   PC_MG          *mg = (PC_MG*)pc->data;
257   PC_MG_Levels   **mglevels = mg->levels;
258 
259   PetscFunctionBegin;
260   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
261   *ksp = mglevels[l]->smoothd;
262   PetscFunctionReturn(0);
263 }
264 
265 #undef __FUNCT__
266 #define __FUNCT__ "PCMGGetSmootherUp"
267 /*@
268    PCMGGetSmootherUp - Gets the KSP context to be used as smoother after
269    coarse grid correction (post-smoother).
270 
271    Not Collective, KSP returned is parallel if PC is
272 
273    Input Parameters:
274 +  pc - the multigrid context
275 -  l  - the level (0 is coarsest) to supply
276 
277    Ouput Parameters:
278 .  ksp - the smoother
279 
280    Level: advanced
281 
282 .keywords: MG, multigrid, get, smoother, up, post-smoother, level
283 
284 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
285 @*/
286 PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
287 {
288   PC_MG          *mg = (PC_MG*)pc->data;
289   PC_MG_Levels   **mglevels = mg->levels;
290   PetscErrorCode ierr;
291   const char     *prefix;
292   MPI_Comm       comm;
293 
294   PetscFunctionBegin;
295   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
296   /*
297      This is called only if user wants a different pre-smoother from post.
298      Thus we check if a different one has already been allocated,
299      if not we allocate it.
300   */
301   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
302   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
303     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
304     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
305     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
306     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
307     ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
308     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
309     ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr);
310   }
311   if (ksp) *ksp = mglevels[l]->smoothu;
312   PetscFunctionReturn(0);
313 }
314 
315 #undef __FUNCT__
316 #define __FUNCT__ "PCMGGetSmootherDown"
317 /*@
318    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
319    coarse grid correction (pre-smoother).
320 
321    Not Collective, KSP returned is parallel if PC is
322 
323    Input Parameters:
324 +  pc - the multigrid context
325 -  l  - the level (0 is coarsest) to supply
326 
327    Ouput Parameters:
328 .  ksp - the smoother
329 
330    Level: advanced
331 
332 .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
333 
334 .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
335 @*/
336 PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
337 {
338   PetscErrorCode ierr;
339   PC_MG          *mg = (PC_MG*)pc->data;
340   PC_MG_Levels   **mglevels = mg->levels;
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
344   /* make sure smoother up and down are different */
345   if (l) {
346     ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
347   }
348   *ksp = mglevels[l]->smoothd;
349   PetscFunctionReturn(0);
350 }
351 
352 #undef __FUNCT__
353 #define __FUNCT__ "PCMGSetCyclesOnLevel"
354 /*@
355    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
356 
357    Logically Collective on PC
358 
359    Input Parameters:
360 +  pc - the multigrid context
361 .  l  - the level (0 is coarsest) this is to be used for
362 -  n  - the number of cycles
363 
364    Level: advanced
365 
366 .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
367 
368 .seealso: PCMGSetCycles()
369 @*/
370 PetscErrorCode  PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
371 {
372   PC_MG          *mg = (PC_MG*)pc->data;
373   PC_MG_Levels   **mglevels = mg->levels;
374 
375   PetscFunctionBegin;
376   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
377   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
378   PetscValidLogicalCollectiveInt(pc,l,2);
379   PetscValidLogicalCollectiveInt(pc,c,3);
380   mglevels[l]->cycles  = c;
381   PetscFunctionReturn(0);
382 }
383 
384 #undef __FUNCT__
385 #define __FUNCT__ "PCMGSetRhs"
386 /*@
387    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
388    on a particular level.
389 
390    Logically Collective on PC and Vec
391 
392    Input Parameters:
393 +  pc - the multigrid context
394 .  l  - the level (0 is coarsest) this is to be used for
395 -  c  - the space
396 
397    Level: advanced
398 
399    Notes: If this is not provided PETSc will automatically generate one.
400 
401           You do not need to keep a reference to this vector if you do
402           not need it PCDestroy() will properly free it.
403 
404 .keywords: MG, multigrid, set, right-hand-side, rhs, level
405 
406 .seealso: PCMGSetX(), PCMGSetR()
407 @*/
408 PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
409 {
410   PetscErrorCode ierr;
411   PC_MG          *mg = (PC_MG*)pc->data;
412   PC_MG_Levels   **mglevels = mg->levels;
413 
414   PetscFunctionBegin;
415   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
416   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
417   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
418   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
419   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
420   mglevels[l]->b  = c;
421   PetscFunctionReturn(0);
422 }
423 
424 #undef __FUNCT__
425 #define __FUNCT__ "PCMGSetX"
426 /*@
427    PCMGSetX - Sets the vector space to be used to store the solution on a
428    particular level.
429 
430    Logically Collective on PC and Vec
431 
432    Input Parameters:
433 +  pc - the multigrid context
434 .  l - the level (0 is coarsest) this is to be used for
435 -  c - the space
436 
437    Level: advanced
438 
439    Notes: If this is not provided PETSc will automatically generate one.
440 
441           You do not need to keep a reference to this vector if you do
442           not need it PCDestroy() will properly free it.
443 
444 .keywords: MG, multigrid, set, solution, level
445 
446 .seealso: PCMGSetRhs(), PCMGSetR()
447 @*/
448 PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
449 {
450   PetscErrorCode ierr;
451   PC_MG          *mg = (PC_MG*)pc->data;
452   PC_MG_Levels   **mglevels = mg->levels;
453 
454   PetscFunctionBegin;
455   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
456   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
457   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
458   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
459   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
460   mglevels[l]->x  = c;
461   PetscFunctionReturn(0);
462 }
463 
464 #undef __FUNCT__
465 #define __FUNCT__ "PCMGSetR"
466 /*@
467    PCMGSetR - Sets the vector space to be used to store the residual on a
468    particular level.
469 
470    Logically Collective on PC and Vec
471 
472    Input Parameters:
473 +  pc - the multigrid context
474 .  l - the level (0 is coarsest) this is to be used for
475 -  c - the space
476 
477    Level: advanced
478 
479    Notes: If this is not provided PETSc will automatically generate one.
480 
481           You do not need to keep a reference to this vector if you do
482           not need it PCDestroy() will properly free it.
483 
484 .keywords: MG, multigrid, set, residual, level
485 @*/
486 PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
487 {
488   PetscErrorCode ierr;
489   PC_MG          *mg = (PC_MG*)pc->data;
490   PC_MG_Levels   **mglevels = mg->levels;
491 
492   PetscFunctionBegin;
493   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
494   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
495   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
496   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
497   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
498   mglevels[l]->r  = c;
499   PetscFunctionReturn(0);
500 }
501