xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1 /*
2     Defines the multigrid preconditioner interface.
3 */
4 #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
5 
6 
7 /*
8        MGMCycle_Private - Given an MG structure created with MGCreate() runs
9                   one multiplicative cycle down through the levels and
10                   back up.
11 
12     Input Parameter:
13 .   mg - structure created with  MGCreate().
14 */
15 #undef __FUNCT__
16 #define __FUNCT__ "MGMCycle_Private"
17 PetscErrorCode MGMCycle_Private(MG *mglevels,PetscTruth *converged)
18 {
19   MG          mg = *mglevels,mgc;
20   PetscErrorCode ierr;
21   int         cycles = mg->cycles;
22   PetscScalar zero = 0.0;
23 
24   PetscFunctionBegin;
25   if (converged) *converged = PETSC_FALSE;
26 
27   if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
28   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);
29   if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
30   if (mg->level) {  /* not the coarsest grid */
31     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
32 
33     /* if on finest level and have convergence criteria set */
34     if (mg->level == mg->levels-1 && mg->ttol) {
35       PetscReal rnorm;
36       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
37       if (rnorm <= mg->ttol) {
38         *converged = PETSC_TRUE;
39         if (rnorm < mg->atol) {
40           PetscLogInfo(0,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",rnorm,mg->atol);
41         } else {
42           PetscLogInfo(0,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",rnorm,mg->ttol);
43         }
44         PetscFunctionReturn(0);
45       }
46     }
47 
48     mgc = *(mglevels - 1);
49     ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr);
50     ierr = VecSet(&zero,mgc->x);CHKERRQ(ierr);
51     while (cycles--) {
52       ierr = MGMCycle_Private(mglevels-1,converged);CHKERRQ(ierr);
53     }
54     ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr);
55     if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
56     ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr);
57     if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
58   }
59   PetscFunctionReturn(0);
60 }
61 
62 /*
63        MGCreate_Private - Creates a MG structure for use with the
64                multigrid code. Level 0 is the coarsest. (But the
65                finest level is stored first in the array).
66 
67 */
68 #undef __FUNCT__
69 #define __FUNCT__ "MGCreate_Private"
70 static PetscErrorCode MGCreate_Private(MPI_Comm comm,int levels,PC pc,MPI_Comm *comms,MG **result)
71 {
72   MG   *mg;
73   PetscErrorCode ierr;
74   int  i,size;
75   char *prefix;
76   PC   ipc;
77 
78   PetscFunctionBegin;
79   ierr = PetscMalloc(levels*sizeof(MG),&mg);CHKERRQ(ierr);
80   PetscLogObjectMemory(pc,levels*(sizeof(MG)+sizeof(struct _MG)));
81 
82   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
83 
84   for (i=0; i<levels; i++) {
85     ierr = PetscNew(struct _MG,&mg[i]);CHKERRQ(ierr);
86     ierr = PetscMemzero(mg[i],sizeof(struct _MG));CHKERRQ(ierr);
87     mg[i]->level  = i;
88     mg[i]->levels = levels;
89     mg[i]->cycles = 1;
90 
91     if (comms) comm = comms[i];
92     ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr);
93     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr);
94     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
95 
96     /* do special stuff for coarse grid */
97     if (!i && levels > 1) {
98       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
99 
100       /* coarse solve is (redundant) LU by default */
101       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
102       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
103       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
104       if (size > 1) {
105         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
106         ierr = PCRedundantGetPC(ipc,&ipc);CHKERRQ(ierr);
107       }
108       ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
109 
110     } else {
111       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,"mg_levels_");CHKERRQ(ierr);
112     }
113     PetscLogObjectParent(pc,mg[i]->smoothd);
114     mg[i]->smoothu         = mg[i]->smoothd;
115     mg[i]->default_smoothu = 10000;
116     mg[i]->default_smoothd = 10000;
117     mg[i]->rtol = 0.0;
118     mg[i]->atol = 0.0;
119     mg[i]->dtol = 0.0;
120     mg[i]->ttol = 0.0;
121     mg[i]->eventsetup = 0;
122     mg[i]->eventsolve = 0;
123   }
124   *result = mg;
125   PetscFunctionReturn(0);
126 }
127 
128 #undef __FUNCT__
129 #define __FUNCT__ "PCDestroy_MG"
130 static PetscErrorCode PCDestroy_MG(PC pc)
131 {
132   MG  *mg = (MG*)pc->data;
133   PetscErrorCode ierr;
134   int i,n = mg[0]->levels;
135 
136   PetscFunctionBegin;
137   for (i=0; i<n; i++) {
138     if (mg[i]->smoothd != mg[i]->smoothu) {
139       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
140     }
141     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
142     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
143   }
144   ierr = PetscFree(mg);CHKERRQ(ierr);
145   PetscFunctionReturn(0);
146 }
147 
148 
149 
150 EXTERN PetscErrorCode MGACycle_Private(MG*);
151 EXTERN PetscErrorCode MGFCycle_Private(MG*);
152 EXTERN PetscErrorCode MGKCycle_Private(MG*);
153 
154 /*
155    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
156              or full cycle of multigrid.
157 
158   Note:
159   A simple wrapper which calls MGMCycle(),MGACycle(), or MGFCycle().
160 */
161 #undef __FUNCT__
162 #define __FUNCT__ "PCApply_MG"
163 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
164 {
165   MG          *mg = (MG*)pc->data;
166   PetscScalar zero = 0.0;
167   PetscErrorCode ierr;
168   int         levels = mg[0]->levels;
169 
170   PetscFunctionBegin;
171   mg[levels-1]->b = b;
172   mg[levels-1]->x = x;
173   if (mg[0]->am == MGMULTIPLICATIVE) {
174     ierr = VecSet(&zero,x);CHKERRQ(ierr);
175     ierr = MGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
176   }
177   else if (mg[0]->am == MGADDITIVE) {
178     ierr = MGACycle_Private(mg);CHKERRQ(ierr);
179   }
180   else if (mg[0]->am == MGKASKADE) {
181     ierr = MGKCycle_Private(mg);CHKERRQ(ierr);
182   }
183   else {
184     ierr = MGFCycle_Private(mg);CHKERRQ(ierr);
185   }
186   PetscFunctionReturn(0);
187 }
188 
189 #undef __FUNCT__
190 #define __FUNCT__ "PCApplyRichardson_MG"
191 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal atol, PetscReal dtol,int its)
192 {
193   MG         *mg = (MG*)pc->data;
194   PetscErrorCode ierr;
195   int levels = mg[0]->levels;
196   PetscTruth converged = PETSC_FALSE;
197 
198   PetscFunctionBegin;
199   mg[levels-1]->b    = b;
200   mg[levels-1]->x    = x;
201 
202   mg[levels-1]->rtol = rtol;
203   mg[levels-1]->atol = atol;
204   mg[levels-1]->dtol = dtol;
205   if (rtol) {
206     /* compute initial residual norm for relative convergence test */
207     PetscReal rnorm;
208     ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
209     ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
210     mg[levels-1]->ttol = PetscMax(rtol*rnorm,atol);
211   } else if (atol) {
212     mg[levels-1]->ttol = atol;
213   } else {
214     mg[levels-1]->ttol = 0.0;
215   }
216 
217   while (its-- && !converged) {
218     ierr = MGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr);
219   }
220   PetscFunctionReturn(0);
221 }
222 
223 #undef __FUNCT__
224 #define __FUNCT__ "PCSetFromOptions_MG"
225 static PetscErrorCode PCSetFromOptions_MG(PC pc)
226 {
227   PetscErrorCode ierr;
228   int indx,m,levels = 1;
229   PetscTruth flg;
230   const char *type[] = {"additive","multiplicative","full","cascade","kascade"};
231 
232   PetscFunctionBegin;
233 
234   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
235     if (!pc->data) {
236       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","MGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
237       ierr = MGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
238     }
239     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
240     if (flg) {
241       ierr = MGSetCycles(pc,m);CHKERRQ(ierr);
242     }
243     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
244     if (flg) {
245       ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
246     }
247     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
248     if (flg) {
249       ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
250     }
251     ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr);
252     if (flg) {
253       MGType mg = (MGType) 0;
254       switch (indx) {
255       case 0:
256         mg = MGADDITIVE;
257         break;
258       case 1:
259         mg = MGMULTIPLICATIVE;
260         break;
261       case 2:
262         mg = MGFULL;
263         break;
264       case 3:
265         mg = MGKASKADE;
266         break;
267       case 4:
268         mg = MGKASKADE;
269         break;
270       }
271       ierr = MGSetType(pc,mg);CHKERRQ(ierr);
272     }
273     ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr);
274     if (flg) {
275       MG   *mg = (MG*)pc->data;
276       int  i;
277       char eventname[128];
278       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
279       levels = mg[0]->levels;
280       for (i=0; i<levels; i++) {
281         sprintf(eventname,"MSetup Level %d",i);
282         ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr);
283         sprintf(eventname,"MGSolve Level %d",i);
284         ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr);
285       }
286     }
287   ierr = PetscOptionsTail();CHKERRQ(ierr);
288   PetscFunctionReturn(0);
289 }
290 
291 #undef __FUNCT__
292 #define __FUNCT__ "PCView_MG"
293 static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
294 {
295   MG         *mg = (MG*)pc->data;
296   PetscErrorCode ierr;
297   int levels = mg[0]->levels,i;
298   const char *cstring;
299   PetscTruth iascii;
300 
301   PetscFunctionBegin;
302   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
303   if (iascii) {
304     if (mg[0]->am == MGMULTIPLICATIVE) cstring = "multiplicative";
305     else if (mg[0]->am == MGADDITIVE)  cstring = "additive";
306     else if (mg[0]->am == MGFULL)      cstring = "full";
307     else if (mg[0]->am == MGKASKADE)   cstring = "Kaskade";
308     else cstring = "unknown";
309     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%d cycles=%d, pre-smooths=%d, post-smooths=%d\n",
310                       cstring,levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
311     for (i=0; i<levels; i++) {
312       ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %d -------------------------------\n",i);CHKERRQ(ierr);
313       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
314       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
315       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
316       if (mg[i]->smoothd == mg[i]->smoothu) {
317         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
318       } else {
319         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %d -------------------------------\n",i);CHKERRQ(ierr);
320         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
321         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
322         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
323       }
324     }
325   } else {
326     SETERRQ1(1,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
327   }
328   PetscFunctionReturn(0);
329 }
330 
331 /*
332     Calls setup for the KSP on each level
333 */
334 #undef __FUNCT__
335 #define __FUNCT__ "PCSetUp_MG"
336 static PetscErrorCode PCSetUp_MG(PC pc)
337 {
338   MG          *mg = (MG*)pc->data;
339   PetscErrorCode ierr;
340   int i,n = mg[0]->levels;
341   PC          cpc;
342   PetscTruth  preonly,lu,redundant,monitor = PETSC_FALSE,dump;
343   PetscViewer ascii;
344   MPI_Comm    comm;
345 
346   PetscFunctionBegin;
347   if (!pc->setupcalled) {
348     ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr);
349 
350     for (i=1; i<n; i++) {
351       if (mg[i]->smoothd) {
352         if (monitor) {
353           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
354           ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
355           ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
356           ierr = KSPSetMonitor(mg[i]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
357         }
358         ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
359       }
360     }
361     for (i=0; i<n; i++) {
362       if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
363         if (monitor) {
364           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
365           ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
366           ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
367           ierr = KSPSetMonitor(mg[i]->smoothu,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
368         }
369         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
370       }
371     }
372   }
373 
374   for (i=1; i<n; i++) {
375     if (mg[i]->smoothd) {
376       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
377       if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
378       ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
379       if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
380     }
381   }
382   for (i=0; i<n; i++) {
383     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
384       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
385       if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
386       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
387       if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
388     }
389   }
390 
391   /*
392       If coarse solver is not direct method then DO NOT USE preonly
393   */
394   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
395   if (preonly) {
396     ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
397     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
398     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
399     if (!lu && !redundant) {
400       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
401     }
402   }
403 
404   if (!pc->setupcalled) {
405     if (monitor) {
406       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
407       ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
408       ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr);
409       ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
410     }
411     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
412   }
413 
414   if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
415   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
416   if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
417 
418   /*
419      Dump the interpolation/restriction matrices to matlab plus the
420    Jacobian/stiffness on each level. This allows Matlab users to
421    easily check if the Galerkin condition A_c = R A_f R^T is satisfied */
422   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
423   if (dump) {
424     for (i=1; i<n; i++) {
425       ierr = MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
426     }
427     for (i=0; i<n; i++) {
428       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
429       ierr = MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
430     }
431   }
432 
433   PetscFunctionReturn(0);
434 }
435 
436 /* -------------------------------------------------------------------------------------*/
437 
438 #undef __FUNCT__
439 #define __FUNCT__ "MGSetLevels"
440 /*@C
441    MGSetLevels - Sets the number of levels to use with MG.
442    Must be called before any other MG routine.
443 
444    Collective on PC
445 
446    Input Parameters:
447 +  pc - the preconditioner context
448 .  levels - the number of levels
449 -  comms - optional communicators for each level; this is to allow solving the coarser problems
450            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
451 
452    Level: intermediate
453 
454    Notes:
455      If the number of levels is one then the multigrid uses the -mg_levels prefix
456   for setting the level options rather than the -mg_coarse prefix.
457 
458 .keywords: MG, set, levels, multigrid
459 
460 .seealso: MGSetType(), MGGetLevels()
461 @*/
462 PetscErrorCode MGSetLevels(PC pc,int levels,MPI_Comm *comms)
463 {
464   PetscErrorCode ierr;
465   MG  *mg;
466 
467   PetscFunctionBegin;
468   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
469 
470   if (pc->data) {
471     SETERRQ(1,"Number levels already set for MG\n\
472     make sure that you call MGSetLevels() before KSPSetFromOptions()");
473   }
474   ierr                     = MGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
475   mg[0]->am                = MGMULTIPLICATIVE;
476   pc->data                 = (void*)mg;
477   pc->ops->applyrichardson = PCApplyRichardson_MG;
478   PetscFunctionReturn(0);
479 }
480 
481 #undef __FUNCT__
482 #define __FUNCT__ "MGGetLevels"
483 /*@
484    MGGetLevels - Gets the number of levels to use with MG.
485 
486    Not Collective
487 
488    Input Parameter:
489 .  pc - the preconditioner context
490 
491    Output parameter:
492 .  levels - the number of levels
493 
494    Level: advanced
495 
496 .keywords: MG, get, levels, multigrid
497 
498 .seealso: MGSetLevels()
499 @*/
500 PetscErrorCode MGGetLevels(PC pc,int *levels)
501 {
502   MG  *mg;
503 
504   PetscFunctionBegin;
505   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
506   PetscValidIntPointer(levels,2);
507 
508   mg      = (MG*)pc->data;
509   *levels = mg[0]->levels;
510   PetscFunctionReturn(0);
511 }
512 
513 #undef __FUNCT__
514 #define __FUNCT__ "MGSetType"
515 /*@
516    MGSetType - Determines the form of multigrid to use:
517    multiplicative, additive, full, or the Kaskade algorithm.
518 
519    Collective on PC
520 
521    Input Parameters:
522 +  pc - the preconditioner context
523 -  form - multigrid form, one of MGMULTIPLICATIVE, MGADDITIVE,
524    MGFULL, MGKASKADE
525 
526    Options Database Key:
527 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
528    additive, full, kaskade
529 
530    Level: advanced
531 
532 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
533 
534 .seealso: MGSetLevels()
535 @*/
536 PetscErrorCode MGSetType(PC pc,MGType form)
537 {
538   MG *mg;
539 
540   PetscFunctionBegin;
541   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
542   mg = (MG*)pc->data;
543 
544   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
545   mg[0]->am = form;
546   if (form == MGMULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
547   else pc->ops->applyrichardson = 0;
548   PetscFunctionReturn(0);
549 }
550 
551 #undef __FUNCT__
552 #define __FUNCT__ "MGSetCycles"
553 /*@
554    MGSetCycles - Sets the type cycles to use.  Use MGSetCyclesOnLevel() for more
555    complicated cycling.
556 
557    Collective on PC
558 
559    Input Parameters:
560 +  mg - the multigrid context
561 -  n - the number of cycles
562 
563    Options Database Key:
564 $  -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle.
565 
566    Level: advanced
567 
568 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
569 
570 .seealso: MGSetCyclesOnLevel()
571 @*/
572 PetscErrorCode MGSetCycles(PC pc,int n)
573 {
574   MG  *mg;
575   int i,levels;
576 
577   PetscFunctionBegin;
578   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
579   mg     = (MG*)pc->data;
580   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
581   levels = mg[0]->levels;
582 
583   for (i=0; i<levels; i++) {
584     mg[i]->cycles  = n;
585   }
586   PetscFunctionReturn(0);
587 }
588 
589 #undef __FUNCT__
590 #define __FUNCT__ "MGCheck"
591 /*@
592    MGCheck - Checks that all components of the MG structure have
593    been set.
594 
595    Collective on PC
596 
597    Input Parameters:
598 .  mg - the MG structure
599 
600    Level: advanced
601 
602 .keywords: MG, check, set, multigrid
603 @*/
604 PetscErrorCode MGCheck(PC pc)
605 {
606   MG  *mg;
607   int i,n,count = 0;
608 
609   PetscFunctionBegin;
610   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
611   mg = (MG*)pc->data;
612 
613   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
614 
615   n = mg[0]->levels;
616 
617   for (i=1; i<n; i++) {
618     if (!mg[i]->restrct) {
619       (*PetscErrorPrintf)("No restrict set level %d \n",n-i); count++;
620     }
621     if (!mg[i]->interpolate) {
622       (*PetscErrorPrintf)("No interpolate set level %d \n",n-i); count++;
623     }
624     if (!mg[i]->residual) {
625       (*PetscErrorPrintf)("No residual set level %d \n",n-i); count++;
626     }
627     if (!mg[i]->smoothu) {
628       (*PetscErrorPrintf)("No smoothup set level %d \n",n-i); count++;
629     }
630     if (!mg[i]->smoothd) {
631       (*PetscErrorPrintf)("No smoothdown set level %d \n",n-i); count++;
632     }
633     if (!mg[i]->r) {
634       (*PetscErrorPrintf)("No r set level %d \n",n-i); count++;
635     }
636     if (!mg[i-1]->x) {
637       (*PetscErrorPrintf)("No x set level %d \n",n-i); count++;
638     }
639     if (!mg[i-1]->b) {
640       (*PetscErrorPrintf)("No b set level %d \n",n-i); count++;
641     }
642   }
643   PetscFunctionReturn(count);
644 }
645 
646 
647 #undef __FUNCT__
648 #define __FUNCT__ "MGSetNumberSmoothDown"
649 /*@
650    MGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
651    use on all levels. Use MGGetSmootherDown() to set different
652    pre-smoothing steps on different levels.
653 
654    Collective on PC
655 
656    Input Parameters:
657 +  mg - the multigrid context
658 -  n - the number of smoothing steps
659 
660    Options Database Key:
661 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
662 
663    Level: advanced
664 
665 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
666 
667 .seealso: MGSetNumberSmoothUp()
668 @*/
669 PetscErrorCode MGSetNumberSmoothDown(PC pc,int n)
670 {
671   MG  *mg;
672   PetscErrorCode ierr;
673   int i,levels;
674 
675   PetscFunctionBegin;
676   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
677   mg     = (MG*)pc->data;
678   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
679   levels = mg[0]->levels;
680 
681   for (i=0; i<levels; i++) {
682     /* make sure smoother up and down are different */
683     ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
684     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
685     mg[i]->default_smoothd = n;
686   }
687   PetscFunctionReturn(0);
688 }
689 
690 #undef __FUNCT__
691 #define __FUNCT__ "MGSetNumberSmoothUp"
692 /*@
693    MGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
694    on all levels. Use MGGetSmootherUp() to set different numbers of
695    post-smoothing steps on different levels.
696 
697    Collective on PC
698 
699    Input Parameters:
700 +  mg - the multigrid context
701 -  n - the number of smoothing steps
702 
703    Options Database Key:
704 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
705 
706    Level: advanced
707 
708    Note: this does not set a value on the coarsest grid, since we assume that
709     there is no seperate smooth up on the coarsest grid. If you want to have a
710     seperate smooth up on the coarsest grid then call MGGetSmoothUp(pc,0,&ksp);
711 
712 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
713 
714 .seealso: MGSetNumberSmoothDown()
715 @*/
716 PetscErrorCode  MGSetNumberSmoothUp(PC pc,int n)
717 {
718   MG  *mg;
719   PetscErrorCode ierr;
720   int i,levels;
721 
722   PetscFunctionBegin;
723   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
724   mg     = (MG*)pc->data;
725   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
726   levels = mg[0]->levels;
727 
728   for (i=1; i<levels; i++) {
729     /* make sure smoother up and down are different */
730     ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
731     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
732     mg[i]->default_smoothu = n;
733   }
734   PetscFunctionReturn(0);
735 }
736 
737 /* ----------------------------------------------------------------------------------------*/
738 
739 /*MC
740    PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional
741     information about the coarser grid matrices and restriction/interpolation operators.
742 
743    Options Database Keys:
744 +  -pc_mg_levels <nlevels> - number of levels including finest
745 .  -pc_mg_cycles 1 or 2 - for V or W-cycle
746 .  -pc_mg_smoothup <n> - number of smoothing steps before interpolation
747 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
748 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
749 .  -pc_mg_log - log information about time spent on each level of the solver
750 .  -pc_mg_monitor - print information on the multigrid convergence
751 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
752                         to the Socket viewer for reading from Matlab.
753 
754    Notes:
755 
756    Level: intermediate
757 
758    Concepts: multigrid
759 
760 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
761            MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(),
762            MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(),
763            MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(),
764            MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR()
765 M*/
766 
767 EXTERN_C_BEGIN
768 #undef __FUNCT__
769 #define __FUNCT__ "PCCreate_MG"
770 PetscErrorCode PCCreate_MG(PC pc)
771 {
772   PetscFunctionBegin;
773   pc->ops->apply          = PCApply_MG;
774   pc->ops->setup          = PCSetUp_MG;
775   pc->ops->destroy        = PCDestroy_MG;
776   pc->ops->setfromoptions = PCSetFromOptions_MG;
777   pc->ops->view           = PCView_MG;
778 
779   pc->data                = (void*)0;
780   PetscFunctionReturn(0);
781 }
782 EXTERN_C_END
783