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