xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 4d91e141dc6054a2aad51cdf1cb23e4e8696c81e)
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_ERR_ORDER,"Number levels already set for MG\n  make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
152   }
153   if (mg->levels) SETERRQ(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_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_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 = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr);
467       ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
468       ierr = MatDestroy(p);CHKERRQ(ierr);
469     }
470 
471     if (!mg->galerkin) {
472       /* each coarse level gets its DM */
473       for (i=n-2; i>-1; i--) {
474         ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
475       }
476     }
477 
478     for (i=n-2; i>-1; i--) {
479       ierr = DMDestroy(dms[i]);CHKERRQ(ierr);
480     }
481     ierr = PetscFree(dms);CHKERRQ(ierr);
482   }
483 
484   if (mg->galerkin) {
485     Mat B;
486     mg->galerkinused = PETSC_TRUE;
487     /* currently only handle case where mat and pmat are the same on coarser levels */
488     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
489     if (!pc->setupcalled) {
490       for (i=n-2; i>-1; i--) {
491         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
492         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
493 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
494         dB   = B;
495       }
496       ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);
497     } else {
498       for (i=n-2; i>-1; i--) {
499         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
500         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
501         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
502         dB   = B;
503       }
504     }
505   }
506 
507   if (!pc->setupcalled) {
508     ierr = PetscOptionsGetTruth(0,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
509 
510     for (i=0; i<n; i++) {
511       if (monitor) {
512         ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr);
513         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
514         ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
515       }
516       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
517     }
518     for (i=1; i<n; i++) {
519       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
520         if (monitor) {
521           ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr);
522           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
523           ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
524         }
525         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
526       }
527     }
528     for (i=1; i<n; i++) {
529       if (!mglevels[i]->residual) {
530         Mat mat;
531         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
532         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
533       }
534       if (mglevels[i]->restrct && !mglevels[i]->interpolate) {
535         ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr);
536       }
537       if (!mglevels[i]->restrct && mglevels[i]->interpolate) {
538         ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr);
539       }
540 #if defined(PETSC_USE_DEBUG)
541       if (!mglevels[i]->restrct || !mglevels[i]->interpolate) {
542         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
543       }
544 #endif
545     }
546     for (i=0; i<n-1; i++) {
547       if (!mglevels[i]->b) {
548         Vec *vec;
549         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
550         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
551         ierr = VecDestroy(*vec);CHKERRQ(ierr);
552         ierr = PetscFree(vec);CHKERRQ(ierr);
553       }
554       if (!mglevels[i]->r && i) {
555         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
556         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
557         ierr = VecDestroy(tvec);CHKERRQ(ierr);
558       }
559       if (!mglevels[i]->x) {
560         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
561         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
562         ierr = VecDestroy(tvec);CHKERRQ(ierr);
563       }
564     }
565     if (n != 1 && !mglevels[n-1]->r) {
566       /* PCMGSetR() on the finest level if user did not supply it */
567       Vec *vec;
568       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
569       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
570       ierr = VecDestroy(*vec);CHKERRQ(ierr);
571       ierr = PetscFree(vec);CHKERRQ(ierr);
572     }
573   }
574 
575 
576   for (i=1; i<n; i++) {
577     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
578       /* if doing only down then initial guess is zero */
579       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
580     }
581     if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
582     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
583     if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
584   }
585   for (i=1; i<n; i++) {
586     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
587       Mat          downmat,downpmat;
588       MatStructure matflag;
589       PetscTruth   opsset;
590 
591       /* check if operators have been set for up, if not use down operators to set them */
592       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
593       if (!opsset) {
594         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
595         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
596       }
597 
598       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
599       if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
600       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
601       if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
602     }
603   }
604 
605   /*
606       If coarse solver is not direct method then DO NOT USE preonly
607   */
608   ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
609   if (preonly) {
610     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
611     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
612     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
613     if (!lu && !redundant && !cholesky) {
614       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
615     }
616   }
617 
618   if (!pc->setupcalled) {
619     if (monitor) {
620       ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr);
621       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
622       ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
623     }
624     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
625   }
626 
627   if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
628   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
629   if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
630 
631   /*
632      Dump the interpolation/restriction matrices plus the
633    Jacobian/stiffness on each level. This allows Matlab users to
634    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
635 
636    Only support one or the other at the same time.
637   */
638 #if defined(PETSC_USE_SOCKET_VIEWER)
639   ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
640   if (dump) {
641     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
642   }
643   dump = PETSC_FALSE;
644 #endif
645   ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
646   if (dump) {
647     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
648   }
649 
650   if (viewer) {
651     for (i=1; i<n; i++) {
652       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
653     }
654     for (i=0; i<n; i++) {
655       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
656       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
657     }
658   }
659   PetscFunctionReturn(0);
660 }
661 
662 /* -------------------------------------------------------------------------------------*/
663 
664 #undef __FUNCT__
665 #define __FUNCT__ "PCMGGetLevels"
666 /*@
667    PCMGGetLevels - Gets the number of levels to use with MG.
668 
669    Not Collective
670 
671    Input Parameter:
672 .  pc - the preconditioner context
673 
674    Output parameter:
675 .  levels - the number of levels
676 
677    Level: advanced
678 
679 .keywords: MG, get, levels, multigrid
680 
681 .seealso: PCMGSetLevels()
682 @*/
683 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
684 {
685   PC_MG *mg = (PC_MG*)pc->data;
686 
687   PetscFunctionBegin;
688   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
689   PetscValidIntPointer(levels,2);
690   *levels = mg->nlevels;
691   PetscFunctionReturn(0);
692 }
693 
694 #undef __FUNCT__
695 #define __FUNCT__ "PCMGSetType"
696 /*@
697    PCMGSetType - Determines the form of multigrid to use:
698    multiplicative, additive, full, or the Kaskade algorithm.
699 
700    Collective on PC
701 
702    Input Parameters:
703 +  pc - the preconditioner context
704 -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
705    PC_MG_FULL, PC_MG_KASKADE
706 
707    Options Database Key:
708 .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
709    additive, full, kaskade
710 
711    Level: advanced
712 
713 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
714 
715 .seealso: PCMGSetLevels()
716 @*/
717 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
718 {
719   PC_MG                   *mg = (PC_MG*)pc->data;
720 
721   PetscFunctionBegin;
722   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
723   mg->am = form;
724   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
725   else pc->ops->applyrichardson = 0;
726   PetscFunctionReturn(0);
727 }
728 
729 #undef __FUNCT__
730 #define __FUNCT__ "PCMGSetCycleType"
731 /*@
732    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
733    complicated cycling.
734 
735    Collective on PC
736 
737    Input Parameters:
738 +  pc - the multigrid context
739 -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
740 
741    Options Database Key:
742 $  -pc_mg_cycle_type v or w
743 
744    Level: advanced
745 
746 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
747 
748 .seealso: PCMGSetCycleTypeOnLevel()
749 @*/
750 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n)
751 {
752   PC_MG        *mg = (PC_MG*)pc->data;
753   PC_MG_Levels **mglevels = mg->levels;
754   PetscInt     i,levels;
755 
756   PetscFunctionBegin;
757   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
758   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
759   levels = mglevels[0]->levels;
760 
761   for (i=0; i<levels; i++) {
762     mglevels[i]->cycles  = n;
763   }
764   PetscFunctionReturn(0);
765 }
766 
767 #undef __FUNCT__
768 #define __FUNCT__ "PCMGMultiplicativeSetCycles"
769 /*@
770    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
771          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
772 
773    Collective on PC
774 
775    Input Parameters:
776 +  pc - the multigrid context
777 -  n - number of cycles (default is 1)
778 
779    Options Database Key:
780 $  -pc_mg_multiplicative_cycles n
781 
782    Level: advanced
783 
784    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
785 
786 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
787 
788 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
789 @*/
790 PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
791 {
792   PC_MG        *mg = (PC_MG*)pc->data;
793   PC_MG_Levels **mglevels = mg->levels;
794   PetscInt     i,levels;
795 
796   PetscFunctionBegin;
797   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
798   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
799   levels = mglevels[0]->levels;
800 
801   for (i=0; i<levels; i++) {
802     mg->cyclesperpcapply  = n;
803   }
804   PetscFunctionReturn(0);
805 }
806 
807 #undef __FUNCT__
808 #define __FUNCT__ "PCMGSetGalerkin"
809 /*@
810    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
811       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
812 
813    Collective on PC
814 
815    Input Parameters:
816 .  pc - the multigrid context
817 
818    Options Database Key:
819 $  -pc_mg_galerkin
820 
821    Level: intermediate
822 
823 .keywords: MG, set, Galerkin
824 
825 .seealso: PCMGGetGalerkin()
826 
827 @*/
828 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
829 {
830   PC_MG        *mg = (PC_MG*)pc->data;
831 
832   PetscFunctionBegin;
833   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
834   mg->galerkin = PETSC_TRUE;
835   PetscFunctionReturn(0);
836 }
837 
838 #undef __FUNCT__
839 #define __FUNCT__ "PCMGGetGalerkin"
840 /*@
841    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
842       A_i-1 = r_i * A_i * r_i^t
843 
844    Not Collective
845 
846    Input Parameter:
847 .  pc - the multigrid context
848 
849    Output Parameter:
850 .  gelerkin - PETSC_TRUE or PETSC_FALSE
851 
852    Options Database Key:
853 $  -pc_mg_galerkin
854 
855    Level: intermediate
856 
857 .keywords: MG, set, Galerkin
858 
859 .seealso: PCMGSetGalerkin()
860 
861 @*/
862 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
863 {
864   PC_MG        *mg = (PC_MG*)pc->data;
865 
866   PetscFunctionBegin;
867   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
868   *galerkin = mg->galerkin;
869   PetscFunctionReturn(0);
870 }
871 
872 #undef __FUNCT__
873 #define __FUNCT__ "PCMGSetNumberSmoothDown"
874 /*@
875    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
876    use on all levels. Use PCMGGetSmootherDown() to set different
877    pre-smoothing steps on different levels.
878 
879    Collective on PC
880 
881    Input Parameters:
882 +  mg - the multigrid context
883 -  n - the number of smoothing steps
884 
885    Options Database Key:
886 .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
887 
888    Level: advanced
889 
890 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
891 
892 .seealso: PCMGSetNumberSmoothUp()
893 @*/
894 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
895 {
896   PC_MG          *mg = (PC_MG*)pc->data;
897   PC_MG_Levels   **mglevels = mg->levels;
898   PetscErrorCode ierr;
899   PetscInt       i,levels;
900 
901   PetscFunctionBegin;
902   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
903   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
904   levels = mglevels[0]->levels;
905 
906   for (i=1; i<levels; i++) {
907     /* make sure smoother up and down are different */
908     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
909     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
910     mg->default_smoothd = n;
911   }
912   PetscFunctionReturn(0);
913 }
914 
915 #undef __FUNCT__
916 #define __FUNCT__ "PCMGSetNumberSmoothUp"
917 /*@
918    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
919    on all levels. Use PCMGGetSmootherUp() to set different numbers of
920    post-smoothing steps on different levels.
921 
922    Collective on PC
923 
924    Input Parameters:
925 +  mg - the multigrid context
926 -  n - the number of smoothing steps
927 
928    Options Database Key:
929 .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
930 
931    Level: advanced
932 
933    Note: this does not set a value on the coarsest grid, since we assume that
934     there is no separate smooth up on the coarsest grid.
935 
936 .keywords: MG, smooth, up, post-smoothing, steps, multigrid
937 
938 .seealso: PCMGSetNumberSmoothDown()
939 @*/
940 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
941 {
942   PC_MG          *mg = (PC_MG*)pc->data;
943   PC_MG_Levels   **mglevels = mg->levels;
944   PetscErrorCode ierr;
945   PetscInt       i,levels;
946 
947   PetscFunctionBegin;
948   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
949   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
950   levels = mglevels[0]->levels;
951 
952   for (i=1; i<levels; i++) {
953     /* make sure smoother up and down are different */
954     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
955     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
956     mg->default_smoothu = n;
957   }
958   PetscFunctionReturn(0);
959 }
960 
961 /* ----------------------------------------------------------------------------------------*/
962 
963 /*MC
964    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
965     information about the coarser grid matrices and restriction/interpolation operators.
966 
967    Options Database Keys:
968 +  -pc_mg_levels <nlevels> - number of levels including finest
969 .  -pc_mg_cycles v or w
970 .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
971 .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
972 .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
973 .  -pc_mg_log - log information about time spent on each level of the solver
974 .  -pc_mg_monitor - print information on the multigrid convergence
975 .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
976 -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
977                         to the Socket viewer for reading from Matlab.
978 
979    Notes:
980 
981    Level: intermediate
982 
983    Concepts: multigrid/multilevel
984 
985 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
986            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
987            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
988            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
989            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
990 M*/
991 
992 EXTERN_C_BEGIN
993 #undef __FUNCT__
994 #define __FUNCT__ "PCCreate_MG"
995 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
996 {
997   PC_MG          *mg;
998   PetscErrorCode ierr;
999 
1000   PetscFunctionBegin;
1001   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1002   pc->data    = (void*)mg;
1003   mg->nlevels = -1;
1004 
1005   pc->ops->apply          = PCApply_MG;
1006   pc->ops->setup          = PCSetUp_MG;
1007   pc->ops->destroy        = PCDestroy_MG;
1008   pc->ops->setfromoptions = PCSetFromOptions_MG;
1009   pc->ops->view           = PCView_MG;
1010   PetscFunctionReturn(0);
1011 }
1012 EXTERN_C_END
1013