xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision fa59e80594ae4f5b96fb1a0d88f82abc4bbaec69)
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       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
488       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
489         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
490       }
491     }
492   }
493   PetscFunctionReturn(0);
494 }
495 
496 #include <petsc-private/dmimpl.h>
497 #include <petsc-private/kspimpl.h>
498 
499 /*
500     Calls setup for the KSP on each level
501 */
502 #undef __FUNCT__
503 #define __FUNCT__ "PCSetUp_MG"
504 PetscErrorCode PCSetUp_MG(PC pc)
505 {
506   PC_MG                   *mg = (PC_MG*)pc->data;
507   PC_MG_Levels            **mglevels = mg->levels;
508   PetscErrorCode          ierr;
509   PetscInt                i,n = mglevels[0]->levels;
510   PC                      cpc;
511   PetscBool               preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset;
512   Mat                     dA,dB;
513   MatStructure            uflag;
514   Vec                     tvec;
515   DM                      *dms;
516   PetscViewer             viewer = 0;
517 
518   PetscFunctionBegin;
519   /* FIX: Move this to PCSetFromOptions_MG? */
520   if (mg->usedmfornumberoflevels) {
521     PetscInt levels;
522     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
523     levels++;
524     if (levels > n) { /* the problem is now being solved on a finer grid */
525       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
526       n    = levels;
527       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
528       mglevels =  mg->levels;
529     }
530   }
531   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
532 
533 
534   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
535   /* so use those from global PC */
536   /* Is this what we always want? What if user wants to keep old one? */
537   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
538   if (opsset) {
539     Mat mmat;
540     ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr);
541     if (mmat == pc->pmat) opsset = PETSC_FALSE;
542   }
543   if (!opsset) {
544     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);
545     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
546   }
547 
548   /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */
549   if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
550     /* construct the interpolation from the DMs */
551     Mat p;
552     Vec rscale;
553     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
554     dms[n-1] = pc->dm;
555     for (i=n-2; i>-1; i--) {
556       DMKSP kdm;
557       ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);
558       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
559       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
560       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
561       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
562        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
563       kdm->computerhs = PETSC_NULL;
564       kdm->rhsctx = PETSC_NULL;
565       if (!mglevels[i+1]->interpolate) {
566 	ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
567 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
568 	if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
569         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
570         ierr = MatDestroy(&p);CHKERRQ(ierr);
571       }
572     }
573 
574     for (i=n-2; i>-1; i--) {
575       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
576     }
577     ierr = PetscFree(dms);CHKERRQ(ierr);
578   }
579 
580   if (pc->dm && !pc->setupcalled) {
581     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
582     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
583     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
584   }
585 
586   if (mg->galerkin == 1) {
587     Mat B;
588     /* currently only handle case where mat and pmat are the same on coarser levels */
589     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
590     if (!pc->setupcalled) {
591       for (i=n-2; i>-1; i--) {
592         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
593         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
594 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
595         dB   = B;
596       }
597       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
598     } else {
599       for (i=n-2; i>-1; i--) {
600         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
601         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
602         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
603         dB   = B;
604       }
605     }
606   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
607     /* need to restrict Jacobian location to coarser meshes for evaluation */
608     for (i=n-2;i>-1; i--) {
609       Mat R;
610       Vec rscale;
611       if (!mglevels[i]->smoothd->dm->x) {
612         Vec *vecs;
613         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr);
614         mglevels[i]->smoothd->dm->x = vecs[0];
615         ierr = PetscFree(vecs);CHKERRQ(ierr);
616       }
617       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
618       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
619       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
620       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
621     }
622   }
623   if (!mg->galerkin && pc->dm) {
624     for (i=n-2;i>=0; i--) {
625       DM dmfine,dmcoarse;
626       Mat Restrict,Inject;
627       Vec rscale;
628       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
629       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
630       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
631       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
632       Inject = PETSC_NULL;      /* Callback should create it if it needs Injection */
633       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
634     }
635   }
636 
637   if (!pc->setupcalled) {
638     for (i=0; i<n; i++) {
639       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
640     }
641     for (i=1; i<n; i++) {
642       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
643         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
644       }
645     }
646     for (i=1; i<n; i++) {
647       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
648       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
649     }
650     for (i=0; i<n-1; i++) {
651       if (!mglevels[i]->b) {
652         Vec *vec;
653         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
654         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
655         ierr = VecDestroy(vec);CHKERRQ(ierr);
656         ierr = PetscFree(vec);CHKERRQ(ierr);
657       }
658       if (!mglevels[i]->r && i) {
659         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
660         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
661         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
662       }
663       if (!mglevels[i]->x) {
664         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
665         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
666         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
667       }
668     }
669     if (n != 1 && !mglevels[n-1]->r) {
670       /* PCMGSetR() on the finest level if user did not supply it */
671       Vec *vec;
672       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
673       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
674       ierr = VecDestroy(vec);CHKERRQ(ierr);
675       ierr = PetscFree(vec);CHKERRQ(ierr);
676     }
677   }
678 
679   if (pc->dm) {
680     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
681     for (i=0; i<n-1; i++){
682       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
683     }
684   }
685 
686   for (i=1; i<n; i++) {
687     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
688       /* if doing only down then initial guess is zero */
689       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
690     }
691     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
692     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
693     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
694     if (!mglevels[i]->residual) {
695       Mat mat;
696       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
697       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
698     }
699   }
700   for (i=1; i<n; i++) {
701     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
702       Mat          downmat,downpmat;
703       MatStructure matflag;
704 
705       /* check if operators have been set for up, if not use down operators to set them */
706       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
707       if (!opsset) {
708         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
709         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
710       }
711 
712       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
713       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
714       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
715       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
716     }
717   }
718 
719   /*
720       If coarse solver is not direct method then DO NOT USE preonly
721   */
722   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
723   if (preonly) {
724     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
725     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
726     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
727     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
728     if (!lu && !redundant && !cholesky && !svd) {
729       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
730     }
731   }
732 
733   if (!pc->setupcalled) {
734     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
735   }
736 
737   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
738   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
739   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
740 
741   /*
742      Dump the interpolation/restriction matrices plus the
743    Jacobian/stiffness on each level. This allows MATLAB users to
744    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
745 
746    Only support one or the other at the same time.
747   */
748 #if defined(PETSC_USE_SOCKET_VIEWER)
749   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
750   if (dump) {
751     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
752   }
753   dump = PETSC_FALSE;
754 #endif
755   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
756   if (dump) {
757 
758     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
759   }
760 
761   if (viewer) {
762     for (i=1; i<n; i++) {
763       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
764     }
765     for (i=0; i<n; i++) {
766       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
767       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
768     }
769   }
770   PetscFunctionReturn(0);
771 }
772 
773 /* -------------------------------------------------------------------------------------*/
774 
775 #undef __FUNCT__
776 #define __FUNCT__ "PCMGGetLevels"
777 /*@
778    PCMGGetLevels - Gets the number of levels to use with MG.
779 
780    Not Collective
781 
782    Input Parameter:
783 .  pc - the preconditioner context
784 
785    Output parameter:
786 .  levels - the number of levels
787 
788    Level: advanced
789 
790 .keywords: MG, get, levels, multigrid
791 
792 .seealso: PCMGSetLevels()
793 @*/
794 PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
795 {
796   PC_MG *mg = (PC_MG*)pc->data;
797 
798   PetscFunctionBegin;
799   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
800   PetscValidIntPointer(levels,2);
801   *levels = mg->nlevels;
802   PetscFunctionReturn(0);
803 }
804 
805 #undef __FUNCT__
806 #define __FUNCT__ "PCMGSetType"
807 /*@
808    PCMGSetType - Determines the form of multigrid to use:
809    multiplicative, additive, full, or the Kaskade algorithm.
810 
811    Logically Collective on PC
812 
813    Input Parameters:
814 +  pc - the preconditioner context
815 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
816    PC_MG_FULL, PC_MG_KASKADE
817 
818    Options Database Key:
819 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
820    additive, full, kaskade
821 
822    Level: advanced
823 
824 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
825 
826 .seealso: PCMGSetLevels()
827 @*/
828 PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
829 {
830   PC_MG                   *mg = (PC_MG*)pc->data;
831 
832   PetscFunctionBegin;
833   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
834   PetscValidLogicalCollectiveEnum(pc,form,2);
835   mg->am = form;
836   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
837   else pc->ops->applyrichardson = 0;
838   PetscFunctionReturn(0);
839 }
840 
841 #undef __FUNCT__
842 #define __FUNCT__ "PCMGSetCycleType"
843 /*@
844    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
845    complicated cycling.
846 
847    Logically Collective on PC
848 
849    Input Parameters:
850 +  pc - the multigrid context
851 -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
852 
853    Options Database Key:
854 $  -pc_mg_cycle_type v or w
855 
856    Level: advanced
857 
858 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
859 
860 .seealso: PCMGSetCycleTypeOnLevel()
861 @*/
862 PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
863 {
864   PC_MG        *mg = (PC_MG*)pc->data;
865   PC_MG_Levels **mglevels = mg->levels;
866   PetscInt     i,levels;
867 
868   PetscFunctionBegin;
869   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
870   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
871   PetscValidLogicalCollectiveInt(pc,n,2);
872   levels = mglevels[0]->levels;
873 
874   for (i=0; i<levels; i++) {
875     mglevels[i]->cycles  = n;
876   }
877   PetscFunctionReturn(0);
878 }
879 
880 #undef __FUNCT__
881 #define __FUNCT__ "PCMGMultiplicativeSetCycles"
882 /*@
883    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
884          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
885 
886    Logically Collective on PC
887 
888    Input Parameters:
889 +  pc - the multigrid context
890 -  n - number of cycles (default is 1)
891 
892    Options Database Key:
893 $  -pc_mg_multiplicative_cycles n
894 
895    Level: advanced
896 
897    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
898 
899 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
900 
901 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
902 @*/
903 PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
904 {
905   PC_MG        *mg = (PC_MG*)pc->data;
906   PC_MG_Levels **mglevels = mg->levels;
907   PetscInt     i,levels;
908 
909   PetscFunctionBegin;
910   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
911   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
912   PetscValidLogicalCollectiveInt(pc,n,2);
913   levels = mglevels[0]->levels;
914 
915   for (i=0; i<levels; i++) {
916     mg->cyclesperpcapply  = n;
917   }
918   PetscFunctionReturn(0);
919 }
920 
921 #undef __FUNCT__
922 #define __FUNCT__ "PCMGSetGalerkin"
923 /*@
924    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
925       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
926 
927    Logically Collective on PC
928 
929    Input Parameters:
930 +  pc - the multigrid context
931 -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
932 
933    Options Database Key:
934 $  -pc_mg_galerkin
935 
936    Level: intermediate
937 
938 .keywords: MG, set, Galerkin
939 
940 .seealso: PCMGGetGalerkin()
941 
942 @*/
943 PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
944 {
945   PC_MG        *mg = (PC_MG*)pc->data;
946 
947   PetscFunctionBegin;
948   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
949   mg->galerkin = use ? 1 : 0;
950   PetscFunctionReturn(0);
951 }
952 
953 #undef __FUNCT__
954 #define __FUNCT__ "PCMGGetGalerkin"
955 /*@
956    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
957       A_i-1 = r_i * A_i * r_i^t
958 
959    Not Collective
960 
961    Input Parameter:
962 .  pc - the multigrid context
963 
964    Output Parameter:
965 .  gelerkin - PETSC_TRUE or PETSC_FALSE
966 
967    Options Database Key:
968 $  -pc_mg_galerkin
969 
970    Level: intermediate
971 
972 .keywords: MG, set, Galerkin
973 
974 .seealso: PCMGSetGalerkin()
975 
976 @*/
977 PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
978 {
979   PC_MG        *mg = (PC_MG*)pc->data;
980 
981   PetscFunctionBegin;
982   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
983   *galerkin = (PetscBool)mg->galerkin;
984   PetscFunctionReturn(0);
985 }
986 
987 #undef __FUNCT__
988 #define __FUNCT__ "PCMGSetNumberSmoothDown"
989 /*@
990    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
991    use on all levels. Use PCMGGetSmootherDown() to set different
992    pre-smoothing steps on different levels.
993 
994    Logically Collective on PC
995 
996    Input Parameters:
997 +  mg - the multigrid context
998 -  n - the number of smoothing steps
999 
1000    Options Database Key:
1001 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
1002 
1003    Level: advanced
1004 
1005 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
1006 
1007 .seealso: PCMGSetNumberSmoothUp()
1008 @*/
1009 PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
1010 {
1011   PC_MG          *mg = (PC_MG*)pc->data;
1012   PC_MG_Levels   **mglevels = mg->levels;
1013   PetscErrorCode ierr;
1014   PetscInt       i,levels;
1015 
1016   PetscFunctionBegin;
1017   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1018   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1019   PetscValidLogicalCollectiveInt(pc,n,2);
1020   levels = mglevels[0]->levels;
1021 
1022   for (i=1; i<levels; i++) {
1023     /* make sure smoother up and down are different */
1024     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1025     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1026     mg->default_smoothd = n;
1027   }
1028   PetscFunctionReturn(0);
1029 }
1030 
1031 #undef __FUNCT__
1032 #define __FUNCT__ "PCMGSetNumberSmoothUp"
1033 /*@
1034    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
1035    on all levels. Use PCMGGetSmootherUp() to set different numbers of
1036    post-smoothing steps on different levels.
1037 
1038    Logically Collective on PC
1039 
1040    Input Parameters:
1041 +  mg - the multigrid context
1042 -  n - the number of smoothing steps
1043 
1044    Options Database Key:
1045 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
1046 
1047    Level: advanced
1048 
1049    Note: this does not set a value on the coarsest grid, since we assume that
1050     there is no separate smooth up on the coarsest grid.
1051 
1052 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1053 
1054 .seealso: PCMGSetNumberSmoothDown()
1055 @*/
1056 PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
1057 {
1058   PC_MG          *mg = (PC_MG*)pc->data;
1059   PC_MG_Levels   **mglevels = mg->levels;
1060   PetscErrorCode ierr;
1061   PetscInt       i,levels;
1062 
1063   PetscFunctionBegin;
1064   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1065   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1066   PetscValidLogicalCollectiveInt(pc,n,2);
1067   levels = mglevels[0]->levels;
1068 
1069   for (i=1; i<levels; i++) {
1070     /* make sure smoother up and down are different */
1071     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1072     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1073     mg->default_smoothu = n;
1074   }
1075   PetscFunctionReturn(0);
1076 }
1077 
1078 /* ----------------------------------------------------------------------------------------*/
1079 
1080 /*MC
1081    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1082     information about the coarser grid matrices and restriction/interpolation operators.
1083 
1084    Options Database Keys:
1085 +  -pc_mg_levels <nlevels> - number of levels including finest
1086 .  -pc_mg_cycles v or w
1087 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
1088 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
1089 .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1090 .  -pc_mg_log - log information about time spent on each level of the solver
1091 .  -pc_mg_monitor - print information on the multigrid convergence
1092 .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1093 .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1094 .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1095                         to the Socket viewer for reading from MATLAB.
1096 -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1097                         to the binary output file called binaryoutput
1098 
1099    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
1100 
1101        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
1102 
1103    Level: intermediate
1104 
1105    Concepts: multigrid/multilevel
1106 
1107 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1108            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1109            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1110            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1111            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1112 M*/
1113 
1114 EXTERN_C_BEGIN
1115 #undef __FUNCT__
1116 #define __FUNCT__ "PCCreate_MG"
1117 PetscErrorCode  PCCreate_MG(PC pc)
1118 {
1119   PC_MG          *mg;
1120   PetscErrorCode ierr;
1121 
1122   PetscFunctionBegin;
1123   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1124   pc->data    = (void*)mg;
1125   mg->nlevels = -1;
1126 
1127   pc->ops->apply          = PCApply_MG;
1128   pc->ops->setup          = PCSetUp_MG;
1129   pc->ops->reset          = PCReset_MG;
1130   pc->ops->destroy        = PCDestroy_MG;
1131   pc->ops->setfromoptions = PCSetFromOptions_MG;
1132   pc->ops->view           = PCView_MG;
1133   PetscFunctionReturn(0);
1134 }
1135 EXTERN_C_END
1136