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