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