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