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