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