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