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