xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision e39fd77f588fe8a9f63ea3dd3ce9052e17ee1991)
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 (n != 1 && !mg[n-1]->r) {
472       /* PCMGSetR() on the finest level if user did not supply it */
473       Vec *vec;
474       ierr = KSPGetVecs(mg[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
475       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
476       ierr = PetscFree(vec);CHKERRQ(ierr);
477     }
478   }
479 
480 
481   for (i=1; i<n; i++) {
482     if (mg[i]->smoothu == mg[i]->smoothd) {
483       /* if doing only down then initial guess is zero */
484       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
485     }
486     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
487     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
488     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
489   }
490   for (i=1; i<n; i++) {
491     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
492       Mat          downmat,downpmat;
493       MatStructure matflag;
494       PetscTruth   opsset;
495 
496       /* check if operators have been set for up, if not use down operators to set them */
497       ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
498       if (!opsset) {
499         ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
500         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
501       }
502 
503       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
504       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
505       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
506       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
507     }
508   }
509 
510   /*
511       If coarse solver is not direct method then DO NOT USE preonly
512   */
513   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
514   if (preonly) {
515     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
516     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
517     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
518     if (!lu && !redundant && !cholesky) {
519       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
520     }
521   }
522 
523   if (!pc->setupcalled) {
524     if (monitor) {
525       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
526       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
527       ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
528     }
529     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
530   }
531 
532   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
533   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
534   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
535 
536   /*
537      Dump the interpolation/restriction matrices plus the
538    Jacobian/stiffness on each level. This allows Matlab users to
539    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
540 
541    Only support one or the other at the same time.
542   */
543 #if defined(PETSC_USE_SOCKET_VIEWER)
544   ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
545   if (dump) {
546     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
547   }
548 #endif
549   ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr);
550   if (dump) {
551     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
552   }
553 
554   if (viewer) {
555     for (i=1; i<n; i++) {
556       ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr);
557     }
558     for (i=0; i<n; i++) {
559       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
560       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
561     }
562   }
563   PetscFunctionReturn(0);
564 }
565 
566 /* -------------------------------------------------------------------------------------*/
567 
568 #undef __FUNCT__
569 #define __FUNCT__ "PCMGSetLevels"
570 /*@C
571    PCMGSetLevels - Sets the number of levels to use with MG.
572    Must be called before any other MG routine.
573 
574    Collective on PC
575 
576    Input Parameters:
577 +  pc - the preconditioner context
578 .  levels - the number of levels
579 -  comms - optional communicators for each level; this is to allow solving the coarser problems
580            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
581 
582    Level: intermediate
583 
584    Notes:
585      If the number of levels is one then the multigrid uses the -mg_levels prefix
586   for setting the level options rather than the -mg_coarse prefix.
587 
588 .keywords: MG, set, levels, multigrid
589 
590 .seealso: PCMGSetType(), PCMGGetLevels()
591 @*/
592 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
593 {
594   PetscErrorCode ierr;
595   PC_MG          **mg=0;
596 
597   PetscFunctionBegin;
598   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
599 
600   if (pc->data) {
601     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
602     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
603   }
604   ierr                     = PCMGCreate_Private(((PetscObject)pc)->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
605   mg[0]->am                = PC_MG_MULTIPLICATIVE;
606   pc->data                 = (void*)mg;
607   pc->ops->applyrichardson = PCApplyRichardson_MG;
608   PetscFunctionReturn(0);
609 }
610 
611 #undef __FUNCT__
612 #define __FUNCT__ "PCMGGetLevels"
613 /*@
614    PCMGGetLevels - Gets the number of levels to use with MG.
615 
616    Not Collective
617 
618    Input Parameter:
619 .  pc - the preconditioner context
620 
621    Output parameter:
622 .  levels - the number of levels
623 
624    Level: advanced
625 
626 .keywords: MG, get, levels, multigrid
627 
628 .seealso: PCMGSetLevels()
629 @*/
630 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
631 {
632   PC_MG  **mg;
633 
634   PetscFunctionBegin;
635   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
636   PetscValidIntPointer(levels,2);
637 
638   mg      = (PC_MG**)pc->data;
639   *levels = mg[0]->levels;
640   PetscFunctionReturn(0);
641 }
642 
643 #undef __FUNCT__
644 #define __FUNCT__ "PCMGSetType"
645 /*@
646    PCMGSetType - Determines the form of multigrid to use:
647    multiplicative, additive, full, or the Kaskade algorithm.
648 
649    Collective on PC
650 
651    Input Parameters:
652 +  pc - the preconditioner context
653 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
654    PC_MG_FULL, PC_MG_KASKADE
655 
656    Options Database Key:
657 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
658    additive, full, kaskade
659 
660    Level: advanced
661 
662 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
663 
664 .seealso: PCMGSetLevels()
665 @*/
666 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
667 {
668   PC_MG **mg;
669 
670   PetscFunctionBegin;
671   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
672   mg = (PC_MG**)pc->data;
673 
674   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
675   mg[0]->am = form;
676   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
677   else pc->ops->applyrichardson = 0;
678   PetscFunctionReturn(0);
679 }
680 
681 #undef __FUNCT__
682 #define __FUNCT__ "PCMGSetCycleType"
683 /*@
684    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
685    complicated cycling.
686 
687    Collective on PC
688 
689    Input Parameters:
690 +  pc - the multigrid context
691 -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
692 
693    Options Database Key:
694 $  -pc_mg_cycle_type v or w
695 
696    Level: advanced
697 
698 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
699 
700 .seealso: PCMGSetCycleTypeOnLevel()
701 @*/
702 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n)
703 {
704   PC_MG    **mg;
705   PetscInt i,levels;
706 
707   PetscFunctionBegin;
708   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
709   mg     = (PC_MG**)pc->data;
710   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
711   levels = mg[0]->levels;
712 
713   for (i=0; i<levels; i++) {
714     mg[i]->cycles  = n;
715   }
716   PetscFunctionReturn(0);
717 }
718 
719 #undef __FUNCT__
720 #define __FUNCT__ "PCMGMultiplicativeSetCycles"
721 /*@
722    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
723          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
724 
725    Collective on PC
726 
727    Input Parameters:
728 +  pc - the multigrid context
729 -  n - number of cycles (default is 1)
730 
731    Options Database Key:
732 $  -pc_mg_multiplicative_cycles n
733 
734    Level: advanced
735 
736    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
737 
738 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
739 
740 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
741 @*/
742 PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
743 {
744   PC_MG    **mg;
745   PetscInt i,levels;
746 
747   PetscFunctionBegin;
748   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
749   mg     = (PC_MG**)pc->data;
750   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
751   levels = mg[0]->levels;
752 
753   for (i=0; i<levels; i++) {
754     mg[i]->cyclesperpcapply  = n;
755   }
756   PetscFunctionReturn(0);
757 }
758 
759 #undef __FUNCT__
760 #define __FUNCT__ "PCMGSetGalerkin"
761 /*@
762    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
763       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
764 
765    Collective on PC
766 
767    Input Parameters:
768 .  pc - the multigrid context
769 
770    Options Database Key:
771 $  -pc_mg_galerkin
772 
773    Level: intermediate
774 
775 .keywords: MG, set, Galerkin
776 
777 .seealso: PCMGGetGalerkin()
778 
779 @*/
780 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
781 {
782   PC_MG    **mg;
783   PetscInt i,levels;
784 
785   PetscFunctionBegin;
786   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
787   mg     = (PC_MG**)pc->data;
788   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
789   levels = mg[0]->levels;
790 
791   for (i=0; i<levels; i++) {
792     mg[i]->galerkin = PETSC_TRUE;
793   }
794   PetscFunctionReturn(0);
795 }
796 
797 #undef __FUNCT__
798 #define __FUNCT__ "PCMGGetGalerkin"
799 /*@
800    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
801       A_i-1 = r_i * A_i * r_i^t
802 
803    Not Collective
804 
805    Input Parameter:
806 .  pc - the multigrid context
807 
808    Output Parameter:
809 .  gelerkin - PETSC_TRUE or PETSC_FALSE
810 
811    Options Database Key:
812 $  -pc_mg_galerkin
813 
814    Level: intermediate
815 
816 .keywords: MG, set, Galerkin
817 
818 .seealso: PCMGSetGalerkin()
819 
820 @*/
821 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
822 {
823   PC_MG    **mg;
824 
825   PetscFunctionBegin;
826   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
827   mg     = (PC_MG**)pc->data;
828   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
829   *galerkin = mg[0]->galerkin;
830   PetscFunctionReturn(0);
831 }
832 
833 #undef __FUNCT__
834 #define __FUNCT__ "PCMGSetNumberSmoothDown"
835 /*@
836    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
837    use on all levels. Use PCMGGetSmootherDown() to set different
838    pre-smoothing steps on different levels.
839 
840    Collective on PC
841 
842    Input Parameters:
843 +  mg - the multigrid context
844 -  n - the number of smoothing steps
845 
846    Options Database Key:
847 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
848 
849    Level: advanced
850 
851 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
852 
853 .seealso: PCMGSetNumberSmoothUp()
854 @*/
855 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
856 {
857   PC_MG          **mg;
858   PetscErrorCode ierr;
859   PetscInt       i,levels;
860 
861   PetscFunctionBegin;
862   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
863   mg     = (PC_MG**)pc->data;
864   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
865   levels = mg[0]->levels;
866 
867   for (i=1; i<levels; i++) {
868     /* make sure smoother up and down are different */
869     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
870     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
871     mg[i]->default_smoothd = n;
872   }
873   PetscFunctionReturn(0);
874 }
875 
876 #undef __FUNCT__
877 #define __FUNCT__ "PCMGSetNumberSmoothUp"
878 /*@
879    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
880    on all levels. Use PCMGGetSmootherUp() to set different numbers of
881    post-smoothing steps on different levels.
882 
883    Collective on PC
884 
885    Input Parameters:
886 +  mg - the multigrid context
887 -  n - the number of smoothing steps
888 
889    Options Database Key:
890 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
891 
892    Level: advanced
893 
894    Note: this does not set a value on the coarsest grid, since we assume that
895     there is no separate smooth up on the coarsest grid.
896 
897 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
898 
899 .seealso: PCMGSetNumberSmoothDown()
900 @*/
901 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
902 {
903   PC_MG          **mg;
904   PetscErrorCode ierr;
905   PetscInt       i,levels;
906 
907   PetscFunctionBegin;
908   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
909   mg     = (PC_MG**)pc->data;
910   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
911   levels = mg[0]->levels;
912 
913   for (i=1; i<levels; i++) {
914     /* make sure smoother up and down are different */
915     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
916     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
917     mg[i]->default_smoothu = n;
918   }
919   PetscFunctionReturn(0);
920 }
921 
922 /* ----------------------------------------------------------------------------------------*/
923 
924 /*MC
925    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
926     information about the coarser grid matrices and restriction/interpolation operators.
927 
928    Options Database Keys:
929 +  -pc_mg_levels <nlevels> - number of levels including finest
930 .  -pc_mg_cycles v or w
931 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
932 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
933 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
934 .  -pc_mg_log - log information about time spent on each level of the solver
935 .  -pc_mg_monitor - print information on the multigrid convergence
936 .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
937 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
938                         to the Socket viewer for reading from Matlab.
939 
940    Notes:
941 
942    Level: intermediate
943 
944    Concepts: multigrid/multilevel
945 
946 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
947            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
948            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
949            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
950            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
951 M*/
952 
953 EXTERN_C_BEGIN
954 #undef __FUNCT__
955 #define __FUNCT__ "PCCreate_MG"
956 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
957 {
958   PetscFunctionBegin;
959   pc->ops->apply          = PCApply_MG;
960   pc->ops->setup          = PCSetUp_MG;
961   pc->ops->destroy        = PCDestroy_MG;
962   pc->ops->setfromoptions = PCSetFromOptions_MG;
963   pc->ops->view           = PCView_MG;
964 
965   pc->data                = (void*)0;
966   PetscFunctionReturn(0);
967 }
968 EXTERN_C_END
969