xref: /petsc/src/ksp/pc/impls/mg/mgadapt.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
1f3b08a26SMatthew G. Knepley #include <petsc/private/pcmgimpl.h>       /*I "petscksp.h" I*/
2f3b08a26SMatthew G. Knepley #include <petscdm.h>
3f3b08a26SMatthew G. Knepley 
4557cf195SMatthew G. Knepley static PetscErrorCode xfunc(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
5557cf195SMatthew G. Knepley {
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);
9557cf195SMatthew G. Knepley   return 0;
10557cf195SMatthew G. Knepley }
11557cf195SMatthew G. Knepley static PetscErrorCode yfunc(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
12557cf195SMatthew G. Knepley {
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);
16557cf195SMatthew G. Knepley   return 0;
17557cf195SMatthew G. Knepley }
18557cf195SMatthew G. Knepley static PetscErrorCode zfunc(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
19557cf195SMatthew G. Knepley {
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);
23557cf195SMatthew G. Knepley   return 0;
24557cf195SMatthew G. Knepley }
25557cf195SMatthew G. Knepley static PetscErrorCode xsin(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
26557cf195SMatthew G. Knepley {
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]);
30557cf195SMatthew G. Knepley   return 0;
31557cf195SMatthew G. Knepley }
32557cf195SMatthew G. Knepley static PetscErrorCode ysin(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
33557cf195SMatthew G. Knepley {
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]);
37557cf195SMatthew G. Knepley   return 0;
38557cf195SMatthew G. Knepley }
39557cf195SMatthew G. Knepley static PetscErrorCode zsin(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
40557cf195SMatthew G. Knepley {
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]);
44557cf195SMatthew G. Knepley   return 0;
45557cf195SMatthew G. Knepley }
46557cf195SMatthew G. Knepley 
47557cf195SMatthew G. Knepley PetscErrorCode DMSetBasisFunction_Internal(PetscInt Nf, PetscBool usePoly, PetscInt dir, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *))
48557cf195SMatthew G. Knepley {
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) {
55557cf195SMatthew G. Knepley       case 0: funcs[f] = xfunc;break;
56557cf195SMatthew G. Knepley       case 1: funcs[f] = yfunc;break;
57557cf195SMatthew G. Knepley       case 2: funcs[f] = zfunc;break;
58*63a3b9bcSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No function for direction %" PetscInt_FMT, dir);
59557cf195SMatthew G. Knepley       }
60557cf195SMatthew G. Knepley     } else {
61557cf195SMatthew G. Knepley       switch (dir) {
62557cf195SMatthew G. Knepley       case 0: funcs[f] = xsin;break;
63557cf195SMatthew G. Knepley       case 1: funcs[f] = ysin;break;
64557cf195SMatthew G. Knepley       case 2: funcs[f] = zsin;break;
65*63a3b9bcSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No function for direction %" PetscInt_FMT, dir);
66557cf195SMatthew G. Knepley       }
67557cf195SMatthew G. Knepley     }
68557cf195SMatthew G. Knepley   }
69557cf195SMatthew G. Knepley   PetscFunctionReturn(0);
70557cf195SMatthew G. Knepley }
71557cf195SMatthew G. Knepley 
72f3b08a26SMatthew G. Knepley static PetscErrorCode PCMGCreateCoarseSpaceDefault_Private(PC pc, PetscInt level, PCMGCoarseSpaceType cstype, DM dm, KSP ksp, PetscInt Nc, const Vec initialGuess[], Vec **coarseSpace)
73f3b08a26SMatthew G. Knepley {
74f3b08a26SMatthew G. Knepley   PetscBool         poly = cstype == PCMG_POLYNOMIAL ? PETSC_TRUE : PETSC_FALSE;
75f3b08a26SMatthew G. Knepley   PetscErrorCode (**funcs)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar*,void*);
76f3b08a26SMatthew G. Knepley   void            **ctxs;
77f3b08a26SMatthew G. Knepley   PetscInt          dim, d, Nf, f, k;
78f3b08a26SMatthew G. Knepley 
79f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dim));
819566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
82*63a3b9bcSJacob 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);
839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &funcs, Nf, &ctxs));
849566063dSJacob Faibussowitsch   if (!*coarseSpace) PetscCall(PetscCalloc1(Nc, coarseSpace));
85f3b08a26SMatthew G. Knepley   for (k = 0; k < Nc/dim; ++k) {
86f3b08a26SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ctxs[f] = &k;}
87f3b08a26SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
889566063dSJacob Faibussowitsch       if (!(*coarseSpace)[k*dim+d]) PetscCall(DMCreateGlobalVector(dm, &(*coarseSpace)[k*dim+d]));
899566063dSJacob Faibussowitsch       PetscCall(DMSetBasisFunction_Internal(Nf, poly, d, funcs));
909566063dSJacob Faibussowitsch       PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, (*coarseSpace)[k*dim+d]));
91f3b08a26SMatthew G. Knepley     }
92f3b08a26SMatthew G. Knepley   }
939566063dSJacob Faibussowitsch   PetscCall(PetscFree2(funcs, ctxs));
94f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
95f3b08a26SMatthew G. Knepley }
96f3b08a26SMatthew G. Knepley 
97f3b08a26SMatthew G. Knepley static PetscErrorCode PCMGCreateCoarseSpace_Polynomial(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, const Vec initialGuess[], Vec **coarseSpace)
98f3b08a26SMatthew G. Knepley {
99f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1009566063dSJacob Faibussowitsch   PetscCall(PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_POLYNOMIAL, dm, ksp, Nc, initialGuess, coarseSpace));
101f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
102f3b08a26SMatthew G. Knepley }
103f3b08a26SMatthew G. Knepley 
104f3b08a26SMatthew G. Knepley PetscErrorCode PCMGCreateCoarseSpace_Harmonic(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, const Vec initialGuess[], Vec **coarseSpace)
105f3b08a26SMatthew G. Knepley {
106f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1079566063dSJacob Faibussowitsch   PetscCall(PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_HARMONIC, dm, ksp, Nc, initialGuess, coarseSpace));
108f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
109f3b08a26SMatthew G. Knepley }
110f3b08a26SMatthew G. Knepley 
111f3b08a26SMatthew G. Knepley /*
112f3b08a26SMatthew G. Knepley   PCMGComputeCoarseSpace_Internal - Compute vectors on level l that must be accurately interpolated.
113f3b08a26SMatthew G. Knepley 
114f3b08a26SMatthew G. Knepley   Input Parameters:
115f3b08a26SMatthew G. Knepley + pc     - The PCMG
116f3b08a26SMatthew G. Knepley . l      - The level
117f3b08a26SMatthew G. Knepley . Nc     - The size of the space (number of vectors)
118f3b08a26SMatthew G. Knepley - cspace - The space from level l-1, or NULL
119f3b08a26SMatthew G. Knepley 
120f3b08a26SMatthew G. Knepley   Output Parameter:
121f3b08a26SMatthew G. Knepley . space  - The space which must be accurately interpolated.
122f3b08a26SMatthew G. Knepley 
123f3b08a26SMatthew G. Knepley   Level: developer
124f3b08a26SMatthew G. Knepley 
125f3b08a26SMatthew G. Knepley   Note: This space is normally used to adapt the interpolator.
126f3b08a26SMatthew G. Knepley 
127f3b08a26SMatthew G. Knepley .seealso: PCMGAdaptInterpolator_Private()
128f3b08a26SMatthew G. Knepley */
129f3b08a26SMatthew G. Knepley PetscErrorCode PCMGComputeCoarseSpace_Internal(PC pc, PetscInt l, PCMGCoarseSpaceType cstype, PetscInt Nc, const Vec cspace[], Vec *space[])
130f3b08a26SMatthew G. Knepley {
131f3b08a26SMatthew G. Knepley   PetscErrorCode (*coarseConstructor)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec*[]);
132f3b08a26SMatthew G. Knepley   DM             dm;
133f3b08a26SMatthew G. Knepley   KSP            smooth;
134f3b08a26SMatthew G. Knepley 
135f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
136f3b08a26SMatthew G. Knepley   switch (cstype) {
137f3b08a26SMatthew G. Knepley   case PCMG_POLYNOMIAL: coarseConstructor = &PCMGCreateCoarseSpace_Polynomial;break;
138f3b08a26SMatthew G. Knepley   case PCMG_HARMONIC:   coarseConstructor = &PCMGCreateCoarseSpace_Harmonic;break;
139f3b08a26SMatthew G. Knepley   case PCMG_EIGENVECTOR:
1409566063dSJacob Faibussowitsch     if (l > 0) PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_MEV", &coarseConstructor));
1419566063dSJacob Faibussowitsch     else       PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_EV", &coarseConstructor));
142f3b08a26SMatthew G. Knepley     break;
143f3b08a26SMatthew G. Knepley   case PCMG_GENERALIZED_EIGENVECTOR:
1449566063dSJacob Faibussowitsch     if (l > 0) PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_MGEV", &coarseConstructor));
1459566063dSJacob Faibussowitsch     else       PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_GEV", &coarseConstructor));
146f3b08a26SMatthew G. Knepley     break;
147*63a3b9bcSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle coarse space type %d", cstype);
148f3b08a26SMatthew G. Knepley   }
1499566063dSJacob Faibussowitsch   PetscCall(PCMGGetSmoother(pc, l, &smooth));
1509566063dSJacob Faibussowitsch   PetscCall(KSPGetDM(smooth, &dm));
1519566063dSJacob Faibussowitsch   PetscCall((*coarseConstructor)(pc, l, dm, smooth, Nc, cspace, space));
152f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
153f3b08a26SMatthew G. Knepley }
154f3b08a26SMatthew G. Knepley 
155f3b08a26SMatthew G. Knepley /*
156f3b08a26SMatthew G. Knepley   PCMGAdaptInterpolator_Internal - Adapt interpolator from level l-1 to level 1
157f3b08a26SMatthew G. Knepley 
158f3b08a26SMatthew G. Knepley   Input Parameters:
159f3b08a26SMatthew G. Knepley + pc      - The PCMG
160f3b08a26SMatthew G. Knepley . l       - The level l
161f3b08a26SMatthew G. Knepley . csmooth - The (coarse) smoother for level l-1
162f3b08a26SMatthew G. Knepley . fsmooth - The (fine) smoother for level l
163f3b08a26SMatthew G. Knepley . Nc      - The size of the subspace used for adaptation
164f3b08a26SMatthew G. Knepley . cspace  - The (coarse) vectors in the subspace for level l-1
165f3b08a26SMatthew G. Knepley - fspace  - The (fine) vectors in the subspace for level l
166f3b08a26SMatthew G. Knepley 
167f3b08a26SMatthew G. Knepley   Level: developer
168f3b08a26SMatthew G. Knepley 
169f3b08a26SMatthew G. Knepley   Note: This routine resets the interpolation and restriction for level l.
170f3b08a26SMatthew G. Knepley 
171f3b08a26SMatthew G. Knepley .seealso: PCMGComputeCoarseSpace_Private()
172f3b08a26SMatthew G. Knepley */
173f3b08a26SMatthew G. Knepley PetscErrorCode PCMGAdaptInterpolator_Internal(PC pc, PetscInt l, KSP csmooth, KSP fsmooth, PetscInt Nc, Vec cspace[], Vec fspace[])
174f3b08a26SMatthew G. Knepley {
175f3b08a26SMatthew G. Knepley   PC_MG         *mg = (PC_MG *) pc->data;
176f3b08a26SMatthew G. Knepley   DM             dm, cdm;
177f3b08a26SMatthew G. Knepley   Mat            Interp, InterpAdapt;
178f3b08a26SMatthew G. Knepley 
179f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
180f3b08a26SMatthew G. Knepley   /* There is no interpolator for the coarse level */
181f3b08a26SMatthew G. Knepley   if (!l) PetscFunctionReturn(0);
1829566063dSJacob Faibussowitsch   PetscCall(KSPGetDM(csmooth, &cdm));
1839566063dSJacob Faibussowitsch   PetscCall(KSPGetDM(fsmooth, &dm));
1849566063dSJacob Faibussowitsch   PetscCall(PCMGGetInterpolation(pc, l, &Interp));
185f3b08a26SMatthew G. Knepley 
1869566063dSJacob Faibussowitsch   PetscCall(DMAdaptInterpolator(cdm, dm, Interp, fsmooth, Nc, fspace, cspace, &InterpAdapt, pc));
1879566063dSJacob Faibussowitsch   if (mg->mespMonitor) PetscCall(DMCheckInterpolator(dm, InterpAdapt, Nc, cspace, fspace, 0.5/* PETSC_SMALL */));
1889566063dSJacob Faibussowitsch   PetscCall(PCMGSetInterpolation(pc, l, InterpAdapt));
1899566063dSJacob Faibussowitsch   PetscCall(PCMGSetRestriction(pc, l, InterpAdapt));
1909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&InterpAdapt));
191f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
192f3b08a26SMatthew G. Knepley }
193f3b08a26SMatthew G. Knepley 
194f3b08a26SMatthew G. Knepley /*
195f3b08a26SMatthew G. Knepley   PCMGRecomputeLevelOperators_Internal - Recomputes Galerkin coarse operator when interpolation is adapted
196f3b08a26SMatthew G. Knepley 
197f3b08a26SMatthew G. Knepley   Input Parameters:
198f3b08a26SMatthew G. Knepley + pc - The PCMG
199f3b08a26SMatthew G. Knepley - l  - The level l
200f3b08a26SMatthew G. Knepley 
201f3b08a26SMatthew G. Knepley   Level: developer
202f3b08a26SMatthew G. Knepley 
203f3b08a26SMatthew G. Knepley   Note: This routine recomputes the Galerkin triple product for the operator on level l.
204f3b08a26SMatthew G. Knepley */
205f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRecomputeLevelOperators_Internal(PC pc, PetscInt l)
206f3b08a26SMatthew G. Knepley {
207f3b08a26SMatthew G. Knepley   Mat              fA, fB;                   /* The system and preconditioning operators on level l+1 */
208f3b08a26SMatthew G. Knepley   Mat              A,  B;                    /* The system and preconditioning operators on level l */
209f3b08a26SMatthew G. Knepley   Mat              Interp, Restrc;           /* The interpolation operator from level l to l+1, and restriction operator from level l+1 to l */
210f3b08a26SMatthew G. Knepley   KSP              smooth, fsmooth;          /* The smoothers on levels l and l+1 */
211f3b08a26SMatthew G. Knepley   PCMGGalerkinType galerkin;                 /* The Galerkin projection flag */
212f3b08a26SMatthew G. Knepley   MatReuse         reuse = MAT_REUSE_MATRIX; /* The matrices are always assumed to be present already */
213f3b08a26SMatthew G. Knepley   PetscBool        doA   = PETSC_FALSE;      /* Updates the system operator */
214f3b08a26SMatthew G. Knepley   PetscBool        doB   = PETSC_FALSE;      /* Updates the preconditioning operator (A == B, then update B) */
215f3b08a26SMatthew G. Knepley   PetscInt         n;                        /* The number of multigrid levels */
216f3b08a26SMatthew G. Knepley 
217f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
2189566063dSJacob Faibussowitsch   PetscCall(PCMGGetGalerkin(pc, &galerkin));
219f3b08a26SMatthew G. Knepley   if (galerkin >= PC_MG_GALERKIN_NONE) PetscFunctionReturn(0);
2209566063dSJacob Faibussowitsch   PetscCall(PCMGGetLevels(pc, &n));
221f3b08a26SMatthew G. Knepley   /* Do not recompute operator for the finest grid */
222f3b08a26SMatthew G. Knepley   if (l == n-1) PetscFunctionReturn(0);
2239566063dSJacob Faibussowitsch   PetscCall(PCMGGetSmoother(pc, l,   &smooth));
2249566063dSJacob Faibussowitsch   PetscCall(KSPGetOperators(smooth, &A, &B));
2259566063dSJacob Faibussowitsch   PetscCall(PCMGGetSmoother(pc, l+1, &fsmooth));
2269566063dSJacob Faibussowitsch   PetscCall(KSPGetOperators(fsmooth, &fA, &fB));
2279566063dSJacob Faibussowitsch   PetscCall(PCMGGetInterpolation(pc, l+1, &Interp));
2289566063dSJacob Faibussowitsch   PetscCall(PCMGGetRestriction(pc, l+1, &Restrc));
229f3b08a26SMatthew G. Knepley   if ((galerkin == PC_MG_GALERKIN_PMAT) ||  (galerkin == PC_MG_GALERKIN_BOTH))                doB = PETSC_TRUE;
230f3b08a26SMatthew G. Knepley   if ((galerkin == PC_MG_GALERKIN_MAT)  || ((galerkin == PC_MG_GALERKIN_BOTH) && (fA != fB))) doA = PETSC_TRUE;
2319566063dSJacob Faibussowitsch   if (doA) PetscCall(MatGalerkin(Restrc, fA, Interp, reuse, 1.0, &A));
2329566063dSJacob Faibussowitsch   if (doB) PetscCall(MatGalerkin(Restrc, fB, Interp, reuse, 1.0, &B));
233f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
234f3b08a26SMatthew G. Knepley }
235