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