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