1*f3b08a26SMatthew G. Knepley #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 2*f3b08a26SMatthew G. Knepley #include <petscdm.h> 3*f3b08a26SMatthew G. Knepley 4*f3b08a26SMatthew G. Knepley static PetscErrorCode PCMGCreateCoarseSpaceDefault_Private(PC pc, PetscInt level, PCMGCoarseSpaceType cstype, DM dm, KSP ksp, PetscInt Nc, const Vec initialGuess[], Vec **coarseSpace) 5*f3b08a26SMatthew G. Knepley { 6*f3b08a26SMatthew G. Knepley PetscBool poly = cstype == PCMG_POLYNOMIAL ? PETSC_TRUE : PETSC_FALSE; 7*f3b08a26SMatthew G. Knepley PetscErrorCode (**funcs)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar*,void*); 8*f3b08a26SMatthew G. Knepley void **ctxs; 9*f3b08a26SMatthew G. Knepley PetscInt dim, d, Nf, f, k; 10*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 11*f3b08a26SMatthew G. Knepley 12*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 13*f3b08a26SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 14*f3b08a26SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 15*f3b08a26SMatthew 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); 16*f3b08a26SMatthew G. Knepley ierr = PetscMalloc2(Nf, &funcs, Nf, &ctxs);CHKERRQ(ierr); 17*f3b08a26SMatthew G. Knepley if (!*coarseSpace) {ierr = PetscCalloc1(Nc, coarseSpace);CHKERRQ(ierr);} 18*f3b08a26SMatthew G. Knepley for (k = 0; k < Nc/dim; ++k) { 19*f3b08a26SMatthew G. Knepley for (f = 0; f < Nf; ++f) {ctxs[f] = &k;} 20*f3b08a26SMatthew G. Knepley for (d = 0; d < dim; ++d) { 21*f3b08a26SMatthew G. Knepley if (!(*coarseSpace)[k*dim+d]) {ierr = DMCreateGlobalVector(dm, &(*coarseSpace)[k*dim+d]);CHKERRQ(ierr);} 22*f3b08a26SMatthew G. Knepley ierr = DMSetBasisFunction_Internal(Nf, poly, d, funcs);CHKERRQ(ierr); 23*f3b08a26SMatthew G. Knepley ierr = DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, (*coarseSpace)[k*dim+d]);CHKERRQ(ierr); 24*f3b08a26SMatthew G. Knepley } 25*f3b08a26SMatthew G. Knepley } 26*f3b08a26SMatthew G. Knepley ierr = PetscFree2(funcs, ctxs);CHKERRQ(ierr); 27*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 28*f3b08a26SMatthew G. Knepley } 29*f3b08a26SMatthew G. Knepley 30*f3b08a26SMatthew G. Knepley static PetscErrorCode PCMGCreateCoarseSpace_Polynomial(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, const Vec initialGuess[], Vec **coarseSpace) 31*f3b08a26SMatthew G. Knepley { 32*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 33*f3b08a26SMatthew G. Knepley 34*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 35*f3b08a26SMatthew G. Knepley ierr = PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_POLYNOMIAL, dm, ksp, Nc, initialGuess, coarseSpace);CHKERRQ(ierr); 36*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 37*f3b08a26SMatthew G. Knepley } 38*f3b08a26SMatthew G. Knepley 39*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGCreateCoarseSpace_Harmonic(PC pc, PetscInt level, DM dm, KSP ksp, PetscInt Nc, const Vec initialGuess[], Vec **coarseSpace) 40*f3b08a26SMatthew G. Knepley { 41*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 42*f3b08a26SMatthew G. Knepley 43*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 44*f3b08a26SMatthew G. Knepley ierr = PCMGCreateCoarseSpaceDefault_Private(pc, level, PCMG_HARMONIC, dm, ksp, Nc, initialGuess, coarseSpace);CHKERRQ(ierr); 45*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 46*f3b08a26SMatthew G. Knepley } 47*f3b08a26SMatthew G. Knepley 48*f3b08a26SMatthew G. Knepley /* 49*f3b08a26SMatthew G. Knepley PCMGComputeCoarseSpace_Internal - Compute vectors on level l that must be accurately interpolated. 50*f3b08a26SMatthew G. Knepley 51*f3b08a26SMatthew G. Knepley Input Parameters: 52*f3b08a26SMatthew G. Knepley + pc - The PCMG 53*f3b08a26SMatthew G. Knepley . l - The level 54*f3b08a26SMatthew G. Knepley . Nc - The size of the space (number of vectors) 55*f3b08a26SMatthew G. Knepley - cspace - The space from level l-1, or NULL 56*f3b08a26SMatthew G. Knepley 57*f3b08a26SMatthew G. Knepley Output Parameter: 58*f3b08a26SMatthew G. Knepley . space - The space which must be accurately interpolated. 59*f3b08a26SMatthew G. Knepley 60*f3b08a26SMatthew G. Knepley Level: developer 61*f3b08a26SMatthew G. Knepley 62*f3b08a26SMatthew G. Knepley Note: This space is normally used to adapt the interpolator. 63*f3b08a26SMatthew G. Knepley 64*f3b08a26SMatthew G. Knepley .seealso: PCMGAdaptInterpolator_Private() 65*f3b08a26SMatthew G. Knepley */ 66*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGComputeCoarseSpace_Internal(PC pc, PetscInt l, PCMGCoarseSpaceType cstype, PetscInt Nc, const Vec cspace[], Vec *space[]) 67*f3b08a26SMatthew G. Knepley { 68*f3b08a26SMatthew G. Knepley PetscErrorCode (*coarseConstructor)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec*[]); 69*f3b08a26SMatthew G. Knepley DM dm; 70*f3b08a26SMatthew G. Knepley KSP smooth; 71*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 72*f3b08a26SMatthew G. Knepley 73*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 74*f3b08a26SMatthew G. Knepley switch (cstype) { 75*f3b08a26SMatthew G. Knepley case PCMG_POLYNOMIAL: coarseConstructor = &PCMGCreateCoarseSpace_Polynomial;break; 76*f3b08a26SMatthew G. Knepley case PCMG_HARMONIC: coarseConstructor = &PCMGCreateCoarseSpace_Harmonic;break; 77*f3b08a26SMatthew G. Knepley case PCMG_EIGENVECTOR: 78*f3b08a26SMatthew G. Knepley if (l > 0) {ierr = PCMGGetCoarseSpaceConstructor("BAMG_MEV", &coarseConstructor);CHKERRQ(ierr);} 79*f3b08a26SMatthew G. Knepley else {ierr = PCMGGetCoarseSpaceConstructor("BAMG_EV", &coarseConstructor);CHKERRQ(ierr);} 80*f3b08a26SMatthew G. Knepley break; 81*f3b08a26SMatthew G. Knepley case PCMG_GENERALIZED_EIGENVECTOR: 82*f3b08a26SMatthew G. Knepley if (l > 0) {ierr = PCMGGetCoarseSpaceConstructor("BAMG_MGEV", &coarseConstructor);CHKERRQ(ierr);} 83*f3b08a26SMatthew G. Knepley else {ierr = PCMGGetCoarseSpaceConstructor("BAMG_GEV", &coarseConstructor);CHKERRQ(ierr);} 84*f3b08a26SMatthew G. Knepley break; 85*f3b08a26SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle coarse space type %D", cstype); 86*f3b08a26SMatthew G. Knepley } 87*f3b08a26SMatthew G. Knepley ierr = PCMGGetSmoother(pc, l, &smooth);CHKERRQ(ierr); 88*f3b08a26SMatthew G. Knepley ierr = KSPGetDM(smooth, &dm);CHKERRQ(ierr); 89*f3b08a26SMatthew G. Knepley ierr = (*coarseConstructor)(pc, l, dm, smooth, Nc, cspace, space);CHKERRQ(ierr); 90*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 91*f3b08a26SMatthew G. Knepley } 92*f3b08a26SMatthew G. Knepley 93*f3b08a26SMatthew G. Knepley /* 94*f3b08a26SMatthew G. Knepley PCMGAdaptInterpolator_Internal - Adapt interpolator from level l-1 to level 1 95*f3b08a26SMatthew G. Knepley 96*f3b08a26SMatthew G. Knepley Input Parameters: 97*f3b08a26SMatthew G. Knepley + pc - The PCMG 98*f3b08a26SMatthew G. Knepley . l - The level l 99*f3b08a26SMatthew G. Knepley . csmooth - The (coarse) smoother for level l-1 100*f3b08a26SMatthew G. Knepley . fsmooth - The (fine) smoother for level l 101*f3b08a26SMatthew G. Knepley . Nc - The size of the subspace used for adaptation 102*f3b08a26SMatthew G. Knepley . cspace - The (coarse) vectors in the subspace for level l-1 103*f3b08a26SMatthew G. Knepley - fspace - The (fine) vectors in the subspace for level l 104*f3b08a26SMatthew G. Knepley 105*f3b08a26SMatthew G. Knepley Level: developer 106*f3b08a26SMatthew G. Knepley 107*f3b08a26SMatthew G. Knepley Note: This routine resets the interpolation and restriction for level l. 108*f3b08a26SMatthew G. Knepley 109*f3b08a26SMatthew G. Knepley .seealso: PCMGComputeCoarseSpace_Private() 110*f3b08a26SMatthew G. Knepley */ 111*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGAdaptInterpolator_Internal(PC pc, PetscInt l, KSP csmooth, KSP fsmooth, PetscInt Nc, Vec cspace[], Vec fspace[]) 112*f3b08a26SMatthew G. Knepley { 113*f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 114*f3b08a26SMatthew G. Knepley DM dm, cdm; 115*f3b08a26SMatthew G. Knepley Mat Interp, InterpAdapt; 116*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 117*f3b08a26SMatthew G. Knepley 118*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 119*f3b08a26SMatthew G. Knepley /* There is no interpolator for the coarse level */ 120*f3b08a26SMatthew G. Knepley if (!l) PetscFunctionReturn(0); 121*f3b08a26SMatthew G. Knepley ierr = KSPGetDM(csmooth, &cdm);CHKERRQ(ierr); 122*f3b08a26SMatthew G. Knepley ierr = KSPGetDM(fsmooth, &dm);CHKERRQ(ierr); 123*f3b08a26SMatthew G. Knepley ierr = PCMGGetInterpolation(pc, l, &Interp);CHKERRQ(ierr); 124*f3b08a26SMatthew G. Knepley 125*f3b08a26SMatthew G. Knepley ierr = DMAdaptInterpolator(cdm, dm, Interp, fsmooth, Nc, fspace, cspace, &InterpAdapt, pc);CHKERRQ(ierr); 126*f3b08a26SMatthew G. Knepley if (mg->mespMonitor) {ierr = DMCheckInterpolator(dm, InterpAdapt, Nc, cspace, fspace, 0.5/*PETSC_SMALL*/);CHKERRQ(ierr);} 127*f3b08a26SMatthew G. Knepley ierr = PCMGSetInterpolation(pc, l, InterpAdapt);CHKERRQ(ierr); 128*f3b08a26SMatthew G. Knepley ierr = PCMGSetRestriction(pc, l, InterpAdapt);CHKERRQ(ierr); 129*f3b08a26SMatthew G. Knepley ierr = MatDestroy(&InterpAdapt);CHKERRQ(ierr); 130*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 131*f3b08a26SMatthew G. Knepley } 132*f3b08a26SMatthew G. Knepley 133*f3b08a26SMatthew G. Knepley /* 134*f3b08a26SMatthew G. Knepley PCMGRecomputeLevelOperators_Internal - Recomputes Galerkin coarse operator when interpolation is adapted 135*f3b08a26SMatthew G. Knepley 136*f3b08a26SMatthew G. Knepley Input Parameters: 137*f3b08a26SMatthew G. Knepley + pc - The PCMG 138*f3b08a26SMatthew G. Knepley - l - The level l 139*f3b08a26SMatthew G. Knepley 140*f3b08a26SMatthew G. Knepley Level: developer 141*f3b08a26SMatthew G. Knepley 142*f3b08a26SMatthew G. Knepley Note: This routine recomputes the Galerkin triple product for the operator on level l. 143*f3b08a26SMatthew G. Knepley */ 144*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRecomputeLevelOperators_Internal(PC pc, PetscInt l) 145*f3b08a26SMatthew G. Knepley { 146*f3b08a26SMatthew G. Knepley Mat fA, fB; /* The system and preconditioning operators on level l+1 */ 147*f3b08a26SMatthew G. Knepley Mat A, B; /* The system and preconditioning operators on level l */ 148*f3b08a26SMatthew G. Knepley Mat Interp, Restrc; /* The interpolation operator from level l to l+1, and restriction operator from level l+1 to l */ 149*f3b08a26SMatthew G. Knepley KSP smooth, fsmooth; /* The smoothers on levels l and l+1 */ 150*f3b08a26SMatthew G. Knepley PCMGGalerkinType galerkin; /* The Galerkin projection flag */ 151*f3b08a26SMatthew G. Knepley MatReuse reuse = MAT_REUSE_MATRIX; /* The matrices are always assumed to be present already */ 152*f3b08a26SMatthew G. Knepley PetscBool doA = PETSC_FALSE; /* Updates the system operator */ 153*f3b08a26SMatthew G. Knepley PetscBool doB = PETSC_FALSE; /* Updates the preconditioning operator (A == B, then update B) */ 154*f3b08a26SMatthew G. Knepley PetscInt n; /* The number of multigrid levels */ 155*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 156*f3b08a26SMatthew G. Knepley 157*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 158*f3b08a26SMatthew G. Knepley ierr = PCMGGetGalerkin(pc, &galerkin);CHKERRQ(ierr); 159*f3b08a26SMatthew G. Knepley if (galerkin >= PC_MG_GALERKIN_NONE) PetscFunctionReturn(0); 160*f3b08a26SMatthew G. Knepley ierr = PCMGGetLevels(pc, &n);CHKERRQ(ierr); 161*f3b08a26SMatthew G. Knepley /* Do not recompute operator for the finest grid */ 162*f3b08a26SMatthew G. Knepley if (l == n-1) PetscFunctionReturn(0); 163*f3b08a26SMatthew G. Knepley ierr = PCMGGetSmoother(pc, l, &smooth);CHKERRQ(ierr); 164*f3b08a26SMatthew G. Knepley ierr = KSPGetOperators(smooth, &A, &B);CHKERRQ(ierr); 165*f3b08a26SMatthew G. Knepley ierr = PCMGGetSmoother(pc, l+1, &fsmooth);CHKERRQ(ierr); 166*f3b08a26SMatthew G. Knepley ierr = KSPGetOperators(fsmooth, &fA, &fB);CHKERRQ(ierr); 167*f3b08a26SMatthew G. Knepley ierr = PCMGGetInterpolation(pc, l+1, &Interp);CHKERRQ(ierr); 168*f3b08a26SMatthew G. Knepley ierr = PCMGGetRestriction(pc, l+1, &Restrc);CHKERRQ(ierr); 169*f3b08a26SMatthew G. Knepley if ((galerkin == PC_MG_GALERKIN_PMAT) || (galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 170*f3b08a26SMatthew G. Knepley if ((galerkin == PC_MG_GALERKIN_MAT) || ((galerkin == PC_MG_GALERKIN_BOTH) && (fA != fB))) doA = PETSC_TRUE; 171*f3b08a26SMatthew G. Knepley if (doA) {ierr = MatGalerkin(Restrc, fA, Interp, reuse, 1.0, &A);CHKERRQ(ierr);} 172*f3b08a26SMatthew G. Knepley if (doB) {ierr = MatGalerkin(Restrc, fB, Interp, reuse, 1.0, &B);CHKERRQ(ierr);} 173*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 174*f3b08a26SMatthew G. Knepley } 175