xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 8fdcdf54d034d7b460de9a85cdb3b4202ec80101)
1 #define PETSCKSP_DLL
2 
3 /*
4     Defines the multigrid preconditioner interface.
5 */
6 #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
7 
8 
9 #undef __FUNCT__
10 #define __FUNCT__ "PCMGMCycle_Private"
11 PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG **mglevels,PetscTruth *converged)
12 {
13   PC_MG          *mg = *mglevels,*mgc;
14   PetscErrorCode ierr;
15   PetscInt       cycles = (PetscInt) mg->cycles;
16 
17   PetscFunctionBegin;
18   if (converged) *converged = PETSC_FALSE;
19 
20   if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
21   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);  /* pre-smooth */
22   if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
23   if (mg->level) {  /* not the coarsest grid */
24     if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
25     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
26     if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
27 
28     /* if on finest level and have convergence criteria set */
29     if (mg->level == mg->levels-1 && mg->ttol) {
30       PetscReal rnorm;
31       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
32       if (rnorm <= mg->ttol) {
33         *converged = PETSC_TRUE;
34         if (rnorm < mg->abstol) {
35           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
36         } else {
37           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
38         }
39         PetscFunctionReturn(0);
40       }
41     }
42 
43     mgc = *(mglevels - 1);
44     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
45     ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr);
46     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
47     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
48     while (cycles--) {
49       ierr = PCMGMCycle_Private(pc,mglevels-1,converged);CHKERRQ(ierr);
50     }
51     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
52     ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr);
53     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
54     if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
55     ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr);    /* post smooth */
56     if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 /*
62        PCMGCreate_Private - Creates a PC_MG structure for use with the
63                multigrid code. Level 0 is the coarsest. (But the
64                finest level is stored first in the array).
65 
66 */
67 #undef __FUNCT__
68 #define __FUNCT__ "PCMGCreate_Private"
69 static PetscErrorCode PCMGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,PC_MG ***result)
70 {
71   PC_MG          **mg;
72   PetscErrorCode ierr;
73   PetscInt       i;
74   PetscMPIInt    size;
75   const char     *prefix;
76   PC             ipc;
77 
78   PetscFunctionBegin;
79   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mg);CHKERRQ(ierr);
80   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
81 
82   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
83 
84   for (i=0; i<levels; i++) {
85     ierr = PetscNewLog(pc,PC_MG,&mg[i]);CHKERRQ(ierr);
86     mg[i]->level           = i;
87     mg[i]->levels          = levels;
88     mg[i]->cycles          = PC_MG_CYCLE_V;
89     mg[i]->galerkin        = PETSC_FALSE;
90     mg[i]->galerkinused    = PETSC_FALSE;
91     mg[i]->default_smoothu = 1;
92     mg[i]->default_smoothd = 1;
93 
94     if (comms) comm = comms[i];
95     ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr);
96     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr);
97     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
98 
99     /* do special stuff for coarse grid */
100     if (!i && levels > 1) {
101       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
102 
103       /* coarse solve is (redundant) LU by default */
104       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
105       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
106       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
107       if (size > 1) {
108         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
109       } else {
110         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
111       }
112 
113     } else {
114       char tprefix[128];
115       sprintf(tprefix,"mg_levels_%d_",(int)i);
116       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr);
117     }
118     ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr);
119     mg[i]->smoothu             = mg[i]->smoothd;
120     mg[i]->rtol                = 0.0;
121     mg[i]->abstol              = 0.0;
122     mg[i]->dtol                = 0.0;
123     mg[i]->ttol                = 0.0;
124     mg[i]->eventsmoothsetup    = 0;
125     mg[i]->eventsmoothsolve    = 0;
126     mg[i]->eventresidual       = 0;
127     mg[i]->eventinterprestrict = 0;
128     mg[i]->cyclesperpcapply    = 1;
129   }
130   *result = mg;
131   PetscFunctionReturn(0);
132 }
133 
134 #undef __FUNCT__
135 #define __FUNCT__ "PCDestroy_MG"
136 static PetscErrorCode PCDestroy_MG(PC pc)
137 {
138   PC_MG          **mg = (PC_MG**)pc->data;
139   PetscErrorCode ierr;
140   PetscInt       i,n;
141 
142   PetscFunctionBegin;
143   if (!mg) PetscFunctionReturn(0);
144   n = mg[0]->levels;
145   for (i=0; i<n-1; i++) {
146     if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);}
147     if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);}
148     if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);}
149     if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);}
150     if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);}
151   }
152 
153   for (i=0; i<n; i++) {
154     if (mg[i]->smoothd != mg[i]->smoothu) {
155       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
156     }
157     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
158     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
159   }
160   ierr = PetscFree(mg);CHKERRQ(ierr);
161   PetscFunctionReturn(0);
162 }
163 
164 
165 
166 EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**);
167 EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG**);
168 EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**);
169 
170 /*
171    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
172              or full cycle of multigrid.
173 
174   Note:
175   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
176 */
177 #undef __FUNCT__
178 #define __FUNCT__ "PCApply_MG"
179 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
180 {
181   PC_MG          **mg = (PC_MG**)pc->data;
182   PetscErrorCode ierr;
183   PetscInt       levels = mg[0]->levels,i;
184 
185   PetscFunctionBegin;
186   mg[levels-1]->b = b;
187   mg[levels-1]->x = x;
188   if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) {
189     Vec tvec;
190     ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr);
191     ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr);
192     ierr = VecDestroy(tvec);CHKERRQ(ierr);
193   }
194   if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
195     ierr = VecSet(x,0.0);CHKERRQ(ierr);
196     for (i=0; i<mg[0]->cyclesperpcapply; i++) {
197       ierr = PCMGMCycle_Private(pc,mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
198     }
199   }
200   else if (mg[0]->am == PC_MG_ADDITIVE) {
201     ierr = PCMGACycle_Private(mg);CHKERRQ(ierr);
202   }
203   else if (mg[0]->am == PC_MG_KASKADE) {
204     ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr);
205   }
206   else {
207     ierr = PCMGFCycle_Private(pc,mg);CHKERRQ(ierr);
208   }
209   PetscFunctionReturn(0);
210 }
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "PCApplyRichardson_MG"
214 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
215 {
216   PC_MG          **mg = (PC_MG**)pc->data;
217   PetscErrorCode ierr;
218   PetscInt       levels = mg[0]->levels;
219   PetscTruth     converged = PETSC_FALSE;
220 
221   PetscFunctionBegin;
222   mg[levels-1]->b    = b;
223   mg[levels-1]->x    = x;
224 
225   mg[levels-1]->rtol = rtol;
226   mg[levels-1]->abstol = abstol;
227   mg[levels-1]->dtol = dtol;
228   if (rtol) {
229     /* compute initial residual norm for relative convergence test */
230     PetscReal rnorm;
231     ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
232     ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
233     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
234   } else if (abstol) {
235     mg[levels-1]->ttol = abstol;
236   } else {
237     mg[levels-1]->ttol = 0.0;
238   }
239 
240   while (its-- && !converged) {
241     ierr = PCMGMCycle_Private(pc,mg+levels-1,&converged);CHKERRQ(ierr);
242   }
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "PCSetFromOptions_MG"
248 PetscErrorCode PCSetFromOptions_MG(PC pc)
249 {
250   PetscErrorCode ierr;
251   PetscInt       m,levels = 1,cycles;
252   PetscTruth     flg;
253   PC_MG          **mg = (PC_MG**)pc->data;
254   PCMGType       mgtype = PC_MG_ADDITIVE;
255   PCMGCycleType  mgctype;
256 
257   PetscFunctionBegin;
258   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
259     if (!pc->data) {
260       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
261       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
262       mg = (PC_MG**)pc->data;
263     }
264     mgctype = (PCMGCycleType) mg[0]->cycles;
265     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
266     if (flg) {
267       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
268     };
269     ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr);
270     if (flg) {
271       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
272     }
273     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
274     if (flg) {
275       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
276     }
277     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
278     if (flg) {
279       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
280     }
281     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
282     if (flg) {
283       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
284     }
285     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
286       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
287       if (flg) {
288 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
289       }
290     }
291     ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr);
292     if (flg) {
293       PetscInt i;
294       char     eventname[128];
295       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
296       levels = mg[0]->levels;
297       for (i=0; i<levels; i++) {
298         sprintf(eventname,"MGSetup Level %d",(int)i);
299         ierr = PetscLogEventRegister(&mg[i]->eventsmoothsetup,eventname,((PetscObject)pc)->cookie);CHKERRQ(ierr);
300         sprintf(eventname,"MGSmooth Level %d",(int)i);
301         ierr = PetscLogEventRegister(&mg[i]->eventsmoothsolve,eventname,((PetscObject)pc)->cookie);CHKERRQ(ierr);
302         if (i) {
303           sprintf(eventname,"MGResid Level %d",(int)i);
304           ierr = PetscLogEventRegister(&mg[i]->eventresidual,eventname,((PetscObject)pc)->cookie);CHKERRQ(ierr);
305           sprintf(eventname,"MGInterp Level %d",(int)i);
306           ierr = PetscLogEventRegister(&mg[i]->eventinterprestrict,eventname,((PetscObject)pc)->cookie);CHKERRQ(ierr);
307         }
308       }
309     }
310   ierr = PetscOptionsTail();CHKERRQ(ierr);
311   PetscFunctionReturn(0);
312 }
313 
314 const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
315 const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
316 
317 #undef __FUNCT__
318 #define __FUNCT__ "PCView_MG"
319 static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
320 {
321   PC_MG          **mg = (PC_MG**)pc->data;
322   PetscErrorCode ierr;
323   PetscInt       levels = mg[0]->levels,i;
324   PetscTruth     iascii;
325 
326   PetscFunctionBegin;
327   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
328   if (iascii) {
329     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s, pre-smooths=%D, post-smooths=%D\n",
330 				  PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w",
331                                   mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
332     if (mg[0]->galerkin) {
333       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
334     }
335     for (i=0; i<levels; i++) {
336       if (!i) {
337         ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
338       } else {
339         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
340       }
341       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
342       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
343       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
344       if (i && mg[i]->smoothd == mg[i]->smoothu) {
345         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
346       } else if (i){
347         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
348         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
349         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
350         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
351       }
352     }
353   } else {
354     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
355   }
356   PetscFunctionReturn(0);
357 }
358 
359 /*
360     Calls setup for the KSP on each level
361 */
362 #undef __FUNCT__
363 #define __FUNCT__ "PCSetUp_MG"
364 static PetscErrorCode PCSetUp_MG(PC pc)
365 {
366   PC_MG                   **mg = (PC_MG**)pc->data;
367   PetscErrorCode          ierr;
368   PetscInt                i,n = mg[0]->levels;
369   PC                      cpc,mpc;
370   PetscTruth              preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset;
371   PetscViewerASCIIMonitor ascii;
372   PetscViewer             viewer = PETSC_NULL;
373   MPI_Comm                comm;
374   Mat                     dA,dB;
375   MatStructure            uflag;
376   Vec                     tvec;
377 
378   PetscFunctionBegin;
379 
380   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
381   /* so use those from global PC */
382   /* Is this what we always want? What if user wants to keep old one? */
383   ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
384   ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
385   ierr = KSPGetPC(mg[n-1]->smoothd,&mpc);CHKERRQ(ierr);
386   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2))) {
387     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);
388     ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
389   }
390 
391   if (mg[0]->galerkin) {
392     Mat B;
393     mg[0]->galerkinused = PETSC_TRUE;
394     /* currently only handle case where mat and pmat are the same on coarser levels */
395     ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
396     if (!pc->setupcalled) {
397       for (i=n-2; i>-1; i--) {
398         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
399         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
400 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
401         dB   = B;
402       }
403       ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);
404     } else {
405       for (i=n-2; i>-1; i--) {
406         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
407         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
408         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
409         dB   = B;
410       }
411     }
412   }
413 
414   if (!pc->setupcalled) {
415     ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr);
416 
417     for (i=0; i<n; i++) {
418       if (monitor) {
419         ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
420         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
421         ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
422       }
423       ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
424     }
425     for (i=1; i<n; i++) {
426       if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) {
427         if (monitor) {
428           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
429           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
430           ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
431         }
432         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
433       }
434     }
435     for (i=1; i<n; i++) {
436       if (!mg[i]->residual) {
437         Mat mat;
438         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
439         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
440       }
441       if (mg[i]->restrct && !mg[i]->interpolate) {
442         ierr = PCMGSetInterpolation(pc,i,mg[i]->restrct);CHKERRQ(ierr);
443       }
444       if (!mg[i]->restrct && mg[i]->interpolate) {
445         ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr);
446       }
447 #if defined(PETSC_USE_DEBUG)
448       if (!mg[i]->restrct || !mg[i]->interpolate) {
449         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
450       }
451 #endif
452     }
453     for (i=0; i<n-1; i++) {
454       if (!mg[i]->b) {
455         Vec *vec;
456         ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
457         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
458         ierr = PetscFree(vec);CHKERRQ(ierr);
459       }
460       if (!mg[i]->r && i) {
461         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
462         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
463         ierr = VecDestroy(tvec);CHKERRQ(ierr);
464       }
465       if (!mg[i]->x) {
466         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
467         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
468         ierr = VecDestroy(tvec);CHKERRQ(ierr);
469       }
470     }
471     if (!mg[n-1]->r) {
472       Vec *vec;
473       ierr = KSPGetVecs(mg[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
474       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
475       ierr = PetscFree(vec);CHKERRQ(ierr);
476     }
477   }
478 
479 
480   for (i=1; i<n; i++) {
481     if (mg[i]->smoothu == mg[i]->smoothd) {
482       /* if doing only down then initial guess is zero */
483       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
484     }
485     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
486     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
487     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
488   }
489   for (i=1; i<n; i++) {
490     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
491       Mat          downmat,downpmat;
492       MatStructure matflag;
493       PetscTruth   opsset;
494 
495       /* check if operators have been set for up, if not use down operators to set them */
496       ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
497       if (!opsset) {
498         ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
499         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
500       }
501 
502       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
503       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
504       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
505       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
506     }
507   }
508 
509   /*
510       If coarse solver is not direct method then DO NOT USE preonly
511   */
512   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
513   if (preonly) {
514     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
515     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
516     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
517     if (!lu && !redundant && !cholesky) {
518       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
519     }
520   }
521 
522   if (!pc->setupcalled) {
523     if (monitor) {
524       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
525       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
526       ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
527     }
528     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
529   }
530 
531   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
532   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
533   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
534 
535   /*
536      Dump the interpolation/restriction matrices plus the
537    Jacobian/stiffness on each level. This allows Matlab users to
538    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
539 
540    Only support one or the other at the same time.
541   */
542 #if defined(PETSC_USE_SOCKET_VIEWER)
543   ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
544   if (dump) {
545     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
546   }
547 #endif
548   ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr);
549   if (dump) {
550     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
551   }
552 
553   if (viewer) {
554     for (i=1; i<n; i++) {
555       ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr);
556     }
557     for (i=0; i<n; i++) {
558       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
559       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
560     }
561   }
562   PetscFunctionReturn(0);
563 }
564 
565 /* -------------------------------------------------------------------------------------*/
566 
567 #undef __FUNCT__
568 #define __FUNCT__ "PCMGSetLevels"
569 /*@C
570    PCMGSetLevels - Sets the number of levels to use with MG.
571    Must be called before any other MG routine.
572 
573    Collective on PC
574 
575    Input Parameters:
576 +  pc - the preconditioner context
577 .  levels - the number of levels
578 -  comms - optional communicators for each level; this is to allow solving the coarser problems
579            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
580 
581    Level: intermediate
582 
583    Notes:
584      If the number of levels is one then the multigrid uses the -mg_levels prefix
585   for setting the level options rather than the -mg_coarse prefix.
586 
587 .keywords: MG, set, levels, multigrid
588 
589 .seealso: PCMGSetType(), PCMGGetLevels()
590 @*/
591 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
592 {
593   PetscErrorCode ierr;
594   PC_MG          **mg=0;
595 
596   PetscFunctionBegin;
597   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
598 
599   if (pc->data) {
600     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
601     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
602   }
603   ierr                     = PCMGCreate_Private(((PetscObject)pc)->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
604   mg[0]->am                = PC_MG_MULTIPLICATIVE;
605   pc->data                 = (void*)mg;
606   pc->ops->applyrichardson = PCApplyRichardson_MG;
607   PetscFunctionReturn(0);
608 }
609 
610 #undef __FUNCT__
611 #define __FUNCT__ "PCMGGetLevels"
612 /*@
613    PCMGGetLevels - Gets the number of levels to use with MG.
614 
615    Not Collective
616 
617    Input Parameter:
618 .  pc - the preconditioner context
619 
620    Output parameter:
621 .  levels - the number of levels
622 
623    Level: advanced
624 
625 .keywords: MG, get, levels, multigrid
626 
627 .seealso: PCMGSetLevels()
628 @*/
629 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
630 {
631   PC_MG  **mg;
632 
633   PetscFunctionBegin;
634   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
635   PetscValidIntPointer(levels,2);
636 
637   mg      = (PC_MG**)pc->data;
638   *levels = mg[0]->levels;
639   PetscFunctionReturn(0);
640 }
641 
642 #undef __FUNCT__
643 #define __FUNCT__ "PCMGSetType"
644 /*@
645    PCMGSetType - Determines the form of multigrid to use:
646    multiplicative, additive, full, or the Kaskade algorithm.
647 
648    Collective on PC
649 
650    Input Parameters:
651 +  pc - the preconditioner context
652 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
653    PC_MG_FULL, PC_MG_KASKADE
654 
655    Options Database Key:
656 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
657    additive, full, kaskade
658 
659    Level: advanced
660 
661 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
662 
663 .seealso: PCMGSetLevels()
664 @*/
665 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
666 {
667   PC_MG **mg;
668 
669   PetscFunctionBegin;
670   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
671   mg = (PC_MG**)pc->data;
672 
673   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
674   mg[0]->am = form;
675   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
676   else pc->ops->applyrichardson = 0;
677   PetscFunctionReturn(0);
678 }
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "PCMGSetCycleType"
682 /*@
683    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
684    complicated cycling.
685 
686    Collective on PC
687 
688    Input Parameters:
689 +  pc - the multigrid context
690 -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
691 
692    Options Database Key:
693 $  -pc_mg_cycle_type v or w
694 
695    Level: advanced
696 
697 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
698 
699 .seealso: PCMGSetCycleTypeOnLevel()
700 @*/
701 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n)
702 {
703   PC_MG    **mg;
704   PetscInt i,levels;
705 
706   PetscFunctionBegin;
707   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
708   mg     = (PC_MG**)pc->data;
709   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
710   levels = mg[0]->levels;
711 
712   for (i=0; i<levels; i++) {
713     mg[i]->cycles  = n;
714   }
715   PetscFunctionReturn(0);
716 }
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "PCMGMultiplicativeSetCycles"
720 /*@
721    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
722          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
723 
724    Collective on PC
725 
726    Input Parameters:
727 +  pc - the multigrid context
728 -  n - number of cycles (default is 1)
729 
730    Options Database Key:
731 $  -pc_mg_multiplicative_cycles n
732 
733    Level: advanced
734 
735    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
736 
737 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
738 
739 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
740 @*/
741 PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
742 {
743   PC_MG    **mg;
744   PetscInt i,levels;
745 
746   PetscFunctionBegin;
747   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
748   mg     = (PC_MG**)pc->data;
749   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
750   levels = mg[0]->levels;
751 
752   for (i=0; i<levels; i++) {
753     mg[i]->cyclesperpcapply  = n;
754   }
755   PetscFunctionReturn(0);
756 }
757 
758 #undef __FUNCT__
759 #define __FUNCT__ "PCMGSetGalerkin"
760 /*@
761    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
762       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
763 
764    Collective on PC
765 
766    Input Parameters:
767 .  pc - the multigrid context
768 
769    Options Database Key:
770 $  -pc_mg_galerkin
771 
772    Level: intermediate
773 
774 .keywords: MG, set, Galerkin
775 
776 .seealso: PCMGGetGalerkin()
777 
778 @*/
779 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
780 {
781   PC_MG    **mg;
782   PetscInt i,levels;
783 
784   PetscFunctionBegin;
785   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
786   mg     = (PC_MG**)pc->data;
787   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
788   levels = mg[0]->levels;
789 
790   for (i=0; i<levels; i++) {
791     mg[i]->galerkin = PETSC_TRUE;
792   }
793   PetscFunctionReturn(0);
794 }
795 
796 #undef __FUNCT__
797 #define __FUNCT__ "PCMGGetGalerkin"
798 /*@
799    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
800       A_i-1 = r_i * A_i * r_i^t
801 
802    Not Collective
803 
804    Input Parameter:
805 .  pc - the multigrid context
806 
807    Output Parameter:
808 .  gelerkin - PETSC_TRUE or PETSC_FALSE
809 
810    Options Database Key:
811 $  -pc_mg_galerkin
812 
813    Level: intermediate
814 
815 .keywords: MG, set, Galerkin
816 
817 .seealso: PCMGSetGalerkin()
818 
819 @*/
820 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
821 {
822   PC_MG    **mg;
823 
824   PetscFunctionBegin;
825   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
826   mg     = (PC_MG**)pc->data;
827   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
828   *galerkin = mg[0]->galerkin;
829   PetscFunctionReturn(0);
830 }
831 
832 #undef __FUNCT__
833 #define __FUNCT__ "PCMGSetNumberSmoothDown"
834 /*@
835    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
836    use on all levels. Use PCMGGetSmootherDown() to set different
837    pre-smoothing steps on different levels.
838 
839    Collective on PC
840 
841    Input Parameters:
842 +  mg - the multigrid context
843 -  n - the number of smoothing steps
844 
845    Options Database Key:
846 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
847 
848    Level: advanced
849 
850 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
851 
852 .seealso: PCMGSetNumberSmoothUp()
853 @*/
854 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
855 {
856   PC_MG          **mg;
857   PetscErrorCode ierr;
858   PetscInt       i,levels;
859 
860   PetscFunctionBegin;
861   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
862   mg     = (PC_MG**)pc->data;
863   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
864   levels = mg[0]->levels;
865 
866   for (i=1; i<levels; i++) {
867     /* make sure smoother up and down are different */
868     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
869     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
870     mg[i]->default_smoothd = n;
871   }
872   PetscFunctionReturn(0);
873 }
874 
875 #undef __FUNCT__
876 #define __FUNCT__ "PCMGSetNumberSmoothUp"
877 /*@
878    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
879    on all levels. Use PCMGGetSmootherUp() to set different numbers of
880    post-smoothing steps on different levels.
881 
882    Collective on PC
883 
884    Input Parameters:
885 +  mg - the multigrid context
886 -  n - the number of smoothing steps
887 
888    Options Database Key:
889 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
890 
891    Level: advanced
892 
893    Note: this does not set a value on the coarsest grid, since we assume that
894     there is no separate smooth up on the coarsest grid.
895 
896 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
897 
898 .seealso: PCMGSetNumberSmoothDown()
899 @*/
900 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
901 {
902   PC_MG          **mg;
903   PetscErrorCode ierr;
904   PetscInt       i,levels;
905 
906   PetscFunctionBegin;
907   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
908   mg     = (PC_MG**)pc->data;
909   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
910   levels = mg[0]->levels;
911 
912   for (i=1; i<levels; i++) {
913     /* make sure smoother up and down are different */
914     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
915     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
916     mg[i]->default_smoothu = n;
917   }
918   PetscFunctionReturn(0);
919 }
920 
921 /* ----------------------------------------------------------------------------------------*/
922 
923 /*MC
924    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
925     information about the coarser grid matrices and restriction/interpolation operators.
926 
927    Options Database Keys:
928 +  -pc_mg_levels <nlevels> - number of levels including finest
929 .  -pc_mg_cycles v or w
930 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
931 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
932 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
933 .  -pc_mg_log - log information about time spent on each level of the solver
934 .  -pc_mg_monitor - print information on the multigrid convergence
935 .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
936 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
937                         to the Socket viewer for reading from Matlab.
938 
939    Notes:
940 
941    Level: intermediate
942 
943    Concepts: multigrid/multilevel
944 
945 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
946            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
947            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
948            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
949            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
950 M*/
951 
952 EXTERN_C_BEGIN
953 #undef __FUNCT__
954 #define __FUNCT__ "PCCreate_MG"
955 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
956 {
957   PetscFunctionBegin;
958   pc->ops->apply          = PCApply_MG;
959   pc->ops->setup          = PCSetUp_MG;
960   pc->ops->destroy        = PCDestroy_MG;
961   pc->ops->setfromoptions = PCSetFromOptions_MG;
962   pc->ops->view           = PCView_MG;
963 
964   pc->data                = (void*)0;
965   PetscFunctionReturn(0);
966 }
967 EXTERN_C_END
968