xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 18ebb3a5280c02fe10d3ba9cebbb96a615bebc44)
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   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
190   if (mg->nlevels == levels) PetscFunctionReturn(0);
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         m,levels = 1,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   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
360   if (!mg->levels) {
361     ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
362     if (!flg && 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   }
369   mglevels = mg->levels;
370 
371   mgctype = (PCMGCycleType) mglevels[0]->cycles;
372   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
373   if (flg) {
374     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
375   }
376   gtype = mg->galerkin;
377   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
378   if (flg) {
379     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
380   }
381   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr);
382   if (flg) {
383     ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
384   }
385   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr);
386   if (flg) {
387     ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
388   }
389   mgtype = mg->am;
390   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
391   if (flg) {
392     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
393   }
394   if (mg->am == PC_MG_MULTIPLICATIVE) {
395     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
396     if (flg) {
397       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
398     }
399   }
400   flg  = PETSC_FALSE;
401   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
402   if (flg) {
403     PetscInt i;
404     char     eventname[128];
405     if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
406     levels = mglevels[0]->levels;
407     for (i=0; i<levels; i++) {
408       sprintf(eventname,"MGSetup Level %d",(int)i);
409       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
410       sprintf(eventname,"MGSmooth Level %d",(int)i);
411       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
412       if (i) {
413         sprintf(eventname,"MGResid Level %d",(int)i);
414         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
415         sprintf(eventname,"MGInterp Level %d",(int)i);
416         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
417       }
418     }
419 
420 #if defined(PETSC_USE_LOG)
421     {
422       const char    *sname = "MG Apply";
423       PetscStageLog stageLog;
424       PetscInt      st;
425 
426       PetscFunctionBegin;
427       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
428       for (st = 0; st < stageLog->numStages; ++st) {
429         PetscBool same;
430 
431         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
432         if (same) mg->stageApply = st;
433       }
434       if (!mg->stageApply) {
435         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
436       }
437     }
438 #endif
439   }
440   ierr = PetscOptionsTail();CHKERRQ(ierr);
441   PetscFunctionReturn(0);
442 }
443 
444 const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
445 const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
446 const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};
447 
448 #include <petscdraw.h>
449 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
450 {
451   PC_MG          *mg        = (PC_MG*)pc->data;
452   PC_MG_Levels   **mglevels = mg->levels;
453   PetscErrorCode ierr;
454   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
455   PetscBool      iascii,isbinary,isdraw;
456 
457   PetscFunctionBegin;
458   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
459   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
460   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
461   if (iascii) {
462     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
463     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
464     if (mg->am == PC_MG_MULTIPLICATIVE) {
465       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
466     }
467     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
468       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
469     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
470       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
471     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
472       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
473     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
474       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
475     } else {
476       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
477     }
478     if (mg->view){
479       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
480     }
481     for (i=0; i<levels; i++) {
482       if (!i) {
483         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
484       } else {
485         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
486       }
487       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
488       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
489       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
490       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
491         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
492       } else if (i) {
493         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
494         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
495         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
496         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
497       }
498     }
499   } else if (isbinary) {
500     for (i=levels-1; i>=0; i--) {
501       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
502       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
503         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
504       }
505     }
506   } else if (isdraw) {
507     PetscDraw draw;
508     PetscReal x,w,y,bottom,th;
509     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
510     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
511     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
512     bottom = y - th;
513     for (i=levels-1; i>=0; i--) {
514       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
515         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
516         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
517         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
518       } else {
519         w    = 0.5*PetscMin(1.0-x,x);
520         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
521         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
522         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
523         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
524         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
525         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
526       }
527       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
528       bottom -= th;
529     }
530   }
531   PetscFunctionReturn(0);
532 }
533 
534 #include <petsc/private/dmimpl.h>
535 #include <petsc/private/kspimpl.h>
536 
537 /*
538     Calls setup for the KSP on each level
539 */
540 PetscErrorCode PCSetUp_MG(PC pc)
541 {
542   PC_MG          *mg        = (PC_MG*)pc->data;
543   PC_MG_Levels   **mglevels = mg->levels;
544   PetscErrorCode ierr;
545   PetscInt       i,n = mglevels[0]->levels;
546   PC             cpc;
547   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
548   Mat            dA,dB;
549   Vec            tvec;
550   DM             *dms;
551   PetscViewer    viewer = 0;
552   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
553 
554   PetscFunctionBegin;
555   /* FIX: Move this to PCSetFromOptions_MG? */
556   if (mg->usedmfornumberoflevels) {
557     PetscInt levels;
558     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
559     levels++;
560     if (levels > n) { /* the problem is now being solved on a finer grid */
561       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
562       n        = levels;
563       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
564       mglevels =  mg->levels;
565     }
566   }
567   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
568 
569 
570   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
571   /* so use those from global PC */
572   /* Is this what we always want? What if user wants to keep old one? */
573   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
574   if (opsset) {
575     Mat mmat;
576     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
577     if (mmat == pc->pmat) opsset = PETSC_FALSE;
578   }
579 
580   if (!opsset) {
581     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
582     if(use_amat){
583       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);
584       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
585     }
586     else {
587       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);
588       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
589     }
590   }
591 
592   for (i=n-1; i>0; i--) {
593     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
594       missinginterpolate = PETSC_TRUE;
595       continue;
596     }
597   }
598 
599   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
600   if (dA == dB) dAeqdB = PETSC_TRUE;
601   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
602     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
603   }
604 
605 
606   /*
607    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
608    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
609   */
610   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
611 	/* construct the interpolation from the DMs */
612     Mat p;
613     Vec rscale;
614     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
615     dms[n-1] = pc->dm;
616     /* Separately create them so we do not get DMKSP interference between levels */
617     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
618 	/*
619 	   Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
620 	   Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
621 	   But it is safe to use -dm_mat_type.
622 
623 	   The mat type should not be hardcoded like this, we need to find a better way.
624     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
625     */
626     for (i=n-2; i>-1; i--) {
627       DMKSP     kdm;
628       PetscBool dmhasrestrict;
629       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
630       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
631       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
632       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
633        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
634       kdm->ops->computerhs = NULL;
635       kdm->rhsctx          = NULL;
636       if (!mglevels[i+1]->interpolate) {
637         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
638         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
639         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
640         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
641         ierr = MatDestroy(&p);CHKERRQ(ierr);
642       }
643       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
644       if (dmhasrestrict && !mglevels[i+1]->restrct){
645         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
646         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
647         ierr = MatDestroy(&p);CHKERRQ(ierr);
648       }
649     }
650 
651     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
652     ierr = PetscFree(dms);CHKERRQ(ierr);
653   }
654 
655   if (pc->dm && !pc->setupcalled) {
656     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
657     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
658     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
659   }
660 
661   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
662     Mat       A,B;
663     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
664     MatReuse  reuse = MAT_INITIAL_MATRIX;
665 
666     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
667     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
668     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
669     for (i=n-2; i>-1; i--) {
670       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");
671       if (!mglevels[i+1]->interpolate) {
672         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
673       }
674       if (!mglevels[i+1]->restrct) {
675         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
676       }
677       if (reuse == MAT_REUSE_MATRIX) {
678         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
679       }
680       if (doA) {
681         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
682       }
683       if (doB) {
684         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
685       }
686       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
687       if (!doA && dAeqdB) {
688         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
689         A = B;
690       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
691         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
692         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
693       }
694       if (!doB && dAeqdB) {
695         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
696         B = A;
697       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
698         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
699         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
700       }
701       if (reuse == MAT_INITIAL_MATRIX) {
702         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
703         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
704         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
705       }
706       dA = A;
707       dB = B;
708     }
709   }
710   if (needRestricts && pc->dm && pc->dm->x) {
711     /* need to restrict Jacobian location to coarser meshes for evaluation */
712     for (i=n-2; i>-1; i--) {
713       Mat R;
714       Vec rscale;
715       if (!mglevels[i]->smoothd->dm->x) {
716         Vec *vecs;
717         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
718         mglevels[i]->smoothd->dm->x = vecs[0];
719         ierr = PetscFree(vecs);CHKERRQ(ierr);
720       }
721       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
722       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
723       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
724       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
725     }
726   }
727   if (needRestricts && pc->dm) {
728     for (i=n-2; i>=0; i--) {
729       DM  dmfine,dmcoarse;
730       Mat Restrict,Inject;
731       Vec rscale;
732       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
733       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
734       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
735       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
736       Inject = NULL;      /* Callback should create it if it needs Injection */
737       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
738     }
739   }
740 
741   if (!pc->setupcalled) {
742     for (i=0; i<n; i++) {
743       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
744     }
745     for (i=1; i<n; i++) {
746       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
747         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
748       }
749     }
750     /* insure that if either interpolation or restriction is set the other other one is set */
751     for (i=1; i<n; i++) {
752       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
753       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
754     }
755     for (i=0; i<n-1; i++) {
756       if (!mglevels[i]->b) {
757         Vec *vec;
758         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
759         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
760         ierr = VecDestroy(vec);CHKERRQ(ierr);
761         ierr = PetscFree(vec);CHKERRQ(ierr);
762       }
763       if (!mglevels[i]->r && i) {
764         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
765         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
766         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
767       }
768       if (!mglevels[i]->x) {
769         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
770         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
771         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
772       }
773     }
774     if (n != 1 && !mglevels[n-1]->r) {
775       /* PCMGSetR() on the finest level if user did not supply it */
776       Vec *vec;
777       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
778       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
779       ierr = VecDestroy(vec);CHKERRQ(ierr);
780       ierr = PetscFree(vec);CHKERRQ(ierr);
781     }
782   }
783 
784   if (pc->dm) {
785     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
786     for (i=0; i<n-1; i++) {
787       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
788     }
789   }
790 
791   for (i=1; i<n; i++) {
792     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
793       /* if doing only down then initial guess is zero */
794       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
795     }
796     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
797     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
798     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
799       pc->failedreason = PC_SUBPC_ERROR;
800     }
801     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
802     if (!mglevels[i]->residual) {
803       Mat mat;
804       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
805       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
806     }
807   }
808   for (i=1; i<n; i++) {
809     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
810       Mat downmat,downpmat;
811 
812       /* check if operators have been set for up, if not use down operators to set them */
813       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
814       if (!opsset) {
815         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
816         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
817       }
818 
819       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
820       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
821       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
822       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
823         pc->failedreason = PC_SUBPC_ERROR;
824       }
825       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
826     }
827   }
828 
829   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
830   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
831   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
832     pc->failedreason = PC_SUBPC_ERROR;
833   }
834   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
835 
836   /*
837      Dump the interpolation/restriction matrices plus the
838    Jacobian/stiffness on each level. This allows MATLAB users to
839    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
840 
841    Only support one or the other at the same time.
842   */
843 #if defined(PETSC_USE_SOCKET_VIEWER)
844   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
845   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
846   dump = PETSC_FALSE;
847 #endif
848   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
849   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
850 
851   if (viewer) {
852     for (i=1; i<n; i++) {
853       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
854     }
855     for (i=0; i<n; i++) {
856       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
857       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
858     }
859   }
860   PetscFunctionReturn(0);
861 }
862 
863 /* -------------------------------------------------------------------------------------*/
864 
865 /*@
866    PCMGGetLevels - Gets the number of levels to use with MG.
867 
868    Not Collective
869 
870    Input Parameter:
871 .  pc - the preconditioner context
872 
873    Output parameter:
874 .  levels - the number of levels
875 
876    Level: advanced
877 
878 .keywords: MG, get, levels, multigrid
879 
880 .seealso: PCMGSetLevels()
881 @*/
882 PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
883 {
884   PC_MG *mg = (PC_MG*)pc->data;
885 
886   PetscFunctionBegin;
887   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
888   PetscValidIntPointer(levels,2);
889   *levels = mg->nlevels;
890   PetscFunctionReturn(0);
891 }
892 
893 /*@
894    PCMGSetType - Determines the form of multigrid to use:
895    multiplicative, additive, full, or the Kaskade algorithm.
896 
897    Logically Collective on PC
898 
899    Input Parameters:
900 +  pc - the preconditioner context
901 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
902    PC_MG_FULL, PC_MG_KASKADE
903 
904    Options Database Key:
905 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
906    additive, full, kaskade
907 
908    Level: advanced
909 
910 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
911 
912 .seealso: PCMGSetLevels()
913 @*/
914 PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
915 {
916   PC_MG *mg = (PC_MG*)pc->data;
917 
918   PetscFunctionBegin;
919   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
920   PetscValidLogicalCollectiveEnum(pc,form,2);
921   mg->am = form;
922   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
923   else pc->ops->applyrichardson = NULL;
924   PetscFunctionReturn(0);
925 }
926 
927 /*@
928    PCMGGetType - Determines the form of multigrid to use:
929    multiplicative, additive, full, or the Kaskade algorithm.
930 
931    Logically Collective on PC
932 
933    Input Parameter:
934 .  pc - the preconditioner context
935 
936    Output Parameter:
937 .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
938 
939 
940    Level: advanced
941 
942 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
943 
944 .seealso: PCMGSetLevels()
945 @*/
946 PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
947 {
948   PC_MG *mg = (PC_MG*)pc->data;
949 
950   PetscFunctionBegin;
951   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
952   *type = mg->am;
953   PetscFunctionReturn(0);
954 }
955 
956 /*@
957    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
958    complicated cycling.
959 
960    Logically Collective on PC
961 
962    Input Parameters:
963 +  pc - the multigrid context
964 -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
965 
966    Options Database Key:
967 .  -pc_mg_cycle_type <v,w> - provide the cycle desired
968 
969    Level: advanced
970 
971 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
972 
973 .seealso: PCMGSetCycleTypeOnLevel()
974 @*/
975 PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
976 {
977   PC_MG        *mg        = (PC_MG*)pc->data;
978   PC_MG_Levels **mglevels = mg->levels;
979   PetscInt     i,levels;
980 
981   PetscFunctionBegin;
982   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
983   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
984   PetscValidLogicalCollectiveEnum(pc,n,2);
985   levels = mglevels[0]->levels;
986 
987   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
988   PetscFunctionReturn(0);
989 }
990 
991 /*@
992    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
993          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
994 
995    Logically Collective on PC
996 
997    Input Parameters:
998 +  pc - the multigrid context
999 -  n - number of cycles (default is 1)
1000 
1001    Options Database Key:
1002 .  -pc_mg_multiplicative_cycles n
1003 
1004    Level: advanced
1005 
1006    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
1007 
1008 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
1009 
1010 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1011 @*/
1012 PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1013 {
1014   PC_MG        *mg        = (PC_MG*)pc->data;
1015 
1016   PetscFunctionBegin;
1017   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1018   PetscValidLogicalCollectiveInt(pc,n,2);
1019   mg->cyclesperpcapply = n;
1020   PetscFunctionReturn(0);
1021 }
1022 
1023 PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1024 {
1025   PC_MG *mg = (PC_MG*)pc->data;
1026 
1027   PetscFunctionBegin;
1028   mg->galerkin = use;
1029   PetscFunctionReturn(0);
1030 }
1031 
1032 /*@
1033    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
1034       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1035 
1036    Logically Collective on PC
1037 
1038    Input Parameters:
1039 +  pc - the multigrid context
1040 -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1041 
1042    Options Database Key:
1043 .  -pc_mg_galerkin <both,pmat,mat,none>
1044 
1045    Level: intermediate
1046 
1047    Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
1048      use the PCMG construction of the coarser grids.
1049 
1050 .keywords: MG, set, Galerkin
1051 
1052 .seealso: PCMGGetGalerkin(), PCMGGalerkinType
1053 
1054 @*/
1055 PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1056 {
1057   PetscErrorCode ierr;
1058 
1059   PetscFunctionBegin;
1060   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1061   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1062   PetscFunctionReturn(0);
1063 }
1064 
1065 /*@
1066    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1067       A_i-1 = r_i * A_i * p_i
1068 
1069    Not Collective
1070 
1071    Input Parameter:
1072 .  pc - the multigrid context
1073 
1074    Output Parameter:
1075 .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL
1076 
1077    Level: intermediate
1078 
1079 .keywords: MG, set, Galerkin
1080 
1081 .seealso: PCMGSetGalerkin(), PCMGGalerkinType
1082 
1083 @*/
1084 PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
1085 {
1086   PC_MG *mg = (PC_MG*)pc->data;
1087 
1088   PetscFunctionBegin;
1089   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1090   *galerkin = mg->galerkin;
1091   PetscFunctionReturn(0);
1092 }
1093 
1094 /*@
1095    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
1096    use on all levels. Use PCMGGetSmootherDown() to set different
1097    pre-smoothing steps on different levels.
1098 
1099    Logically Collective on PC
1100 
1101    Input Parameters:
1102 +  mg - the multigrid context
1103 -  n - the number of smoothing steps
1104 
1105    Options Database Key:
1106 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
1107 
1108    Level: advanced
1109     If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the
1110    up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same
1111    number of smoothing steps for pre and post smoothing.
1112 
1113 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
1114 
1115 .seealso: PCMGSetNumberSmoothUp(), PCMGSetNumberSmooth()
1116 @*/
1117 PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
1118 {
1119   PC_MG          *mg        = (PC_MG*)pc->data;
1120   PC_MG_Levels   **mglevels = mg->levels;
1121   PetscErrorCode ierr;
1122   PetscInt       i,levels;
1123 
1124   PetscFunctionBegin;
1125   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1126   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1127   PetscValidLogicalCollectiveInt(pc,n,2);
1128   levels = mglevels[0]->levels;
1129 
1130   for (i=1; i<levels; i++) {
1131     PetscInt nc;
1132     ierr = KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);CHKERRQ(ierr);
1133     if (nc == n) continue;
1134 
1135     /* make sure smoother up and down are different */
1136     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1137     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1138 
1139     mg->default_smoothd = n;
1140   }
1141   PetscFunctionReturn(0);
1142 }
1143 
1144 /*@
1145    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
1146    on all levels. Use PCMGGetSmootherUp() to set different numbers of
1147    post-smoothing steps on different levels.
1148 
1149    Logically Collective on PC
1150 
1151    Input Parameters:
1152 +  mg - the multigrid context
1153 -  n - the number of smoothing steps
1154 
1155    Options Database Key:
1156 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
1157 
1158    Level: advanced
1159 
1160    Notes: this does not set a value on the coarsest grid, since we assume that
1161     there is no separate smooth up on the coarsest grid.
1162 
1163     If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the
1164    up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same
1165    number of smoothing steps for pre and post smoothing.
1166 
1167 
1168 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1169 
1170 .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmooth()
1171 @*/
1172 PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
1173 {
1174   PC_MG          *mg        = (PC_MG*)pc->data;
1175   PC_MG_Levels   **mglevels = mg->levels;
1176   PetscErrorCode ierr;
1177   PetscInt       i,levels;
1178 
1179   PetscFunctionBegin;
1180   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1181   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1182   PetscValidLogicalCollectiveInt(pc,n,2);
1183   levels = mglevels[0]->levels;
1184 
1185   for (i=1; i<levels; i++) {
1186     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
1187       PetscInt nc;
1188       ierr = KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);CHKERRQ(ierr);
1189       if (nc == n) continue;
1190     }
1191 
1192     /* make sure smoother up and down are different */
1193     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1194     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1195 
1196     mg->default_smoothu = n;
1197   }
1198   PetscFunctionReturn(0);
1199 }
1200 
1201 /*@
1202    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1203    on all levels. Use PCMGSetSmoothUp() and PCMGSetSmoothDown() set different numbers of
1204    pre ad post-smoothing steps
1205 
1206    Logically Collective on PC
1207 
1208    Input Parameters:
1209 +  mg - the multigrid context
1210 -  n - the number of smoothing steps
1211 
1212    Options Database Key:
1213 +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
1214 .  -pc_mg_smooth_down <n> - Sets number of pre-smoothing steps (if setting different pre and post amounts)
1215 -  -pc_mg_smooth_up <n> - Sets number of post-smoothing steps (if setting different pre and post amounts)
1216 
1217    Level: advanced
1218 
1219    Notes: this does not set a value on the coarsest grid, since we assume that
1220     there is no separate smooth up on the coarsest grid.
1221 
1222 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1223 
1224 .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmoothUp()
1225 @*/
1226 PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
1227 {
1228   PC_MG          *mg        = (PC_MG*)pc->data;
1229   PC_MG_Levels   **mglevels = mg->levels;
1230   PetscErrorCode ierr;
1231   PetscInt       i,levels;
1232 
1233   PetscFunctionBegin;
1234   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1235   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1236   PetscValidLogicalCollectiveInt(pc,n,2);
1237   levels = mglevels[0]->levels;
1238 
1239   for (i=1; i<levels; i++) {
1240     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1241     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
1242     mg->default_smoothu = n;
1243     mg->default_smoothd = n;
1244   }
1245   PetscFunctionReturn(0);
1246 }
1247 
1248 /* ----------------------------------------------------------------------------------------*/
1249 
1250 /*MC
1251    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1252     information about the coarser grid matrices and restriction/interpolation operators.
1253 
1254    Options Database Keys:
1255 +  -pc_mg_levels <nlevels> - number of levels including finest
1256 .  -pc_mg_cycle_type <v,w> -
1257 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
1258 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
1259 .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1260 .  -pc_mg_log - log information about time spent on each level of the solver
1261 .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1262 .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1263 .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1264                         to the Socket viewer for reading from MATLAB.
1265 -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1266                         to the binary output file called binaryoutput
1267 
1268    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
1269 
1270        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
1271 
1272        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
1273        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
1274        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
1275        residual is computed at the end of each cycle.
1276 
1277    Level: intermediate
1278 
1279    Concepts: multigrid/multilevel
1280 
1281 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1282            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1283            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1284            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1285            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1286 M*/
1287 
1288 PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1289 {
1290   PC_MG          *mg;
1291   PetscErrorCode ierr;
1292 
1293   PetscFunctionBegin;
1294   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1295   pc->data     = (void*)mg;
1296   mg->nlevels  = -1;
1297   mg->am       = PC_MG_MULTIPLICATIVE;
1298   mg->galerkin = PC_MG_GALERKIN_NONE;
1299 
1300   pc->useAmat = PETSC_TRUE;
1301 
1302   pc->ops->apply          = PCApply_MG;
1303   pc->ops->setup          = PCSetUp_MG;
1304   pc->ops->reset          = PCReset_MG;
1305   pc->ops->destroy        = PCDestroy_MG;
1306   pc->ops->setfromoptions = PCSetFromOptions_MG;
1307   pc->ops->view           = PCView_MG;
1308 
1309   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
1310   PetscFunctionReturn(0);
1311 }
1312