xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 87e657c68951b264cdcd22b0bb0e60ff85785b6e)
1 
2 /*
3     Defines the multigrid preconditioner interface.
4 */
5 #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscpcmg.h" I*/
6 
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "PCMGMCycle_Private"
10 PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
11 {
12   PC_MG          *mg = (PC_MG*)pc->data;
13   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
14   PetscErrorCode ierr;
15   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
16 
17   PetscFunctionBegin;
18 
19   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
20   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
21   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
22   if (mglevels->level) {  /* not the coarsest grid */
23     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
24     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
25     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
26 
27     /* if on finest level and have convergence criteria set */
28     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
29       PetscReal rnorm;
30       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
31       if (rnorm <= mg->ttol) {
32         if (rnorm < mg->abstol) {
33           *reason = PCRICHARDSON_CONVERGED_ATOL;
34           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
35         } else {
36           *reason = PCRICHARDSON_CONVERGED_RTOL;
37           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
38         }
39         PetscFunctionReturn(0);
40       }
41     }
42 
43     mgc = *(mglevelsin - 1);
44     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
45     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
46     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
47     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
48     while (cycles--) {
49       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
50     }
51     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
52     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
53     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
54     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
55     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
56     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "PCApplyRichardson_MG"
63 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)
64 {
65   PC_MG          *mg = (PC_MG*)pc->data;
66   PC_MG_Levels   **mglevels = mg->levels;
67   PetscErrorCode ierr;
68   PetscInt       levels = mglevels[0]->levels,i;
69 
70   PetscFunctionBegin;
71   mglevels[levels-1]->b    = b;
72   mglevels[levels-1]->x    = x;
73 
74   mg->rtol = rtol;
75   mg->abstol = abstol;
76   mg->dtol = dtol;
77   if (rtol) {
78     /* compute initial residual norm for relative convergence test */
79     PetscReal rnorm;
80     if (zeroguess) {
81       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
82     } else {
83       ierr               = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
84       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
85     }
86     mg->ttol = PetscMax(rtol*rnorm,abstol);
87   } else if (abstol) {
88     mg->ttol = abstol;
89   } else {
90     mg->ttol = 0.0;
91   }
92 
93   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
94      stop prematurely do to small residual */
95   for (i=1; i<levels; i++) {
96     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
97     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
98       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
99     }
100   }
101 
102   *reason = (PCRichardsonConvergedReason)0;
103   for (i=0; i<its; i++) {
104     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
105     if (*reason) break;
106   }
107   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
108   *outits = i;
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "PCReset_MG"
114 PetscErrorCode PCReset_MG(PC pc)
115 {
116   PC_MG          *mg = (PC_MG*)pc->data;
117   PC_MG_Levels   **mglevels = mg->levels;
118   PetscErrorCode ierr;
119   PetscInt       i,n;
120 
121   PetscFunctionBegin;
122   if (mglevels) {
123     n = mglevels[0]->levels;
124     for (i=0; i<n-1; i++) {
125       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
126       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
127       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
128       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
129       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
130       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
131     }
132 
133     for (i=0; i<n; i++) {
134       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
135       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
136 	ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
137       }
138       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
139     }
140   }
141   PetscFunctionReturn(0);
142 }
143 
144 #undef __FUNCT__
145 #define __FUNCT__ "PCMGSetLevels"
146 /*@C
147    PCMGSetLevels - Sets the number of levels to use with MG.
148    Must be called before any other MG routine.
149 
150    Logically Collective on PC
151 
152    Input Parameters:
153 +  pc - the preconditioner context
154 .  levels - the number of levels
155 -  comms - optional communicators for each level; this is to allow solving the coarser problems
156            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
157 
158    Level: intermediate
159 
160    Notes:
161      If the number of levels is one then the multigrid uses the -mg_levels prefix
162   for setting the level options rather than the -mg_coarse prefix.
163 
164 .keywords: MG, set, levels, multigrid
165 
166 .seealso: PCMGSetType(), PCMGGetLevels()
167 @*/
168 PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
169 {
170   PetscErrorCode ierr;
171   PC_MG          *mg = (PC_MG*)pc->data;
172   MPI_Comm       comm = ((PetscObject)pc)->comm;
173   PC_MG_Levels   **mglevels = mg->levels;
174   PetscInt       i;
175   PetscMPIInt    size;
176   const char     *prefix;
177   PC             ipc;
178   PetscInt       n;
179 
180   PetscFunctionBegin;
181   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
182   PetscValidLogicalCollectiveInt(pc,levels,2);
183   if (mg->nlevels == levels) PetscFunctionReturn(0);
184   if (mglevels) {
185     /* changing the number of levels so free up the previous stuff */
186     ierr = PCReset_MG(pc);CHKERRQ(ierr);
187     n = mglevels[0]->levels;
188     for (i=0; i<n; i++) {
189       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
190 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
191       }
192       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
193       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
194     }
195     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
196   }
197 
198   mg->nlevels      = levels;
199   mg->galerkin     = PETSC_FALSE;
200   mg->galerkinused = PETSC_FALSE;
201 
202   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr);
203   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
204 
205   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
206 
207   for (i=0; i<levels; i++) {
208     ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr);
209     mglevels[i]->level           = i;
210     mglevels[i]->levels          = levels;
211     mglevels[i]->cycles          = PC_MG_CYCLE_V;
212     mg->default_smoothu = 1;
213     mg->default_smoothd = 1;
214     mglevels[i]->eventsmoothsetup    = 0;
215     mglevels[i]->eventsmoothsolve    = 0;
216     mglevels[i]->eventresidual       = 0;
217     mglevels[i]->eventinterprestrict = 0;
218 
219     if (comms) comm = comms[i];
220     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
221     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
222     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
223     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
224 
225     /* do special stuff for coarse grid */
226     if (!i && levels > 1) {
227       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
228 
229       /* coarse solve is (redundant) LU by default */
230       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
231       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
232       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
233       if (size > 1) {
234         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
235       } else {
236         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
237       }
238 
239     } else {
240       char tprefix[128];
241       sprintf(tprefix,"mg_levels_%d_",(int)i);
242       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
243     }
244     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
245     mglevels[i]->smoothu    = mglevels[i]->smoothd;
246     mg->rtol                = 0.0;
247     mg->abstol              = 0.0;
248     mg->dtol                = 0.0;
249     mg->ttol                = 0.0;
250     mg->cyclesperpcapply    = 1;
251   }
252   mg->am          = PC_MG_MULTIPLICATIVE;
253   mg->levels      = mglevels;
254   pc->ops->applyrichardson = PCApplyRichardson_MG;
255   PetscFunctionReturn(0);
256 }
257 
258 
259 #undef __FUNCT__
260 #define __FUNCT__ "PCDestroy_MG"
261 PetscErrorCode PCDestroy_MG(PC pc)
262 {
263   PetscErrorCode ierr;
264   PC_MG          *mg = (PC_MG*)pc->data;
265   PC_MG_Levels   **mglevels = mg->levels;
266   PetscInt       i,n;
267 
268   PetscFunctionBegin;
269   ierr = PCReset_MG(pc);CHKERRQ(ierr);
270   if (mglevels) {
271     n = mglevels[0]->levels;
272     for (i=0; i<n; i++) {
273       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
274 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
275       }
276       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
277       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
278     }
279     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
280   }
281   ierr = PetscFree(pc->data);CHKERRQ(ierr);
282   PetscFunctionReturn(0);
283 }
284 
285 
286 
287 extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
288 extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
289 extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
290 
291 /*
292    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
293              or full cycle of multigrid.
294 
295   Note:
296   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
297 */
298 #undef __FUNCT__
299 #define __FUNCT__ "PCApply_MG"
300 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
301 {
302   PC_MG          *mg = (PC_MG*)pc->data;
303   PC_MG_Levels   **mglevels = mg->levels;
304   PetscErrorCode ierr;
305   PetscInt       levels = mglevels[0]->levels,i;
306 
307   PetscFunctionBegin;
308 
309   /* When the DM is supplying the matrix then it will not exist until here */
310   for (i=0; i<levels; i++) {
311     if (!mglevels[i]->A) {
312       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
313       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
314     }
315   }
316 
317   mglevels[levels-1]->b = b;
318   mglevels[levels-1]->x = x;
319   if (mg->am == PC_MG_MULTIPLICATIVE) {
320     ierr = VecSet(x,0.0);CHKERRQ(ierr);
321     for (i=0; i<mg->cyclesperpcapply; i++) {
322       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
323     }
324   }
325   else if (mg->am == PC_MG_ADDITIVE) {
326     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
327   }
328   else if (mg->am == PC_MG_KASKADE) {
329     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
330   }
331   else {
332     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
333   }
334   PetscFunctionReturn(0);
335 }
336 
337 
338 #undef __FUNCT__
339 #define __FUNCT__ "PCSetFromOptions_MG"
340 PetscErrorCode PCSetFromOptions_MG(PC pc)
341 {
342   PetscErrorCode ierr;
343   PetscInt       m,levels = 1,cycles;
344   PetscBool      flg;
345   PC_MG          *mg = (PC_MG*)pc->data;
346   PC_MG_Levels   **mglevels = mg->levels;
347   PCMGType       mgtype;
348   PCMGCycleType  mgctype;
349 
350   PetscFunctionBegin;
351   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
352     if (!mg->levels) {
353       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
354       if (!flg && pc->dm) {
355         ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
356         levels++;
357         mg->usedmfornumberoflevels = PETSC_TRUE;
358       }
359       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
360     }
361     mglevels = mg->levels;
362 
363     mgctype = (PCMGCycleType) mglevels[0]->cycles;
364     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
365     if (flg) {
366       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
367     };
368     flg  = PETSC_FALSE;
369     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
370     if (flg) {
371       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
372     }
373     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
374     if (flg) {
375       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
376     }
377     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
378     if (flg) {
379       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
380     }
381     mgtype = mg->am;
382     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
383     if (flg) {
384       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
385     }
386     if (mg->am == PC_MG_MULTIPLICATIVE) {
387       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
388       if (flg) {
389 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
390       }
391     }
392     flg  = PETSC_FALSE;
393     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
394     if (flg) {
395       PetscInt i;
396       char     eventname[128];
397       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
398       levels = mglevels[0]->levels;
399       for (i=0; i<levels; i++) {
400         sprintf(eventname,"MGSetup Level %d",(int)i);
401         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
402         sprintf(eventname,"MGSmooth Level %d",(int)i);
403         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
404         if (i) {
405           sprintf(eventname,"MGResid Level %d",(int)i);
406           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
407           sprintf(eventname,"MGInterp Level %d",(int)i);
408           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
409         }
410       }
411     }
412   ierr = PetscOptionsTail();CHKERRQ(ierr);
413   PetscFunctionReturn(0);
414 }
415 
416 const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
417 const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
418 
419 #undef __FUNCT__
420 #define __FUNCT__ "PCView_MG"
421 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
422 {
423   PC_MG          *mg = (PC_MG*)pc->data;
424   PC_MG_Levels   **mglevels = mg->levels;
425   PetscErrorCode ierr;
426   PetscInt       levels = mglevels[0]->levels,i;
427   PetscBool      iascii;
428 
429   PetscFunctionBegin;
430   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
431   if (iascii) {
432     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,(mglevels[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr);
433     if (mg->am == PC_MG_MULTIPLICATIVE) {
434       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
435     }
436     if (mg->galerkin) {
437       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
438     } else {
439       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
440     }
441     for (i=0; i<levels; i++) {
442       if (!i) {
443         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
444       } else {
445         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
446       }
447       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
448       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
449       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
450       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
451         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
452       } else if (i){
453         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr);
454         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
455         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
456         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
457       }
458     }
459   } else {
460     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
461   }
462   PetscFunctionReturn(0);
463 }
464 
465 #include <private/dmimpl.h>
466 #include <private/kspimpl.h>
467 
468 /*
469     Calls setup for the KSP on each level
470 */
471 #undef __FUNCT__
472 #define __FUNCT__ "PCSetUp_MG"
473 PetscErrorCode PCSetUp_MG(PC pc)
474 {
475   PC_MG                   *mg = (PC_MG*)pc->data;
476   PC_MG_Levels            **mglevels = mg->levels;
477   PetscErrorCode          ierr;
478   PetscInt                i,n = mglevels[0]->levels;
479   PC                      cpc,mpc;
480   PetscBool               preonly,lu,redundant,cholesky,svd,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset;
481   PetscViewerASCIIMonitor ascii;
482   PetscViewer             viewer = PETSC_NULL;
483   MPI_Comm                comm;
484   Mat                     dA,dB;
485   MatStructure            uflag;
486   Vec                     tvec;
487   DM                      *dms;
488 
489   PetscFunctionBegin;
490   if (mg->usedmfornumberoflevels) {
491     PetscInt levels;
492     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
493     levels++;
494     if (levels > n) { /* the problem is now being solved on a finer grid */
495       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
496       n    = levels;
497       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
498       mglevels =  mg->levels;
499     }
500   }
501 
502 
503   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
504   /* so use those from global PC */
505   /* Is this what we always want? What if user wants to keep old one? */
506   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
507   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
508   ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr);
509   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) {
510     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);
511     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
512   }
513 
514   if (pc->dm && !pc->setupcalled) {
515     /* construct the interpolation from the DMs */
516     Mat p;
517     Vec rscale;
518     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
519     dms[n-1] = pc->dm;
520     for (i=n-2; i>-1; i--) {
521       ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr);
522       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
523       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
524       ierr = DMSetFunction(dms[i],0);
525       ierr = DMSetInitialGuess(dms[i],0);
526       if (!mglevels[i+1]->interpolate) {
527 	ierr = DMGetInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
528 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
529 	ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);
530         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
531         ierr = MatDestroy(&p);CHKERRQ(ierr);
532       }
533     }
534 
535     for (i=n-2; i>-1; i--) {
536       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
537     }
538     ierr = PetscFree(dms);CHKERRQ(ierr);
539 
540     /* finest smoother also gets DM but it is not active */
541     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
542     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
543   }
544 
545   if (mg->galerkin) {
546     Mat B;
547     mg->galerkinused = PETSC_TRUE;
548     /* currently only handle case where mat and pmat are the same on coarser levels */
549     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
550     if (!pc->setupcalled) {
551       for (i=n-2; i>-1; i--) {
552         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
553         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
554 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
555         dB   = B;
556       }
557       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
558     } else {
559       for (i=n-2; i>-1; i--) {
560         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
561         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
562         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
563         dB   = B;
564       }
565     }
566   } else if (pc->dm && pc->dm->x) {
567     /* need to restrict Jacobian location to coarser meshes for evaluation */
568     for (i=n-2;i>-1; i--) {
569       if (!mglevels[i]->smoothd->dm->x) {
570         Vec *vecs;
571         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr);
572         mglevels[i]->smoothd->dm->x = vecs[0];
573         ierr = PetscFree(vecs);CHKERRQ(ierr);
574       }
575       ierr = MatRestrict(mglevels[i+1]->interpolate,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
576       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,mglevels[i+1]->rscale);CHKERRQ(ierr);
577     }
578   }
579 
580   if (!pc->setupcalled) {
581     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
582 
583     for (i=0; i<n; i++) {
584       if (monitor) {
585         ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr);
586         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
587         ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
588       }
589       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
590     }
591     for (i=1; i<n; i++) {
592       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
593         if (monitor) {
594           ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr);
595           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
596           ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
597         }
598         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
599       }
600     }
601     for (i=1; i<n; i++) {
602       if (mglevels[i]->restrct && !mglevels[i]->interpolate) {
603         ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr);
604       }
605       if (!mglevels[i]->restrct && mglevels[i]->interpolate) {
606         ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr);
607       }
608 #if defined(PETSC_USE_DEBUG)
609       if (!mglevels[i]->restrct || !mglevels[i]->interpolate) {
610         SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
611       }
612 #endif
613     }
614     for (i=0; i<n-1; i++) {
615       if (!mglevels[i]->b) {
616         Vec *vec;
617         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
618         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
619         ierr = VecDestroy(vec);CHKERRQ(ierr);
620         ierr = PetscFree(vec);CHKERRQ(ierr);
621       }
622       if (!mglevels[i]->r && i) {
623         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
624         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
625         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
626       }
627       if (!mglevels[i]->x) {
628         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
629         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
630         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
631       }
632     }
633     if (n != 1 && !mglevels[n-1]->r) {
634       /* PCMGSetR() on the finest level if user did not supply it */
635       Vec *vec;
636       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
637       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
638       ierr = VecDestroy(vec);CHKERRQ(ierr);
639       ierr = PetscFree(vec);CHKERRQ(ierr);
640     }
641   }
642 
643 
644   for (i=1; i<n; i++) {
645     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
646       /* if doing only down then initial guess is zero */
647       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
648     }
649     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
650     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
651     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
652     if (!mglevels[i]->residual) {
653       Mat mat;
654       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
655       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
656     }
657   }
658   for (i=1; i<n; i++) {
659     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
660       Mat          downmat,downpmat;
661       MatStructure matflag;
662       PetscBool    opsset;
663 
664       /* check if operators have been set for up, if not use down operators to set them */
665       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
666       if (!opsset) {
667         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
668         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
669       }
670 
671       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
672       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
673       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
674       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
675     }
676   }
677 
678   /*
679       If coarse solver is not direct method then DO NOT USE preonly
680   */
681   ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
682   if (preonly) {
683     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
684     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
685     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
686     ierr = PetscTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
687     if (!lu && !redundant && !cholesky && !svd) {
688       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
689     }
690   }
691 
692   if (!pc->setupcalled) {
693     if (monitor) {
694       ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr);
695       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
696       ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
697     }
698     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
699   }
700 
701   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
702   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
703   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
704 
705   /*
706      Dump the interpolation/restriction matrices plus the
707    Jacobian/stiffness on each level. This allows MATLAB users to
708    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
709 
710    Only support one or the other at the same time.
711   */
712 #if defined(PETSC_USE_SOCKET_VIEWER)
713   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
714   if (dump) {
715     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
716   }
717   dump = PETSC_FALSE;
718 #endif
719   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
720   if (dump) {
721     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
722   }
723 
724   if (viewer) {
725     for (i=1; i<n; i++) {
726       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
727     }
728     for (i=0; i<n; i++) {
729       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
730       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
731     }
732   }
733   PetscFunctionReturn(0);
734 }
735 
736 /* -------------------------------------------------------------------------------------*/
737 
738 #undef __FUNCT__
739 #define __FUNCT__ "PCMGGetLevels"
740 /*@
741    PCMGGetLevels - Gets the number of levels to use with MG.
742 
743    Not Collective
744 
745    Input Parameter:
746 .  pc - the preconditioner context
747 
748    Output parameter:
749 .  levels - the number of levels
750 
751    Level: advanced
752 
753 .keywords: MG, get, levels, multigrid
754 
755 .seealso: PCMGSetLevels()
756 @*/
757 PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
758 {
759   PC_MG *mg = (PC_MG*)pc->data;
760 
761   PetscFunctionBegin;
762   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
763   PetscValidIntPointer(levels,2);
764   *levels = mg->nlevels;
765   PetscFunctionReturn(0);
766 }
767 
768 #undef __FUNCT__
769 #define __FUNCT__ "PCMGSetType"
770 /*@
771    PCMGSetType - Determines the form of multigrid to use:
772    multiplicative, additive, full, or the Kaskade algorithm.
773 
774    Logically Collective on PC
775 
776    Input Parameters:
777 +  pc - the preconditioner context
778 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
779    PC_MG_FULL, PC_MG_KASKADE
780 
781    Options Database Key:
782 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
783    additive, full, kaskade
784 
785    Level: advanced
786 
787 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
788 
789 .seealso: PCMGSetLevels()
790 @*/
791 PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
792 {
793   PC_MG                   *mg = (PC_MG*)pc->data;
794 
795   PetscFunctionBegin;
796   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
797   PetscValidLogicalCollectiveEnum(pc,form,2);
798   mg->am = form;
799   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
800   else pc->ops->applyrichardson = 0;
801   PetscFunctionReturn(0);
802 }
803 
804 #undef __FUNCT__
805 #define __FUNCT__ "PCMGSetCycleType"
806 /*@
807    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
808    complicated cycling.
809 
810    Logically Collective on PC
811 
812    Input Parameters:
813 +  pc - the multigrid context
814 -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
815 
816    Options Database Key:
817 $  -pc_mg_cycle_type v or w
818 
819    Level: advanced
820 
821 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
822 
823 .seealso: PCMGSetCycleTypeOnLevel()
824 @*/
825 PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
826 {
827   PC_MG        *mg = (PC_MG*)pc->data;
828   PC_MG_Levels **mglevels = mg->levels;
829   PetscInt     i,levels;
830 
831   PetscFunctionBegin;
832   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
833   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
834   PetscValidLogicalCollectiveInt(pc,n,2);
835   levels = mglevels[0]->levels;
836 
837   for (i=0; i<levels; i++) {
838     mglevels[i]->cycles  = n;
839   }
840   PetscFunctionReturn(0);
841 }
842 
843 #undef __FUNCT__
844 #define __FUNCT__ "PCMGMultiplicativeSetCycles"
845 /*@
846    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
847          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
848 
849    Logically Collective on PC
850 
851    Input Parameters:
852 +  pc - the multigrid context
853 -  n - number of cycles (default is 1)
854 
855    Options Database Key:
856 $  -pc_mg_multiplicative_cycles n
857 
858    Level: advanced
859 
860    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
861 
862 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
863 
864 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
865 @*/
866 PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
867 {
868   PC_MG        *mg = (PC_MG*)pc->data;
869   PC_MG_Levels **mglevels = mg->levels;
870   PetscInt     i,levels;
871 
872   PetscFunctionBegin;
873   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
874   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
875   PetscValidLogicalCollectiveInt(pc,n,2);
876   levels = mglevels[0]->levels;
877 
878   for (i=0; i<levels; i++) {
879     mg->cyclesperpcapply  = n;
880   }
881   PetscFunctionReturn(0);
882 }
883 
884 #undef __FUNCT__
885 #define __FUNCT__ "PCMGSetGalerkin"
886 /*@
887    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
888       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
889 
890    Logically Collective on PC
891 
892    Input Parameters:
893 .  pc - the multigrid context
894 
895    Options Database Key:
896 $  -pc_mg_galerkin
897 
898    Level: intermediate
899 
900 .keywords: MG, set, Galerkin
901 
902 .seealso: PCMGGetGalerkin()
903 
904 @*/
905 PetscErrorCode  PCMGSetGalerkin(PC pc)
906 {
907   PC_MG        *mg = (PC_MG*)pc->data;
908 
909   PetscFunctionBegin;
910   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
911   mg->galerkin = PETSC_TRUE;
912   PetscFunctionReturn(0);
913 }
914 
915 #undef __FUNCT__
916 #define __FUNCT__ "PCMGGetGalerkin"
917 /*@
918    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
919       A_i-1 = r_i * A_i * r_i^t
920 
921    Not Collective
922 
923    Input Parameter:
924 .  pc - the multigrid context
925 
926    Output Parameter:
927 .  gelerkin - PETSC_TRUE or PETSC_FALSE
928 
929    Options Database Key:
930 $  -pc_mg_galerkin
931 
932    Level: intermediate
933 
934 .keywords: MG, set, Galerkin
935 
936 .seealso: PCMGSetGalerkin()
937 
938 @*/
939 PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
940 {
941   PC_MG        *mg = (PC_MG*)pc->data;
942 
943   PetscFunctionBegin;
944   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
945   *galerkin = mg->galerkin;
946   PetscFunctionReturn(0);
947 }
948 
949 #undef __FUNCT__
950 #define __FUNCT__ "PCMGSetNumberSmoothDown"
951 /*@
952    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
953    use on all levels. Use PCMGGetSmootherDown() to set different
954    pre-smoothing steps on different levels.
955 
956    Logically Collective on PC
957 
958    Input Parameters:
959 +  mg - the multigrid context
960 -  n - the number of smoothing steps
961 
962    Options Database Key:
963 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
964 
965    Level: advanced
966 
967 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
968 
969 .seealso: PCMGSetNumberSmoothUp()
970 @*/
971 PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
972 {
973   PC_MG          *mg = (PC_MG*)pc->data;
974   PC_MG_Levels   **mglevels = mg->levels;
975   PetscErrorCode ierr;
976   PetscInt       i,levels;
977 
978   PetscFunctionBegin;
979   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
980   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
981   PetscValidLogicalCollectiveInt(pc,n,2);
982   levels = mglevels[0]->levels;
983 
984   for (i=1; i<levels; i++) {
985     /* make sure smoother up and down are different */
986     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
987     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
988     mg->default_smoothd = n;
989   }
990   PetscFunctionReturn(0);
991 }
992 
993 #undef __FUNCT__
994 #define __FUNCT__ "PCMGSetNumberSmoothUp"
995 /*@
996    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
997    on all levels. Use PCMGGetSmootherUp() to set different numbers of
998    post-smoothing steps on different levels.
999 
1000    Logically Collective on PC
1001 
1002    Input Parameters:
1003 +  mg - the multigrid context
1004 -  n - the number of smoothing steps
1005 
1006    Options Database Key:
1007 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
1008 
1009    Level: advanced
1010 
1011    Note: this does not set a value on the coarsest grid, since we assume that
1012     there is no separate smooth up on the coarsest grid.
1013 
1014 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1015 
1016 .seealso: PCMGSetNumberSmoothDown()
1017 @*/
1018 PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
1019 {
1020   PC_MG          *mg = (PC_MG*)pc->data;
1021   PC_MG_Levels   **mglevels = mg->levels;
1022   PetscErrorCode ierr;
1023   PetscInt       i,levels;
1024 
1025   PetscFunctionBegin;
1026   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1027   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1028   PetscValidLogicalCollectiveInt(pc,n,2);
1029   levels = mglevels[0]->levels;
1030 
1031   for (i=1; i<levels; i++) {
1032     /* make sure smoother up and down are different */
1033     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1034     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1035     mg->default_smoothu = n;
1036   }
1037   PetscFunctionReturn(0);
1038 }
1039 
1040 /* ----------------------------------------------------------------------------------------*/
1041 
1042 /*MC
1043    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1044     information about the coarser grid matrices and restriction/interpolation operators.
1045 
1046    Options Database Keys:
1047 +  -pc_mg_levels <nlevels> - number of levels including finest
1048 .  -pc_mg_cycles v or w
1049 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
1050 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
1051 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
1052 .  -pc_mg_log - log information about time spent on each level of the solver
1053 .  -pc_mg_monitor - print information on the multigrid convergence
1054 .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
1055 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1056                         to the Socket viewer for reading from MATLAB.
1057 
1058    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
1059 
1060    Level: intermediate
1061 
1062    Concepts: multigrid/multilevel
1063 
1064 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC
1065            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1066            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1067            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1068            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1069 M*/
1070 
1071 EXTERN_C_BEGIN
1072 #undef __FUNCT__
1073 #define __FUNCT__ "PCCreate_MG"
1074 PetscErrorCode  PCCreate_MG(PC pc)
1075 {
1076   PC_MG          *mg;
1077   PetscErrorCode ierr;
1078 
1079   PetscFunctionBegin;
1080   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1081   pc->data    = (void*)mg;
1082   mg->nlevels = -1;
1083 
1084   pc->ops->apply          = PCApply_MG;
1085   pc->ops->setup          = PCSetUp_MG;
1086   pc->ops->reset          = PCReset_MG;
1087   pc->ops->destroy        = PCDestroy_MG;
1088   pc->ops->setfromoptions = PCSetFromOptions_MG;
1089   pc->ops->view           = PCView_MG;
1090   PetscFunctionReturn(0);
1091 }
1092 EXTERN_C_END
1093