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