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