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_GDSW: 162 coarseConstructor = &PCMGGDSWCreateCoarseSpace_Private; 163 break; 164 case PCMG_ADAPT_NONE: 165 break; 166 default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot handle coarse space type %d", cstype); 167 } 168 if (coarseConstructor) { 169 PetscCall(PCMGGetSmoother(pc, l, &smooth)); 170 PetscCall(KSPGetDM(smooth, &dm)); 171 PetscCall((*coarseConstructor)(pc, l, dm, smooth, Nc, cspace, space)); 172 } 173 PetscFunctionReturn(0); 174 } 175 176 /* 177 PCMGAdaptInterpolator_Internal - Adapt interpolator from level l-1 to level l 178 179 Input Parameters: 180 + pc - The PCMG 181 . l - The level l 182 . csmooth - The (coarse) smoother for level l-1 183 . fsmooth - The (fine) smoother for level l 184 . cspace - The (coarse) vectors in the subspace for level l-1 185 - fspace - The (fine) vectors in the subspace for level l 186 187 Level: developer 188 189 Note: This routine resets the interpolation and restriction for level l. 190 191 .seealso: `PCMGComputeCoarseSpace_Private()` 192 */ 193 PetscErrorCode PCMGAdaptInterpolator_Internal(PC pc, PetscInt l, KSP csmooth, KSP fsmooth, Mat cspace, Mat fspace) 194 { 195 PC_MG *mg = (PC_MG *) pc->data; 196 DM dm, cdm; 197 Mat Interp, InterpAdapt; 198 199 PetscFunctionBegin; 200 /* There is no interpolator for the coarse level */ 201 if (!l) PetscFunctionReturn(0); 202 PetscCall(KSPGetDM(csmooth, &cdm)); 203 PetscCall(KSPGetDM(fsmooth, &dm)); 204 PetscCall(PCMGGetInterpolation(pc, l, &Interp)); 205 if (Interp == fspace && !cspace) PetscFunctionReturn(0); 206 PetscCall(DMAdaptInterpolator(cdm, dm, Interp, fsmooth, fspace, cspace, &InterpAdapt, pc)); 207 if (mg->mespMonitor) PetscCall(DMCheckInterpolator(dm, InterpAdapt, cspace, fspace, 0.5/* PETSC_SMALL */)); 208 PetscCall(PCMGSetInterpolation(pc, l, InterpAdapt)); 209 PetscCall(PCMGSetRestriction(pc, l, InterpAdapt)); /* MATT: Remove????? */ 210 PetscCall(MatDestroy(&InterpAdapt)); 211 PetscFunctionReturn(0); 212 } 213 214 /* 215 PCMGRecomputeLevelOperators_Internal - Recomputes Galerkin coarse operator when interpolation is adapted 216 217 Input Parameters: 218 + pc - The PCMG 219 - l - The level l 220 221 Level: developer 222 223 Note: This routine recomputes the Galerkin triple product for the operator on level l. 224 */ 225 PetscErrorCode PCMGRecomputeLevelOperators_Internal(PC pc, PetscInt l) 226 { 227 Mat fA, fB; /* The system and preconditioning operators on level l+1 */ 228 Mat A, B; /* The system and preconditioning operators on level l */ 229 Mat Interp, Restrc; /* The interpolation operator from level l to l+1, and restriction operator from level l+1 to l */ 230 KSP smooth, fsmooth; /* The smoothers on levels l and l+1 */ 231 PCMGGalerkinType galerkin; /* The Galerkin projection flag */ 232 MatReuse reuse = MAT_REUSE_MATRIX; /* The matrices are always assumed to be present already */ 233 PetscBool doA = PETSC_FALSE; /* Updates the system operator */ 234 PetscBool doB = PETSC_FALSE; /* Updates the preconditioning operator (A == B, then update B) */ 235 PetscInt n; /* The number of multigrid levels */ 236 237 PetscFunctionBegin; 238 PetscCall(PCMGGetGalerkin(pc, &galerkin)); 239 if (galerkin >= PC_MG_GALERKIN_NONE) PetscFunctionReturn(0); 240 PetscCall(PCMGGetLevels(pc, &n)); 241 /* Do not recompute operator for the finest grid */ 242 if (l == n-1) PetscFunctionReturn(0); 243 PetscCall(PCMGGetSmoother(pc, l, &smooth)); 244 PetscCall(KSPGetOperators(smooth, &A, &B)); 245 PetscCall(PCMGGetSmoother(pc, l+1, &fsmooth)); 246 PetscCall(KSPGetOperators(fsmooth, &fA, &fB)); 247 PetscCall(PCMGGetInterpolation(pc, l+1, &Interp)); 248 PetscCall(PCMGGetRestriction(pc, l+1, &Restrc)); 249 if ((galerkin == PC_MG_GALERKIN_PMAT) || (galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 250 if ((galerkin == PC_MG_GALERKIN_MAT) || ((galerkin == PC_MG_GALERKIN_BOTH) && (fA != fB))) doA = PETSC_TRUE; 251 if (doA) PetscCall(MatGalerkin(Restrc, fA, Interp, reuse, 1.0, &A)); 252 if (doB) PetscCall(MatGalerkin(Restrc, fB, Interp, reuse, 1.0, &B)); 253 PetscFunctionReturn(0); 254 } 255