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