xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 7b6bb2c608b6fc6714ef38fda02c2dbb91c82665)
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,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     if (!lu && !redundant && !cholesky) {
664       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
665     }
666   }
667 
668   if (!pc->setupcalled) {
669     if (monitor) {
670       ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr);
671       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
672       ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
673     }
674     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
675   }
676 
677   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
678   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
679   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
680 
681   /*
682      Dump the interpolation/restriction matrices plus the
683    Jacobian/stiffness on each level. This allows MATLAB users to
684    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
685 
686    Only support one or the other at the same time.
687   */
688 #if defined(PETSC_USE_SOCKET_VIEWER)
689   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
690   if (dump) {
691     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
692   }
693   dump = PETSC_FALSE;
694 #endif
695   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
696   if (dump) {
697     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
698   }
699 
700   if (viewer) {
701     for (i=1; i<n; i++) {
702       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
703     }
704     for (i=0; i<n; i++) {
705       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
706       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
707     }
708   }
709   PetscFunctionReturn(0);
710 }
711 
712 /* -------------------------------------------------------------------------------------*/
713 
714 #undef __FUNCT__
715 #define __FUNCT__ "PCMGGetLevels"
716 /*@
717    PCMGGetLevels - Gets the number of levels to use with MG.
718 
719    Not Collective
720 
721    Input Parameter:
722 .  pc - the preconditioner context
723 
724    Output parameter:
725 .  levels - the number of levels
726 
727    Level: advanced
728 
729 .keywords: MG, get, levels, multigrid
730 
731 .seealso: PCMGSetLevels()
732 @*/
733 PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
734 {
735   PC_MG *mg = (PC_MG*)pc->data;
736 
737   PetscFunctionBegin;
738   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
739   PetscValidIntPointer(levels,2);
740   *levels = mg->nlevels;
741   PetscFunctionReturn(0);
742 }
743 
744 #undef __FUNCT__
745 #define __FUNCT__ "PCMGSetType"
746 /*@
747    PCMGSetType - Determines the form of multigrid to use:
748    multiplicative, additive, full, or the Kaskade algorithm.
749 
750    Logically Collective on PC
751 
752    Input Parameters:
753 +  pc - the preconditioner context
754 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
755    PC_MG_FULL, PC_MG_KASKADE
756 
757    Options Database Key:
758 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
759    additive, full, kaskade
760 
761    Level: advanced
762 
763 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
764 
765 .seealso: PCMGSetLevels()
766 @*/
767 PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
768 {
769   PC_MG                   *mg = (PC_MG*)pc->data;
770 
771   PetscFunctionBegin;
772   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
773   PetscValidLogicalCollectiveEnum(pc,form,2);
774   mg->am = form;
775   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
776   else pc->ops->applyrichardson = 0;
777   PetscFunctionReturn(0);
778 }
779 
780 #undef __FUNCT__
781 #define __FUNCT__ "PCMGSetCycleType"
782 /*@
783    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
784    complicated cycling.
785 
786    Logically Collective on PC
787 
788    Input Parameters:
789 +  pc - the multigrid context
790 -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
791 
792    Options Database Key:
793 $  -pc_mg_cycle_type v or w
794 
795    Level: advanced
796 
797 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
798 
799 .seealso: PCMGSetCycleTypeOnLevel()
800 @*/
801 PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
802 {
803   PC_MG        *mg = (PC_MG*)pc->data;
804   PC_MG_Levels **mglevels = mg->levels;
805   PetscInt     i,levels;
806 
807   PetscFunctionBegin;
808   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
809   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
810   PetscValidLogicalCollectiveInt(pc,n,2);
811   levels = mglevels[0]->levels;
812 
813   for (i=0; i<levels; i++) {
814     mglevels[i]->cycles  = n;
815   }
816   PetscFunctionReturn(0);
817 }
818 
819 #undef __FUNCT__
820 #define __FUNCT__ "PCMGMultiplicativeSetCycles"
821 /*@
822    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
823          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
824 
825    Logically Collective on PC
826 
827    Input Parameters:
828 +  pc - the multigrid context
829 -  n - number of cycles (default is 1)
830 
831    Options Database Key:
832 $  -pc_mg_multiplicative_cycles n
833 
834    Level: advanced
835 
836    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
837 
838 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
839 
840 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
841 @*/
842 PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
843 {
844   PC_MG        *mg = (PC_MG*)pc->data;
845   PC_MG_Levels **mglevels = mg->levels;
846   PetscInt     i,levels;
847 
848   PetscFunctionBegin;
849   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
850   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
851   PetscValidLogicalCollectiveInt(pc,n,2);
852   levels = mglevels[0]->levels;
853 
854   for (i=0; i<levels; i++) {
855     mg->cyclesperpcapply  = n;
856   }
857   PetscFunctionReturn(0);
858 }
859 
860 #undef __FUNCT__
861 #define __FUNCT__ "PCMGSetGalerkin"
862 /*@
863    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
864       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
865 
866    Logically Collective on PC
867 
868    Input Parameters:
869 .  pc - the multigrid context
870 
871    Options Database Key:
872 $  -pc_mg_galerkin
873 
874    Level: intermediate
875 
876 .keywords: MG, set, Galerkin
877 
878 .seealso: PCMGGetGalerkin()
879 
880 @*/
881 PetscErrorCode  PCMGSetGalerkin(PC pc)
882 {
883   PC_MG        *mg = (PC_MG*)pc->data;
884 
885   PetscFunctionBegin;
886   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
887   mg->galerkin = PETSC_TRUE;
888   PetscFunctionReturn(0);
889 }
890 
891 #undef __FUNCT__
892 #define __FUNCT__ "PCMGGetGalerkin"
893 /*@
894    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
895       A_i-1 = r_i * A_i * r_i^t
896 
897    Not Collective
898 
899    Input Parameter:
900 .  pc - the multigrid context
901 
902    Output Parameter:
903 .  gelerkin - PETSC_TRUE or PETSC_FALSE
904 
905    Options Database Key:
906 $  -pc_mg_galerkin
907 
908    Level: intermediate
909 
910 .keywords: MG, set, Galerkin
911 
912 .seealso: PCMGSetGalerkin()
913 
914 @*/
915 PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
916 {
917   PC_MG        *mg = (PC_MG*)pc->data;
918 
919   PetscFunctionBegin;
920   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
921   *galerkin = mg->galerkin;
922   PetscFunctionReturn(0);
923 }
924 
925 #undef __FUNCT__
926 #define __FUNCT__ "PCMGSetNumberSmoothDown"
927 /*@
928    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
929    use on all levels. Use PCMGGetSmootherDown() to set different
930    pre-smoothing steps on different levels.
931 
932    Logically Collective on PC
933 
934    Input Parameters:
935 +  mg - the multigrid context
936 -  n - the number of smoothing steps
937 
938    Options Database Key:
939 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
940 
941    Level: advanced
942 
943 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
944 
945 .seealso: PCMGSetNumberSmoothUp()
946 @*/
947 PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
948 {
949   PC_MG          *mg = (PC_MG*)pc->data;
950   PC_MG_Levels   **mglevels = mg->levels;
951   PetscErrorCode ierr;
952   PetscInt       i,levels;
953 
954   PetscFunctionBegin;
955   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
956   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
957   PetscValidLogicalCollectiveInt(pc,n,2);
958   levels = mglevels[0]->levels;
959 
960   for (i=1; i<levels; i++) {
961     /* make sure smoother up and down are different */
962     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
963     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
964     mg->default_smoothd = n;
965   }
966   PetscFunctionReturn(0);
967 }
968 
969 #undef __FUNCT__
970 #define __FUNCT__ "PCMGSetNumberSmoothUp"
971 /*@
972    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
973    on all levels. Use PCMGGetSmootherUp() to set different numbers of
974    post-smoothing steps on different levels.
975 
976    Logically Collective on PC
977 
978    Input Parameters:
979 +  mg - the multigrid context
980 -  n - the number of smoothing steps
981 
982    Options Database Key:
983 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
984 
985    Level: advanced
986 
987    Note: this does not set a value on the coarsest grid, since we assume that
988     there is no separate smooth up on the coarsest grid.
989 
990 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
991 
992 .seealso: PCMGSetNumberSmoothDown()
993 @*/
994 PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
995 {
996   PC_MG          *mg = (PC_MG*)pc->data;
997   PC_MG_Levels   **mglevels = mg->levels;
998   PetscErrorCode ierr;
999   PetscInt       i,levels;
1000 
1001   PetscFunctionBegin;
1002   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1003   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1004   PetscValidLogicalCollectiveInt(pc,n,2);
1005   levels = mglevels[0]->levels;
1006 
1007   for (i=1; i<levels; i++) {
1008     /* make sure smoother up and down are different */
1009     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1010     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1011     mg->default_smoothu = n;
1012   }
1013   PetscFunctionReturn(0);
1014 }
1015 
1016 /* ----------------------------------------------------------------------------------------*/
1017 
1018 /*MC
1019    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1020     information about the coarser grid matrices and restriction/interpolation operators.
1021 
1022    Options Database Keys:
1023 +  -pc_mg_levels <nlevels> - number of levels including finest
1024 .  -pc_mg_cycles v or w
1025 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
1026 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
1027 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
1028 .  -pc_mg_log - log information about time spent on each level of the solver
1029 .  -pc_mg_monitor - print information on the multigrid convergence
1030 .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
1031 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1032                         to the Socket viewer for reading from MATLAB.
1033 
1034    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
1035 
1036    Level: intermediate
1037 
1038    Concepts: multigrid/multilevel
1039 
1040 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC
1041            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1042            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1043            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1044            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1045 M*/
1046 
1047 EXTERN_C_BEGIN
1048 #undef __FUNCT__
1049 #define __FUNCT__ "PCCreate_MG"
1050 PetscErrorCode  PCCreate_MG(PC pc)
1051 {
1052   PC_MG          *mg;
1053   PetscErrorCode ierr;
1054 
1055   PetscFunctionBegin;
1056   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1057   pc->data    = (void*)mg;
1058   mg->nlevels = -1;
1059 
1060   pc->ops->apply          = PCApply_MG;
1061   pc->ops->setup          = PCSetUp_MG;
1062   pc->ops->reset          = PCReset_MG;
1063   pc->ops->destroy        = PCDestroy_MG;
1064   pc->ops->setfromoptions = PCSetFromOptions_MG;
1065   pc->ops->view           = PCView_MG;
1066   PetscFunctionReturn(0);
1067 }
1068 EXTERN_C_END
1069