xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision c2be2410037ee5fc98eb6412cbeb90529affecc5)
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 /*
10        MGMCycle_Private - Given an MG structure created with MGCreate() runs
11                   one multiplicative cycle down through the levels and
12                   back up.
13 
14     Input Parameter:
15 .   mg - structure created with  MGCreate().
16 */
17 #undef __FUNCT__
18 #define __FUNCT__ "MGMCycle_Private"
19 PetscErrorCode MGMCycle_Private(MG *mglevels,PetscTruth *converged)
20 {
21   MG             mg = *mglevels,mgc;
22   PetscErrorCode ierr;
23   PetscInt       cycles = mg->cycles;
24   PetscScalar    zero = 0.0;
25 
26   PetscFunctionBegin;
27   if (converged) *converged = PETSC_FALSE;
28 
29   if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
30   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);
31   if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
32   if (mg->level) {  /* not the coarsest grid */
33     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
34 
35     /* if on finest level and have convergence criteria set */
36     if (mg->level == mg->levels-1 && mg->ttol) {
37       PetscReal rnorm;
38       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
39       if (rnorm <= mg->ttol) {
40         *converged = PETSC_TRUE;
41         if (rnorm < mg->abstol) {
42           ierr = PetscLogInfo((0,"MGMCycle_Private:Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",rnorm,mg->abstol));CHKERRQ(ierr);
43         } else {
44           ierr = PetscLogInfo((0,"MGMCycle_Private:Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",rnorm,mg->ttol));CHKERRQ(ierr);
45         }
46         PetscFunctionReturn(0);
47       }
48     }
49 
50     mgc = *(mglevels - 1);
51     ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr);
52     ierr = VecSet(&zero,mgc->x);CHKERRQ(ierr);
53     while (cycles--) {
54       ierr = MGMCycle_Private(mglevels-1,converged);CHKERRQ(ierr);
55     }
56     ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr);
57     if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
58     ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr);
59     if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);}
60   }
61   PetscFunctionReturn(0);
62 }
63 
64 /*
65        MGCreate_Private - Creates a MG structure for use with the
66                multigrid code. Level 0 is the coarsest. (But the
67                finest level is stored first in the array).
68 
69 */
70 #undef __FUNCT__
71 #define __FUNCT__ "MGCreate_Private"
72 static PetscErrorCode MGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,MG **result)
73 {
74   MG             *mg;
75   PetscErrorCode ierr;
76   PetscInt       i;
77   PetscMPIInt    size;
78   char           *prefix;
79   PC             ipc;
80 
81   PetscFunctionBegin;
82   ierr = PetscMalloc(levels*sizeof(MG),&mg);CHKERRQ(ierr);
83   ierr = PetscLogObjectMemory(pc,levels*(sizeof(MG)+sizeof(struct _MG)));CHKERRQ(ierr);
84 
85   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
86 
87   for (i=0; i<levels; i++) {
88     ierr = PetscNew(struct _MG,&mg[i]);CHKERRQ(ierr);
89     mg[i]->level           = i;
90     mg[i]->levels          = levels;
91     mg[i]->cycles          = 1;
92     mg[i]->galerkin        = PETSC_FALSE;
93     mg[i]->galerkinused    = PETSC_FALSE;
94     mg[i]->default_smoothu = 1;
95     mg[i]->default_smoothd = 1;
96 
97     if (comms) comm = comms[i];
98     ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr);
99     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr);
100     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
101 
102     /* do special stuff for coarse grid */
103     if (!i && levels > 1) {
104       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
105 
106       /* coarse solve is (redundant) LU by default */
107       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
108       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
109       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
110       if (size > 1) {
111         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
112         ierr = PCRedundantGetPC(ipc,&ipc);CHKERRQ(ierr);
113       }
114       ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
115 
116     } else {
117       char tprefix[128];
118       sprintf(tprefix,"mg_levels_%d_",(int)i);
119       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr);
120     }
121     ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr);
122     mg[i]->smoothu         = mg[i]->smoothd;
123     mg[i]->rtol = 0.0;
124     mg[i]->abstol = 0.0;
125     mg[i]->dtol = 0.0;
126     mg[i]->ttol = 0.0;
127     mg[i]->eventsetup = 0;
128     mg[i]->eventsolve = 0;
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   MG             *mg = (MG*)pc->data;
139   PetscErrorCode ierr;
140   PetscInt       i,n = mg[0]->levels;
141 
142   PetscFunctionBegin;
143   if (mg[0]->galerkinused) {
144     Mat B;
145     for (i=0; i<n-1; i++) {
146       ierr = KSPGetOperators(mg[i]->smoothd,0,&B,0);CHKERRQ(ierr);
147       ierr = MatDestroy(B);CHKERRQ(ierr);
148     }
149   }
150 
151   for (i=0; i<n; i++) {
152     if (mg[i]->smoothd != mg[i]->smoothu) {
153       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
154     }
155     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
156     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
157   }
158   ierr = PetscFree(mg);CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 
163 
164 EXTERN PetscErrorCode MGACycle_Private(MG*);
165 EXTERN PetscErrorCode MGFCycle_Private(MG*);
166 EXTERN PetscErrorCode MGKCycle_Private(MG*);
167 
168 /*
169    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
170              or full cycle of multigrid.
171 
172   Note:
173   A simple wrapper which calls MGMCycle(),MGACycle(), or MGFCycle().
174 */
175 #undef __FUNCT__
176 #define __FUNCT__ "PCApply_MG"
177 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
178 {
179   MG             *mg = (MG*)pc->data;
180   PetscScalar    zero = 0.0;
181   PetscErrorCode ierr;
182   PetscInt       levels = mg[0]->levels;
183 
184   PetscFunctionBegin;
185   mg[levels-1]->b = b;
186   mg[levels-1]->x = x;
187   if (mg[0]->am == MGMULTIPLICATIVE) {
188     ierr = VecSet(&zero,x);CHKERRQ(ierr);
189     ierr = MGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
190   }
191   else if (mg[0]->am == MGADDITIVE) {
192     ierr = MGACycle_Private(mg);CHKERRQ(ierr);
193   }
194   else if (mg[0]->am == MGKASKADE) {
195     ierr = MGKCycle_Private(mg);CHKERRQ(ierr);
196   }
197   else {
198     ierr = MGFCycle_Private(mg);CHKERRQ(ierr);
199   }
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "PCApplyRichardson_MG"
205 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
206 {
207   MG             *mg = (MG*)pc->data;
208   PetscErrorCode ierr;
209   PetscInt       levels = mg[0]->levels;
210   PetscTruth     converged = PETSC_FALSE;
211 
212   PetscFunctionBegin;
213   mg[levels-1]->b    = b;
214   mg[levels-1]->x    = x;
215 
216   mg[levels-1]->rtol = rtol;
217   mg[levels-1]->abstol = abstol;
218   mg[levels-1]->dtol = dtol;
219   if (rtol) {
220     /* compute initial residual norm for relative convergence test */
221     PetscReal rnorm;
222     ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
223     ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
224     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
225   } else if (abstol) {
226     mg[levels-1]->ttol = abstol;
227   } else {
228     mg[levels-1]->ttol = 0.0;
229   }
230 
231   while (its-- && !converged) {
232     ierr = MGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr);
233   }
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "PCSetFromOptions_MG"
239 static PetscErrorCode PCSetFromOptions_MG(PC pc)
240 {
241   PetscErrorCode ierr;
242   PetscInt       indx,m,levels = 1;
243   PetscTruth     flg;
244   const char     *type[] = {"additive","multiplicative","full","cascade","kascade"};
245 
246   PetscFunctionBegin;
247 
248   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
249     if (!pc->data) {
250       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","MGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
251       ierr = MGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
252     }
253     ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr);
254     if (flg) {
255       ierr = MGSetCycles(pc,m);CHKERRQ(ierr);
256     }
257     ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","MGSetGalerkin",&flg);CHKERRQ(ierr);
258     if (flg) {
259       ierr = MGSetGalerkin(pc);CHKERRQ(ierr);
260     }
261     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
262     if (flg) {
263       ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
264     }
265     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
266     if (flg) {
267       ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
268     }
269     ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr);
270     if (flg) {
271       MGType mg = (MGType) 0;
272       switch (indx) {
273       case 0:
274         mg = MGADDITIVE;
275         break;
276       case 1:
277         mg = MGMULTIPLICATIVE;
278         break;
279       case 2:
280         mg = MGFULL;
281         break;
282       case 3:
283         mg = MGKASKADE;
284         break;
285       case 4:
286         mg = MGKASKADE;
287         break;
288       }
289       ierr = MGSetType(pc,mg);CHKERRQ(ierr);
290     }
291     ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr);
292     if (flg) {
293       MG   *mg = (MG*)pc->data;
294       PetscInt i;
295       char eventname[128];
296       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
297       levels = mg[0]->levels;
298       for (i=0; i<levels; i++) {
299         sprintf(eventname,"MSetup Level %d",(int)i);
300         ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr);
301         sprintf(eventname,"MGSolve Level %d to 0",(int)i);
302         ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr);
303       }
304     }
305   ierr = PetscOptionsTail();CHKERRQ(ierr);
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "PCView_MG"
311 static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
312 {
313   MG             *mg = (MG*)pc->data;
314   PetscErrorCode ierr;
315   PetscInt       levels = mg[0]->levels,i;
316   const char     *cstring;
317   PetscTruth     iascii;
318 
319   PetscFunctionBegin;
320   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
321   if (iascii) {
322     if (mg[0]->am == MGMULTIPLICATIVE) cstring = "multiplicative";
323     else if (mg[0]->am == MGADDITIVE)  cstring = "additive";
324     else if (mg[0]->am == MGFULL)      cstring = "full";
325     else if (mg[0]->am == MGKASKADE)   cstring = "Kaskade";
326     else cstring = "unknown";
327     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%D, pre-smooths=%D, post-smooths=%D\n",
328                       cstring,levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
329     if (mg[0]->galerkin) {
330       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
331     }
332     for (i=0; i<levels; i++) {
333       if (!i) {
334         ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
335       } else {
336         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
337       }
338       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
339       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
340       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
341       if (i && mg[i]->smoothd == mg[i]->smoothu) {
342         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
343       } else if (i){
344         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
345         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
346         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
347         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
348       }
349     }
350   } else {
351     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
352   }
353   PetscFunctionReturn(0);
354 }
355 
356 /*
357     Calls setup for the KSP on each level
358 */
359 #undef __FUNCT__
360 #define __FUNCT__ "PCSetUp_MG"
361 static PetscErrorCode PCSetUp_MG(PC pc)
362 {
363   MG             *mg = (MG*)pc->data;
364   PetscErrorCode ierr;
365   PetscInt       i,n = mg[0]->levels;
366   PC             cpc;
367   PetscTruth     preonly,lu,redundant,monitor = PETSC_FALSE,dump;
368   PetscViewer    ascii;
369   MPI_Comm       comm;
370   Mat            dA,dB;
371   MatStructure   uflag;
372 
373   PetscFunctionBegin;
374   if (!pc->setupcalled) {
375     ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr);
376 
377     for (i=0; i<n; i++) {
378       if (monitor) {
379         ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
380         ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
381         ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
382         ierr = KSPSetMonitor(mg[i]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
383       }
384       ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
385     }
386     for (i=1; i<n; i++) {
387       if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
388         if (monitor) {
389           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
390           ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
391           ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr);
392           ierr = KSPSetMonitor(mg[i]->smoothu,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
393         }
394         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
395       }
396     }
397   }
398 
399   /* If user did not provide fine grid operators, use those from PC */
400   /* BUG BUG BUG This will work ONLY the first time called: hence if the user changes
401      the PC matrices between solves PCMG will continue to use first set provided */
402   ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
403   if (!dA  && !dB) {
404     ierr = PetscLogInfo((pc,"PCSetUp_MG: Using outer operators to define finest grid operator \n  because MGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n"));
405     ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,uflag);CHKERRQ(ierr);
406   }
407 
408   if (mg[0]->galerkin) {
409     mg[0]->galerkinused = PETSC_TRUE;
410     /* currently only handle case where mat and pmat are the same on coarser levels */
411     Mat B;
412     ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
413     if (!pc->setupcalled) {
414       for (i=n-2; i>-1; i--) {
415         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
416         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
417         dB   = B;
418       }
419     } else {
420       for (i=n-2; i>-1; i--) {
421         ierr = KSPGetOperators(mg[i]->smoothd,0,&B,0);CHKERRQ(ierr);
422         ierr = MatPtAP(dB,mg[i]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
423         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
424         dB   = B;
425       }
426     }
427   }
428 
429   for (i=1; i<n; i++) {
430     if (mg[i]->smoothu == mg[i]->smoothd) {
431       /* if doing only down then initial guess is zero */
432       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
433     }
434     if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
435     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
436     if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
437   }
438   for (i=1; i<n; i++) {
439     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
440       PC           uppc,downpc;
441       Mat          downmat,downpmat,upmat,uppmat;
442       MatStructure matflag;
443 
444       /* check if operators have been set for up, if not use down operators to set them */
445       ierr = KSPGetPC(mg[i]->smoothu,&uppc);CHKERRQ(ierr);
446       ierr = PCGetOperators(uppc,&upmat,&uppmat,PETSC_NULL);CHKERRQ(ierr);
447       if (!upmat) {
448         ierr = KSPGetPC(mg[i]->smoothd,&downpc);CHKERRQ(ierr);
449         ierr = PCGetOperators(downpc,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
450         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
451       }
452 
453       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
454       if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
455       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
456       if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
457     }
458   }
459 
460   /*
461       If coarse solver is not direct method then DO NOT USE preonly
462   */
463   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
464   if (preonly) {
465     ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
466     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
467     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
468     if (!lu && !redundant) {
469       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
470     }
471   }
472 
473   if (!pc->setupcalled) {
474     if (monitor) {
475       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
476       ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr);
477       ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr);
478       ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(PetscErrorCode(*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
479     }
480     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
481   }
482 
483   if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
484   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
485   if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);}
486 
487   /*
488      Dump the interpolation/restriction matrices to matlab plus the
489    Jacobian/stiffness on each level. This allows Matlab users to
490    easily check if the Galerkin condition A_c = R A_f R^T is satisfied */
491   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
492   if (dump) {
493     for (i=1; i<n; i++) {
494       ierr = MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
495     }
496     for (i=0; i<n; i++) {
497       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
498       ierr = MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr);
499     }
500   }
501   ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr);
502   if (dump) {
503     for (i=1; i<n; i++) {
504       ierr = MatView(mg[i]->restrct,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr);
505     }
506     for (i=0; i<n; i++) {
507       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
508       ierr = MatView(pc->mat,PETSC_VIEWER_BINARY_(pc->comm));CHKERRQ(ierr);
509     }
510   }
511   PetscFunctionReturn(0);
512 }
513 
514 /* -------------------------------------------------------------------------------------*/
515 
516 #undef __FUNCT__
517 #define __FUNCT__ "MGSetLevels"
518 /*@C
519    MGSetLevels - Sets the number of levels to use with MG.
520    Must be called before any other MG routine.
521 
522    Collective on PC
523 
524    Input Parameters:
525 +  pc - the preconditioner context
526 .  levels - the number of levels
527 -  comms - optional communicators for each level; this is to allow solving the coarser problems
528            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
529 
530    Level: intermediate
531 
532    Notes:
533      If the number of levels is one then the multigrid uses the -mg_levels prefix
534   for setting the level options rather than the -mg_coarse prefix.
535 
536 .keywords: MG, set, levels, multigrid
537 
538 .seealso: MGSetType(), MGGetLevels()
539 @*/
540 PetscErrorCode PETSCKSP_DLLEXPORT MGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
541 {
542   PetscErrorCode ierr;
543   MG             *mg;
544 
545   PetscFunctionBegin;
546   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
547 
548   if (pc->data) {
549     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
550     make sure that you call MGSetLevels() before KSPSetFromOptions()");
551   }
552   ierr                     = MGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
553   mg[0]->am                = MGMULTIPLICATIVE;
554   pc->data                 = (void*)mg;
555   pc->ops->applyrichardson = PCApplyRichardson_MG;
556   PetscFunctionReturn(0);
557 }
558 
559 #undef __FUNCT__
560 #define __FUNCT__ "MGGetLevels"
561 /*@
562    MGGetLevels - Gets the number of levels to use with MG.
563 
564    Not Collective
565 
566    Input Parameter:
567 .  pc - the preconditioner context
568 
569    Output parameter:
570 .  levels - the number of levels
571 
572    Level: advanced
573 
574 .keywords: MG, get, levels, multigrid
575 
576 .seealso: MGSetLevels()
577 @*/
578 PetscErrorCode PETSCKSP_DLLEXPORT MGGetLevels(PC pc,PetscInt *levels)
579 {
580   MG  *mg;
581 
582   PetscFunctionBegin;
583   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
584   PetscValidIntPointer(levels,2);
585 
586   mg      = (MG*)pc->data;
587   *levels = mg[0]->levels;
588   PetscFunctionReturn(0);
589 }
590 
591 #undef __FUNCT__
592 #define __FUNCT__ "MGSetType"
593 /*@
594    MGSetType - Determines the form of multigrid to use:
595    multiplicative, additive, full, or the Kaskade algorithm.
596 
597    Collective on PC
598 
599    Input Parameters:
600 +  pc - the preconditioner context
601 -  form - multigrid form, one of MGMULTIPLICATIVE, MGADDITIVE,
602    MGFULL, MGKASKADE
603 
604    Options Database Key:
605 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
606    additive, full, kaskade
607 
608    Level: advanced
609 
610 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
611 
612 .seealso: MGSetLevels()
613 @*/
614 PetscErrorCode PETSCKSP_DLLEXPORT MGSetType(PC pc,MGType form)
615 {
616   MG *mg;
617 
618   PetscFunctionBegin;
619   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
620   mg = (MG*)pc->data;
621 
622   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
623   mg[0]->am = form;
624   if (form == MGMULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
625   else pc->ops->applyrichardson = 0;
626   PetscFunctionReturn(0);
627 }
628 
629 #undef __FUNCT__
630 #define __FUNCT__ "MGSetCycles"
631 /*@
632    MGSetCycles - Sets the type cycles to use.  Use MGSetCyclesOnLevel() for more
633    complicated cycling.
634 
635    Collective on PC
636 
637    Input Parameters:
638 +  pc - the multigrid context
639 -  n - the number of cycles
640 
641    Options Database Key:
642 $  -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle.
643 
644    Level: advanced
645 
646 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
647 
648 .seealso: MGSetCyclesOnLevel()
649 @*/
650 PetscErrorCode PETSCKSP_DLLEXPORT MGSetCycles(PC pc,PetscInt n)
651 {
652   MG       *mg;
653   PetscInt i,levels;
654 
655   PetscFunctionBegin;
656   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
657   mg     = (MG*)pc->data;
658   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
659   levels = mg[0]->levels;
660 
661   for (i=0; i<levels; i++) {
662     mg[i]->cycles  = n;
663   }
664   PetscFunctionReturn(0);
665 }
666 
667 #undef __FUNCT__
668 #define __FUNCT__ "MGSetGalerkin"
669 /*@
670    MGSetGalerkin - Causes the coarser grid matrices to be computed from the
671       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
672 
673    Collective on PC
674 
675    Input Parameters:
676 +  pc - the multigrid context
677 -  n - the number of cycles
678 
679    Options Database Key:
680 $  -pc_mg_galerkin
681 
682    Level: intermediate
683 
684 .keywords: MG, set, Galerkin
685 
686 @*/
687 PetscErrorCode PETSCKSP_DLLEXPORT MGSetGalerkin(PC pc)
688 {
689   MG       *mg;
690   PetscInt i,levels;
691 
692   PetscFunctionBegin;
693   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
694   mg     = (MG*)pc->data;
695   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
696   levels = mg[0]->levels;
697 
698   for (i=0; i<levels; i++) {
699     mg[i]->galerkin = PETSC_TRUE;
700   }
701   PetscFunctionReturn(0);
702 }
703 
704 #undef __FUNCT__
705 #define __FUNCT__ "MGCheck"
706 /*@
707    MGCheck - Checks that all components of the MG structure have
708    been set.
709 
710    Collective on PC
711 
712    Input Parameters:
713 .  mg - the MG structure
714 
715    Level: advanced
716 
717 .keywords: MG, check, set, multigrid
718 @*/
719 PetscErrorCode PETSCKSP_DLLEXPORT MGCheck(PC pc)
720 {
721   MG       *mg;
722   PetscInt i,n,count = 0;
723 
724   PetscFunctionBegin;
725   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
726   mg = (MG*)pc->data;
727 
728   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
729 
730   n = mg[0]->levels;
731 
732   for (i=1; i<n; i++) {
733     if (!mg[i]->restrct) {
734       (*PetscErrorPrintf)("No restrict set level %D \n",n-i); count++;
735     }
736     if (!mg[i]->interpolate) {
737       (*PetscErrorPrintf)("No interpolate set level %D \n",n-i); count++;
738     }
739     if (!mg[i]->r) {
740       (*PetscErrorPrintf)("No r set level %D \n",n-i); count++;
741     }
742     if (!mg[i-1]->x) {
743       (*PetscErrorPrintf)("No x set level %D \n",n-i); count++;
744     }
745     if (!mg[i-1]->b) {
746       (*PetscErrorPrintf)("No b set level %D \n",n-i); count++;
747     }
748   }
749   PetscFunctionReturn(count);
750 }
751 
752 
753 #undef __FUNCT__
754 #define __FUNCT__ "MGSetNumberSmoothDown"
755 /*@
756    MGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
757    use on all levels. Use MGGetSmootherDown() to set different
758    pre-smoothing steps on different levels.
759 
760    Collective on PC
761 
762    Input Parameters:
763 +  mg - the multigrid context
764 -  n - the number of smoothing steps
765 
766    Options Database Key:
767 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
768 
769    Level: advanced
770 
771 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
772 
773 .seealso: MGSetNumberSmoothUp()
774 @*/
775 PetscErrorCode PETSCKSP_DLLEXPORT MGSetNumberSmoothDown(PC pc,PetscInt n)
776 {
777   MG             *mg;
778   PetscErrorCode ierr;
779   PetscInt       i,levels;
780 
781   PetscFunctionBegin;
782   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
783   mg     = (MG*)pc->data;
784   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
785   levels = mg[0]->levels;
786 
787   for (i=0; i<levels; i++) {
788     /* make sure smoother up and down are different */
789     ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
790     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
791     mg[i]->default_smoothd = n;
792   }
793   PetscFunctionReturn(0);
794 }
795 
796 #undef __FUNCT__
797 #define __FUNCT__ "MGSetNumberSmoothUp"
798 /*@
799    MGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
800    on all levels. Use MGGetSmootherUp() to set different numbers of
801    post-smoothing steps on different levels.
802 
803    Collective on PC
804 
805    Input Parameters:
806 +  mg - the multigrid context
807 -  n - the number of smoothing steps
808 
809    Options Database Key:
810 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
811 
812    Level: advanced
813 
814    Note: this does not set a value on the coarsest grid, since we assume that
815     there is no seperate smooth up on the coarsest grid.
816 
817 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
818 
819 .seealso: MGSetNumberSmoothDown()
820 @*/
821 PetscErrorCode PETSCKSP_DLLEXPORT MGSetNumberSmoothUp(PC pc,PetscInt n)
822 {
823   MG             *mg;
824   PetscErrorCode ierr;
825   PetscInt       i,levels;
826 
827   PetscFunctionBegin;
828   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
829   mg     = (MG*)pc->data;
830   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
831   levels = mg[0]->levels;
832 
833   for (i=1; i<levels; i++) {
834     /* make sure smoother up and down are different */
835     ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
836     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
837     mg[i]->default_smoothu = n;
838   }
839   PetscFunctionReturn(0);
840 }
841 
842 /* ----------------------------------------------------------------------------------------*/
843 
844 /*MC
845    PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional
846     information about the coarser grid matrices and restriction/interpolation operators.
847 
848    Options Database Keys:
849 +  -pc_mg_levels <nlevels> - number of levels including finest
850 .  -pc_mg_cycles 1 or 2 - for V or W-cycle
851 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
852 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
853 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
854 .  -pc_mg_log - log information about time spent on each level of the solver
855 .  -pc_mg_monitor - print information on the multigrid convergence
856 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
857                         to the Socket viewer for reading from Matlab.
858 
859    Notes:
860 
861    Level: intermediate
862 
863    Concepts: multigrid
864 
865 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
866            MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(),
867            MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(),
868            MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(),
869            MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR()
870 M*/
871 
872 EXTERN_C_BEGIN
873 #undef __FUNCT__
874 #define __FUNCT__ "PCCreate_MG"
875 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
876 {
877   PetscFunctionBegin;
878   pc->ops->apply          = PCApply_MG;
879   pc->ops->setup          = PCSetUp_MG;
880   pc->ops->destroy        = PCDestroy_MG;
881   pc->ops->setfromoptions = PCSetFromOptions_MG;
882   pc->ops->view           = PCView_MG;
883 
884   pc->data                = (void*)0;
885   PetscFunctionReturn(0);
886 }
887 EXTERN_C_END
888