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