xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 6e111a19f6677190c8cb13236301fcb65e0e3d3b)
1 
2 /*
3     Defines the multigrid preconditioner interface.
4 */
5 #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscpcmg.h" I*/
6 
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "PCMGMCycle_Private"
10 PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
11 {
12   PC_MG          *mg = (PC_MG*)pc->data;
13   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
14   PetscErrorCode ierr;
15   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
16 
17   PetscFunctionBegin;
18   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
19   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
20   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
21   if (mglevels->level) {  /* not the coarsest grid */
22     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
23     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
24     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
25 
26     /* if on finest level and have convergence criteria set */
27     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
28       PetscReal rnorm;
29       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
30       if (rnorm <= mg->ttol) {
31         if (rnorm < mg->abstol) {
32           *reason = PCRICHARDSON_CONVERGED_ATOL;
33           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
34         } else {
35           *reason = PCRICHARDSON_CONVERGED_RTOL;
36           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
37         }
38         PetscFunctionReturn(0);
39       }
40     }
41 
42     mgc = *(mglevelsin - 1);
43     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
44     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
45     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
46     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
47     while (cycles--) {
48       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
49     }
50     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
51     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
52     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
53     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
54     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
55     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
56   }
57   PetscFunctionReturn(0);
58 }
59 
60 #undef __FUNCT__
61 #define __FUNCT__ "PCApplyRichardson_MG"
62 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool  zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
63 {
64   PC_MG          *mg = (PC_MG*)pc->data;
65   PC_MG_Levels   **mglevels = mg->levels;
66   PetscErrorCode ierr;
67   PetscInt       levels = mglevels[0]->levels,i;
68 
69   PetscFunctionBegin;
70   mglevels[levels-1]->b    = b;
71   mglevels[levels-1]->x    = x;
72 
73   mg->rtol = rtol;
74   mg->abstol = abstol;
75   mg->dtol = dtol;
76   if (rtol) {
77     /* compute initial residual norm for relative convergence test */
78     PetscReal rnorm;
79     if (zeroguess) {
80       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
81     } else {
82       ierr               = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
83       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
84     }
85     mg->ttol = PetscMax(rtol*rnorm,abstol);
86   } else if (abstol) {
87     mg->ttol = abstol;
88   } else {
89     mg->ttol = 0.0;
90   }
91 
92   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
93      stop prematurely due to small residual */
94   for (i=1; i<levels; i++) {
95     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
96     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
97       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
98     }
99   }
100 
101   *reason = (PCRichardsonConvergedReason)0;
102   for (i=0; i<its; i++) {
103     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
104     if (*reason) break;
105   }
106   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
107   *outits = i;
108   PetscFunctionReturn(0);
109 }
110 
111 #undef __FUNCT__
112 #define __FUNCT__ "PCReset_MG"
113 PetscErrorCode PCReset_MG(PC pc)
114 {
115   PC_MG          *mg = (PC_MG*)pc->data;
116   PC_MG_Levels   **mglevels = mg->levels;
117   PetscErrorCode ierr;
118   PetscInt       i,n;
119 
120   PetscFunctionBegin;
121   if (mglevels) {
122     n = mglevels[0]->levels;
123     for (i=0; i<n-1; i++) {
124       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
125       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
126       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
127       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
128       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
129       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
130     }
131 
132     for (i=0; i<n; i++) {
133       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
134       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
135         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
136       }
137       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
138     }
139   }
140   PetscFunctionReturn(0);
141 }
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "PCMGSetLevels"
145 /*@C
146    PCMGSetLevels - Sets the number of levels to use with MG.
147    Must be called before any other MG routine.
148 
149    Logically Collective on PC
150 
151    Input Parameters:
152 +  pc - the preconditioner context
153 .  levels - the number of levels
154 -  comms - optional communicators for each level; this is to allow solving the coarser problems
155            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
156 
157    Level: intermediate
158 
159    Notes:
160      If the number of levels is one then the multigrid uses the -mg_levels prefix
161   for setting the level options rather than the -mg_coarse prefix.
162 
163 .keywords: MG, set, levels, multigrid
164 
165 .seealso: PCMGSetType(), PCMGGetLevels()
166 @*/
167 PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
168 {
169   PetscErrorCode ierr;
170   PC_MG          *mg = (PC_MG*)pc->data;
171   MPI_Comm       comm = ((PetscObject)pc)->comm;
172   PC_MG_Levels   **mglevels = mg->levels;
173   PetscInt       i;
174   PetscMPIInt    size;
175   const char     *prefix;
176   PC             ipc;
177   PetscInt       n;
178 
179   PetscFunctionBegin;
180   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
181   PetscValidLogicalCollectiveInt(pc,levels,2);
182   if (mg->nlevels == levels) PetscFunctionReturn(0);
183   if (mglevels) {
184     /* changing the number of levels so free up the previous stuff */
185     ierr = PCReset_MG(pc);CHKERRQ(ierr);
186     n = mglevels[0]->levels;
187     for (i=0; i<n; i++) {
188       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
189         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
190       }
191       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
192       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
193     }
194     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
195   }
196 
197   mg->nlevels      = levels;
198 
199   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr);
200   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
201 
202   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
203 
204   mg->stageApply = 0;
205   for (i=0; i<levels; i++) {
206     ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr);
207     mglevels[i]->level           = i;
208     mglevels[i]->levels          = levels;
209     mglevels[i]->cycles          = PC_MG_CYCLE_V;
210     mg->default_smoothu = 2;
211     mg->default_smoothd = 2;
212     mglevels[i]->eventsmoothsetup    = 0;
213     mglevels[i]->eventsmoothsolve    = 0;
214     mglevels[i]->eventresidual       = 0;
215     mglevels[i]->eventinterprestrict = 0;
216 
217     if (comms) comm = comms[i];
218     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
219     ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
220     ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
221     ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
222     ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
223     ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
224     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
225     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i?mg->default_smoothd:1);CHKERRQ(ierr);
226     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
227 
228     /* do special stuff for coarse grid */
229     if (!i && levels > 1) {
230       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
231 
232       /* coarse solve is (redundant) LU by default */
233       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
234       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
235       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
236       if (size > 1) {
237         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
238       } else {
239         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
240       }
241 
242     } else {
243       char tprefix[128];
244       sprintf(tprefix,"mg_levels_%d_",(int)i);
245       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
246     }
247     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
248     mglevels[i]->smoothu    = mglevels[i]->smoothd;
249     mg->rtol                = 0.0;
250     mg->abstol              = 0.0;
251     mg->dtol                = 0.0;
252     mg->ttol                = 0.0;
253     mg->cyclesperpcapply    = 1;
254   }
255   mg->am          = PC_MG_MULTIPLICATIVE;
256   mg->levels      = mglevels;
257   pc->ops->applyrichardson = PCApplyRichardson_MG;
258   PetscFunctionReturn(0);
259 }
260 
261 
262 #undef __FUNCT__
263 #define __FUNCT__ "PCDestroy_MG"
264 PetscErrorCode PCDestroy_MG(PC pc)
265 {
266   PetscErrorCode ierr;
267   PC_MG          *mg = (PC_MG*)pc->data;
268   PC_MG_Levels   **mglevels = mg->levels;
269   PetscInt       i,n;
270 
271   PetscFunctionBegin;
272   ierr = PCReset_MG(pc);CHKERRQ(ierr);
273   if (mglevels) {
274     n = mglevels[0]->levels;
275     for (i=0; i<n; i++) {
276       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
277         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
278       }
279       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
280       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
281     }
282     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
283   }
284   ierr = PetscFree(pc->data);CHKERRQ(ierr);
285   PetscFunctionReturn(0);
286 }
287 
288 
289 
290 extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
291 extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
292 extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
293 
294 /*
295    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
296              or full cycle of multigrid.
297 
298   Note:
299   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
300 */
301 #undef __FUNCT__
302 #define __FUNCT__ "PCApply_MG"
303 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
304 {
305   PC_MG          *mg = (PC_MG*)pc->data;
306   PC_MG_Levels   **mglevels = mg->levels;
307   PetscErrorCode ierr;
308   PetscInt       levels = mglevels[0]->levels,i;
309 
310   PetscFunctionBegin;
311   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
312   /* When the DM is supplying the matrix then it will not exist until here */
313   for (i=0; i<levels; i++) {
314     if (!mglevels[i]->A) {
315       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
316       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
317     }
318   }
319 
320   mglevels[levels-1]->b = b;
321   mglevels[levels-1]->x = x;
322   if (mg->am == PC_MG_MULTIPLICATIVE) {
323     ierr = VecSet(x,0.0);CHKERRQ(ierr);
324     for (i=0; i<mg->cyclesperpcapply; i++) {
325       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
326     }
327   }
328   else if (mg->am == PC_MG_ADDITIVE) {
329     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
330   }
331   else if (mg->am == PC_MG_KASKADE) {
332     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
333   }
334   else {
335     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
336   }
337   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
338   PetscFunctionReturn(0);
339 }
340 
341 
342 #undef __FUNCT__
343 #define __FUNCT__ "PCSetFromOptions_MG"
344 PetscErrorCode PCSetFromOptions_MG(PC pc)
345 {
346   PetscErrorCode ierr;
347   PetscInt       m,levels = 1,cycles;
348   PetscBool      flg,set;
349   PC_MG          *mg = (PC_MG*)pc->data;
350   PC_MG_Levels   **mglevels = mg->levels;
351   PCMGType       mgtype;
352   PCMGCycleType  mgctype;
353 
354   PetscFunctionBegin;
355   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
356     if (!mg->levels) {
357       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
358       if (!flg && pc->dm) {
359         ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
360         levels++;
361         mg->usedmfornumberoflevels = PETSC_TRUE;
362       }
363       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
364     }
365     mglevels = mg->levels;
366 
367     mgctype = (PCMGCycleType) mglevels[0]->cycles;
368     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
369     if (flg) {
370       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
371     };
372     flg  = PETSC_FALSE;
373     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
374     if (set) {
375       ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
376     }
377     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
378     if (flg) {
379       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
380     }
381     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
382     if (flg) {
383       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
384     }
385     mgtype = mg->am;
386     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
387     if (flg) {
388       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
389     }
390     if (mg->am == PC_MG_MULTIPLICATIVE) {
391       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
392       if (flg) {
393         ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
394       }
395     }
396     flg  = PETSC_FALSE;
397     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
398     if (flg) {
399       PetscInt i;
400       char     eventname[128];
401       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
402       levels = mglevels[0]->levels;
403       for (i=0; i<levels; i++) {
404         sprintf(eventname,"MGSetup Level %d",(int)i);
405         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
406         sprintf(eventname,"MGSmooth Level %d",(int)i);
407         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
408         if (i) {
409           sprintf(eventname,"MGResid Level %d",(int)i);
410           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
411           sprintf(eventname,"MGInterp Level %d",(int)i);
412           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
413         }
414       }
415 
416 #if defined(PETSC_USE_LOG)
417       {
418         const char   *sname = "MG Apply";
419         PetscStageLog stageLog;
420         PetscInt      st;
421 
422         PetscFunctionBegin;
423         ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
424         for (st = 0; st < stageLog->numStages; ++st) {
425           PetscBool same;
426 
427           ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
428           if (same) {mg->stageApply = st;}
429         }
430         if (!mg->stageApply) {
431           ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
432         }
433       }
434 #endif
435     }
436   ierr = PetscOptionsTail();CHKERRQ(ierr);
437   PetscFunctionReturn(0);
438 }
439 
440 const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
441 const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
442 
443 #undef __FUNCT__
444 #define __FUNCT__ "PCView_MG"
445 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
446 {
447   PC_MG          *mg = (PC_MG*)pc->data;
448   PC_MG_Levels   **mglevels = mg->levels;
449   PetscErrorCode ierr;
450   PetscInt       levels = mglevels[0]->levels,i;
451   PetscBool      iascii,isbinary,isdraw;
452 
453   PetscFunctionBegin;
454   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
455   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
456   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
457   if (iascii) {
458     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,(mglevels[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr);
459     if (mg->am == PC_MG_MULTIPLICATIVE) {
460       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
461     }
462     if (mg->galerkin) {
463       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
464     } else {
465       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
466     }
467     for (i=0; i<levels; i++) {
468       if (!i) {
469         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
470       } else {
471         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
472       }
473       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
474       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
475       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
476       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
477         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
478       } else if (i) {
479         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
480         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
481         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
482         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
483       }
484     }
485   } else if (isbinary) {
486     for (i=levels-1; i>=0; i--) {
487       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
488       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
489         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
490       }
491     }
492   } else if (isdraw) {
493     PetscDraw draw;
494     PetscReal x,w,y,bottom,th;
495     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
496     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
497     ierr = PetscDrawStringGetSize(draw,PETSC_NULL,&th);CHKERRQ(ierr);
498     bottom = y - th;
499     for (i=levels-1; i>=0; i--) {
500       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
501         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
502         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
503         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
504       } else {
505         w = 0.5*PetscMin(1.0-x,x);
506         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
507         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
508         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
509         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
510         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
511         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
512       }
513       ierr = PetscDrawGetBoundingBox(draw,PETSC_NULL,&bottom,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
514       bottom -= th;
515     }
516   }
517   PetscFunctionReturn(0);
518 }
519 
520 #include <petsc-private/dmimpl.h>
521 #include <petsc-private/kspimpl.h>
522 
523 /*
524     Calls setup for the KSP on each level
525 */
526 #undef __FUNCT__
527 #define __FUNCT__ "PCSetUp_MG"
528 PetscErrorCode PCSetUp_MG(PC pc)
529 {
530   PC_MG                   *mg = (PC_MG*)pc->data;
531   PC_MG_Levels            **mglevels = mg->levels;
532   PetscErrorCode          ierr;
533   PetscInt                i,n = mglevels[0]->levels;
534   PC                      cpc;
535   PetscBool               preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset;
536   Mat                     dA,dB;
537   MatStructure            uflag;
538   Vec                     tvec;
539   DM                      *dms;
540   PetscViewer             viewer = 0;
541 
542   PetscFunctionBegin;
543   /* FIX: Move this to PCSetFromOptions_MG? */
544   if (mg->usedmfornumberoflevels) {
545     PetscInt levels;
546     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
547     levels++;
548     if (levels > n) { /* the problem is now being solved on a finer grid */
549       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
550       n    = levels;
551       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
552       mglevels =  mg->levels;
553     }
554   }
555   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
556 
557 
558   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
559   /* so use those from global PC */
560   /* Is this what we always want? What if user wants to keep old one? */
561   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
562   if (opsset) {
563     Mat mmat;
564     ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr);
565     if (mmat == pc->pmat) opsset = PETSC_FALSE;
566   }
567   if (!opsset) {
568     ierr = PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr);
569     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
570   }
571 
572   /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */
573   if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
574     /* construct the interpolation from the DMs */
575     Mat p;
576     Vec rscale;
577     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
578     dms[n-1] = pc->dm;
579     for (i=n-2; i>-1; i--) {
580       DMKSP kdm;
581       ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);
582       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
583       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
584       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
585       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
586        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
587       kdm->ops->computerhs = PETSC_NULL;
588       kdm->rhsctx          = PETSC_NULL;
589       if (!mglevels[i+1]->interpolate) {
590         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
591         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
592         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
593         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
594         ierr = MatDestroy(&p);CHKERRQ(ierr);
595       }
596     }
597 
598     for (i=n-2; i>-1; i--) {
599       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
600     }
601     ierr = PetscFree(dms);CHKERRQ(ierr);
602   }
603 
604   if (pc->dm && !pc->setupcalled) {
605     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
606     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
607     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
608   }
609 
610   if (mg->galerkin == 1) {
611     Mat B;
612     /* currently only handle case where mat and pmat are the same on coarser levels */
613     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
614     if (!pc->setupcalled) {
615       for (i=n-2; i>-1; i--) {
616         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
617         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
618         if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
619         dB   = B;
620       }
621       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
622     } else {
623       for (i=n-2; i>-1; i--) {
624         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
625         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
626         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
627         dB   = B;
628       }
629     }
630   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
631     /* need to restrict Jacobian location to coarser meshes for evaluation */
632     for (i=n-2;i>-1; i--) {
633       Mat R;
634       Vec rscale;
635       if (!mglevels[i]->smoothd->dm->x) {
636         Vec *vecs;
637         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr);
638         mglevels[i]->smoothd->dm->x = vecs[0];
639         ierr = PetscFree(vecs);CHKERRQ(ierr);
640       }
641       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
642       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
643       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
644       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
645     }
646   }
647   if (!mg->galerkin && pc->dm) {
648     for (i=n-2;i>=0; i--) {
649       DM dmfine,dmcoarse;
650       Mat Restrict,Inject;
651       Vec rscale;
652       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
653       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
654       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
655       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
656       Inject = PETSC_NULL;      /* Callback should create it if it needs Injection */
657       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
658     }
659   }
660 
661   if (!pc->setupcalled) {
662     for (i=0; i<n; i++) {
663       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
664     }
665     for (i=1; i<n; i++) {
666       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
667         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
668       }
669     }
670     for (i=1; i<n; i++) {
671       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
672       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
673     }
674     for (i=0; i<n-1; i++) {
675       if (!mglevels[i]->b) {
676         Vec *vec;
677         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
678         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
679         ierr = VecDestroy(vec);CHKERRQ(ierr);
680         ierr = PetscFree(vec);CHKERRQ(ierr);
681       }
682       if (!mglevels[i]->r && i) {
683         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
684         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
685         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
686       }
687       if (!mglevels[i]->x) {
688         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
689         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
690         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
691       }
692     }
693     if (n != 1 && !mglevels[n-1]->r) {
694       /* PCMGSetR() on the finest level if user did not supply it */
695       Vec *vec;
696       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
697       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
698       ierr = VecDestroy(vec);CHKERRQ(ierr);
699       ierr = PetscFree(vec);CHKERRQ(ierr);
700     }
701   }
702 
703   if (pc->dm) {
704     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
705     for (i=0; i<n-1; i++) {
706       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
707     }
708   }
709 
710   for (i=1; i<n; i++) {
711     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
712       /* if doing only down then initial guess is zero */
713       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
714     }
715     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
716     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
717     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
718     if (!mglevels[i]->residual) {
719       Mat mat;
720       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
721       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
722     }
723   }
724   for (i=1; i<n; i++) {
725     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
726       Mat          downmat,downpmat;
727       MatStructure matflag;
728 
729       /* check if operators have been set for up, if not use down operators to set them */
730       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
731       if (!opsset) {
732         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
733         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
734       }
735 
736       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
737       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
738       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
739       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
740     }
741   }
742 
743   /*
744       If coarse solver is not direct method then DO NOT USE preonly
745   */
746   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
747   if (preonly) {
748     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
749     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
750     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
751     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
752     if (!lu && !redundant && !cholesky && !svd) {
753       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
754     }
755   }
756 
757   if (!pc->setupcalled) {
758     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
759   }
760 
761   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
762   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
763   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
764 
765   /*
766      Dump the interpolation/restriction matrices plus the
767    Jacobian/stiffness on each level. This allows MATLAB users to
768    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
769 
770    Only support one or the other at the same time.
771   */
772 #if defined(PETSC_USE_SOCKET_VIEWER)
773   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
774   if (dump) {
775     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
776   }
777   dump = PETSC_FALSE;
778 #endif
779   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
780   if (dump) {
781 
782     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
783   }
784 
785   if (viewer) {
786     for (i=1; i<n; i++) {
787       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
788     }
789     for (i=0; i<n; i++) {
790       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
791       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
792     }
793   }
794   PetscFunctionReturn(0);
795 }
796 
797 /* -------------------------------------------------------------------------------------*/
798 
799 #undef __FUNCT__
800 #define __FUNCT__ "PCMGGetLevels"
801 /*@
802    PCMGGetLevels - Gets the number of levels to use with MG.
803 
804    Not Collective
805 
806    Input Parameter:
807 .  pc - the preconditioner context
808 
809    Output parameter:
810 .  levels - the number of levels
811 
812    Level: advanced
813 
814 .keywords: MG, get, levels, multigrid
815 
816 .seealso: PCMGSetLevels()
817 @*/
818 PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
819 {
820   PC_MG *mg = (PC_MG*)pc->data;
821 
822   PetscFunctionBegin;
823   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
824   PetscValidIntPointer(levels,2);
825   *levels = mg->nlevels;
826   PetscFunctionReturn(0);
827 }
828 
829 #undef __FUNCT__
830 #define __FUNCT__ "PCMGSetType"
831 /*@
832    PCMGSetType - Determines the form of multigrid to use:
833    multiplicative, additive, full, or the Kaskade algorithm.
834 
835    Logically Collective on PC
836 
837    Input Parameters:
838 +  pc - the preconditioner context
839 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
840    PC_MG_FULL, PC_MG_KASKADE
841 
842    Options Database Key:
843 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
844    additive, full, kaskade
845 
846    Level: advanced
847 
848 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
849 
850 .seealso: PCMGSetLevels()
851 @*/
852 PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
853 {
854   PC_MG                   *mg = (PC_MG*)pc->data;
855 
856   PetscFunctionBegin;
857   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
858   PetscValidLogicalCollectiveEnum(pc,form,2);
859   mg->am = form;
860   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
861   else pc->ops->applyrichardson = 0;
862   PetscFunctionReturn(0);
863 }
864 
865 #undef __FUNCT__
866 #define __FUNCT__ "PCMGSetCycleType"
867 /*@
868    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
869    complicated cycling.
870 
871    Logically Collective on PC
872 
873    Input Parameters:
874 +  pc - the multigrid context
875 -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
876 
877    Options Database Key:
878 $  -pc_mg_cycle_type v or w
879 
880    Level: advanced
881 
882 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
883 
884 .seealso: PCMGSetCycleTypeOnLevel()
885 @*/
886 PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
887 {
888   PC_MG        *mg = (PC_MG*)pc->data;
889   PC_MG_Levels **mglevels = mg->levels;
890   PetscInt     i,levels;
891 
892   PetscFunctionBegin;
893   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
894   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
895   PetscValidLogicalCollectiveInt(pc,n,2);
896   levels = mglevels[0]->levels;
897 
898   for (i=0; i<levels; i++) {
899     mglevels[i]->cycles  = n;
900   }
901   PetscFunctionReturn(0);
902 }
903 
904 #undef __FUNCT__
905 #define __FUNCT__ "PCMGMultiplicativeSetCycles"
906 /*@
907    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
908          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
909 
910    Logically Collective on PC
911 
912    Input Parameters:
913 +  pc - the multigrid context
914 -  n - number of cycles (default is 1)
915 
916    Options Database Key:
917 $  -pc_mg_multiplicative_cycles n
918 
919    Level: advanced
920 
921    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
922 
923 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
924 
925 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
926 @*/
927 PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
928 {
929   PC_MG        *mg = (PC_MG*)pc->data;
930   PC_MG_Levels **mglevels = mg->levels;
931   PetscInt     i,levels;
932 
933   PetscFunctionBegin;
934   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
935   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
936   PetscValidLogicalCollectiveInt(pc,n,2);
937   levels = mglevels[0]->levels;
938 
939   for (i=0; i<levels; i++) {
940     mg->cyclesperpcapply  = n;
941   }
942   PetscFunctionReturn(0);
943 }
944 
945 #undef __FUNCT__
946 #define __FUNCT__ "PCMGSetGalerkin"
947 /*@
948    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
949       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
950 
951    Logically Collective on PC
952 
953    Input Parameters:
954 +  pc - the multigrid context
955 -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
956 
957    Options Database Key:
958 $  -pc_mg_galerkin
959 
960    Level: intermediate
961 
962 .keywords: MG, set, Galerkin
963 
964 .seealso: PCMGGetGalerkin()
965 
966 @*/
967 PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
968 {
969   PC_MG        *mg = (PC_MG*)pc->data;
970 
971   PetscFunctionBegin;
972   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
973   mg->galerkin = use ? 1 : 0;
974   PetscFunctionReturn(0);
975 }
976 
977 #undef __FUNCT__
978 #define __FUNCT__ "PCMGGetGalerkin"
979 /*@
980    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
981       A_i-1 = r_i * A_i * r_i^t
982 
983    Not Collective
984 
985    Input Parameter:
986 .  pc - the multigrid context
987 
988    Output Parameter:
989 .  gelerkin - PETSC_TRUE or PETSC_FALSE
990 
991    Options Database Key:
992 $  -pc_mg_galerkin
993 
994    Level: intermediate
995 
996 .keywords: MG, set, Galerkin
997 
998 .seealso: PCMGSetGalerkin()
999 
1000 @*/
1001 PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
1002 {
1003   PC_MG        *mg = (PC_MG*)pc->data;
1004 
1005   PetscFunctionBegin;
1006   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1007   *galerkin = (PetscBool)mg->galerkin;
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 #undef __FUNCT__
1012 #define __FUNCT__ "PCMGSetNumberSmoothDown"
1013 /*@
1014    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
1015    use on all levels. Use PCMGGetSmootherDown() to set different
1016    pre-smoothing steps on different levels.
1017 
1018    Logically Collective on PC
1019 
1020    Input Parameters:
1021 +  mg - the multigrid context
1022 -  n - the number of smoothing steps
1023 
1024    Options Database Key:
1025 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
1026 
1027    Level: advanced
1028 
1029 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
1030 
1031 .seealso: PCMGSetNumberSmoothUp()
1032 @*/
1033 PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
1034 {
1035   PC_MG          *mg = (PC_MG*)pc->data;
1036   PC_MG_Levels   **mglevels = mg->levels;
1037   PetscErrorCode ierr;
1038   PetscInt       i,levels;
1039 
1040   PetscFunctionBegin;
1041   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1042   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1043   PetscValidLogicalCollectiveInt(pc,n,2);
1044   levels = mglevels[0]->levels;
1045 
1046   for (i=1; i<levels; i++) {
1047     /* make sure smoother up and down are different */
1048     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1049     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1050     mg->default_smoothd = n;
1051   }
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 #undef __FUNCT__
1056 #define __FUNCT__ "PCMGSetNumberSmoothUp"
1057 /*@
1058    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
1059    on all levels. Use PCMGGetSmootherUp() to set different numbers of
1060    post-smoothing steps on different levels.
1061 
1062    Logically Collective on PC
1063 
1064    Input Parameters:
1065 +  mg - the multigrid context
1066 -  n - the number of smoothing steps
1067 
1068    Options Database Key:
1069 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
1070 
1071    Level: advanced
1072 
1073    Note: this does not set a value on the coarsest grid, since we assume that
1074     there is no separate smooth up on the coarsest grid.
1075 
1076 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1077 
1078 .seealso: PCMGSetNumberSmoothDown()
1079 @*/
1080 PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
1081 {
1082   PC_MG          *mg = (PC_MG*)pc->data;
1083   PC_MG_Levels   **mglevels = mg->levels;
1084   PetscErrorCode ierr;
1085   PetscInt       i,levels;
1086 
1087   PetscFunctionBegin;
1088   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1089   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1090   PetscValidLogicalCollectiveInt(pc,n,2);
1091   levels = mglevels[0]->levels;
1092 
1093   for (i=1; i<levels; i++) {
1094     /* make sure smoother up and down are different */
1095     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1096     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1097     mg->default_smoothu = n;
1098   }
1099   PetscFunctionReturn(0);
1100 }
1101 
1102 /* ----------------------------------------------------------------------------------------*/
1103 
1104 /*MC
1105    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1106     information about the coarser grid matrices and restriction/interpolation operators.
1107 
1108    Options Database Keys:
1109 +  -pc_mg_levels <nlevels> - number of levels including finest
1110 .  -pc_mg_cycles v or w
1111 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
1112 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
1113 .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1114 .  -pc_mg_log - log information about time spent on each level of the solver
1115 .  -pc_mg_monitor - print information on the multigrid convergence
1116 .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1117 .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1118 .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1119                         to the Socket viewer for reading from MATLAB.
1120 -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1121                         to the binary output file called binaryoutput
1122 
1123    Notes: By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES
1124 
1125        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
1126 
1127    Level: intermediate
1128 
1129    Concepts: multigrid/multilevel
1130 
1131 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1132            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1133            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1134            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1135            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1136 M*/
1137 
1138 EXTERN_C_BEGIN
1139 #undef __FUNCT__
1140 #define __FUNCT__ "PCCreate_MG"
1141 PetscErrorCode  PCCreate_MG(PC pc)
1142 {
1143   PC_MG          *mg;
1144   PetscErrorCode ierr;
1145 
1146   PetscFunctionBegin;
1147   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1148   pc->data    = (void*)mg;
1149   mg->nlevels = -1;
1150 
1151   pc->ops->apply          = PCApply_MG;
1152   pc->ops->setup          = PCSetUp_MG;
1153   pc->ops->reset          = PCReset_MG;
1154   pc->ops->destroy        = PCDestroy_MG;
1155   pc->ops->setfromoptions = PCSetFromOptions_MG;
1156   pc->ops->view           = PCView_MG;
1157   PetscFunctionReturn(0);
1158 }
1159 EXTERN_C_END
1160