xref: /petsc/src/ksp/pc/impls/mg/mgfunc.c (revision dba47a550923b04c7c4ebbb735eb62a1b3e4e9ae)
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 .keywords:  multigrid, set, interpolate, level
121 
122 .seealso: MGSetRestriction()
123 @*/
124 PetscErrorCode PETSCKSP_DLLEXPORT MGSetInterpolate(PC pc,PetscInt l,Mat mat)
125 {
126   MG *mg = (MG*)pc->data;
127 
128   PetscFunctionBegin;
129   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
130   mg[l]->interpolate = mat;
131   PetscFunctionReturn(0);
132 }
133 
134 #undef __FUNCT__
135 #define __FUNCT__ "MGSetRestriction"
136 /*@
137    MGSetRestriction - Sets the function to be used to restrict vector
138    from level l to l-1.
139 
140    Collective on PC and Mat
141 
142    Input Parameters:
143 +  pc - the multigrid context
144 .  mat - the restriction matrix
145 -  l - the level (0 is coarsest) to supply
146 
147    Level: advanced
148 
149    Notes:
150           Usually this is the same matrix used also to set the interpolation
151     for the same level.
152 
153           One can pass in the interpolation matrix or its transpose; PETSc figures
154     out from the matrix size which one it is.
155 
156 .keywords: MG, set, multigrid, restriction, level
157 
158 .seealso: MGSetInterpolate()
159 @*/
160 PetscErrorCode PETSCKSP_DLLEXPORT MGSetRestriction(PC pc,PetscInt l,Mat mat)
161 {
162   MG *mg = (MG*)pc->data;
163 
164   PetscFunctionBegin;
165   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
166   mg[l]->restrct  = mat;
167   PetscFunctionReturn(0);
168 }
169 
170 #undef __FUNCT__
171 #define __FUNCT__ "MGGetSmoother"
172 /*@C
173    MGGetSmoother - Gets the KSP context to be used as smoother for
174    both pre- and post-smoothing.  Call both MGGetSmootherUp() and
175    MGGetSmootherDown() to use different functions for pre- and
176    post-smoothing.
177 
178    Not Collective, KSP returned is parallel if PC is
179 
180    Input Parameters:
181 +  pc - the multigrid context
182 -  l - the level (0 is coarsest) to supply
183 
184    Ouput Parameters:
185 .  ksp - the smoother
186 
187    Level: advanced
188 
189 .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother
190 
191 .seealso: MGGetSmootherUp(), MGGetSmootherDown()
192 @*/
193 PetscErrorCode PETSCKSP_DLLEXPORT MGGetSmoother(PC pc,PetscInt l,KSP *ksp)
194 {
195   MG *mg = (MG*)pc->data;
196 
197   PetscFunctionBegin;
198   *ksp = mg[l]->smoothd;
199   PetscFunctionReturn(0);
200 }
201 
202 #undef __FUNCT__
203 #define __FUNCT__ "MGGetSmootherUp"
204 /*@C
205    MGGetSmootherUp - Gets the KSP context to be used as smoother after
206    coarse grid correction (post-smoother).
207 
208    Not Collective, KSP returned is parallel if PC is
209 
210    Input Parameters:
211 +  pc - the multigrid context
212 -  l  - the level (0 is coarsest) to supply
213 
214    Ouput Parameters:
215 .  ksp - the smoother
216 
217    Level: advanced
218 
219 .keywords: MG, multigrid, get, smoother, up, post-smoother, level
220 
221 .seealso: MGGetSmootherUp(), MGGetSmootherDown()
222 @*/
223 PetscErrorCode PETSCKSP_DLLEXPORT MGGetSmootherUp(PC pc,PetscInt l,KSP *ksp)
224 {
225   MG             *mg = (MG*)pc->data;
226   PetscErrorCode ierr;
227   char           *prefix;
228   MPI_Comm       comm;
229 
230   PetscFunctionBegin;
231   /*
232      This is called only if user wants a different pre-smoother from post.
233      Thus we check if a different one has already been allocated,
234      if not we allocate it.
235   */
236   if (!l) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid");
237   if (mg[l]->smoothu == mg[l]->smoothd) {
238     ierr = PetscObjectGetComm((PetscObject)mg[l]->smoothd,&comm);CHKERRQ(ierr);
239     ierr = KSPGetOptionsPrefix(mg[l]->smoothd,&prefix);CHKERRQ(ierr);
240     ierr = KSPCreate(comm,&mg[l]->smoothu);CHKERRQ(ierr);
241     ierr = KSPSetTolerances(mg[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
242     ierr = KSPSetOptionsPrefix(mg[l]->smoothu,prefix);CHKERRQ(ierr);
243     ierr = PetscLogObjectParent(pc,mg[l]->smoothu);CHKERRQ(ierr);
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 PETSCKSP_DLLEXPORT MGGetSmootherDown(PC pc,PetscInt 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 PETSCKSP_DLLEXPORT MGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt 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 PETSCKSP_DLLEXPORT MGSetRhs(PC pc,PetscInt 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 PETSCKSP_DLLEXPORT MGSetX(PC pc,PetscInt 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 PETSCKSP_DLLEXPORT MGSetR(PC pc,PetscInt 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