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