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