1
0
mirror of https://github.com/paboyle/Grid.git synced 2024-11-10 07:55:35 +00:00

WilsonMG: Huge refactor into something that could be considered an algorithm

This commit is contained in:
Daniel Richtmann 2018-03-19 10:59:42 +01:00
parent 1cfed3de7c
commit 0f6009a29f
No known key found for this signature in database
GPG Key ID: B33C490AF0772057

View File

@ -69,55 +69,67 @@ public:
}
};
class myclass : Serializable {
// clang-format off
struct MultigridParams : Serializable {
public:
// clang-format off
GRID_SERIALIZABLE_CLASS_MEMBERS(myclass,
int, domaindecompose,
int, domainsize,
int, coarsegrids,
int, order,
int, Ls,
double, mq,
double, lo,
double, hi,
int, steps);
// clang-format on
myclass(){};
GRID_SERIALIZABLE_CLASS_MEMBERS(MultigridParams,
int, nLevels,
std::vector<std::vector<int>>, blockSizes);
MultigridParams(){};
};
myclass params;
MultigridParams mgParams;
// clang-format on
template<int nbasis> struct CoarseGrids {
struct LevelInfo {
public:
std::vector<std::vector<int>> LattSizes;
std::vector<std::vector<int>> Seeds;
std::vector<GridCartesian *> Grids;
std::vector<GridParallelRNG> PRNGs;
CoarseGrids(std::vector<std::vector<int>> const &blockSizes, int coarsegrids) {
LevelInfo(GridCartesian *FineGrid, MultigridParams const &Params) {
auto nCoarseLevels = Params.blockSizes.size();
assert(blockSizes.size() == coarsegrids);
assert(nCoarseLevels == Params.nLevels - 1);
std::cout << GridLogMessage << "Constructing " << coarsegrids << " CoarseGrids" << std::endl;
// set up values for finest grid
Grids.push_back(FineGrid);
Seeds.push_back({1, 2, 3, 4});
PRNGs.push_back(GridParallelRNG(Grids.back()));
PRNGs.back().SeedFixedIntegers(Seeds.back());
for(int cl = 0; cl < coarsegrids; ++cl) { // may be a bit ugly and slow but not perf critical
// need to differentiate between first and other coarse levels in size calculation
LattSizes.push_back({cl == 0 ? GridDefaultLatt() : LattSizes[cl - 1]});
Seeds.push_back(std::vector<int>(LattSizes[cl].size()));
// set up values for coarser grids
for(int level = 1; level < Params.nLevels; ++level) {
auto Nd = Grids[level - 1]->_ndimension;
auto tmp = Grids[level - 1]->_fdimensions;
assert(tmp.size() == Nd);
for(int d = 0; d < LattSizes[cl].size(); ++d) {
LattSizes[cl][d] = LattSizes[cl][d] / blockSizes[cl][d];
Seeds[cl][d] = (cl + 1) * LattSizes[cl].size() + d + 1;
// calculation unimportant, just to get. e.g., {5, 6, 7, 8} for first coarse level and so on
Seeds.push_back(std::vector<int>(Nd));
for(int d = 0; d < Nd; ++d) {
tmp[d] /= Params.blockSizes[level - 1][d];
Seeds[level][d] = (level)*Nd + d + 1;
}
Grids.push_back(SpaceTimeGrid::makeFourDimGrid(LattSizes[cl], GridDefaultSimd(Nd, vComplex::Nsimd()), GridDefaultMpi()));
PRNGs.push_back(GridParallelRNG(Grids[cl]));
Grids.push_back(SpaceTimeGrid::makeFourDimGrid(tmp, GridDefaultSimd(Nd, vComplex::Nsimd()), GridDefaultMpi()));
PRNGs.push_back(GridParallelRNG(Grids[level]));
PRNGs[cl].SeedFixedIntegers(Seeds[cl]);
PRNGs[level].SeedFixedIntegers(Seeds[level]);
}
std::cout << GridLogMessage << "cl = " << cl << ": LattSize = " << LattSizes[cl] << std::endl;
std::cout << GridLogMessage << "cl = " << cl << ": Seeds = " << Seeds[cl] << std::endl;
std::cout << GridLogMessage << "Constructed " << Params.nLevels << " levels" << std::endl;
// The construction above corresponds to the finest level having level == 0
// (simply because it's not as ugly to implement), but we need it the
// other way round (i.e., the coarsest level to have level == 0) for the MG
// Preconditioner -> reverse the vectors
std::reverse(Seeds.begin(), Seeds.end());
std::reverse(Grids.begin(), Grids.end());
std::reverse(PRNGs.begin(), PRNGs.end());
for(int level = 0; level < Params.nLevels; ++level) {
std::cout << GridLogMessage << "level = " << level << ":" << std::endl;
Grids[level]->show_decomposition();
}
}
};
@ -221,116 +233,177 @@ template<class Field> void testLinearOperator(LinearOperatorBase<Field> &LinOp,
}
}
// template < class Fobj, class CComplex, int coarseSpins, int nbasis, class Matrix >
// class MultiGridPreconditioner : public LinearFunction< Lattice< Fobj > > {
template<class Fobj, class CComplex, int nbasis, class Matrix> class MultiGridPreconditioner : public LinearFunction<Lattice<Fobj>> {
template<class Fobj, class CoarseScalar, int nCoarseSpins, int nBasis, int level, class Matrix>
class MultiGridPreconditioner : public LinearFunction<Lattice<Fobj>> {
public:
typedef Aggregation<Fobj, CComplex, nbasis> Aggregates;
typedef CoarsenedMatrix<Fobj, CComplex, nbasis> CoarseOperator;
/////////////////////////////////////////////
// Type Definitions
/////////////////////////////////////////////
typedef typename Aggregation<Fobj, CComplex, nbasis>::siteVector siteVector;
typedef typename Aggregation<Fobj, CComplex, nbasis>::CoarseScalar CoarseScalar;
typedef typename Aggregation<Fobj, CComplex, nbasis>::CoarseVector CoarseVector;
typedef typename Aggregation<Fobj, CComplex, nbasis>::CoarseMatrix CoarseMatrix;
typedef typename Aggregation<Fobj, CComplex, nbasis>::FineField FineField;
typedef LinearOperatorBase<FineField> FineOperator;
typedef Aggregation<Fobj, CoarseScalar, nBasis> Aggregates;
typedef CoarsenedMatrix<Fobj, CoarseScalar, nBasis> CoarseMatrix;
typedef typename Aggregates::CoarseVector CoarseVector;
typedef typename Aggregates::siteVector CoarseSiteVector;
typedef Matrix FineMatrix;
typedef typename Aggregates::FineField FineVector;
typedef MultiGridPreconditioner<CoarseSiteVector, CoarseScalar, nCoarseSpins, nBasis, level - 1, CoarseMatrix> NextPreconditionerLevel;
Aggregates & _Aggregates;
CoarseOperator &_CoarseOperator;
Matrix & _FineMatrix;
FineOperator & _FineOperator;
Matrix & _SmootherMatrix;
FineOperator & _SmootherOperator;
/////////////////////////////////////////////
// Member Data
/////////////////////////////////////////////
// Constructor
MultiGridPreconditioner(Aggregates & Agg,
CoarseOperator &Coarse,
FineOperator & Fine,
Matrix & FineMatrix,
FineOperator & Smooth,
Matrix & SmootherMatrix)
: _Aggregates(Agg)
, _CoarseOperator(Coarse)
, _FineOperator(Fine)
, _FineMatrix(FineMatrix)
, _SmootherOperator(Smooth)
, _SmootherMatrix(SmootherMatrix) {}
LevelInfo & _LevelInfo;
FineMatrix & _FineMatrix;
FineMatrix & _SmootherMatrix;
Aggregates _Aggregates;
CoarseMatrix _CoarseMatrix;
std::unique_ptr<NextPreconditionerLevel> _NextPreconditionerLevel;
void operator()(const FineField &in, FineField &out) {
/////////////////////////////////////////////
// Member Functions
/////////////////////////////////////////////
CoarseVector coarseSrc(_CoarseOperator.Grid());
CoarseVector coarseTmp(_CoarseOperator.Grid());
CoarseVector coarseSol(_CoarseOperator.Grid());
coarseSol = zero;
GeneralisedMinimalResidual<CoarseVector> coarseGMRES(5.0e-2, 100, 25, false);
GeneralisedMinimalResidual<FineField> fineGMRES(5.0e-2, 100, 25, false);
HermitianLinearOperator<CoarseOperator, CoarseVector> coarseHermOp(_CoarseOperator);
MdagMLinearOperator<CoarseOperator, CoarseVector> coarseMdagMOp(_CoarseOperator);
MdagMLinearOperator<Matrix, FineField> fineMdagMOp(_SmootherMatrix);
FineField fineTmp1(in._grid);
FineField fineTmp2(in._grid);
RealD Ni = norm2(in);
// no pre smoothing for now
auto preSmootherNorm = 0;
auto preSmootherResidual = 0;
RealD r;
// Project to coarse grid, solve, project back to fine grid
_Aggregates.ProjectToSubspace(coarseSrc, in);
coarseGMRES(coarseMdagMOp, coarseSrc, coarseSol);
_Aggregates.PromoteFromSubspace(coarseSol, out);
// Recompute error
_FineOperator.Op(out, fineTmp1);
fineTmp1 = in - fineTmp1;
r = norm2(fineTmp1);
auto coarseResidual = std::sqrt(r / Ni);
// Apply smoother, use GMRES for the moment
fineGMRES(fineMdagMOp, in, out);
// Recompute error
_FineOperator.Op(out, fineTmp1);
fineTmp1 = in - fineTmp1;
r = norm2(fineTmp1);
auto postSmootherResidual = std::sqrt(r / Ni);
std::cout << GridLogIterative << "Input norm = " << Ni << " Pre-Smoother norm " << preSmootherNorm
<< " Pre-Smoother residual = " << preSmootherResidual << " Coarse residual = " << coarseResidual
<< " Post-Smoother residual = " << postSmootherResidual << std::endl;
MultiGridPreconditioner(LevelInfo &LvlInfo, FineMatrix &FineMat, FineMatrix &SmootherMat)
: _LevelInfo(LvlInfo)
, _FineMatrix(FineMat)
, _SmootherMatrix(SmootherMat)
, _Aggregates(_LevelInfo.Grids[level - 1], _LevelInfo.Grids[level], 0)
, _CoarseMatrix(*_LevelInfo.Grids[level - 1]) {
_NextPreconditionerLevel = std::unique_ptr<NextPreconditionerLevel>(
new NextPreconditionerLevel(_LevelInfo, _CoarseMatrix, _CoarseMatrix));
}
void runChecks(CoarseGrids<nbasis> &cGrids, int whichCoarseGrid) {
void setup() {
/////////////////////////////////////////////
// Some stuff we need for the checks below //
/////////////////////////////////////////////
auto tolerance = 1e-13; // TODO: this obviously depends on the precision we use, current value is for double
Gamma g5(Gamma::Algebra::Gamma5);
MdagMLinearOperator<FineMatrix, FineVector> fineMdagMOp(_FineMatrix);
std::vector<CoarseVector> cTmps(4, _CoarseOperator.Grid());
std::vector<FineField> fTmps(2, _Aggregates.subspace[0]._grid); // atm only for one coarser grid
_Aggregates.CreateSubspace(_LevelInfo.PRNGs[level], fineMdagMOp /*, nb */); // NOTE: Don't specify nb to see the orthogonalization check
// need to construct an operator, since _CoarseOperator is not a LinearOperator but only a matrix (the name is a bit misleading)
MdagMLinearOperator<CoarseOperator, CoarseVector> MdagMOp(_CoarseOperator);
// TestVectorAnalyzer<FineVector, nbasis> fineTVA;
// fineTVA(fineMdagMOp, _Aggregates.subspace);
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "MG correctness check: 0 == (1 - P R) v" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
static_assert((nBasis & 0x1) == 0, "MG Preconditioner only supports an even number of basis vectors");
int nb = nBasis / 2;
for(
int n = 0; n < nb;
n++) { // TODO: to get this to work for more than two levels, I would need to either implement coarse spins or have a template specialization of this class also for the finest level
_Aggregates.subspace[n + nb] = g5 * _Aggregates.subspace[n];
}
_CoarseMatrix.CoarsenOperator(_LevelInfo.Grids[level], fineMdagMOp, _Aggregates);
_NextPreconditionerLevel->setup();
}
virtual void operator()(Lattice<Fobj> const &in, Lattice<Fobj> &out) {
// TODO: implement a W-cycle and a toggle to switch between the cycling strategies
vCycle(in, out);
// kCycle(in, out);
}
void vCycle(Lattice<Fobj> const &in, Lattice<Fobj> &out) {
RealD inputNorm = norm2(in);
CoarseVector coarseSrc(_LevelInfo.Grids[level - 1]);
CoarseVector coarseSol(_LevelInfo.Grids[level - 1]);
coarseSol = zero;
FineVector fineTmp(in._grid);
TrivialPrecon<FineVector> fineTrivialPreconditioner;
FlexibleGeneralisedMinimalResidual<FineVector> fineFGMRES(1.0e-14, 1, fineTrivialPreconditioner, 1, false);
MdagMLinearOperator<FineMatrix, FineVector> fineMdagMOp(_FineMatrix);
MdagMLinearOperator<FineMatrix, FineVector> fineSmootherMdagMOp(_SmootherMatrix);
_Aggregates.ProjectToSubspace(coarseSrc, in);
(*_NextPreconditionerLevel)(coarseSrc, coarseSol);
_Aggregates.PromoteFromSubspace(coarseSol, out);
fineMdagMOp.Op(out, fineTmp);
fineTmp = in - fineTmp;
auto r = norm2(fineTmp);
auto residualAfterCoarseGridCorrection = std::sqrt(r / inputNorm);
fineFGMRES(fineSmootherMdagMOp, in, out);
fineMdagMOp.Op(out, fineTmp);
fineTmp = in - fineTmp;
r = norm2(fineTmp);
auto residualAfterPostSmoother = std::sqrt(r / inputNorm);
std::cout << GridLogMG << " Level " << level << ": Input norm = " << std::sqrt(inputNorm)
<< " Coarse residual = " << residualAfterCoarseGridCorrection << " Post-Smoother residual = " << residualAfterPostSmoother
<< std::endl;
}
void kCycle(Lattice<Fobj> const &in, Lattice<Fobj> &out) {
RealD inputNorm = norm2(in);
CoarseVector coarseSrc(_LevelInfo.Grids[level - 1]);
CoarseVector coarseSol(_LevelInfo.Grids[level - 1]);
coarseSol = zero;
FineVector fineTmp(in._grid);
TrivialPrecon<FineVector> fineTrivialPreconditioner;
FlexibleGeneralisedMinimalResidual<FineVector> fineFGMRES(1.0e-14, 1, fineTrivialPreconditioner, 1, false);
FlexibleGeneralisedMinimalResidual<CoarseVector> coarseFGMRES(1.0e-14, 1, *_NextPreconditionerLevel, 1, false);
MdagMLinearOperator<FineMatrix, FineVector> fineMdagMOp(_FineMatrix);
MdagMLinearOperator<FineMatrix, FineVector> fineSmootherMdagMOp(_SmootherMatrix);
MdagMLinearOperator<CoarseMatrix, CoarseVector> coarseMdagMOp(_CoarseMatrix);
_Aggregates.ProjectToSubspace(coarseSrc, in);
coarseFGMRES(coarseMdagMOp, coarseSrc, coarseSol);
_Aggregates.PromoteFromSubspace(coarseSol, out);
fineMdagMOp.Op(out, fineTmp);
fineTmp = in - fineTmp;
auto r = norm2(fineTmp);
auto residualAfterCoarseGridCorrection = std::sqrt(r / inputNorm);
fineFGMRES(fineSmootherMdagMOp, in, out);
fineMdagMOp.Op(out, fineTmp);
fineTmp = in - fineTmp;
r = norm2(fineTmp);
auto residualAfterPostSmoother = std::sqrt(r / inputNorm);
std::cout << GridLogMG << " Level " << level << ": Input norm = " << std::sqrt(inputNorm)
<< " Coarse residual = " << residualAfterCoarseGridCorrection << " Post-Smoother residual = " << residualAfterPostSmoother
<< std::endl;
}
void runChecks() {
auto tolerance = 1e-13; // TODO: this obviously depends on the precision we use, current value is for double
auto coarseLevel = level - 1;
std::vector<FineVector> fineTmps(2, _LevelInfo.Grids[level]);
std::vector<CoarseVector> coarseTmps(4, _LevelInfo.Grids[level - 1]);
MdagMLinearOperator<FineMatrix, FineVector> fineMdagMOp(_FineMatrix);
MdagMLinearOperator<CoarseMatrix, CoarseVector> coarseMdagMOp(_CoarseMatrix);
std::cout << GridLogMG << " Level " << level << ": **************************************************" << std::endl;
std::cout << GridLogMG << " Level " << level << ": MG correctness check: 0 == (1 - P R) v" << std::endl;
std::cout << GridLogMG << " Level " << level << ": **************************************************" << std::endl;
for(auto i = 0; i < _Aggregates.subspace.size(); ++i) {
_Aggregates.ProjectToSubspace(cTmps[0], _Aggregates.subspace[i]); // R v_i
_Aggregates.PromoteFromSubspace(cTmps[0], fTmps[0]); // P R v_i
_Aggregates.ProjectToSubspace(coarseTmps[0], _Aggregates.subspace[i]); // R v_i
_Aggregates.PromoteFromSubspace(coarseTmps[0], fineTmps[0]); // P R v_i
fTmps[1] = _Aggregates.subspace[i] - fTmps[0]; // v_i - P R v_i
auto deviation = std::sqrt(norm2(fTmps[1]) / norm2(_Aggregates.subspace[i]));
fineTmps[1] = _Aggregates.subspace[i] - fineTmps[0]; // v_i - P R v_i
auto deviation = std::sqrt(norm2(fineTmps[1]) / norm2(_Aggregates.subspace[i]));
std::cout << GridLogMessage << "Vector " << i << ": norm2(v_i) = " << norm2(_Aggregates.subspace[i])
<< " | norm2(R v_i) = " << norm2(cTmps[0]) << " | norm2(P R v_i) = " << norm2(fTmps[0])
std::cout << GridLogMG << " Level " << level << ": Vector " << i << ": norm2(v_i) = " << norm2(_Aggregates.subspace[i])
<< " | norm2(R v_i) = " << norm2(coarseTmps[0]) << " | norm2(P R v_i) = " << norm2(fineTmps[0])
<< " | relative deviation = " << deviation;
if(deviation > tolerance) {
@ -341,44 +414,20 @@ public:
}
}
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "MG correctness check: 0 == (1 - R P) v_c" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMG << " Level " << level << ": **************************************************" << std::endl;
std::cout << GridLogMG << " Level " << level << ": MG correctness check: 0 == (1 - R P) v_c" << std::endl;
std::cout << GridLogMG << " Level " << level << ": **************************************************" << std::endl;
random(cGrids.PRNGs[whichCoarseGrid], cTmps[0]);
random(_LevelInfo.PRNGs[coarseLevel], coarseTmps[0]);
_Aggregates.PromoteFromSubspace(cTmps[0], fTmps[0]); // P v_c
_Aggregates.ProjectToSubspace(cTmps[1], fTmps[0]); // R P v_c
_Aggregates.PromoteFromSubspace(coarseTmps[0], fineTmps[0]); // P v_c
_Aggregates.ProjectToSubspace(coarseTmps[1], fineTmps[0]); // R P v_c
cTmps[2] = cTmps[0] - cTmps[1]; // v_c - R P v_c
auto deviation = std::sqrt(norm2(cTmps[2]) / norm2(cTmps[0]));
coarseTmps[2] = coarseTmps[0] - coarseTmps[1]; // v_c - R P v_c
auto deviation = std::sqrt(norm2(coarseTmps[2]) / norm2(coarseTmps[0]));
std::cout << GridLogMessage << "norm2(v_c) = " << norm2(cTmps[0]) << " | norm2(R P v_c) = " << norm2(cTmps[1])
<< " | norm2(P v_c) = " << norm2(fTmps[0]) << " | relative deviation = " << deviation;
if(deviation > tolerance) {
std::cout << " > " << tolerance << " -> check failed" << std::endl;
// abort();
} else {
std::cout << " < " << tolerance << " -> check passed" << std::endl;
}
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "MG correctness check: 0 == (R D P - D_c) v_c" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
random(cGrids.PRNGs[whichCoarseGrid], cTmps[0]);
_Aggregates.PromoteFromSubspace(cTmps[0], fTmps[0]); // P v_c
_FineOperator.Op(fTmps[0], fTmps[1]); // D P v_c
_Aggregates.ProjectToSubspace(cTmps[1], fTmps[1]); // R D P v_c
MdagMOp.Op(cTmps[0], cTmps[2]); // D_c v_c
cTmps[3] = cTmps[1] - cTmps[2]; // R D P v_c - D_c v_c
deviation = std::sqrt(norm2(cTmps[3]) / norm2(cTmps[1]));
std::cout << GridLogMessage << "norm2(R D P v_c) = " << norm2(cTmps[1]) << " | norm2(D_c v_c) = " << norm2(cTmps[2])
std::cout << GridLogMG << " Level " << level << ": norm2(v_c) = " << norm2(coarseTmps[0])
<< " | norm2(R P v_c) = " << norm2(coarseTmps[1]) << " | norm2(P v_c) = " << norm2(fineTmps[0])
<< " | relative deviation = " << deviation;
if(deviation > tolerance) {
@ -388,54 +437,115 @@ public:
std::cout << " < " << tolerance << " -> check passed" << std::endl;
}
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "MG correctness check: 0 == |(Im(v_c^dag D_c^dag D_c v_c)|" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMG << " Level " << level << ": **************************************************" << std::endl;
std::cout << GridLogMG << " Level " << level << ": MG correctness check: 0 == (R D P - D_c) v_c" << std::endl;
std::cout << GridLogMG << " Level " << level << ": **************************************************" << std::endl;
random(cGrids.PRNGs[whichCoarseGrid], cTmps[0]);
random(_LevelInfo.PRNGs[coarseLevel], coarseTmps[0]);
MdagMOp.Op(cTmps[0], cTmps[1]); // D_c v_c
MdagMOp.AdjOp(cTmps[1], cTmps[2]); // D_c^dag D_c v_c
_Aggregates.PromoteFromSubspace(coarseTmps[0], fineTmps[0]); // P v_c
fineMdagMOp.Op(fineTmps[0], fineTmps[1]); // D P v_c
_Aggregates.ProjectToSubspace(coarseTmps[1], fineTmps[1]); // R D P v_c
auto dot = innerProduct(cTmps[0], cTmps[2]); //v_c^dag D_c^dag D_c v_c
coarseMdagMOp.Op(coarseTmps[0], coarseTmps[2]); // D_c v_c
coarseTmps[3] = coarseTmps[1] - coarseTmps[2]; // R D P v_c - D_c v_c
deviation = std::sqrt(norm2(coarseTmps[3]) / norm2(coarseTmps[1]));
std::cout << GridLogMG << " Level " << level << ": norm2(R D P v_c) = " << norm2(coarseTmps[1])
<< " | norm2(D_c v_c) = " << norm2(coarseTmps[2]) << " | relative deviation = " << deviation;
if(deviation > tolerance) {
std::cout << " > " << tolerance << " -> check failed" << std::endl;
// abort();
} else {
std::cout << " < " << tolerance << " -> check passed" << std::endl;
}
std::cout << GridLogMG << " Level " << level << ": **************************************************" << std::endl;
std::cout << GridLogMG << " Level " << level << ": MG correctness check: 0 == |(Im(v_c^dag D_c^dag D_c v_c)|" << std::endl;
std::cout << GridLogMG << " Level " << level << ": **************************************************" << std::endl;
random(_LevelInfo.PRNGs[coarseLevel], coarseTmps[0]);
coarseMdagMOp.Op(coarseTmps[0], coarseTmps[1]); // D_c v_c
coarseMdagMOp.AdjOp(coarseTmps[1], coarseTmps[2]); // D_c^dag D_c v_c
auto dot = innerProduct(coarseTmps[0], coarseTmps[2]); //v_c^dag D_c^dag D_c v_c
deviation = abs(imag(dot)) / abs(real(dot));
std::cout << GridLogMessage << "Re(v_c^dag D_c^dag D_c v_c) = " << real(dot) << " | Im(v_c^dag D_c^dag D_c v_c) = " << imag(dot)
<< " | relative deviation = " << deviation;
std::cout << GridLogMG << " Level " << level << ": Re(v_c^dag D_c^dag D_c v_c) = " << real(dot)
<< " | Im(v_c^dag D_c^dag D_c v_c) = " << imag(dot) << " | relative deviation = " << deviation;
if(deviation > tolerance) {
std::cout << " > " << tolerance << " -> check failed" << std::endl;
// abort();
} else {
std::cout << " < " << tolerance << " -> check passed" << std::endl;
std::cout << " < " << tolerance << " -> check passed"
<< std::endl; // TODO: this check will work only when I got Mdag in CoarsenedMatrix to work
}
_NextPreconditionerLevel->runChecks();
}
};
// Specialize the coarsest level, this corresponds to counting downwards with level: coarsest = 0, finest = N
template<class Fobj, class CoarseScalar, int nCoarseSpins, int nbasis, class Matrix>
class MultiGridPreconditioner<Fobj, CoarseScalar, nCoarseSpins, nbasis, 0, Matrix> : public LinearFunction<Lattice<Fobj>> {
public:
/////////////////////////////////////////////
// Type Definitions
/////////////////////////////////////////////
typedef Matrix FineMatrix;
typedef Lattice<Fobj> FineVector;
/////////////////////////////////////////////
// Member Data
/////////////////////////////////////////////
LevelInfo & _LevelInfo;
FineMatrix &_FineMatrix;
FineMatrix &_SmootherMatrix;
/////////////////////////////////////////////
// Member Functions
/////////////////////////////////////////////
MultiGridPreconditioner(LevelInfo &LvlInfo, FineMatrix &FineMat, FineMatrix &SmootherMat)
: _LevelInfo(LvlInfo), _FineMatrix(FineMat), _SmootherMatrix(SmootherMat) {}
void setup() {}
virtual void operator()(Lattice<Fobj> const &in, Lattice<Fobj> &out) {
TrivialPrecon<FineVector> fineTrivialPreconditioner;
FlexibleGeneralisedMinimalResidual<FineVector> fineFGMRES(1.0e-14, 1, fineTrivialPreconditioner, 1, false);
MdagMLinearOperator<FineMatrix, FineVector> fineMdagMOp(_FineMatrix);
fineFGMRES(fineMdagMOp, in, out);
}
void runChecks() {}
};
template<class Fobj, class CoarseScalar, int nCoarseSpins, int nbasis, class Matrix>
using FourLevelMGPreconditioner = MultiGridPreconditioner<Fobj, CoarseScalar, nCoarseSpins, nbasis, 4 - 1, Matrix>;
template<class Fobj, class CoarseScalar, int nCoarseSpins, int nbasis, class Matrix>
using ThreeLevelMGPreconditioner = MultiGridPreconditioner<Fobj, CoarseScalar, nCoarseSpins, nbasis, 3 - 1, Matrix>;
template<class Fobj, class CoarseScalar, int nCoarseSpins, int nbasis, class Matrix>
using TwoLevelMGPreconditioner = MultiGridPreconditioner<Fobj, CoarseScalar, nCoarseSpins, nbasis, 2 - 1, Matrix>;
template<class Fobj, class CoarseScalar, int nCoarseSpins, int nbasis, int nlevel, class Matrix>
using NLevelMGPreconditioner = MultiGridPreconditioner<Fobj, CoarseScalar, nCoarseSpins, nbasis, nlevel - 1, Matrix>;
int main(int argc, char **argv) {
Grid_init(&argc, &argv);
params.domainsize = 1;
params.coarsegrids = 1;
params.domaindecompose = 0;
params.order = 30;
params.Ls = 1;
params.mq = -0.5;
params.lo = 0.5;
params.hi = 70.0;
params.steps = 1;
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Params: " << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << params << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Set up some fine level stuff: " << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
GridCartesian *FGrid = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd, vComplex::Nsimd()), GridDefaultMpi());
GridRedBlackCartesian *FrbGrid = SpaceTimeGrid::makeFourDimRedBlackGrid(FGrid);
@ -451,312 +561,54 @@ int main(int argc, char **argv) {
LatticeGaugeField Umu(FGrid); SU3::HotConfiguration(fPRNG, Umu);
// clang-format on
RealD mass = params.mq;
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Set up some coarser levels stuff: " << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
const int nbasis = 20; // fix the number of test vector to the same
// number on every level for now
//////////////////////////////////////////
// toggle to run two/three level method
//////////////////////////////////////////
// two-level algorithm
std::vector<std::vector<int>> blockSizes({{2, 2, 2, 2}});
CoarseGrids<nbasis> coarseGrids(blockSizes, 1);
// // three-level algorithm
// std::vector<std::vector<int>> blockSizes({{2, 2, 2, 2}, {2, 2, 1, 1}});
// CoarseGrids<nbasis> coarseGrids(blockSizes, 2);
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Some typedefs" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
// typedefs for transition from fine to first coarsened grid
typedef vSpinColourVector FineSiteVector;
typedef vTComplex CoarseSiteScalar;
typedef Aggregation<FineSiteVector, CoarseSiteScalar, nbasis> Subspace;
typedef CoarsenedMatrix<FineSiteVector, CoarseSiteScalar, nbasis> CoarseOperator;
typedef CoarseOperator::CoarseVector CoarseVector;
typedef CoarseOperator::siteVector CoarseSiteVector;
typedef TestVectorAnalyzer<LatticeFermion, nbasis> FineTVA;
typedef MultiGridPreconditioner<FineSiteVector, CoarseSiteScalar, nbasis, WilsonFermionR> FineMGPreconditioner;
typedef TrivialPrecon<LatticeFermion> FineTrivialPreconditioner;
// typedefs for transition from a coarse to the next coarser grid (some defs remain the same)
typedef Aggregation<CoarseSiteVector, CoarseSiteScalar, nbasis> SubSubSpace;
typedef CoarsenedMatrix<CoarseSiteVector, CoarseSiteScalar, nbasis> CoarseCoarseOperator;
typedef CoarseCoarseOperator::CoarseVector CoarseCoarseVector;
typedef CoarseCoarseOperator::siteVector CoarseCoarseSiteVector;
typedef TestVectorAnalyzer<CoarseVector, nbasis> CoarseTVA;
typedef MultiGridPreconditioner<CoarseSiteVector, CoarseSiteScalar, nbasis, CoarseOperator> CoarseMGPreconditioner;
typedef TrivialPrecon<CoarseVector> CoarseTrivialPreconditioner;
static_assert(std::is_same<CoarseVector, CoarseCoarseVector>::value, "CoarseVector and CoarseCoarseVector must be of the same type");
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Building the wilson operator on the fine grid" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
RealD mass = 0.5;
const int nbasis = 20;
WilsonFermionR Dw(Umu, *FGrid, *FrbGrid, mass);
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Setting up linear operators" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
// mgParams.blockSizes = {{2, 2, 2, 2}, {2, 2, 1, 1}, {1, 1, 2, 1}};
// mgParams.blockSizes = {{2, 2, 2, 2}, {2, 2, 1, 1}};
mgParams.blockSizes = {{2, 2, 2, 2}};
mgParams.nLevels = mgParams.blockSizes.size() + 1;
std::cout << mgParams << std::endl;
LevelInfo levelInfo(FGrid, mgParams);
MdagMLinearOperator<WilsonFermionR, LatticeFermion> FineMdagMOp(Dw);
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Calling Aggregation class to build subspaces" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
TrivialPrecon<LatticeFermion> TrivialPrecon;
TwoLevelMGPreconditioner<vSpinColourVector, vTComplex, 1, nbasis, WilsonFermionR> TwoLevelMGPrecon(levelInfo, Dw, Dw);
// ThreeLevelMGPreconditioner<vSpinColourVector, vTComplex, 1, nbasis, WilsonFermionR> ThreeLevelMGPrecon(levelInfo, Dw, Dw);
// FourLevelMGPreconditioner<vSpinColourVector, vTComplex, 1, nbasis, WilsonFermionR> FourLevelMGPrecon(levelInfo, Dw, Dw);
// NLevelMGPreconditioner<vSpinColourVector, vTComplex, 1, nbasis, 4, WilsonFermionR> FourLevelMGPrecon(levelInfo, Dw, Dw);
Subspace FineAggregates(coarseGrids.Grids[0], FGrid, 0);
TwoLevelMGPrecon.setup();
TwoLevelMGPrecon.runChecks();
assert((nbasis & 0x1) == 0);
int nb = nbasis / 2;
std::cout << GridLogMessage << " nbasis/2 = " << nb << std::endl;
// ThreeLevelMGPrecon.setup();
// ThreeLevelMGPrecon.runChecks();
FineAggregates.CreateSubspace(fPRNG, FineMdagMOp /*, nb */); // Don't specify nb to see the orthogonalization check
// FourLevelMGPrecon.setup();
// FourLevelMGPrecon.runChecks();
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Test vector analysis after initial creation of subspace" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
FineTVA fineTVA;
fineTVA(FineMdagMOp, FineAggregates.subspace);
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Projecting subspace to definite chirality" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
for(int n = 0; n < nb; n++) {
FineAggregates.subspace[n + nb] = g5 * FineAggregates.subspace[n];
}
auto coarseSites = 1;
for(auto const &elem : coarseGrids.LattSizes[0]) coarseSites *= elem;
std::cout << GridLogMessage << "Norms of MG test vectors after chiral projection (coarse sites = " << coarseSites << ")" << std::endl;
for(int n = 0; n < nbasis; n++) {
std::cout << GridLogMessage << "vec[" << n << "] = " << norm2(FineAggregates.subspace[n]) << std::endl;
}
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Building coarse representation of Dirac operator" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
CoarseOperator Dc(*coarseGrids.Grids[0]);
Dc.CoarsenOperator(FGrid, FineMdagMOp, FineAggregates);
MdagMLinearOperator<CoarseOperator, CoarseVector> CoarseMdagMOp(Dc);
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Test vector analysis after construction of coarse Dirac operator" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
fineTVA(FineMdagMOp, FineAggregates.subspace);
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Testing the linear operators" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
// clang-format off
testLinearOperator(FineMdagMOp, FGrid, "FineMdagMOp"); std::cout << GridLogMessage << std::endl;
testLinearOperator(CoarseMdagMOp, coarseGrids.Grids[0], "CoarseMdagMOp"); std::cout << GridLogMessage << std::endl;
// clang-format on
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Building coarse vectors" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
CoarseVector coarseSource(coarseGrids.Grids[0]);
CoarseVector coarseResult(coarseGrids.Grids[0]);
gaussian(coarseGrids.PRNGs[0], coarseSource);
coarseResult = zero;
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Building some coarse space solvers" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::vector<std::unique_ptr<OperatorFunction<CoarseVector>>> dummyCoarseSolvers;
dummyCoarseSolvers.emplace_back(new GeneralisedMinimalResidual<CoarseVector>(5.0e-2, 100, 8, false));
dummyCoarseSolvers.emplace_back(new MinimalResidual<CoarseVector>(5.0e-2, 100, 0.8, false));
dummyCoarseSolvers.emplace_back(new ConjugateGradient<CoarseVector>(5.0e-2, 100, false));
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Testing some coarse space solvers" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "checking norm of coarse src " << norm2(coarseSource) << std::endl;
for(auto const &solver : dummyCoarseSolvers) {
coarseResult = zero;
(*solver)(CoarseMdagMOp, coarseSource, coarseResult);
}
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Building a multigrid preconditioner" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
FineMGPreconditioner FineMGPrecon(FineAggregates, Dc, FineMdagMOp, Dw, FineMdagMOp, Dw);
FineTrivialPreconditioner FineSimplePrecon;
FineMGPrecon.runChecks(coarseGrids, 0);
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Building krylov subspace solvers w/ & w/o MG Preconditioner" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
// NLevelMGPrecon.setup();
// NLevelMGPrecon.runChecks();
std::vector<std::unique_ptr<OperatorFunction<LatticeFermion>>> solvers;
solvers.emplace_back(new FlexibleGeneralisedMinimalResidual<LatticeFermion>(1.0e-12, 4000000, FineSimplePrecon, 25, false));
solvers.emplace_back(new FlexibleGeneralisedMinimalResidual<LatticeFermion>(1.0e-12, 100, FineMGPrecon, 25, false));
solvers.emplace_back(new PrecGeneralisedConjugateResidual<LatticeFermion>(1.0e-12, 4000000, FineSimplePrecon, 25, 25));
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Testing the (un)?preconditioned solvers" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
solvers.emplace_back(new FlexibleGeneralisedMinimalResidual<LatticeFermion>(1.0e-12, 50000, TrivialPrecon, 1000, false));
solvers.emplace_back(new FlexibleGeneralisedMinimalResidual<LatticeFermion>(1.0e-12, 50000, TwoLevelMGPrecon, 1000, false));
// solvers.emplace_back(new FlexibleGeneralisedMinimalResidual<LatticeFermion>(1.0e-12, 50000, ThreeLevelMGPrecon, 1000, false));
// solvers.emplace_back(new FlexibleGeneralisedMinimalResidual<LatticeFermion>(1.0e-12, 50000, FourLevelMGPrecon, 1000, false));
// solvers.emplace_back(new FlexibleGeneralisedMinimalResidual<LatticeFermion>(1.0e-12, 50000, NLevelMGPrecon, 1000, false));
for(auto const &solver : solvers) {
std::cout << GridLogMessage << "checking norm of fine src " << norm2(src) << std::endl;
std::cout << "Starting with a new solver" << std::endl;
result = zero;
(*solver)(FineMdagMOp, src, result);
std::cout << std::endl;
}
#if 0
if(coarseGrids.LattSizes.size() == 2) {
std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "Some testing for construction of a second coarse level" << std::endl;
std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Calling Aggregation class to build subspaces" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
SubSubSpace CoarseAggregates(coarseGrids.Grids[1], coarseGrids.Grids[0], 0);
CoarseAggregates.CreateSubspace(coarseGrids.PRNGs[0], CoarseMdagMOp);
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Test vector analysis after initial creation of subspace" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// // this doesn't work because this function applies g5 to a vector, which
// // doesn't work for coarse vectors atm -> FIXME
// CoarseTVA coarseTVA;
// coarseTVA(CoarseMdagMOp, CoarseAggregates.subspace);
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Projecting subspace to definite chirality" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// // cannot apply g5 to coarse vectors atm -> FIXME
// for(int n=0;n<nb;n++){
// CoarseAggregates.subspace[n+nb] = g5 * CoarseAggregates.subspace[n];
// std::cout<<GridLogMessage<<n<<" subspace "<<norm2(CoarseAggregates.subspace[n+nb])<<" "<<norm2(CoarseAggregates.subspace[n]) <<std::endl;
// }
auto coarseCoarseSites = 1;
for(auto const &elem : coarseGrids.LattSizes[1]) coarseCoarseSites *= elem;
std::cout << GridLogMessage << "Norms of MG test vectors after chiral projection (coarse coarse sites = " << coarseCoarseSites << ")"
<< std::endl;
for(int n = 0; n < nbasis; n++) {
std::cout << GridLogMessage << "vec[" << n << "] = " << norm2(CoarseAggregates.subspace[n]) << std::endl;
}
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Building coarse coarse representation of Dirac operator" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
CoarseCoarseOperator Dcc(*coarseGrids.Grids[1]);
Dcc.CoarsenOperator(coarseGrids.Grids[0], CoarseMdagMOp, CoarseAggregates);
MdagMLinearOperator<CoarseCoarseOperator, CoarseCoarseVector> CoarseCoarseMdagMOp(Dcc);
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Test vector analysis after construction of coarse Dirac operator" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// // this doesn't work because this function applies g5 to a vector, which
// // doesn't work for coarse vectors atm -> FIXME
// coarseTVA(CoarseMdagMOp, CoarseAggregates.subspace);
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Testing the linear operators" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// clang-format off
testLinearOperator(CoarseMdagMOp, coarseGrids.Grids[0], "CoarseMdagMOp");
testLinearOperator(CoarseCoarseMdagMOp, coarseGrids.Grids[1], "CoarseCoarseMdagMOp");
// clang-format on
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Building coarse coarse vectors" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
CoarseCoarseVector coarseCoarseSource(coarseGrids.Grids[1]);
CoarseCoarseVector coarseCoarseResult(coarseGrids.Grids[1]);
gaussian(coarseGrids.PRNGs[1], coarseCoarseSource);
coarseCoarseResult = zero;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Building some coarse space solvers" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
std::vector<std::unique_ptr<OperatorFunction<CoarseCoarseVector>>> dummyCoarseCoarseSolvers;
dummyCoarseCoarseSolvers.emplace_back(new GeneralisedMinimalResidual<CoarseCoarseVector>(5.0e-2, 100, 8, false));
dummyCoarseCoarseSolvers.emplace_back(new MinimalResidual<CoarseCoarseVector>(5.0e-2, 100, 0.8, false));
dummyCoarseCoarseSolvers.emplace_back(new ConjugateGradient<CoarseCoarseVector>(5.0e-2, 100, false));
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Testing some coarse coarse space solvers" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
std::cout << GridLogMessage << "checking norm of coarse coarse src " << norm2(coarseCoarseSource) << std::endl;
for(auto const &solver : dummyCoarseCoarseSolvers) {
coarseCoarseResult = zero;
(*solver)(CoarseCoarseMdagMOp, coarseCoarseSource, coarseCoarseResult);
}
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Building a multigrid preconditioner" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
CoarseMGPreconditioner CoarseMGPrecon(CoarseAggregates, Dcc, CoarseMdagMOp, Dc, CoarseMdagMOp, Dc);
CoarseTrivialPreconditioner CoarseSimplePrecon;
CoarseMGPrecon.runChecks(coarseGrids, 1);
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Building krylov subspace solvers w/ & w/o MG Preconditioner" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
std::vector<std::unique_ptr<OperatorFunction<CoarseVector>>> solvers;
solvers.emplace_back(new FlexibleGeneralisedMinimalResidual<CoarseVector>(1.0e-12, 4000000, CoarseSimplePrecon, 25, false));
solvers.emplace_back(new FlexibleGeneralisedMinimalResidual<CoarseVector>(1.0e-12, 100, CoarseMGPrecon, 25, false));
solvers.emplace_back(new PrecGeneralisedConjugateResidual<CoarseVector>(1.0e-12, 4000000, CoarseSimplePrecon, 25, 25));
// std::cout << GridLogMessage << "**************************************************" << std::endl;
// std::cout << GridLogMessage << "Testing the (un)?preconditioned solvers" << std::endl;
// std::cout << GridLogMessage << "**************************************************" << std::endl;
for(auto const &solver : solvers) {
std::cout << GridLogMessage << "checking norm of fine src " << norm2(coarseSource) << std::endl;
coarseResult = zero;
(*solver)(CoarseMdagMOp, coarseSource, coarseResult);
std::cout << std::endl;
}
}
#endif
Grid_finalize();
}