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