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