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