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