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