xref: /petsc/src/ksp/pc/impls/mg/mgadapt.c (revision 2b3cbbda7ef6bfb59aeed26278d0c91b4282b9fd)
1 #include <petsc/private/pcmgimpl.h>       /*I "petscksp.h" I*/
2 #include <petscdm.h>
3 
4 static PetscErrorCode xfunc(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
5 {
6   PetscInt k = *((PetscInt *) ctx), c;
7 
8   for (c = 0; c < Nc; ++c) u[c] = PetscPowRealInt(coords[0], k);
9   return 0;
10 }
11 static PetscErrorCode yfunc(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
12 {
13   PetscInt k = *((PetscInt *) ctx), c;
14 
15   for (c = 0; c < Nc; ++c) u[c] = PetscPowRealInt(coords[1], k);
16   return 0;
17 }
18 static PetscErrorCode zfunc(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
19 {
20   PetscInt k = *((PetscInt *) ctx), c;
21 
22   for (c = 0; c < Nc; ++c) u[c] = PetscPowRealInt(coords[2], k);
23   return 0;
24 }
25 static PetscErrorCode xsin(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
26 {
27   PetscInt k = *((PetscInt *) ctx), c;
28 
29   for (c = 0; c < Nc; ++c) u[c] = PetscSinReal(PETSC_PI*(k+1)*coords[0]);
30   return 0;
31 }
32 static PetscErrorCode ysin(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
33 {
34   PetscInt k = *((PetscInt *) ctx), c;
35 
36   for (c = 0; c < Nc; ++c) u[c] = PetscSinReal(PETSC_PI*(k+1)*coords[1]);
37   return 0;
38 }
39 static PetscErrorCode zsin(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx)
40 {
41   PetscInt k = *((PetscInt *) ctx), c;
42 
43   for (c = 0; c < Nc; ++c) u[c] = PetscSinReal(PETSC_PI*(k+1)*coords[2]);
44   return 0;
45 }
46 
47 PetscErrorCode DMSetBasisFunction_Internal(PetscInt Nf, PetscBool usePoly, PetscInt dir, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *))
48 {
49   PetscInt f;
50 
51   PetscFunctionBeginUser;
52   for (f = 0; f < Nf; ++f) {
53     if (usePoly) {
54       switch (dir) {
55       case 0: funcs[f] = xfunc;break;
56       case 1: funcs[f] = yfunc;break;
57       case 2: funcs[f] = zfunc;break;
58       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No function for direction %" PetscInt_FMT, dir);
59       }
60     } else {
61       switch (dir) {
62       case 0: funcs[f] = xsin;break;
63       case 1: funcs[f] = ysin;break;
64       case 2: funcs[f] = zsin;break;
65       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No function for direction %" PetscInt_FMT, dir);
66       }
67     }
68   }
69   PetscFunctionReturn(0);
70 }
71 
72 static PetscErrorCode PCMGCreateCoarseSpaceDefault_Private(PC pc, PetscInt level, PCMGCoarseSpaceType cstype, DM dm, KSP ksp, PetscInt Nc, Mat initialGuess, Mat *coarseSpace)
73 {
74   PetscBool         poly = cstype == PCMG_ADAPT_POLYNOMIAL ? PETSC_TRUE : PETSC_FALSE;
75   PetscErrorCode (**funcs)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar*,void*);
76   void            **ctxs;
77   PetscInt          dim, d, Nf, f, k, m, M;
78   Vec               tmp;
79 
80   PetscFunctionBegin;
81   Nc = Nc < 0 ? 6 : Nc;
82   PetscCall(DMGetCoordinateDim(dm, &dim));
83   PetscCall(DMGetNumFields(dm, &Nf));
84   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);
85   PetscCall(PetscMalloc2(Nf, &funcs, Nf, &ctxs));
86   PetscCall(DMGetGlobalVector(dm, &tmp));
87   PetscCall(VecGetSize(tmp, &M));
88   PetscCall(VecGetLocalSize(tmp, &m));
89   PetscCall(MatCreateDense(PetscObjectComm((PetscObject) pc), m, PETSC_DECIDE, M, Nc, NULL, coarseSpace));
90   PetscCall(DMRestoreGlobalVector(dm, &tmp));
91   for (k = 0; k < Nc/dim; ++k) {
92     for (f = 0; f < Nf; ++f) ctxs[f] = &k;
93     for (d = 0; d < dim; ++d) {
94       PetscCall(MatDenseGetColumnVecWrite(*coarseSpace,k*dim+d,&tmp));
95       PetscCall(DMSetBasisFunction_Internal(Nf, poly, d, funcs));
96       PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, tmp));
97       PetscCall(MatDenseRestoreColumnVecWrite(*coarseSpace,k*dim+d,&tmp));
98     }
99   }
100   PetscCall(PetscFree2(funcs, ctxs));
101   PetscFunctionReturn(0);
102 }
103 
104 static PetscErrorCode PCMGCreateCoarseSpace_Polynomial(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, Mat initialGuess, Mat *coarseSpace)
105 {
106   PetscFunctionBegin;
107   PetscCall(PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_ADAPT_POLYNOMIAL, dm, ksp, Nc, initialGuess, coarseSpace));
108   PetscFunctionReturn(0);
109 }
110 
111 PetscErrorCode PCMGCreateCoarseSpace_Harmonic(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, Mat initialGuess, Mat *coarseSpace)
112 {
113   PetscFunctionBegin;
114   PetscCall(PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_ADAPT_HARMONIC, dm, ksp, Nc, initialGuess, coarseSpace));
115   PetscFunctionReturn(0);
116 }
117 
118 /*
119   PCMGComputeCoarseSpace_Internal - Compute vectors on level l that must be accurately interpolated.
120 
121   Input Parameters:
122 + pc     - The PCMG
123 . l      - The level
124 . Nc     - The number of vectors requested
125 - cspace - The initial guess for the space, or NULL
126 
127   Output Parameter:
128 . space  - The space which must be accurately interpolated.
129 
130   Level: developer
131 
132   Note: This space is normally used to adapt the interpolator. If Nc is negative, an adaptive choice can be made.
133 
134 .seealso: `PCMGAdaptInterpolator_Private()`
135 */
136 PetscErrorCode PCMGComputeCoarseSpace_Internal(PC pc, PetscInt l, PCMGCoarseSpaceType cstype, PetscInt Nc, Mat cspace, Mat *space)
137 {
138   PetscErrorCode (*coarseConstructor)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*) = NULL;
139   DM             dm;
140   KSP            smooth;
141 
142   PetscFunctionBegin;
143   *space = NULL;
144   switch (cstype) {
145   case PCMG_ADAPT_POLYNOMIAL:
146     coarseConstructor = &PCMGCreateCoarseSpace_Polynomial;
147     break;
148   case PCMG_ADAPT_HARMONIC:
149     coarseConstructor = &PCMGCreateCoarseSpace_Harmonic;
150     break;
151   case PCMG_ADAPT_EIGENVECTOR:
152     Nc = Nc < 0 ? 6 : Nc;
153     if (l > 0) PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_MEV", &coarseConstructor));
154     else       PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_EV", &coarseConstructor));
155     break;
156   case PCMG_ADAPT_GENERALIZED_EIGENVECTOR:
157     Nc = Nc < 0 ? 6 : Nc;
158     if (l > 0) PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_MGEV", &coarseConstructor));
159     else       PetscCall(PCMGGetCoarseSpaceConstructor("BAMG_GEV", &coarseConstructor));
160     break;
161   case PCMG_ADAPT_NONE:
162     break;
163   default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot handle coarse space type %d", cstype);
164   }
165   if (coarseConstructor) {
166     PetscCall(PCMGGetSmoother(pc, l, &smooth));
167     PetscCall(KSPGetDM(smooth, &dm));
168     PetscCall((*coarseConstructor)(pc, l, dm, smooth, Nc, cspace, space));
169   }
170   PetscFunctionReturn(0);
171 }
172 
173 /*
174   PCMGAdaptInterpolator_Internal - Adapt interpolator from level l-1 to level l
175 
176   Input Parameters:
177 + pc      - The PCMG
178 . l       - The level l
179 . csmooth - The (coarse) smoother for level l-1
180 . fsmooth - The (fine) smoother for level l
181 . cspace  - The (coarse) vectors in the subspace for level l-1
182 - fspace  - The (fine) vectors in the subspace for level l
183 
184   Level: developer
185 
186   Note: This routine resets the interpolation and restriction for level l.
187 
188 .seealso: `PCMGComputeCoarseSpace_Private()`
189 */
190 PetscErrorCode PCMGAdaptInterpolator_Internal(PC pc, PetscInt l, KSP csmooth, KSP fsmooth, Mat cspace, Mat fspace)
191 {
192   PC_MG         *mg = (PC_MG *) pc->data;
193   DM             dm, cdm;
194   Mat            Interp, InterpAdapt;
195 
196   PetscFunctionBegin;
197   /* There is no interpolator for the coarse level */
198   if (!l) PetscFunctionReturn(0);
199   PetscCall(KSPGetDM(csmooth, &cdm));
200   PetscCall(KSPGetDM(fsmooth, &dm));
201   PetscCall(PCMGGetInterpolation(pc, l, &Interp));
202   if (Interp == fspace && !cspace) PetscFunctionReturn(0);
203   PetscCall(DMAdaptInterpolator(cdm, dm, Interp, fsmooth, fspace, cspace, &InterpAdapt, pc));
204   if (mg->mespMonitor) PetscCall(DMCheckInterpolator(dm, InterpAdapt, cspace, fspace, 0.5/* PETSC_SMALL */));
205   PetscCall(PCMGSetInterpolation(pc, l, InterpAdapt));
206   PetscCall(PCMGSetRestriction(pc, l, InterpAdapt)); /* MATT: Remove????? */
207   PetscCall(MatDestroy(&InterpAdapt));
208   PetscFunctionReturn(0);
209 }
210 
211 /*
212   PCMGRecomputeLevelOperators_Internal - Recomputes Galerkin coarse operator when interpolation is adapted
213 
214   Input Parameters:
215 + pc - The PCMG
216 - l  - The level l
217 
218   Level: developer
219 
220   Note: This routine recomputes the Galerkin triple product for the operator on level l.
221 */
222 PetscErrorCode PCMGRecomputeLevelOperators_Internal(PC pc, PetscInt l)
223 {
224   Mat              fA, fB;                   /* The system and preconditioning operators on level l+1 */
225   Mat              A,  B;                    /* The system and preconditioning operators on level l */
226   Mat              Interp, Restrc;           /* The interpolation operator from level l to l+1, and restriction operator from level l+1 to l */
227   KSP              smooth, fsmooth;          /* The smoothers on levels l and l+1 */
228   PCMGGalerkinType galerkin;                 /* The Galerkin projection flag */
229   MatReuse         reuse = MAT_REUSE_MATRIX; /* The matrices are always assumed to be present already */
230   PetscBool        doA   = PETSC_FALSE;      /* Updates the system operator */
231   PetscBool        doB   = PETSC_FALSE;      /* Updates the preconditioning operator (A == B, then update B) */
232   PetscInt         n;                        /* The number of multigrid levels */
233 
234   PetscFunctionBegin;
235   PetscCall(PCMGGetGalerkin(pc, &galerkin));
236   if (galerkin >= PC_MG_GALERKIN_NONE) PetscFunctionReturn(0);
237   PetscCall(PCMGGetLevels(pc, &n));
238   /* Do not recompute operator for the finest grid */
239   if (l == n-1) PetscFunctionReturn(0);
240   PetscCall(PCMGGetSmoother(pc, l,   &smooth));
241   PetscCall(KSPGetOperators(smooth, &A, &B));
242   PetscCall(PCMGGetSmoother(pc, l+1, &fsmooth));
243   PetscCall(KSPGetOperators(fsmooth, &fA, &fB));
244   PetscCall(PCMGGetInterpolation(pc, l+1, &Interp));
245   PetscCall(PCMGGetRestriction(pc, l+1, &Restrc));
246   if ((galerkin == PC_MG_GALERKIN_PMAT) ||  (galerkin == PC_MG_GALERKIN_BOTH))                doB = PETSC_TRUE;
247   if ((galerkin == PC_MG_GALERKIN_MAT)  || ((galerkin == PC_MG_GALERKIN_BOTH) && (fA != fB))) doA = PETSC_TRUE;
248   if (doA) PetscCall(MatGalerkin(Restrc, fA, Interp, reuse, 1.0, &A));
249   if (doB) PetscCall(MatGalerkin(Restrc, fB, Interp, reuse, 1.0, &B));
250   PetscFunctionReturn(0);
251 }
252