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