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