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