1f3b08a26SMatthew G. Knepley #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/
2f3b08a26SMatthew G. Knepley #include <petscdm.h>
3f3b08a26SMatthew G. Knepley
xfunc(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)4*2a8381b2SBarry Smith static PetscErrorCode xfunc(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
5d71ae5a4SJacob Faibussowitsch {
6557cf195SMatthew G. Knepley PetscInt k = *((PetscInt *)ctx), c;
7557cf195SMatthew G. Knepley
8557cf195SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = PetscPowRealInt(coords[0], k);
93ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
10557cf195SMatthew G. Knepley }
yfunc(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)11*2a8381b2SBarry Smith static PetscErrorCode yfunc(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
12d71ae5a4SJacob Faibussowitsch {
13557cf195SMatthew G. Knepley PetscInt k = *((PetscInt *)ctx), c;
14557cf195SMatthew G. Knepley
15557cf195SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = PetscPowRealInt(coords[1], k);
163ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
17557cf195SMatthew G. Knepley }
zfunc(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)18*2a8381b2SBarry Smith static PetscErrorCode zfunc(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
19d71ae5a4SJacob Faibussowitsch {
20557cf195SMatthew G. Knepley PetscInt k = *((PetscInt *)ctx), c;
21557cf195SMatthew G. Knepley
22557cf195SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = PetscPowRealInt(coords[2], k);
233ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
24557cf195SMatthew G. Knepley }
xsin(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)25*2a8381b2SBarry Smith static PetscErrorCode xsin(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
26d71ae5a4SJacob Faibussowitsch {
27557cf195SMatthew G. Knepley PetscInt k = *((PetscInt *)ctx), c;
28557cf195SMatthew G. Knepley
29557cf195SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = PetscSinReal(PETSC_PI * (k + 1) * coords[0]);
303ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
31557cf195SMatthew G. Knepley }
ysin(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)32*2a8381b2SBarry Smith static PetscErrorCode ysin(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
33d71ae5a4SJacob Faibussowitsch {
34557cf195SMatthew G. Knepley PetscInt k = *((PetscInt *)ctx), c;
35557cf195SMatthew G. Knepley
36557cf195SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = PetscSinReal(PETSC_PI * (k + 1) * coords[1]);
373ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
38557cf195SMatthew G. Knepley }
zsin(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)39*2a8381b2SBarry Smith static PetscErrorCode zsin(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
40d71ae5a4SJacob Faibussowitsch {
41557cf195SMatthew G. Knepley PetscInt k = *((PetscInt *)ctx), c;
42557cf195SMatthew G. Knepley
43557cf195SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = PetscSinReal(PETSC_PI * (k + 1) * coords[2]);
443ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
45557cf195SMatthew G. Knepley }
46557cf195SMatthew G. Knepley
DMSetBasisFunction_Internal(PetscInt Nf,PetscBool usePoly,PetscInt dir,PetscErrorCode (** funcs)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *))47d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetBasisFunction_Internal(PetscInt Nf, PetscBool usePoly, PetscInt dir, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *))
48d71ae5a4SJacob Faibussowitsch {
49557cf195SMatthew G. Knepley PetscInt f;
50557cf195SMatthew G. Knepley
51557cf195SMatthew G. Knepley PetscFunctionBeginUser;
52557cf195SMatthew G. Knepley for (f = 0; f < Nf; ++f) {
53557cf195SMatthew G. Knepley if (usePoly) {
54557cf195SMatthew G. Knepley switch (dir) {
55d71ae5a4SJacob Faibussowitsch case 0:
56d71ae5a4SJacob Faibussowitsch funcs[f] = xfunc;
57d71ae5a4SJacob Faibussowitsch break;
58d71ae5a4SJacob Faibussowitsch case 1:
59d71ae5a4SJacob Faibussowitsch funcs[f] = yfunc;
60d71ae5a4SJacob Faibussowitsch break;
61d71ae5a4SJacob Faibussowitsch case 2:
62d71ae5a4SJacob Faibussowitsch funcs[f] = zfunc;
63d71ae5a4SJacob Faibussowitsch break;
64d71ae5a4SJacob Faibussowitsch default:
65d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No function for direction %" PetscInt_FMT, dir);
66557cf195SMatthew G. Knepley }
67557cf195SMatthew G. Knepley } else {
68557cf195SMatthew G. Knepley switch (dir) {
69d71ae5a4SJacob Faibussowitsch case 0:
70d71ae5a4SJacob Faibussowitsch funcs[f] = xsin;
71d71ae5a4SJacob Faibussowitsch break;
72d71ae5a4SJacob Faibussowitsch case 1:
73d71ae5a4SJacob Faibussowitsch funcs[f] = ysin;
74d71ae5a4SJacob Faibussowitsch break;
75d71ae5a4SJacob Faibussowitsch case 2:
76d71ae5a4SJacob Faibussowitsch funcs[f] = zsin;
77d71ae5a4SJacob Faibussowitsch break;
78d71ae5a4SJacob Faibussowitsch default:
79d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No function for direction %" PetscInt_FMT, dir);
80557cf195SMatthew G. Knepley }
81557cf195SMatthew G. Knepley }
82557cf195SMatthew G. Knepley }
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
84557cf195SMatthew G. Knepley }
85557cf195SMatthew G. Knepley
PCMGCreateCoarseSpaceDefault_Private(PC pc,PetscInt level,PCMGCoarseSpaceType cstype,DM dm,KSP ksp,PetscInt Nc,Mat initialGuess,Mat * coarseSpace)86d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMGCreateCoarseSpaceDefault_Private(PC pc, PetscInt level, PCMGCoarseSpaceType cstype, DM dm, KSP ksp, PetscInt Nc, Mat initialGuess, Mat *coarseSpace)
87d71ae5a4SJacob Faibussowitsch {
882b3cbbdaSStefano Zampini PetscBool poly = cstype == PCMG_ADAPT_POLYNOMIAL ? PETSC_TRUE : PETSC_FALSE;
89f3b08a26SMatthew G. Knepley PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
90f3b08a26SMatthew G. Knepley void **ctxs;
912b3cbbdaSStefano Zampini PetscInt dim, d, Nf, f, k, m, M;
922b3cbbdaSStefano Zampini Vec tmp;
93f3b08a26SMatthew G. Knepley
94f3b08a26SMatthew G. Knepley PetscFunctionBegin;
952b3cbbdaSStefano Zampini Nc = Nc < 0 ? 6 : Nc;
969566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dim));
979566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf));
9863a3b9bcSJacob Faibussowitsch PetscCheck(Nc % dim == 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "The number of coarse vectors %" PetscInt_FMT " must be divisible by the dimension %" PetscInt_FMT, Nc, dim);
999566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nf, &funcs, Nf, &ctxs));
1002b3cbbdaSStefano Zampini PetscCall(DMGetGlobalVector(dm, &tmp));
1012b3cbbdaSStefano Zampini PetscCall(VecGetSize(tmp, &M));
1022b3cbbdaSStefano Zampini PetscCall(VecGetLocalSize(tmp, &m));
1032b3cbbdaSStefano Zampini PetscCall(MatCreateDense(PetscObjectComm((PetscObject)pc), m, PETSC_DECIDE, M, Nc, NULL, coarseSpace));
1042b3cbbdaSStefano Zampini PetscCall(DMRestoreGlobalVector(dm, &tmp));
105f3b08a26SMatthew G. Knepley for (k = 0; k < Nc / dim; ++k) {
1062b3cbbdaSStefano Zampini for (f = 0; f < Nf; ++f) ctxs[f] = &k;
107f3b08a26SMatthew G. Knepley for (d = 0; d < dim; ++d) {
1082b3cbbdaSStefano Zampini PetscCall(MatDenseGetColumnVecWrite(*coarseSpace, k * dim + d, &tmp));
1099566063dSJacob Faibussowitsch PetscCall(DMSetBasisFunction_Internal(Nf, poly, d, funcs));
1102b3cbbdaSStefano Zampini PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, tmp));
1112b3cbbdaSStefano Zampini PetscCall(MatDenseRestoreColumnVecWrite(*coarseSpace, k * dim + d, &tmp));
112f3b08a26SMatthew G. Knepley }
113f3b08a26SMatthew G. Knepley }
1149566063dSJacob Faibussowitsch PetscCall(PetscFree2(funcs, ctxs));
1153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
116f3b08a26SMatthew G. Knepley }
117f3b08a26SMatthew G. Knepley
PCMGCreateCoarseSpace_Polynomial(PC pc,PetscInt level,DM dm,KSP ksp,PetscInt Nc,Mat initialGuess,Mat * coarseSpace)118d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMGCreateCoarseSpace_Polynomial(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, Mat initialGuess, Mat *coarseSpace)
119d71ae5a4SJacob Faibussowitsch {
120f3b08a26SMatthew G. Knepley PetscFunctionBegin;
1212b3cbbdaSStefano Zampini PetscCall(PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_ADAPT_POLYNOMIAL, dm, ksp, Nc, initialGuess, coarseSpace));
1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
123f3b08a26SMatthew G. Knepley }
124f3b08a26SMatthew G. Knepley
PCMGCreateCoarseSpace_Harmonic(PC pc,PetscInt level,DM dm,KSP ksp,PetscInt Nc,Mat initialGuess,Mat * coarseSpace)12566976f2fSJacob Faibussowitsch static PetscErrorCode PCMGCreateCoarseSpace_Harmonic(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, Mat initialGuess, Mat *coarseSpace)
126d71ae5a4SJacob Faibussowitsch {
127f3b08a26SMatthew G. Knepley PetscFunctionBegin;
1282b3cbbdaSStefano Zampini PetscCall(PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_ADAPT_HARMONIC, dm, ksp, Nc, initialGuess, coarseSpace));
1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
130f3b08a26SMatthew G. Knepley }
131f3b08a26SMatthew G. Knepley
132f3b08a26SMatthew G. Knepley /*
133f3b08a26SMatthew G. Knepley PCMGComputeCoarseSpace_Internal - Compute vectors on level l that must be accurately interpolated.
134f3b08a26SMatthew G. Knepley
135f3b08a26SMatthew G. Knepley Input Parameters:
1362fe279fdSBarry Smith + pc - The `PCMG`
137f3b08a26SMatthew G. Knepley . l - The level
1382b3cbbdaSStefano Zampini . Nc - The number of vectors requested
1392fe279fdSBarry Smith - cspace - The initial guess for the space, or `NULL`
140f3b08a26SMatthew G. Knepley
141f3b08a26SMatthew G. Knepley Output Parameter:
142f3b08a26SMatthew G. Knepley . space - The space which must be accurately interpolated.
143f3b08a26SMatthew G. Knepley
144f3b08a26SMatthew G. Knepley Level: developer
145f3b08a26SMatthew G. Knepley
1462b3cbbdaSStefano Zampini Note: This space is normally used to adapt the interpolator. If Nc is negative, an adaptive choice can be made.
147f3b08a26SMatthew G. Knepley
148562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGAdaptInterpolator_Private()`
149f3b08a26SMatthew G. Knepley */
PCMGComputeCoarseSpace_Internal(PC pc,PetscInt l,PCMGCoarseSpaceType cstype,PetscInt Nc,Mat cspace,Mat * space)150d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGComputeCoarseSpace_Internal(PC pc, PetscInt l, PCMGCoarseSpaceType cstype, PetscInt Nc, Mat cspace, Mat *space)
151d71ae5a4SJacob Faibussowitsch {
1522b3cbbdaSStefano Zampini PetscErrorCode (*coarseConstructor)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *) = NULL;
153f3b08a26SMatthew G. Knepley DM dm;
154f3b08a26SMatthew G. Knepley KSP smooth;
155f3b08a26SMatthew G. Knepley
156f3b08a26SMatthew G. Knepley PetscFunctionBegin;
1572b3cbbdaSStefano Zampini *space = NULL;
158f3b08a26SMatthew G. Knepley switch (cstype) {
159d71ae5a4SJacob Faibussowitsch case PCMG_ADAPT_POLYNOMIAL:
160d71ae5a4SJacob Faibussowitsch coarseConstructor = &PCMGCreateCoarseSpace_Polynomial;
161d71ae5a4SJacob Faibussowitsch break;
162d71ae5a4SJacob Faibussowitsch case PCMG_ADAPT_HARMONIC:
163d71ae5a4SJacob Faibussowitsch coarseConstructor = &PCMGCreateCoarseSpace_Harmonic;
164d71ae5a4SJacob Faibussowitsch break;
1652b3cbbdaSStefano Zampini case PCMG_ADAPT_EIGENVECTOR:
1662b3cbbdaSStefano Zampini Nc = Nc < 0 ? 6 : Nc;
1679566063dSJacob Faibussowitsch if (l > 0) PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_MEV", &coarseConstructor));
1689566063dSJacob Faibussowitsch else PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_EV", &coarseConstructor));
169f3b08a26SMatthew G. Knepley break;
1702b3cbbdaSStefano Zampini case PCMG_ADAPT_GENERALIZED_EIGENVECTOR:
1712b3cbbdaSStefano Zampini Nc = Nc < 0 ? 6 : Nc;
1729566063dSJacob Faibussowitsch if (l > 0) PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_MGEV", &coarseConstructor));
1739566063dSJacob Faibussowitsch else PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_GEV", &coarseConstructor));
174f3b08a26SMatthew G. Knepley break;
175d71ae5a4SJacob Faibussowitsch case PCMG_ADAPT_GDSW:
176d71ae5a4SJacob Faibussowitsch coarseConstructor = &PCMGGDSWCreateCoarseSpace_Private;
177d71ae5a4SJacob Faibussowitsch break;
178d71ae5a4SJacob Faibussowitsch case PCMG_ADAPT_NONE:
179d71ae5a4SJacob Faibussowitsch break;
180d71ae5a4SJacob Faibussowitsch default:
181d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot handle coarse space type %d", cstype);
182f3b08a26SMatthew G. Knepley }
1832b3cbbdaSStefano Zampini if (coarseConstructor) {
1849566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, l, &smooth));
1859566063dSJacob Faibussowitsch PetscCall(KSPGetDM(smooth, &dm));
1869566063dSJacob Faibussowitsch PetscCall((*coarseConstructor)(pc, l, dm, smooth, Nc, cspace, space));
1872b3cbbdaSStefano Zampini }
1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
189f3b08a26SMatthew G. Knepley }
190f3b08a26SMatthew G. Knepley
191f3b08a26SMatthew G. Knepley /*
1922b3cbbdaSStefano Zampini PCMGAdaptInterpolator_Internal - Adapt interpolator from level l-1 to level l
193f3b08a26SMatthew G. Knepley
194f3b08a26SMatthew G. Knepley Input Parameters:
195f3b08a26SMatthew G. Knepley + pc - The PCMG
196f3b08a26SMatthew G. Knepley . l - The level l
197f3b08a26SMatthew G. Knepley . csmooth - The (coarse) smoother for level l-1
198f3b08a26SMatthew G. Knepley . fsmooth - The (fine) smoother for level l
199f3b08a26SMatthew G. Knepley . cspace - The (coarse) vectors in the subspace for level l-1
200f3b08a26SMatthew G. Knepley - fspace - The (fine) vectors in the subspace for level l
201f3b08a26SMatthew G. Knepley
202f3b08a26SMatthew G. Knepley Level: developer
203f3b08a26SMatthew G. Knepley
204f3b08a26SMatthew G. Knepley Note: This routine resets the interpolation and restriction for level l.
205f3b08a26SMatthew G. Knepley
206562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGComputeCoarseSpace_Private()`
207f3b08a26SMatthew G. Knepley */
PCMGAdaptInterpolator_Internal(PC pc,PetscInt l,KSP csmooth,KSP fsmooth,Mat cspace,Mat fspace)208d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGAdaptInterpolator_Internal(PC pc, PetscInt l, KSP csmooth, KSP fsmooth, Mat cspace, Mat fspace)
209d71ae5a4SJacob Faibussowitsch {
210f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *)pc->data;
211f3b08a26SMatthew G. Knepley DM dm, cdm;
212f3b08a26SMatthew G. Knepley Mat Interp, InterpAdapt;
213f3b08a26SMatthew G. Knepley
214f3b08a26SMatthew G. Knepley PetscFunctionBegin;
215f3b08a26SMatthew G. Knepley /* There is no interpolator for the coarse level */
2163ba16761SJacob Faibussowitsch if (!l) PetscFunctionReturn(PETSC_SUCCESS);
2179566063dSJacob Faibussowitsch PetscCall(KSPGetDM(csmooth, &cdm));
2189566063dSJacob Faibussowitsch PetscCall(KSPGetDM(fsmooth, &dm));
2199566063dSJacob Faibussowitsch PetscCall(PCMGGetInterpolation(pc, l, &Interp));
2203ba16761SJacob Faibussowitsch if (Interp == fspace && !cspace) PetscFunctionReturn(PETSC_SUCCESS);
2212b3cbbdaSStefano Zampini PetscCall(DMAdaptInterpolator(cdm, dm, Interp, fsmooth, fspace, cspace, &InterpAdapt, pc));
2222b3cbbdaSStefano Zampini if (mg->mespMonitor) PetscCall(DMCheckInterpolator(dm, InterpAdapt, cspace, fspace, 0.5 /* PETSC_SMALL */));
2239566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, l, InterpAdapt));
2242b3cbbdaSStefano Zampini PetscCall(PCMGSetRestriction(pc, l, InterpAdapt)); /* MATT: Remove????? */
2259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&InterpAdapt));
2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
227f3b08a26SMatthew G. Knepley }
228f3b08a26SMatthew G. Knepley
229f3b08a26SMatthew G. Knepley /*
230f3b08a26SMatthew G. Knepley PCMGRecomputeLevelOperators_Internal - Recomputes Galerkin coarse operator when interpolation is adapted
231f3b08a26SMatthew G. Knepley
232f3b08a26SMatthew G. Knepley Note: This routine recomputes the Galerkin triple product for the operator on level l.
233f3b08a26SMatthew G. Knepley */
PCMGRecomputeLevelOperators_Internal(PC pc,PetscInt l)234d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGRecomputeLevelOperators_Internal(PC pc, PetscInt l)
235d71ae5a4SJacob Faibussowitsch {
236f3b08a26SMatthew G. Knepley Mat fA, fB; /* The system and preconditioning operators on level l+1 */
237f3b08a26SMatthew G. Knepley Mat A, B; /* The system and preconditioning operators on level l */
238f3b08a26SMatthew G. Knepley Mat Interp, Restrc; /* The interpolation operator from level l to l+1, and restriction operator from level l+1 to l */
239f3b08a26SMatthew G. Knepley KSP smooth, fsmooth; /* The smoothers on levels l and l+1 */
240f3b08a26SMatthew G. Knepley PCMGGalerkinType galerkin; /* The Galerkin projection flag */
241f3b08a26SMatthew G. Knepley MatReuse reuse = MAT_REUSE_MATRIX; /* The matrices are always assumed to be present already */
242f3b08a26SMatthew G. Knepley PetscBool doA = PETSC_FALSE; /* Updates the system operator */
243f3b08a26SMatthew G. Knepley PetscBool doB = PETSC_FALSE; /* Updates the preconditioning operator (A == B, then update B) */
244f3b08a26SMatthew G. Knepley PetscInt n; /* The number of multigrid levels */
245f3b08a26SMatthew G. Knepley
246f3b08a26SMatthew G. Knepley PetscFunctionBegin;
2479566063dSJacob Faibussowitsch PetscCall(PCMGGetGalerkin(pc, &galerkin));
2483ba16761SJacob Faibussowitsch if (galerkin >= PC_MG_GALERKIN_NONE) PetscFunctionReturn(PETSC_SUCCESS);
2499566063dSJacob Faibussowitsch PetscCall(PCMGGetLevels(pc, &n));
250f3b08a26SMatthew G. Knepley /* Do not recompute operator for the finest grid */
2513ba16761SJacob Faibussowitsch if (l == n - 1) PetscFunctionReturn(PETSC_SUCCESS);
2529566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, l, &smooth));
2539566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(smooth, &A, &B));
2549566063dSJacob Faibussowitsch PetscCall(PCMGGetSmoother(pc, l + 1, &fsmooth));
2559566063dSJacob Faibussowitsch PetscCall(KSPGetOperators(fsmooth, &fA, &fB));
2569566063dSJacob Faibussowitsch PetscCall(PCMGGetInterpolation(pc, l + 1, &Interp));
2579566063dSJacob Faibussowitsch PetscCall(PCMGGetRestriction(pc, l + 1, &Restrc));
258f3b08a26SMatthew G. Knepley if ((galerkin == PC_MG_GALERKIN_PMAT) || (galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
259f3b08a26SMatthew G. Knepley if ((galerkin == PC_MG_GALERKIN_MAT) || ((galerkin == PC_MG_GALERKIN_BOTH) && (fA != fB))) doA = PETSC_TRUE;
2609566063dSJacob Faibussowitsch if (doA) PetscCall(MatGalerkin(Restrc, fA, Interp, reuse, 1.0, &A));
2619566063dSJacob Faibussowitsch if (doB) PetscCall(MatGalerkin(Restrc, fB, Interp, reuse, 1.0, &B));
2623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
263f3b08a26SMatthew G. Knepley }
264