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