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