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