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