xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 84df9cb40eca90ea9b18a456fab7a4ecc7f6c1a4)
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    Notes: calling this will result in a different pre and post smoother so you may need to
283          set options on the pre smoother also
284 
285 .keywords: MG, multigrid, get, smoother, up, post-smoother, level
286 
287 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown()
288 @*/
289 PetscErrorCode  PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
290 {
291   PC_MG          *mg = (PC_MG*)pc->data;
292   PC_MG_Levels   **mglevels = mg->levels;
293   PetscErrorCode ierr;
294   const char     *prefix;
295   MPI_Comm       comm;
296 
297   PetscFunctionBegin;
298   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
299   /*
300      This is called only if user wants a different pre-smoother from post.
301      Thus we check if a different one has already been allocated,
302      if not we allocate it.
303   */
304   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
305   if (mglevels[l]->smoothu == mglevels[l]->smoothd) {
306     ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr);
307     ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr);
308     ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr);
309     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr);
310     ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
311     ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr);
312     ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr);
313   }
314   if (ksp) *ksp = mglevels[l]->smoothu;
315   PetscFunctionReturn(0);
316 }
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "PCMGGetSmootherDown"
320 /*@
321    PCMGGetSmootherDown - Gets the KSP context to be used as smoother before
322    coarse grid correction (pre-smoother).
323 
324    Not Collective, KSP returned is parallel if PC is
325 
326    Input Parameters:
327 +  pc - the multigrid context
328 -  l  - the level (0 is coarsest) to supply
329 
330    Ouput Parameters:
331 .  ksp - the smoother
332 
333    Level: advanced
334 
335    Notes: calling this will result in a different pre and post smoother so you may need to
336          set options on the post smoother also
337 
338 .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
339 
340 .seealso: PCMGGetSmootherUp(), PCMGGetSmoother()
341 @*/
342 PetscErrorCode  PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
343 {
344   PetscErrorCode ierr;
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   /* make sure smoother up and down are different */
351   if (l) {
352     ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
353   }
354   *ksp = mglevels[l]->smoothd;
355   PetscFunctionReturn(0);
356 }
357 
358 #undef __FUNCT__
359 #define __FUNCT__ "PCMGSetCyclesOnLevel"
360 /*@
361    PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level.
362 
363    Logically Collective on PC
364 
365    Input Parameters:
366 +  pc - the multigrid context
367 .  l  - the level (0 is coarsest) this is to be used for
368 -  n  - the number of cycles
369 
370    Level: advanced
371 
372 .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
373 
374 .seealso: PCMGSetCycles()
375 @*/
376 PetscErrorCode  PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
377 {
378   PC_MG          *mg = (PC_MG*)pc->data;
379   PC_MG_Levels   **mglevels = mg->levels;
380 
381   PetscFunctionBegin;
382   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
383   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
384   PetscValidLogicalCollectiveInt(pc,l,2);
385   PetscValidLogicalCollectiveInt(pc,c,3);
386   mglevels[l]->cycles  = c;
387   PetscFunctionReturn(0);
388 }
389 
390 #undef __FUNCT__
391 #define __FUNCT__ "PCMGSetRhs"
392 /*@
393    PCMGSetRhs - Sets the vector space to be used to store the right-hand side
394    on a particular level.
395 
396    Logically Collective on PC and Vec
397 
398    Input Parameters:
399 +  pc - the multigrid context
400 .  l  - the level (0 is coarsest) this is to be used for
401 -  c  - the space
402 
403    Level: advanced
404 
405    Notes: If this is not provided PETSc will automatically generate one.
406 
407           You do not need to keep a reference to this vector if you do
408           not need it PCDestroy() will properly free it.
409 
410 .keywords: MG, multigrid, set, right-hand-side, rhs, level
411 
412 .seealso: PCMGSetX(), PCMGSetR()
413 @*/
414 PetscErrorCode  PCMGSetRhs(PC pc,PetscInt l,Vec c)
415 {
416   PetscErrorCode ierr;
417   PC_MG          *mg = (PC_MG*)pc->data;
418   PC_MG_Levels   **mglevels = mg->levels;
419 
420   PetscFunctionBegin;
421   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
422   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
423   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level");
424   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
425   ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr);
426   mglevels[l]->b  = c;
427   PetscFunctionReturn(0);
428 }
429 
430 #undef __FUNCT__
431 #define __FUNCT__ "PCMGSetX"
432 /*@
433    PCMGSetX - Sets the vector space to be used to store the solution on a
434    particular level.
435 
436    Logically Collective on PC and Vec
437 
438    Input Parameters:
439 +  pc - the multigrid context
440 .  l - the level (0 is coarsest) this is to be used for
441 -  c - the space
442 
443    Level: advanced
444 
445    Notes: If this is not provided PETSc will automatically generate one.
446 
447           You do not need to keep a reference to this vector if you do
448           not need it PCDestroy() will properly free it.
449 
450 .keywords: MG, multigrid, set, solution, level
451 
452 .seealso: PCMGSetRhs(), PCMGSetR()
453 @*/
454 PetscErrorCode  PCMGSetX(PC pc,PetscInt l,Vec c)
455 {
456   PetscErrorCode ierr;
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   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
463   if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level");
464   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
465   ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr);
466   mglevels[l]->x  = c;
467   PetscFunctionReturn(0);
468 }
469 
470 #undef __FUNCT__
471 #define __FUNCT__ "PCMGSetR"
472 /*@
473    PCMGSetR - Sets the vector space to be used to store the residual on a
474    particular level.
475 
476    Logically Collective on PC and Vec
477 
478    Input Parameters:
479 +  pc - the multigrid context
480 .  l - the level (0 is coarsest) this is to be used for
481 -  c - the space
482 
483    Level: advanced
484 
485    Notes: If this is not provided PETSc will automatically generate one.
486 
487           You do not need to keep a reference to this vector if you do
488           not need it PCDestroy() will properly free it.
489 
490 .keywords: MG, multigrid, set, residual, level
491 @*/
492 PetscErrorCode  PCMGSetR(PC pc,PetscInt l,Vec c)
493 {
494   PetscErrorCode ierr;
495   PC_MG          *mg = (PC_MG*)pc->data;
496   PC_MG_Levels   **mglevels = mg->levels;
497 
498   PetscFunctionBegin;
499   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
500   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
501   if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid");
502   ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr);
503   ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr);
504   mglevels[l]->r  = c;
505   PetscFunctionReturn(0);
506 }
507