xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
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,int 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,int 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,int 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,int 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,int 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   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
236 
237   if (mg[l]->smoothu == mg[l]->smoothd) {
238     ierr = PetscObjectGetComm((PetscObject)mg[l]->smoothd,&comm);CHKERRQ(ierr);
239     ierr = KSPCreate(comm,&mg[l]->smoothu);CHKERRQ(ierr);
240     ierr = KSPSetTolerances(mg[l]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
241     ierr = KSPSetOptionsPrefix(mg[l]->smoothu,prefix);CHKERRQ(ierr);
242     ierr = KSPAppendOptionsPrefix(mg[l]->smoothd,"mg_levels_");CHKERRQ(ierr);
243     PetscLogObjectParent(pc,mg[l]->smoothu);
244   }
245   if (ksp) *ksp = mg[l]->smoothu;
246   PetscFunctionReturn(0);
247 }
248 
249 #undef __FUNCT__
250 #define __FUNCT__ "MGGetSmootherDown"
251 /*@C
252    MGGetSmootherDown - Gets the KSP context to be used as smoother before
253    coarse grid correction (pre-smoother).
254 
255    Not Collective, KSP returned is parallel if PC is
256 
257    Input Parameters:
258 +  pc - the multigrid context
259 -  l  - the level (0 is coarsest) to supply
260 
261    Ouput Parameters:
262 .  ksp - the smoother
263 
264    Level: advanced
265 
266 .keywords: MG, multigrid, get, smoother, down, pre-smoother, level
267 
268 .seealso: MGGetSmootherUp(), MGGetSmoother()
269 @*/
270 PetscErrorCode MGGetSmootherDown(PC pc,int l,KSP *ksp)
271 {
272   PetscErrorCode ierr;
273   MG  *mg = (MG*)pc->data;
274 
275   PetscFunctionBegin;
276   /* make sure smoother up and down are different */
277   ierr = MGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr);
278   *ksp = mg[l]->smoothd;
279   PetscFunctionReturn(0);
280 }
281 
282 #undef __FUNCT__
283 #define __FUNCT__ "MGSetCyclesOnLevel"
284 /*@
285    MGSetCyclesOnLevel - Sets the number of cycles to run on this level.
286 
287    Collective on PC
288 
289    Input Parameters:
290 +  pc - the multigrid context
291 .  l  - the level (0 is coarsest) this is to be used for
292 -  n  - the number of cycles
293 
294    Level: advanced
295 
296 .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level
297 
298 .seealso: MGSetCycles()
299 @*/
300 PetscErrorCode MGSetCyclesOnLevel(PC pc,int l,int c)
301 {
302   MG *mg = (MG*)pc->data;
303 
304   PetscFunctionBegin;
305   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
306   mg[l]->cycles  = c;
307   PetscFunctionReturn(0);
308 }
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "MGSetRhs"
312 /*@
313    MGSetRhs - Sets the vector space to be used to store the right-hand side
314    on a particular level.  The user should free this space at the conclusion
315    of multigrid use.
316 
317    Collective on PC and Vec
318 
319    Input Parameters:
320 +  pc - the multigrid context
321 .  l  - the level (0 is coarsest) this is to be used for
322 -  c  - the space
323 
324    Level: advanced
325 
326 .keywords: MG, multigrid, set, right-hand-side, rhs, level
327 
328 .seealso: MGSetX(), MGSetR()
329 @*/
330 PetscErrorCode MGSetRhs(PC pc,int l,Vec c)
331 {
332   MG *mg = (MG*)pc->data;
333 
334   PetscFunctionBegin;
335   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
336   mg[l]->b  = c;
337   PetscFunctionReturn(0);
338 }
339 
340 #undef __FUNCT__
341 #define __FUNCT__ "MGSetX"
342 /*@
343    MGSetX - Sets the vector space to be used to store the solution on a
344    particular level.  The user should free this space at the conclusion
345    of multigrid use.
346 
347    Collective on PC and Vec
348 
349    Input Parameters:
350 +  pc - the multigrid context
351 .  l - the level (0 is coarsest) this is to be used for
352 -  c - the space
353 
354    Level: advanced
355 
356 .keywords: MG, multigrid, set, solution, level
357 
358 .seealso: MGSetRhs(), MGSetR()
359 @*/
360 PetscErrorCode MGSetX(PC pc,int l,Vec c)
361 {
362   MG *mg = (MG*)pc->data;
363 
364   PetscFunctionBegin;
365   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
366   mg[l]->x  = c;
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "MGSetR"
372 /*@
373    MGSetR - Sets the vector space to be used to store the residual on a
374    particular level.  The user should free this space at the conclusion of
375    multigrid use.
376 
377    Collective on PC and Vec
378 
379    Input Parameters:
380 +  pc - the multigrid context
381 .  l - the level (0 is coarsest) this is to be used for
382 -  c - the space
383 
384    Level: advanced
385 
386 .keywords: MG, multigrid, set, residual, level
387 @*/
388 PetscErrorCode MGSetR(PC pc,int l,Vec c)
389 {
390   MG *mg = (MG*)pc->data;
391 
392   PetscFunctionBegin;
393   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
394   mg[l]->r  = c;
395   PetscFunctionReturn(0);
396 }
397 
398 
399 
400 
401 
402 
403 
404 
405