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