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