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