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