xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision bcab024551cbaaaa7f49e357634a3584f91cb6d5)
1 
2 /*
3     Defines the multigrid preconditioner interface.
4 */
5 #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscksp.h" I*/
6 #include <petscdm.h>
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   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
19   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
20   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
21   if (mglevels->level) {  /* not the coarsest grid */
22     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
23     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
24     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
25 
26     /* if on finest level and have convergence criteria set */
27     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
28       PetscReal rnorm;
29       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
30       if (rnorm <= mg->ttol) {
31         if (rnorm < mg->abstol) {
32           *reason = PCRICHARDSON_CONVERGED_ATOL;
33           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr);
34         } else {
35           *reason = PCRICHARDSON_CONVERGED_RTOL;
36           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);CHKERRQ(ierr);
37         }
38         PetscFunctionReturn(0);
39       }
40     }
41 
42     mgc = *(mglevelsin - 1);
43     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
44     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
45     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
46     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
47     while (cycles--) {
48       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
49     }
50     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
51     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
52     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
53     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
54     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
55     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
56   }
57   PetscFunctionReturn(0);
58 }
59 
60 #undef __FUNCT__
61 #define __FUNCT__ "PCApplyRichardson_MG"
62 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)
63 {
64   PC_MG          *mg        = (PC_MG*)pc->data;
65   PC_MG_Levels   **mglevels = mg->levels;
66   PetscErrorCode ierr;
67   PetscInt       levels = mglevels[0]->levels,i;
68 
69   PetscFunctionBegin;
70   /* When the DM is supplying the matrix then it will not exist until here */
71   for (i=0; i<levels; i++) {
72     if (!mglevels[i]->A) {
73       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
74       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
75     }
76   }
77   mglevels[levels-1]->b = b;
78   mglevels[levels-1]->x = x;
79 
80   mg->rtol   = rtol;
81   mg->abstol = abstol;
82   mg->dtol   = dtol;
83   if (rtol) {
84     /* compute initial residual norm for relative convergence test */
85     PetscReal rnorm;
86     if (zeroguess) {
87       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
88     } else {
89       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
90       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
91     }
92     mg->ttol = PetscMax(rtol*rnorm,abstol);
93   } else if (abstol) mg->ttol = abstol;
94   else mg->ttol = 0.0;
95 
96   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
97      stop prematurely due to small residual */
98   for (i=1; i<levels; i++) {
99     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
100     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
101       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
102     }
103   }
104 
105   *reason = (PCRichardsonConvergedReason)0;
106   for (i=0; i<its; i++) {
107     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
108     if (*reason) break;
109   }
110   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
111   *outits = i;
112   PetscFunctionReturn(0);
113 }
114 
115 #undef __FUNCT__
116 #define __FUNCT__ "PCReset_MG"
117 PetscErrorCode PCReset_MG(PC pc)
118 {
119   PC_MG          *mg        = (PC_MG*)pc->data;
120   PC_MG_Levels   **mglevels = mg->levels;
121   PetscErrorCode ierr;
122   PetscInt       i,n;
123 
124   PetscFunctionBegin;
125   if (mglevels) {
126     n = mglevels[0]->levels;
127     for (i=0; i<n-1; i++) {
128       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
129       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
130       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
131       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
132       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
133       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
134     }
135 
136     for (i=0; i<n; i++) {
137       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
138       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
139         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
140       }
141       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
142     }
143   }
144   PetscFunctionReturn(0);
145 }
146 
147 #undef __FUNCT__
148 #define __FUNCT__ "PCMGSetLevels"
149 /*@C
150    PCMGSetLevels - Sets the number of levels to use with MG.
151    Must be called before any other MG routine.
152 
153    Logically Collective on PC
154 
155    Input Parameters:
156 +  pc - the preconditioner context
157 .  levels - the number of levels
158 -  comms - optional communicators for each level; this is to allow solving the coarser problems
159            on smaller sets of processors. Use NULL_OBJECT for default in Fortran
160 
161    Level: intermediate
162 
163    Notes:
164      If the number of levels is one then the multigrid uses the -mg_levels prefix
165   for setting the level options rather than the -mg_coarse prefix.
166 
167 .keywords: MG, set, levels, multigrid
168 
169 .seealso: PCMGSetType(), PCMGGetLevels()
170 @*/
171 PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
172 {
173   PetscErrorCode ierr;
174   PC_MG          *mg        = (PC_MG*)pc->data;
175   MPI_Comm       comm;
176   PC_MG_Levels   **mglevels = mg->levels;
177   PetscInt       i;
178   PetscMPIInt    size;
179   const char     *prefix;
180   PC             ipc;
181   PetscInt       n;
182 
183   PetscFunctionBegin;
184   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
185   PetscValidLogicalCollectiveInt(pc,levels,2);
186   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
187   if (mg->nlevels == levels) PetscFunctionReturn(0);
188   if (mglevels) {
189     /* changing the number of levels so free up the previous stuff */
190     ierr = PCReset_MG(pc);CHKERRQ(ierr);
191     n    = mglevels[0]->levels;
192     for (i=0; i<n; i++) {
193       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
194         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
195       }
196       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
197       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
198     }
199     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
200   }
201 
202   mg->nlevels = levels;
203 
204   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
205   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
206 
207   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
208 
209   mg->stageApply = 0;
210   for (i=0; i<levels; i++) {
211     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
212 
213     mglevels[i]->level               = i;
214     mglevels[i]->levels              = levels;
215     mglevels[i]->cycles              = PC_MG_CYCLE_V;
216     mg->default_smoothu              = 2;
217     mg->default_smoothd              = 2;
218     mglevels[i]->eventsmoothsetup    = 0;
219     mglevels[i]->eventsmoothsolve    = 0;
220     mglevels[i]->eventresidual       = 0;
221     mglevels[i]->eventinterprestrict = 0;
222 
223     if (comms) comm = comms[i];
224     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
225     ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
226     ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
227     ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
228     ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
229     ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
230     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
231     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);CHKERRQ(ierr);
232     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
233 
234     /* do special stuff for coarse grid */
235     if (!i && levels > 1) {
236       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
237 
238       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
239       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
240       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
241       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
242       if (size > 1) {
243         KSP innerksp;
244         PC  innerpc;
245         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
246         ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr);
247         ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr);
248         ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
249       } else {
250         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
251         ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
252       }
253     } else {
254       char tprefix[128];
255       sprintf(tprefix,"mg_levels_%d_",(int)i);
256       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
257     }
258     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
259 
260     mglevels[i]->smoothu = mglevels[i]->smoothd;
261     mg->rtol             = 0.0;
262     mg->abstol           = 0.0;
263     mg->dtol             = 0.0;
264     mg->ttol             = 0.0;
265     mg->cyclesperpcapply = 1;
266   }
267   mg->am                   = PC_MG_MULTIPLICATIVE;
268   mg->levels               = mglevels;
269   pc->ops->applyrichardson = PCApplyRichardson_MG;
270   PetscFunctionReturn(0);
271 }
272 
273 
274 #undef __FUNCT__
275 #define __FUNCT__ "PCDestroy_MG"
276 PetscErrorCode PCDestroy_MG(PC pc)
277 {
278   PetscErrorCode ierr;
279   PC_MG          *mg        = (PC_MG*)pc->data;
280   PC_MG_Levels   **mglevels = mg->levels;
281   PetscInt       i,n;
282 
283   PetscFunctionBegin;
284   ierr = PCReset_MG(pc);CHKERRQ(ierr);
285   if (mglevels) {
286     n = mglevels[0]->levels;
287     for (i=0; i<n; i++) {
288       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
289         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
290       }
291       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
292       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
293     }
294     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
295   }
296   ierr = PetscFree(pc->data);CHKERRQ(ierr);
297   PetscFunctionReturn(0);
298 }
299 
300 
301 
302 extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
303 extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
304 extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
305 
306 /*
307    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
308              or full cycle of multigrid.
309 
310   Note:
311   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
312 */
313 #undef __FUNCT__
314 #define __FUNCT__ "PCApply_MG"
315 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
316 {
317   PC_MG          *mg        = (PC_MG*)pc->data;
318   PC_MG_Levels   **mglevels = mg->levels;
319   PetscErrorCode ierr;
320   PetscInt       levels = mglevels[0]->levels,i;
321 
322   PetscFunctionBegin;
323   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
324   /* When the DM is supplying the matrix then it will not exist until here */
325   for (i=0; i<levels; i++) {
326     if (!mglevels[i]->A) {
327       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
328       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
329     }
330   }
331 
332   mglevels[levels-1]->b = b;
333   mglevels[levels-1]->x = x;
334   if (mg->am == PC_MG_MULTIPLICATIVE) {
335     ierr = VecSet(x,0.0);CHKERRQ(ierr);
336     for (i=0; i<mg->cyclesperpcapply; i++) {
337       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
338     }
339   } else if (mg->am == PC_MG_ADDITIVE) {
340     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
341   } else if (mg->am == PC_MG_KASKADE) {
342     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
343   } else {
344     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
345   }
346   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
347   PetscFunctionReturn(0);
348 }
349 
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "PCSetFromOptions_MG"
353 PetscErrorCode PCSetFromOptions_MG(PC pc)
354 {
355   PetscErrorCode ierr;
356   PetscInt       m,levels = 1,cycles;
357   PetscBool      flg,set;
358   PC_MG          *mg        = (PC_MG*)pc->data;
359   PC_MG_Levels   **mglevels = mg->levels;
360   PCMGType       mgtype;
361   PCMGCycleType  mgctype;
362 
363   PetscFunctionBegin;
364   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
365   if (!mg->levels) {
366     ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
367     if (!flg && pc->dm) {
368       ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
369       levels++;
370       mg->usedmfornumberoflevels = PETSC_TRUE;
371     }
372     ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
373   }
374   mglevels = mg->levels;
375 
376   mgctype = (PCMGCycleType) mglevels[0]->cycles;
377   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
378   if (flg) {
379     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
380   }
381   flg  = PETSC_FALSE;
382   ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
383   if (set) {
384     ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
385   }
386   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr);
387   if (flg) {
388     ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
389   }
390   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr);
391   if (flg) {
392     ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
393   }
394   mgtype = mg->am;
395   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
396   if (flg) {
397     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
398   }
399   if (mg->am == PC_MG_MULTIPLICATIVE) {
400     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
401     if (flg) {
402       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
403     }
404   }
405   flg  = PETSC_FALSE;
406   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
407   if (flg) {
408     PetscInt i;
409     char     eventname[128];
410     if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
411     levels = mglevels[0]->levels;
412     for (i=0; i<levels; i++) {
413       sprintf(eventname,"MGSetup Level %d",(int)i);
414       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
415       sprintf(eventname,"MGSmooth Level %d",(int)i);
416       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
417       if (i) {
418         sprintf(eventname,"MGResid Level %d",(int)i);
419         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
420         sprintf(eventname,"MGInterp Level %d",(int)i);
421         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
422       }
423     }
424 
425 #if defined(PETSC_USE_LOG)
426     {
427       const char    *sname = "MG Apply";
428       PetscStageLog stageLog;
429       PetscInt      st;
430 
431       PetscFunctionBegin;
432       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
433       for (st = 0; st < stageLog->numStages; ++st) {
434         PetscBool same;
435 
436         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
437         if (same) mg->stageApply = st;
438       }
439       if (!mg->stageApply) {
440         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
441       }
442     }
443 #endif
444   }
445   ierr = PetscOptionsTail();CHKERRQ(ierr);
446   PetscFunctionReturn(0);
447 }
448 
449 const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
450 const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
451 
452 #include <petscdraw.h>
453 #undef __FUNCT__
454 #define __FUNCT__ "PCView_MG"
455 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
456 {
457   PC_MG          *mg        = (PC_MG*)pc->data;
458   PC_MG_Levels   **mglevels = mg->levels;
459   PetscErrorCode ierr;
460   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
461   PetscBool      iascii,isbinary,isdraw;
462 
463   PetscFunctionBegin;
464   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
465   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
466   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
467   if (iascii) {
468     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
469     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
470     if (mg->am == PC_MG_MULTIPLICATIVE) {
471       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
472     }
473     if (mg->galerkin) {
474       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
475     } else {
476       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
477     }
478     for (i=0; i<levels; i++) {
479       if (!i) {
480         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
481       } else {
482         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
483       }
484       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
485       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
486       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
487       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
488         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
489       } else if (i) {
490         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
491         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
492         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
493         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
494       }
495     }
496   } else if (isbinary) {
497     for (i=levels-1; i>=0; i--) {
498       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
499       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
500         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
501       }
502     }
503   } else if (isdraw) {
504     PetscDraw draw;
505     PetscReal x,w,y,bottom,th;
506     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
507     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
508     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
509     bottom = y - th;
510     for (i=levels-1; i>=0; i--) {
511       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
512         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
513         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
514         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
515       } else {
516         w    = 0.5*PetscMin(1.0-x,x);
517         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
518         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
519         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
520         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
521         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
522         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
523       }
524       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
525       bottom -= th;
526     }
527   }
528   PetscFunctionReturn(0);
529 }
530 
531 #include <petsc-private/dmimpl.h>
532 #include <petsc-private/kspimpl.h>
533 
534 /*
535     Calls setup for the KSP on each level
536 */
537 #undef __FUNCT__
538 #define __FUNCT__ "PCSetUp_MG"
539 PetscErrorCode PCSetUp_MG(PC pc)
540 {
541   PC_MG          *mg        = (PC_MG*)pc->data;
542   PC_MG_Levels   **mglevels = mg->levels;
543   PetscErrorCode ierr;
544   PetscInt       i,n = mglevels[0]->levels;
545   PC             cpc;
546   PetscBool      preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset,use_amat;
547   Mat            dA,dB;
548   Vec            tvec;
549   DM             *dms;
550   PetscViewer    viewer = 0;
551 
552   PetscFunctionBegin;
553   /* FIX: Move this to PCSetFromOptions_MG? */
554   if (mg->usedmfornumberoflevels) {
555     PetscInt levels;
556     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
557     levels++;
558     if (levels > n) { /* the problem is now being solved on a finer grid */
559       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
560       n        = levels;
561       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
562       mglevels =  mg->levels;
563     }
564   }
565   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
566 
567 
568   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
569   /* so use those from global PC */
570   /* Is this what we always want? What if user wants to keep old one? */
571   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
572   if (opsset) {
573     Mat mmat;
574     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
575     if (mmat == pc->pmat) opsset = PETSC_FALSE;
576   }
577 
578   if (!opsset) {
579     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
580     if(use_amat){
581       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);
582       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
583     }
584     else {
585       ierr = PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr);
586       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
587     }
588   }
589 
590   /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */
591   if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
592     /* construct the interpolation from the DMs */
593     Mat p;
594     Vec rscale;
595     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
596     dms[n-1] = pc->dm;
597     for (i=n-2; i>-1; i--) {
598       DMKSP kdm;
599       ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);
600       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
601       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
602       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
603       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
604        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
605       kdm->ops->computerhs = NULL;
606       kdm->rhsctx          = NULL;
607       if (!mglevels[i+1]->interpolate) {
608         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
609         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
610         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
611         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
612         ierr = MatDestroy(&p);CHKERRQ(ierr);
613       }
614     }
615 
616     for (i=n-2; i>-1; i--) {
617       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
618     }
619     ierr = PetscFree(dms);CHKERRQ(ierr);
620   }
621 
622   if (pc->dm && !pc->setupcalled) {
623     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
624     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
625     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
626   }
627 
628   if (mg->galerkin == 1) {
629     Mat B;
630     /* currently only handle case where mat and pmat are the same on coarser levels */
631     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
632     if (!pc->setupcalled) {
633       for (i=n-2; i>-1; i--) {
634         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct || !mglevels[i+1]->restrct) {
635           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
636         } else {
637           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
638         }
639         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
640         if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
641         dB = B;
642       }
643       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
644     } else {
645       for (i=n-2; i>-1; i--) {
646         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
647         if (mglevels[i+1]->interpolate == mglevels[i+1]->restrct || !mglevels[i+1]->restrct) {
648           ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
649         } else {
650           ierr = MatMatMatMult(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
651         }
652         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B);CHKERRQ(ierr);
653         dB   = B;
654       }
655     }
656   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
657     /* need to restrict Jacobian location to coarser meshes for evaluation */
658     for (i=n-2; i>-1; i--) {
659       Mat R;
660       Vec rscale;
661       if (!mglevels[i]->smoothd->dm->x) {
662         Vec *vecs;
663         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
664 
665         mglevels[i]->smoothd->dm->x = vecs[0];
666 
667         ierr = PetscFree(vecs);CHKERRQ(ierr);
668       }
669       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
670       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
671       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
672       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
673     }
674   }
675   if (!mg->galerkin && pc->dm) {
676     for (i=n-2; i>=0; i--) {
677       DM  dmfine,dmcoarse;
678       Mat Restrict,Inject;
679       Vec rscale;
680       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
681       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
682       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
683       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
684       Inject = NULL;      /* Callback should create it if it needs Injection */
685       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
686     }
687   }
688 
689   if (!pc->setupcalled) {
690     for (i=0; i<n; i++) {
691       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
692     }
693     for (i=1; i<n; i++) {
694       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
695         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
696       }
697     }
698     for (i=1; i<n; i++) {
699       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
700       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
701     }
702     for (i=0; i<n-1; i++) {
703       if (!mglevels[i]->b) {
704         Vec *vec;
705         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
706         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
707         ierr = VecDestroy(vec);CHKERRQ(ierr);
708         ierr = PetscFree(vec);CHKERRQ(ierr);
709       }
710       if (!mglevels[i]->r && i) {
711         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
712         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
713         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
714       }
715       if (!mglevels[i]->x) {
716         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
717         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
718         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
719       }
720     }
721     if (n != 1 && !mglevels[n-1]->r) {
722       /* PCMGSetR() on the finest level if user did not supply it */
723       Vec *vec;
724       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
725       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
726       ierr = VecDestroy(vec);CHKERRQ(ierr);
727       ierr = PetscFree(vec);CHKERRQ(ierr);
728     }
729   }
730 
731   if (pc->dm) {
732     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
733     for (i=0; i<n-1; i++) {
734       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
735     }
736   }
737 
738   for (i=1; i<n; i++) {
739     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
740       /* if doing only down then initial guess is zero */
741       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
742     }
743     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
744     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
745     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
746     if (!mglevels[i]->residual) {
747       Mat mat;
748       ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&mat);CHKERRQ(ierr);
749       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
750     }
751   }
752   for (i=1; i<n; i++) {
753     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
754       Mat          downmat,downpmat;
755 
756       /* check if operators have been set for up, if not use down operators to set them */
757       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
758       if (!opsset) {
759         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
760         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
761       }
762 
763       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
764       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
765       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
766       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
767     }
768   }
769 
770   /*
771       If coarse solver is not direct method then DO NOT USE preonly
772   */
773   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
774   if (preonly) {
775     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
776     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
777     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
778     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
779     if (!lu && !redundant && !cholesky && !svd) {
780       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
781     }
782   }
783 
784   if (!pc->setupcalled) {
785     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
786   }
787 
788   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
789   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
790   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
791 
792   /*
793      Dump the interpolation/restriction matrices plus the
794    Jacobian/stiffness on each level. This allows MATLAB users to
795    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
796 
797    Only support one or the other at the same time.
798   */
799 #if defined(PETSC_USE_SOCKET_VIEWER)
800   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
801   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
802   dump = PETSC_FALSE;
803 #endif
804   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
805   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
806 
807   if (viewer) {
808     for (i=1; i<n; i++) {
809       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
810     }
811     for (i=0; i<n; i++) {
812       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
813       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
814     }
815   }
816   PetscFunctionReturn(0);
817 }
818 
819 /* -------------------------------------------------------------------------------------*/
820 
821 #undef __FUNCT__
822 #define __FUNCT__ "PCMGGetLevels"
823 /*@
824    PCMGGetLevels - Gets the number of levels to use with MG.
825 
826    Not Collective
827 
828    Input Parameter:
829 .  pc - the preconditioner context
830 
831    Output parameter:
832 .  levels - the number of levels
833 
834    Level: advanced
835 
836 .keywords: MG, get, levels, multigrid
837 
838 .seealso: PCMGSetLevels()
839 @*/
840 PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
841 {
842   PC_MG *mg = (PC_MG*)pc->data;
843 
844   PetscFunctionBegin;
845   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
846   PetscValidIntPointer(levels,2);
847   *levels = mg->nlevels;
848   PetscFunctionReturn(0);
849 }
850 
851 #undef __FUNCT__
852 #define __FUNCT__ "PCMGSetType"
853 /*@
854    PCMGSetType - Determines the form of multigrid to use:
855    multiplicative, additive, full, or the Kaskade algorithm.
856 
857    Logically Collective on PC
858 
859    Input Parameters:
860 +  pc - the preconditioner context
861 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
862    PC_MG_FULL, PC_MG_KASKADE
863 
864    Options Database Key:
865 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
866    additive, full, kaskade
867 
868    Level: advanced
869 
870 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
871 
872 .seealso: PCMGSetLevels()
873 @*/
874 PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
875 {
876   PC_MG *mg = (PC_MG*)pc->data;
877 
878   PetscFunctionBegin;
879   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
880   PetscValidLogicalCollectiveEnum(pc,form,2);
881   mg->am = form;
882   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
883   else pc->ops->applyrichardson = NULL;
884   PetscFunctionReturn(0);
885 }
886 
887 /*@
888    PCMGGetType - Determines the form of multigrid to use:
889    multiplicative, additive, full, or the Kaskade algorithm.
890 
891    Logically Collective on PC
892 
893    Input Parameter:
894 .  pc - the preconditioner context
895 
896    Output Parameter:
897 .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
898 
899 
900    Level: advanced
901 
902 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
903 
904 .seealso: PCMGSetLevels()
905 @*/
906 PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
907 {
908   PC_MG *mg = (PC_MG*)pc->data;
909 
910   PetscFunctionBegin;
911   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
912   *type = mg->am;
913   PetscFunctionReturn(0);
914 }
915 
916 #undef __FUNCT__
917 #define __FUNCT__ "PCMGSetCycleType"
918 /*@
919    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
920    complicated cycling.
921 
922    Logically Collective on PC
923 
924    Input Parameters:
925 +  pc - the multigrid context
926 -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
927 
928    Options Database Key:
929 .  -pc_mg_cycle_type <v,w>
930 
931    Level: advanced
932 
933 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
934 
935 .seealso: PCMGSetCycleTypeOnLevel()
936 @*/
937 PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
938 {
939   PC_MG        *mg        = (PC_MG*)pc->data;
940   PC_MG_Levels **mglevels = mg->levels;
941   PetscInt     i,levels;
942 
943   PetscFunctionBegin;
944   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
945   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
946   PetscValidLogicalCollectiveInt(pc,n,2);
947   levels = mglevels[0]->levels;
948 
949   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
950   PetscFunctionReturn(0);
951 }
952 
953 #undef __FUNCT__
954 #define __FUNCT__ "PCMGMultiplicativeSetCycles"
955 /*@
956    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
957          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
958 
959    Logically Collective on PC
960 
961    Input Parameters:
962 +  pc - the multigrid context
963 -  n - number of cycles (default is 1)
964 
965    Options Database Key:
966 .  -pc_mg_multiplicative_cycles n
967 
968    Level: advanced
969 
970    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
971 
972 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
973 
974 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
975 @*/
976 PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
977 {
978   PC_MG        *mg        = (PC_MG*)pc->data;
979   PC_MG_Levels **mglevels = mg->levels;
980   PetscInt     i,levels;
981 
982   PetscFunctionBegin;
983   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
984   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
985   PetscValidLogicalCollectiveInt(pc,n,2);
986   levels = mglevels[0]->levels;
987 
988   for (i=0; i<levels; i++) mg->cyclesperpcapply = n;
989   PetscFunctionReturn(0);
990 }
991 
992 #undef __FUNCT__
993 #define __FUNCT__ "PCMGSetGalerkin"
994 /*@
995    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
996       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
997 
998    Logically Collective on PC
999 
1000    Input Parameters:
1001 +  pc - the multigrid context
1002 -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
1003 
1004    Options Database Key:
1005 .  -pc_mg_galerkin
1006 
1007    Level: intermediate
1008 
1009 .keywords: MG, set, Galerkin
1010 
1011 .seealso: PCMGGetGalerkin()
1012 
1013 @*/
1014 PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
1015 {
1016   PC_MG *mg = (PC_MG*)pc->data;
1017 
1018   PetscFunctionBegin;
1019   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1020   mg->galerkin = use ? 1 : 0;
1021   PetscFunctionReturn(0);
1022 }
1023 
1024 #undef __FUNCT__
1025 #define __FUNCT__ "PCMGGetGalerkin"
1026 /*@
1027    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1028       A_i-1 = r_i * A_i * p_i
1029 
1030    Not Collective
1031 
1032    Input Parameter:
1033 .  pc - the multigrid context
1034 
1035    Output Parameter:
1036 .  galerkin - PETSC_TRUE or PETSC_FALSE
1037 
1038    Options Database Key:
1039 .  -pc_mg_galerkin
1040 
1041    Level: intermediate
1042 
1043 .keywords: MG, set, Galerkin
1044 
1045 .seealso: PCMGSetGalerkin()
1046 
1047 @*/
1048 PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
1049 {
1050   PC_MG *mg = (PC_MG*)pc->data;
1051 
1052   PetscFunctionBegin;
1053   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1054   *galerkin = (PetscBool)mg->galerkin;
1055   PetscFunctionReturn(0);
1056 }
1057 
1058 #undef __FUNCT__
1059 #define __FUNCT__ "PCMGSetNumberSmoothDown"
1060 /*@
1061    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
1062    use on all levels. Use PCMGGetSmootherDown() to set different
1063    pre-smoothing steps on different levels.
1064 
1065    Logically Collective on PC
1066 
1067    Input Parameters:
1068 +  mg - the multigrid context
1069 -  n - the number of smoothing steps
1070 
1071    Options Database Key:
1072 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
1073 
1074    Level: advanced
1075 
1076 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
1077 
1078 .seealso: PCMGSetNumberSmoothUp()
1079 @*/
1080 PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
1081 {
1082   PC_MG          *mg        = (PC_MG*)pc->data;
1083   PC_MG_Levels   **mglevels = mg->levels;
1084   PetscErrorCode ierr;
1085   PetscInt       i,levels;
1086 
1087   PetscFunctionBegin;
1088   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1089   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1090   PetscValidLogicalCollectiveInt(pc,n,2);
1091   levels = mglevels[0]->levels;
1092 
1093   for (i=1; i<levels; i++) {
1094     /* make sure smoother up and down are different */
1095     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1096     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1097 
1098     mg->default_smoothd = n;
1099   }
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 #undef __FUNCT__
1104 #define __FUNCT__ "PCMGSetNumberSmoothUp"
1105 /*@
1106    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
1107    on all levels. Use PCMGGetSmootherUp() to set different numbers of
1108    post-smoothing steps on different levels.
1109 
1110    Logically Collective on PC
1111 
1112    Input Parameters:
1113 +  mg - the multigrid context
1114 -  n - the number of smoothing steps
1115 
1116    Options Database Key:
1117 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
1118 
1119    Level: advanced
1120 
1121    Note: this does not set a value on the coarsest grid, since we assume that
1122     there is no separate smooth up on the coarsest grid.
1123 
1124 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1125 
1126 .seealso: PCMGSetNumberSmoothDown()
1127 @*/
1128 PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
1129 {
1130   PC_MG          *mg        = (PC_MG*)pc->data;
1131   PC_MG_Levels   **mglevels = mg->levels;
1132   PetscErrorCode ierr;
1133   PetscInt       i,levels;
1134 
1135   PetscFunctionBegin;
1136   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1137   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1138   PetscValidLogicalCollectiveInt(pc,n,2);
1139   levels = mglevels[0]->levels;
1140 
1141   for (i=1; i<levels; i++) {
1142     /* make sure smoother up and down are different */
1143     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1144     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1145 
1146     mg->default_smoothu = n;
1147   }
1148   PetscFunctionReturn(0);
1149 }
1150 
1151 /* ----------------------------------------------------------------------------------------*/
1152 
1153 /*MC
1154    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1155     information about the coarser grid matrices and restriction/interpolation operators.
1156 
1157    Options Database Keys:
1158 +  -pc_mg_levels <nlevels> - number of levels including finest
1159 .  -pc_mg_cycles <v,w> -
1160 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
1161 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
1162 .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1163 .  -pc_mg_log - log information about time spent on each level of the solver
1164 .  -pc_mg_monitor - print information on the multigrid convergence
1165 .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1166 .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1167 .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1168                         to the Socket viewer for reading from MATLAB.
1169 -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1170                         to the binary output file called binaryoutput
1171 
1172    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
1173 
1174        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
1175 
1176    Level: intermediate
1177 
1178    Concepts: multigrid/multilevel
1179 
1180 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1181            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1182            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1183            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1184            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1185 M*/
1186 
1187 #undef __FUNCT__
1188 #define __FUNCT__ "PCCreate_MG"
1189 PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1190 {
1191   PC_MG          *mg;
1192   PetscErrorCode ierr;
1193 
1194   PetscFunctionBegin;
1195   ierr        = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1196   pc->data    = (void*)mg;
1197   mg->nlevels = -1;
1198 
1199   pc->useAmat = PETSC_TRUE;
1200 
1201   pc->ops->apply          = PCApply_MG;
1202   pc->ops->setup          = PCSetUp_MG;
1203   pc->ops->reset          = PCReset_MG;
1204   pc->ops->destroy        = PCDestroy_MG;
1205   pc->ops->setfromoptions = PCSetFromOptions_MG;
1206   pc->ops->view           = PCView_MG;
1207   PetscFunctionReturn(0);
1208 }
1209