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