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