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