xref: /petsc/src/ksp/pc/impls/mg/mgadapt.c (revision 2b3cbbda7ef6bfb59aeed26278d0c91b4282b9fd)
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;
5863a3b9bcSJacob 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;
6563a3b9bcSJacob 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 
72*2b3cbbdaSStefano Zampini static PetscErrorCode PCMGCreateCoarseSpaceDefault_Private(PC pc, PetscInt level, PCMGCoarseSpaceType cstype, DM dm, KSP ksp, PetscInt Nc, Mat initialGuess, Mat *coarseSpace)
73f3b08a26SMatthew G. Knepley {
74*2b3cbbdaSStefano Zampini   PetscBool         poly = cstype == PCMG_ADAPT_POLYNOMIAL ? PETSC_TRUE : PETSC_FALSE;
75f3b08a26SMatthew G. Knepley   PetscErrorCode (**funcs)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar*,void*);
76f3b08a26SMatthew G. Knepley   void            **ctxs;
77*2b3cbbdaSStefano Zampini   PetscInt          dim, d, Nf, f, k, m, M;
78*2b3cbbdaSStefano Zampini   Vec               tmp;
79f3b08a26SMatthew G. Knepley 
80f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
81*2b3cbbdaSStefano Zampini   Nc = Nc < 0 ? 6 : Nc;
829566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dim));
839566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
8463a3b9bcSJacob 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);
859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &funcs, Nf, &ctxs));
86*2b3cbbdaSStefano Zampini   PetscCall(DMGetGlobalVector(dm, &tmp));
87*2b3cbbdaSStefano Zampini   PetscCall(VecGetSize(tmp, &M));
88*2b3cbbdaSStefano Zampini   PetscCall(VecGetLocalSize(tmp, &m));
89*2b3cbbdaSStefano Zampini   PetscCall(MatCreateDense(PetscObjectComm((PetscObject) pc), m, PETSC_DECIDE, M, Nc, NULL, coarseSpace));
90*2b3cbbdaSStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &tmp));
91f3b08a26SMatthew G. Knepley   for (k = 0; k < Nc/dim; ++k) {
92*2b3cbbdaSStefano Zampini     for (f = 0; f < Nf; ++f) ctxs[f] = &k;
93f3b08a26SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
94*2b3cbbdaSStefano Zampini       PetscCall(MatDenseGetColumnVecWrite(*coarseSpace,k*dim+d,&tmp));
959566063dSJacob Faibussowitsch       PetscCall(DMSetBasisFunction_Internal(Nf, poly, d, funcs));
96*2b3cbbdaSStefano Zampini       PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, tmp));
97*2b3cbbdaSStefano Zampini       PetscCall(MatDenseRestoreColumnVecWrite(*coarseSpace,k*dim+d,&tmp));
98f3b08a26SMatthew G. Knepley     }
99f3b08a26SMatthew G. Knepley   }
1009566063dSJacob Faibussowitsch   PetscCall(PetscFree2(funcs, ctxs));
101f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
102f3b08a26SMatthew G. Knepley }
103f3b08a26SMatthew G. Knepley 
104*2b3cbbdaSStefano Zampini static PetscErrorCode PCMGCreateCoarseSpace_Polynomial(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, Mat initialGuess, Mat *coarseSpace)
105f3b08a26SMatthew G. Knepley {
106f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
107*2b3cbbdaSStefano Zampini   PetscCall(PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_ADAPT_POLYNOMIAL, dm, ksp, Nc, initialGuess, coarseSpace));
108f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
109f3b08a26SMatthew G. Knepley }
110f3b08a26SMatthew G. Knepley 
111*2b3cbbdaSStefano Zampini PetscErrorCode PCMGCreateCoarseSpace_Harmonic(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, Mat initialGuess, Mat *coarseSpace)
112f3b08a26SMatthew G. Knepley {
113f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
114*2b3cbbdaSStefano Zampini   PetscCall(PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_ADAPT_HARMONIC, dm, ksp, Nc, initialGuess, coarseSpace));
115f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
116f3b08a26SMatthew G. Knepley }
117f3b08a26SMatthew G. Knepley 
118f3b08a26SMatthew G. Knepley /*
119f3b08a26SMatthew G. Knepley   PCMGComputeCoarseSpace_Internal - Compute vectors on level l that must be accurately interpolated.
120f3b08a26SMatthew G. Knepley 
121f3b08a26SMatthew G. Knepley   Input Parameters:
122f3b08a26SMatthew G. Knepley + pc     - The PCMG
123f3b08a26SMatthew G. Knepley . l      - The level
124*2b3cbbdaSStefano Zampini . Nc     - The number of vectors requested
125*2b3cbbdaSStefano Zampini - cspace - The initial guess for the space, or NULL
126f3b08a26SMatthew G. Knepley 
127f3b08a26SMatthew G. Knepley   Output Parameter:
128f3b08a26SMatthew G. Knepley . space  - The space which must be accurately interpolated.
129f3b08a26SMatthew G. Knepley 
130f3b08a26SMatthew G. Knepley   Level: developer
131f3b08a26SMatthew G. Knepley 
132*2b3cbbdaSStefano Zampini   Note: This space is normally used to adapt the interpolator. If Nc is negative, an adaptive choice can be made.
133f3b08a26SMatthew G. Knepley 
134db781477SPatrick Sanan .seealso: `PCMGAdaptInterpolator_Private()`
135f3b08a26SMatthew G. Knepley */
136*2b3cbbdaSStefano Zampini PetscErrorCode PCMGComputeCoarseSpace_Internal(PC pc, PetscInt l, PCMGCoarseSpaceType cstype, PetscInt Nc, Mat cspace, Mat *space)
137f3b08a26SMatthew G. Knepley {
138*2b3cbbdaSStefano Zampini   PetscErrorCode (*coarseConstructor)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*) = NULL;
139f3b08a26SMatthew G. Knepley   DM             dm;
140f3b08a26SMatthew G. Knepley   KSP            smooth;
141f3b08a26SMatthew G. Knepley 
142f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
143*2b3cbbdaSStefano Zampini   *space = NULL;
144f3b08a26SMatthew G. Knepley   switch (cstype) {
145*2b3cbbdaSStefano Zampini   case PCMG_ADAPT_POLYNOMIAL:
146*2b3cbbdaSStefano Zampini     coarseConstructor = &PCMGCreateCoarseSpace_Polynomial;
147*2b3cbbdaSStefano Zampini     break;
148*2b3cbbdaSStefano Zampini   case PCMG_ADAPT_HARMONIC:
149*2b3cbbdaSStefano Zampini     coarseConstructor = &PCMGCreateCoarseSpace_Harmonic;
150*2b3cbbdaSStefano Zampini     break;
151*2b3cbbdaSStefano Zampini   case PCMG_ADAPT_EIGENVECTOR:
152*2b3cbbdaSStefano Zampini     Nc = Nc < 0 ? 6 : Nc;
1539566063dSJacob Faibussowitsch     if (l > 0) PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_MEV", &coarseConstructor));
1549566063dSJacob Faibussowitsch     else       PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_EV", &coarseConstructor));
155f3b08a26SMatthew G. Knepley     break;
156*2b3cbbdaSStefano Zampini   case PCMG_ADAPT_GENERALIZED_EIGENVECTOR:
157*2b3cbbdaSStefano Zampini     Nc = Nc < 0 ? 6 : Nc;
1589566063dSJacob Faibussowitsch     if (l > 0) PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_MGEV", &coarseConstructor));
1599566063dSJacob Faibussowitsch     else       PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_GEV", &coarseConstructor));
160f3b08a26SMatthew G. Knepley     break;
161*2b3cbbdaSStefano Zampini   case PCMG_ADAPT_NONE:
162*2b3cbbdaSStefano Zampini     break;
163*2b3cbbdaSStefano Zampini   default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot handle coarse space type %d", cstype);
164f3b08a26SMatthew G. Knepley   }
165*2b3cbbdaSStefano Zampini   if (coarseConstructor) {
1669566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmoother(pc, l, &smooth));
1679566063dSJacob Faibussowitsch     PetscCall(KSPGetDM(smooth, &dm));
1689566063dSJacob Faibussowitsch     PetscCall((*coarseConstructor)(pc, l, dm, smooth, Nc, cspace, space));
169*2b3cbbdaSStefano Zampini   }
170f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
171f3b08a26SMatthew G. Knepley }
172f3b08a26SMatthew G. Knepley 
173f3b08a26SMatthew G. Knepley /*
174*2b3cbbdaSStefano Zampini   PCMGAdaptInterpolator_Internal - Adapt interpolator from level l-1 to level l
175f3b08a26SMatthew G. Knepley 
176f3b08a26SMatthew G. Knepley   Input Parameters:
177f3b08a26SMatthew G. Knepley + pc      - The PCMG
178f3b08a26SMatthew G. Knepley . l       - The level l
179f3b08a26SMatthew G. Knepley . csmooth - The (coarse) smoother for level l-1
180f3b08a26SMatthew G. Knepley . fsmooth - The (fine) smoother for level l
181f3b08a26SMatthew G. Knepley . cspace  - The (coarse) vectors in the subspace for level l-1
182f3b08a26SMatthew G. Knepley - fspace  - The (fine) vectors in the subspace for level l
183f3b08a26SMatthew G. Knepley 
184f3b08a26SMatthew G. Knepley   Level: developer
185f3b08a26SMatthew G. Knepley 
186f3b08a26SMatthew G. Knepley   Note: This routine resets the interpolation and restriction for level l.
187f3b08a26SMatthew G. Knepley 
188db781477SPatrick Sanan .seealso: `PCMGComputeCoarseSpace_Private()`
189f3b08a26SMatthew G. Knepley */
190*2b3cbbdaSStefano Zampini PetscErrorCode PCMGAdaptInterpolator_Internal(PC pc, PetscInt l, KSP csmooth, KSP fsmooth, Mat cspace, Mat fspace)
191f3b08a26SMatthew G. Knepley {
192f3b08a26SMatthew G. Knepley   PC_MG         *mg = (PC_MG *) pc->data;
193f3b08a26SMatthew G. Knepley   DM             dm, cdm;
194f3b08a26SMatthew G. Knepley   Mat            Interp, InterpAdapt;
195f3b08a26SMatthew G. Knepley 
196f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
197f3b08a26SMatthew G. Knepley   /* There is no interpolator for the coarse level */
198f3b08a26SMatthew G. Knepley   if (!l) PetscFunctionReturn(0);
1999566063dSJacob Faibussowitsch   PetscCall(KSPGetDM(csmooth, &cdm));
2009566063dSJacob Faibussowitsch   PetscCall(KSPGetDM(fsmooth, &dm));
2019566063dSJacob Faibussowitsch   PetscCall(PCMGGetInterpolation(pc, l, &Interp));
202*2b3cbbdaSStefano Zampini   if (Interp == fspace && !cspace) PetscFunctionReturn(0);
203*2b3cbbdaSStefano Zampini   PetscCall(DMAdaptInterpolator(cdm, dm, Interp, fsmooth, fspace, cspace, &InterpAdapt, pc));
204*2b3cbbdaSStefano Zampini   if (mg->mespMonitor) PetscCall(DMCheckInterpolator(dm, InterpAdapt, cspace, fspace, 0.5/* PETSC_SMALL */));
2059566063dSJacob Faibussowitsch   PetscCall(PCMGSetInterpolation(pc, l, InterpAdapt));
206*2b3cbbdaSStefano Zampini   PetscCall(PCMGSetRestriction(pc, l, InterpAdapt)); /* MATT: Remove????? */
2079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&InterpAdapt));
208f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
209f3b08a26SMatthew G. Knepley }
210f3b08a26SMatthew G. Knepley 
211f3b08a26SMatthew G. Knepley /*
212f3b08a26SMatthew G. Knepley   PCMGRecomputeLevelOperators_Internal - Recomputes Galerkin coarse operator when interpolation is adapted
213f3b08a26SMatthew G. Knepley 
214f3b08a26SMatthew G. Knepley   Input Parameters:
215f3b08a26SMatthew G. Knepley + pc - The PCMG
216f3b08a26SMatthew G. Knepley - l  - The level l
217f3b08a26SMatthew G. Knepley 
218f3b08a26SMatthew G. Knepley   Level: developer
219f3b08a26SMatthew G. Knepley 
220f3b08a26SMatthew G. Knepley   Note: This routine recomputes the Galerkin triple product for the operator on level l.
221f3b08a26SMatthew G. Knepley */
222f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRecomputeLevelOperators_Internal(PC pc, PetscInt l)
223f3b08a26SMatthew G. Knepley {
224f3b08a26SMatthew G. Knepley   Mat              fA, fB;                   /* The system and preconditioning operators on level l+1 */
225f3b08a26SMatthew G. Knepley   Mat              A,  B;                    /* The system and preconditioning operators on level l */
226f3b08a26SMatthew G. Knepley   Mat              Interp, Restrc;           /* The interpolation operator from level l to l+1, and restriction operator from level l+1 to l */
227f3b08a26SMatthew G. Knepley   KSP              smooth, fsmooth;          /* The smoothers on levels l and l+1 */
228f3b08a26SMatthew G. Knepley   PCMGGalerkinType galerkin;                 /* The Galerkin projection flag */
229f3b08a26SMatthew G. Knepley   MatReuse         reuse = MAT_REUSE_MATRIX; /* The matrices are always assumed to be present already */
230f3b08a26SMatthew G. Knepley   PetscBool        doA   = PETSC_FALSE;      /* Updates the system operator */
231f3b08a26SMatthew G. Knepley   PetscBool        doB   = PETSC_FALSE;      /* Updates the preconditioning operator (A == B, then update B) */
232f3b08a26SMatthew G. Knepley   PetscInt         n;                        /* The number of multigrid levels */
233f3b08a26SMatthew G. Knepley 
234f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
2359566063dSJacob Faibussowitsch   PetscCall(PCMGGetGalerkin(pc, &galerkin));
236f3b08a26SMatthew G. Knepley   if (galerkin >= PC_MG_GALERKIN_NONE) PetscFunctionReturn(0);
2379566063dSJacob Faibussowitsch   PetscCall(PCMGGetLevels(pc, &n));
238f3b08a26SMatthew G. Knepley   /* Do not recompute operator for the finest grid */
239f3b08a26SMatthew G. Knepley   if (l == n-1) PetscFunctionReturn(0);
2409566063dSJacob Faibussowitsch   PetscCall(PCMGGetSmoother(pc, l,   &smooth));
2419566063dSJacob Faibussowitsch   PetscCall(KSPGetOperators(smooth, &A, &B));
2429566063dSJacob Faibussowitsch   PetscCall(PCMGGetSmoother(pc, l+1, &fsmooth));
2439566063dSJacob Faibussowitsch   PetscCall(KSPGetOperators(fsmooth, &fA, &fB));
2449566063dSJacob Faibussowitsch   PetscCall(PCMGGetInterpolation(pc, l+1, &Interp));
2459566063dSJacob Faibussowitsch   PetscCall(PCMGGetRestriction(pc, l+1, &Restrc));
246f3b08a26SMatthew G. Knepley   if ((galerkin == PC_MG_GALERKIN_PMAT) ||  (galerkin == PC_MG_GALERKIN_BOTH))                doB = PETSC_TRUE;
247f3b08a26SMatthew G. Knepley   if ((galerkin == PC_MG_GALERKIN_MAT)  || ((galerkin == PC_MG_GALERKIN_BOTH) && (fA != fB))) doA = PETSC_TRUE;
2489566063dSJacob Faibussowitsch   if (doA) PetscCall(MatGalerkin(Restrc, fA, Interp, reuse, 1.0, &A));
2499566063dSJacob Faibussowitsch   if (doB) PetscCall(MatGalerkin(Restrc, fB, Interp, reuse, 1.0, &B));
250f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
251f3b08a26SMatthew G. Knepley }
252