xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 8bc8193efbc389280f83b3d41dffa9e2d23e2ace)
1 
2 #include "src/ksp/pc/impls/mg/mgimpl.h"       /*I "petscksp.h" I*/
3                           /*I "petscmg.h"   I*/
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "MGDefaultResidual"
7 /*@C
8    MGDefaultResidual - 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: MGSetResidual()
25 @*/
26 PetscErrorCode MGDefaultResidual(Mat mat,Vec b,Vec x,Vec r)
27 {
28   PetscErrorCode ierr;
29   PetscScalar    mone = -1.0;
30 
31   PetscFunctionBegin;
32   ierr = MatMult(mat,x,r);CHKERRQ(ierr);
33   ierr = VecAYPX(&mone,b,r);CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 /* ---------------------------------------------------------------------------*/
38 
39 #undef __FUNCT__
40 #define __FUNCT__ "MGGetCoarseSolve"
41 /*@C
42    MGGetCoarseSolve - Gets the solver context to be used on the coarse grid.
43 
44    Not Collective
45 
46    Input Parameter:
47 .  pc - the multigrid context
48 
49    Output Parameter:
50 .  ksp - the coarse grid solver context
51 
52    Level: advanced
53 
54 .keywords: MG, multigrid, get, coarse grid
55 @*/
56 PetscErrorCode MGGetCoarseSolve(PC pc,KSP *ksp)
57 {
58   MG *mg = (MG*)pc->data;
59 
60   PetscFunctionBegin;
61   *ksp =  mg[0]->smoothd;
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "MGSetResidual"
67 /*@C
68    MGSetResidual - Sets the function to be used to calculate the residual
69    on the lth level.
70 
71    Collective on PC and Mat
72 
73    Input Parameters:
74 +  pc       - the multigrid context
75 .  l        - the level (0 is coarsest) to supply
76 .  residual - function used to form residual (usually MGDefaultResidual)
77 -  mat      - matrix associated with residual
78 
79    Level: advanced
80 
81 .keywords:  MG, set, multigrid, residual, level
82 
83 .seealso: MGDefaultResidual()
84 @*/
85 PetscErrorCode MGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat)
86 {
87   MG *mg = (MG*)pc->data;
88 
89   PetscFunctionBegin;
90   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
91 
92   mg[l]->residual = residual;
93   mg[l]->A        = mat;
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "MGSetInterpolate"
99 /*@
100    MGSetInterpolate - Sets the function to be used to calculate the
101    interpolation on the lth level.
102 
103    Collective on PC and Mat
104 
105    Input Parameters:
106 +  pc  - the multigrid context
107 .  mat - the interpolation operator
108 -  l   - the level (0 is coarsest) to supply
109 
110    Level: advanced
111 
112    Notes:
113           Usually this is the same matrix used also to set the restriction
114     for the same level.
115 
116           One can pass in the interpolation matrix or its transpose; PETSc figures
117     out from the matrix size which one it is.
118 
119 .keywords:  multigrid, set, interpolate, level
120 
121 .seealso: MGSetRestriction()
122 @*/
123 PetscErrorCode MGSetInterpolate(PC pc,PetscInt l,Mat mat)
124 {
125   MG *mg = (MG*)pc->data;
126 
127   PetscFunctionBegin;
128   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
129   mg[l]->interpolate = mat;
130   PetscFunctionReturn(0);
131 }
132 
133 #undef __FUNCT__
134 #define __FUNCT__ "MGSetRestriction"
135 /*@
136    MGSetRestriction - Sets the function to be used to restrict vector
137    from level l to l-1.
138 
139    Collective on PC and Mat
140 
141    Input Parameters:
142 +  pc - the multigrid context
143 .  mat - the restriction matrix
144 -  l - the level (0 is coarsest) to supply
145 
146    Level: advanced
147 
148    Notes:
149           Usually this is the same matrix used also to set the interpolation
150     for the same level.
151 
152           One can pass in the interpolation matrix or its transpose; PETSc figures
153     out from the matrix size which one it is.
154 
155 .keywords: MG, set, multigrid, restriction, level
156 
157 .seealso: MGSetInterpolate()
158 @*/
159 PetscErrorCode MGSetRestriction(PC pc,PetscInt l,Mat mat)
160 {
161   MG *mg = (MG*)pc->data;
162 
163   PetscFunctionBegin;
164   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
165   mg[l]->restrct  = mat;
166   PetscFunctionReturn(0);
167 }
168 
169 #undef __FUNCT__
170 #define __FUNCT__ "MGGetSmoother"
171 /*@C
172    MGGetSmoother - Gets the KSP context to be used as smoother for
173    both pre- and post-smoothing.  Call both MGGetSmootherUp() and
174    MGGetSmootherDown() to use different functions for pre- and
175    post-smoothing.
176 
177    Not Collective, KSP returned is parallel if PC is
178 
179    Input Parameters:
180 +  pc - the multigrid context
181 -  l - the level (0 is coarsest) to supply
182 
183    Ouput Parameters:
184 .  ksp - the smoother
185 
186    Level: advanced
187 
188 .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
189 
190 .seealso: MGGetSmootherUp(), MGGetSmootherDown()
191 @*/
192 PetscErrorCode MGGetSmoother(PC pc,PetscInt l,KSP *ksp)
193 {
194   MG *mg = (MG*)pc->data;
195 
196   PetscFunctionBegin;
197   *ksp = mg[l]->smoothd;
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "MGGetSmootherUp"
203 /*@C
204    MGGetSmootherUp - Gets the KSP context to be used as smoother after
205    coarse grid correction (post-smoother).
206 
207    Not Collective, KSP returned is parallel if PC is
208 
209    Input Parameters:
210 +  pc - the multigrid context
211 -  l  - the level (0 is coarsest) to supply
212 
213    Ouput Parameters:
214 .  ksp - the smoother
215 
216    Level: advanced
217 
218 .keywords: MG, multigrid, get, smoother, up, post-smoother, level
219 
220 .seealso: MGGetSmootherUp(), MGGetSmootherDown()
221 @*/
222 PetscErrorCode MGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
223 {
224   MG             *mg = (MG*)pc->data;
225   PetscErrorCode ierr;
226   char           *prefix;
227   MPI_Comm       comm;
228 
229   PetscFunctionBegin;
230   /*
231      This is called only if user wants a different pre-smoother from post.
232      Thus we check if a different one has already been allocated,
233      if not we allocate it.
234   */
235 
236   if (mg[l]->smoothu == mg[l]->smoothd) {
237     ierr = PetscObjectGetComm((PetscObject)mg[l]->smoothd,&comm);CHKERRQ(ierr);
238     ierr = KSPGetOptionsPrefix(mg[l]->smoothd,&prefix);CHKERRQ(ierr);
239     ierr = KSPCreate(comm,&mg[l]->smoothu);CHKERRQ(ierr);
240     ierr = KSPSetTolerances(mg[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
241     ierr = KSPSetOptionsPrefix(mg[l]->smoothu,prefix);CHKERRQ(ierr);
242     PetscLogObjectParent(pc,mg[l]->smoothu);
243   }
244   if (ksp) *ksp = mg[l]->smoothu;
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNCT__
249 #define __FUNCT__ "MGGetSmootherDown"
250 /*@C
251    MGGetSmootherDown - Gets the KSP context to be used as smoother before
252    coarse grid correction (pre-smoother).
253 
254    Not Collective, KSP returned is parallel if PC is
255 
256    Input Parameters:
257 +  pc - the multigrid context
258 -  l  - the level (0 is coarsest) to supply
259 
260    Ouput Parameters:
261 .  ksp - the smoother
262 
263    Level: advanced
264 
265 .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
266 
267 .seealso: MGGetSmootherUp(), MGGetSmoother()
268 @*/
269 PetscErrorCode MGGetSmootherDown(PC pc,PetscInt l,KSP *ksp)
270 {
271   PetscErrorCode ierr;
272   MG             *mg = (MG*)pc->data;
273 
274   PetscFunctionBegin;
275   /* make sure smoother up and down are different */
276   ierr = MGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
277   *ksp = mg[l]->smoothd;
278   PetscFunctionReturn(0);
279 }
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "MGSetCyclesOnLevel"
283 /*@
284    MGSetCyclesOnLevel - Sets the number of cycles to run on this level.
285 
286    Collective on PC
287 
288    Input Parameters:
289 +  pc - the multigrid context
290 .  l  - the level (0 is coarsest) this is to be used for
291 -  n  - the number of cycles
292 
293    Level: advanced
294 
295 .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
296 
297 .seealso: MGSetCycles()
298 @*/
299 PetscErrorCode MGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c)
300 {
301   MG *mg = (MG*)pc->data;
302 
303   PetscFunctionBegin;
304   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
305   mg[l]->cycles  = c;
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "MGSetRhs"
311 /*@
312    MGSetRhs - Sets the vector space to be used to store the right-hand side
313    on a particular level.  The user should free this space at the conclusion
314    of multigrid use.
315 
316    Collective on PC and Vec
317 
318    Input Parameters:
319 +  pc - the multigrid context
320 .  l  - the level (0 is coarsest) this is to be used for
321 -  c  - the space
322 
323    Level: advanced
324 
325 .keywords: MG, multigrid, set, right-hand-side, rhs, level
326 
327 .seealso: MGSetX(), MGSetR()
328 @*/
329 PetscErrorCode MGSetRhs(PC pc,PetscInt l,Vec c)
330 {
331   MG *mg = (MG*)pc->data;
332 
333   PetscFunctionBegin;
334   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
335   mg[l]->b  = c;
336   PetscFunctionReturn(0);
337 }
338 
339 #undef __FUNCT__
340 #define __FUNCT__ "MGSetX"
341 /*@
342    MGSetX - Sets the vector space to be used to store the solution on a
343    particular level.  The user should free this space at the conclusion
344    of multigrid use.
345 
346    Collective on PC and Vec
347 
348    Input Parameters:
349 +  pc - the multigrid context
350 .  l - the level (0 is coarsest) this is to be used for
351 -  c - the space
352 
353    Level: advanced
354 
355 .keywords: MG, multigrid, set, solution, level
356 
357 .seealso: MGSetRhs(), MGSetR()
358 @*/
359 PetscErrorCode MGSetX(PC pc,PetscInt l,Vec c)
360 {
361   MG *mg = (MG*)pc->data;
362 
363   PetscFunctionBegin;
364   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
365   mg[l]->x  = c;
366   PetscFunctionReturn(0);
367 }
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "MGSetR"
371 /*@
372    MGSetR - Sets the vector space to be used to store the residual on a
373    particular level.  The user should free this space at the conclusion of
374    multigrid use.
375 
376    Collective on PC and Vec
377 
378    Input Parameters:
379 +  pc - the multigrid context
380 .  l - the level (0 is coarsest) this is to be used for
381 -  c - the space
382 
383    Level: advanced
384 
385 .keywords: MG, multigrid, set, residual, level
386 @*/
387 PetscErrorCode MGSetR(PC pc,PetscInt l,Vec c)
388 {
389   MG *mg = (MG*)pc->data;
390 
391   PetscFunctionBegin;
392   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
393   mg[l]->r  = c;
394   PetscFunctionReturn(0);
395 }
396 
397 
398 
399 
400 
401 
402 
403 
404