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