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