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