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