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