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