1
0
mirror of https://github.com/paboyle/Grid.git synced 2024-09-20 09:15:38 +01:00
Grid/tests/solver/Test_dwf_hdcr.cc

675 lines
25 KiB
C++
Raw Normal View History

/*************************************************************************************
Grid physics library, www.github.com/paboyle/Grid
Source file: ./tests/Test_dwf_hdcr.cc
Copyright (C) 2015
Author: Antonin Portelli <antonin.portelli@me.com>
Author: Peter Boyle <paboyle@ph.ed.ac.uk>
Author: paboyle <paboyle@ph.ed.ac.uk>
This program is free software; you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation; either version 2 of the License, or
(at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License along
with this program; if not, write to the Free Software Foundation, Inc.,
51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
See the full license in the file "LICENSE" in the top level distribution directory
*************************************************************************************/
/* END LEGAL */
2016-07-07 22:31:07 +01:00
#include <Grid/Grid.h>
#include <Grid/algorithms/iterative/PrecGeneralisedConjugateResidual.h>
2015-07-01 13:04:02 +01:00
//#include <algorithms/iterative/PrecConjugateResidual.h>
2015-06-22 12:49:44 +01:00
using namespace std;
using namespace Grid;
2015-11-29 10:59:11 +00:00
class myclass: Serializable {
public:
2015-12-08 13:54:00 +00:00
GRID_SERIALIZABLE_CLASS_MEMBERS(myclass,
int, domaindecompose,
int, domainsize,
int, order,
2016-03-16 09:31:04 +00:00
int, Ls,
double, mq,
double, lo,
double, hi,
int, steps);
myclass(){};
};
RealD InverseApproximation(RealD x){
return 1.0/x;
}
2019-12-28 15:32:35 +00:00
template<class Fobj,class CComplex,int nbasis, class Matrix, class Guesser>
2015-06-22 12:49:44 +01:00
class MultiGridPreconditioner : public LinearFunction< Lattice<Fobj> > {
public:
typedef Aggregation<Fobj,CComplex,nbasis> Aggregates;
typedef CoarsenedMatrix<Fobj,CComplex,nbasis> CoarseOperator;
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;
Aggregates & _Aggregates;
CoarseOperator & _CoarseOperator;
Matrix & _FineMatrix;
2015-06-22 12:49:44 +01:00
FineOperator & _FineOperator;
Matrix & _SmootherMatrix;
FineOperator & _SmootherOperator;
2019-12-28 15:32:35 +00:00
Guesser & _Guess;
double cheby_hi;
double cheby_lo;
int cheby_ord;
myclass _params;
2015-06-22 12:49:44 +01:00
// Constructor
MultiGridPreconditioner(Aggregates &Agg, CoarseOperator &Coarse,
FineOperator &Fine,Matrix &FineMatrix,
2019-12-28 15:32:35 +00:00
FineOperator &Smooth,Matrix &SmootherMatrix,
Guesser &Guess_,
myclass params_)
2015-06-22 12:49:44 +01:00
: _Aggregates(Agg),
_CoarseOperator(Coarse),
_FineOperator(Fine),
_FineMatrix(FineMatrix),
_SmootherOperator(Smooth),
2019-12-28 15:32:35 +00:00
_SmootherMatrix(SmootherMatrix),
_Guess(Guess_),
_params(params_)
2015-06-22 12:49:44 +01:00
{
}
void PowerMethod(const FineField &in) {
2018-01-27 00:04:12 +00:00
FineField p1(in.Grid());
FineField p2(in.Grid());
MdagMLinearOperator<Matrix,FineField> fMdagMOp(_FineMatrix);
p1=in;
2019-12-28 15:32:35 +00:00
for(int i=0;i<50;i++){
RealD absp1=std::sqrt(norm2(p1));
fMdagMOp.HermOp(p1,p2);// this is the G5 herm bit
// _FineOperator.Op(p1,p2);// this is the G5 herm bit
RealD absp2=std::sqrt(norm2(p2));
if(i%10==9)
std::cout<<GridLogMessage << "Power method on mdagm "<<i<<" " << absp2/absp1<<std::endl;
p1=p2*(1.0/std::sqrt(absp2));
}
}
2019-12-28 15:32:35 +00:00
void operator()(const FineField &in, FineField & out ) {
2019-12-09 07:55:45 +00:00
operatorCheby(in,out);
2019-12-28 15:32:35 +00:00
//operatorADEF2(in,out);
}
////////////////////////////////////////////////////////////////////////
// ADEF2: [PTM+Q] in = [1 - Q A] M in + Q in = Min + Q [ in -A Min]
// ADEF1: [MP+Q ] in =M [1 - A Q] in + Q in
////////////////////////////////////////////////////////////////////////
#if 1
void operatorADEF2(const FineField &in, FineField & out) {
CoarseVector Csrc(_CoarseOperator.Grid());
CoarseVector Ctmp(_CoarseOperator.Grid());
CoarseVector Csol(_CoarseOperator.Grid());
2019-12-28 15:32:35 +00:00
ConjugateGradient<CoarseVector> CG(1.0e-3,1000,false);
ConjugateGradient<FineField> fCG(1.0e-3,15,false);
HermitianLinearOperator<CoarseOperator,CoarseVector> HermOp(_CoarseOperator);
MdagMLinearOperator<CoarseOperator,CoarseVector> MdagMOp(_CoarseOperator);
MdagMLinearOperator<Matrix,FineField> fMdagMOp(_FineMatrix);
2018-01-27 00:04:12 +00:00
FineField tmp(in.Grid());
FineField res(in.Grid());
FineField Min(in.Grid());
// Monitor completeness of low mode space
_Aggregates.ProjectToSubspace (Csrc,in);
_Aggregates.PromoteFromSubspace(Csrc,out);
std::cout<<GridLogMessage<<"Coarse Grid Preconditioner\nCompleteness in: "<<std::sqrt(norm2(out)/norm2(in))<<std::endl;
// [PTM+Q] in = [1 - Q A] M in + Q in = Min + Q [ in -A Min]
_FineOperator.Op(in,tmp);// this is the G5 herm bit
fCG(fMdagMOp,tmp,Min); // solves MdagM = g5 M g5M
// Monitor completeness of low mode space
_Aggregates.ProjectToSubspace (Csrc,Min);
_Aggregates.PromoteFromSubspace(Csrc,out);
std::cout<<GridLogMessage<<"Completeness Min: "<<std::sqrt(norm2(out)/norm2(Min))<<std::endl;
_FineOperator.Op(Min,tmp);
tmp = in - tmp; // in - A Min
_Aggregates.ProjectToSubspace (Csrc,tmp);
HermOp.AdjOp(Csrc,Ctmp);// Normal equations
2019-12-28 15:32:35 +00:00
Csol=Zero();
CG(MdagMOp,Ctmp,Csol);
HermOp.Op(Csol,Ctmp);
Ctmp=Ctmp-Csrc;
std::cout<<GridLogMessage<<"coarse space true residual "<<std::sqrt(norm2(Ctmp)/norm2(Csrc))<<std::endl;
_Aggregates.PromoteFromSubspace(Csol,out);
_FineOperator.Op(out,res);
res=res-tmp;
std::cout<<GridLogMessage<<"promoted sol residual "<<std::sqrt(norm2(res)/norm2(tmp))<<std::endl;
_Aggregates.ProjectToSubspace (Csrc,res);
std::cout<<GridLogMessage<<"coarse space proj of residual "<<norm2(Csrc)<<std::endl;
out = out+Min; // additive coarse space correction
// out = Min; // no additive coarse space correction
_FineOperator.Op(out,tmp);
tmp=tmp-in; // tmp is new residual
std::cout<<GridLogMessage<< " Preconditioner in " << norm2(in)<<std::endl;
std::cout<<GridLogMessage<< " Preconditioner out " << norm2(out)<<std::endl;
std::cout<<GridLogMessage<<"preconditioner thinks residual is "<<std::sqrt(norm2(tmp)/norm2(in))<<std::endl;
}
#endif
// ADEF1: [MP+Q ] in =M [1 - A Q] in + Q in
#if 1
void operatorADEF1(const FineField &in, FineField & out) {
CoarseVector Csrc(_CoarseOperator.Grid());
CoarseVector Ctmp(_CoarseOperator.Grid());
2018-01-27 23:46:02 +00:00
CoarseVector Csol(_CoarseOperator.Grid()); Csol=Zero();
ConjugateGradient<CoarseVector> CG(1.0e-10,100000);
ConjugateGradient<FineField> fCG(1.0e-3,1000);
HermitianLinearOperator<CoarseOperator,CoarseVector> HermOp(_CoarseOperator);
MdagMLinearOperator<CoarseOperator,CoarseVector> MdagMOp(_CoarseOperator);
ShiftedMdagMLinearOperator<Matrix,FineField> fMdagMOp(_FineMatrix,0.1);
2018-01-27 00:04:12 +00:00
FineField tmp(in.Grid());
FineField res(in.Grid());
FineField Qin(in.Grid());
// Monitor completeness of low mode space
// _Aggregates.ProjectToSubspace (Csrc,in);
// _Aggregates.PromoteFromSubspace(Csrc,out);
// std::cout<<GridLogMessage<<"Coarse Grid Preconditioner\nCompleteness in: "<<std::sqrt(norm2(out)/norm2(in))<<std::endl;
_Aggregates.ProjectToSubspace (Csrc,in);
HermOp.AdjOp(Csrc,Ctmp);// Normal equations
CG(MdagMOp,Ctmp,Csol);
_Aggregates.PromoteFromSubspace(Csol,Qin);
// Qin=0;
_FineOperator.Op(Qin,tmp);// A Q in
tmp = in - tmp; // in - A Q in
_FineOperator.Op(tmp,res);// this is the G5 herm bit
fCG(fMdagMOp,res,out); // solves MdagM = g5 M g5M
out = out + Qin;
_FineOperator.Op(out,tmp);
tmp=tmp-in; // tmp is new residual
std::cout<<GridLogMessage<<"preconditioner thinks residual is "<<std::sqrt(norm2(tmp)/norm2(in))<<std::endl;
2015-06-22 12:49:44 +01:00
}
#endif
void SmootherTest (const FineField & in){
2018-01-27 00:04:12 +00:00
FineField vec1(in.Grid());
FineField vec2(in.Grid());
RealD lo[3] = { 0.5, 1.0, 2.0};
// MdagMLinearOperator<Matrix,FineField> fMdagMOp(_FineMatrix);
ShiftedMdagMLinearOperator<Matrix,FineField> fMdagMOp(_SmootherMatrix,0.0);
RealD Ni,r;
Ni = norm2(in);
for(int ilo=0;ilo<3;ilo++){
for(int ord=5;ord<50;ord*=2){
2019-12-09 07:55:45 +00:00
std::cout << " lo "<<lo[ilo]<<" order "<<ord<<std::endl;
_SmootherOperator.AdjOp(in,vec1);
Chebyshev<FineField> Cheby (lo[ilo],70.0,ord,InverseApproximation);
Cheby(fMdagMOp,vec1,vec2); // solves MdagM = g5 M g5M
_FineOperator.Op(vec2,vec1);// this is the G5 herm bit
vec1 = in - vec1; // tmp = in - A Min
r=norm2(vec1);
std::cout<<GridLogMessage << "Smoother resid "<<std::sqrt(r/Ni)<<std::endl;
}
}
}
void operatorCheby(const FineField &in, FineField & out) {
CoarseVector Csrc(_CoarseOperator.Grid());
CoarseVector Ctmp(_CoarseOperator.Grid());
2019-12-28 15:32:35 +00:00
CoarseVector Csol(_CoarseOperator.Grid());
ConjugateGradient<CoarseVector> CG(3.0e-2,100000);
HermitianLinearOperator<CoarseOperator,CoarseVector> HermOp(_CoarseOperator);
MdagMLinearOperator<CoarseOperator,CoarseVector> MdagMOp(_CoarseOperator);
// MdagMLinearOperator<Matrix,FineField> fMdagMOp(_FineMatrix);
ShiftedMdagMLinearOperator<Matrix,FineField> fMdagMOp(_SmootherMatrix,0.0);
2018-01-27 00:04:12 +00:00
FineField vec1(in.Grid());
FineField vec2(in.Grid());
2019-12-28 15:32:35 +00:00
Chebyshev<FineField> Cheby (_params.lo,_params.hi,_params.order,InverseApproximation);
Chebyshev<FineField> ChebyAccu(_params.lo,_params.hi,_params.order,InverseApproximation);
// _Aggregates.ProjectToSubspace (Csrc,in);
// _Aggregates.PromoteFromSubspace(Csrc,out);
// std::cout<<GridLogMessage<<"Completeness: "<<std::sqrt(norm2(out)/norm2(in))<<std::endl;
// ofstream fout("smoother");
// Cheby.csv(fout);
// V11 multigrid.
// Use a fixed chebyshev and hope hermiticity helps.
// To make a working smoother for indefinite operator
// must multiply by "Mdag" (ouch loses all low mode content)
// and apply to poly approx of (mdagm)^-1.
// so that we end up with an odd polynomial.
RealD Ni = norm2(in);
std::cout<<GridLogMessage << "Smoother calling Cheby" <<std::endl;
_SmootherOperator.AdjOp(in,vec1);// this is the G5 herm bit
ChebyAccu(fMdagMOp,vec1,out); // solves MdagM = g5 M g5M
std::cout<<GridLogMessage << "Smoother called Cheby" <<std::endl;
// Update with residual for out
_FineOperator.Op(out,vec1);// this is the G5 herm bit
vec1 = in - vec1; // tmp = in - A Min
RealD r = norm2(vec1);
std::cout<<GridLogMessage << "Smoother resid "<<std::sqrt(r/Ni)<< " " << r << " " << Ni <<std::endl;
std::cout<<GridLogMessage << "ProjectToSubspace" <<std::endl;
_Aggregates.ProjectToSubspace (Csrc,vec1);
std::cout<<GridLogMessage << "ProjectToSubspaceDone" <<std::endl;
2019-12-09 07:55:45 +00:00
HermOp.AdjOp(Csrc,Ctmp);// Normal equations // This appears to be zero.
2019-12-28 15:32:35 +00:00
_Guess(Ctmp,Csol);
CG(MdagMOp,Ctmp,Csol);
std::cout<<GridLogMessage << "PromoteFromSubspace" <<std::endl;
_Aggregates.PromoteFromSubspace(Csol,vec1); // Ass^{-1} [in - A Min]_s
2019-12-09 07:55:45 +00:00
// Q = Q[in - A Min]
std::cout<<GridLogMessage << "PromoteFromSubspaceDone" <<std::endl;
out = out+vec1;
// Three preconditioner smoothing -- hermitian if C3 = C1
// Recompute error
_FineOperator.Op(out,vec1);// this is the G5 herm bit
std::cout<<GridLogMessage << "FineOp" <<std::endl;
vec1 = in - vec1; // tmp = in - A Min
r=norm2(vec1);
std::cout<<GridLogMessage << "Coarse resid "<<std::sqrt(r/Ni)<<std::endl;
// Reapply smoother
std::cout<<GridLogMessage << "Smoother calling Cheby" <<std::endl;
_SmootherOperator.Op(vec1,vec2); // this is the G5 herm bit
ChebyAccu(fMdagMOp,vec2,vec1); // solves MdagM = g5 M g5M
std::cout<<GridLogMessage << "Smoother called Cheby" <<std::endl;
out =out+vec1;
vec1 = in - vec1; // tmp = in - A Min
r=norm2(vec1);
std::cout<<GridLogMessage << "Smoother resid "<<std::sqrt(r/Ni)<<std::endl;
}
2015-06-22 12:49:44 +01:00
};
int main (int argc, char ** argv)
{
Grid_init(&argc,&argv);
2019-12-28 15:32:35 +00:00
myclass params;
myclass cparams;
2015-11-29 10:59:11 +00:00
XmlReader RD("params.xml");
read(RD,"params",params);
std::cout<<"Params: Order "<<params.order<<"["<<params.lo<<","<<params.hi<<"]"<< " steps "<<params.steps<<std::endl;
2016-03-16 09:31:04 +00:00
const int Ls=params.Ls;
2015-06-22 12:49:44 +01:00
GridCartesian * UGrid = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd,vComplex::Nsimd()),GridDefaultMpi());
2015-06-22 12:49:44 +01:00
GridRedBlackCartesian * UrbGrid = SpaceTimeGrid::makeFourDimRedBlackGrid(UGrid);
GridCartesian * FGrid = SpaceTimeGrid::makeFiveDimGrid(Ls,UGrid);
GridRedBlackCartesian * FrbGrid = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls,UGrid);
///////////////////////////////////////////////////
// Construct a coarsened grid; utility for this?
///////////////////////////////////////////////////
std::vector<int> block ({2,2,2,2});
const int nbasis= 32;
2018-02-24 22:18:33 +00:00
auto clatt = GridDefaultLatt();
2015-06-22 12:49:44 +01:00
for(int d=0;d<clatt.size();d++){
clatt[d] = clatt[d]/block[d];
2015-06-22 12:49:44 +01:00
}
2019-12-28 15:32:35 +00:00
GridCartesian *Coarse4d = SpaceTimeGrid::makeFourDimGrid(clatt, GridDefaultSimd(Nd,vComplex::Nsimd()),GridDefaultMpi());;
2015-06-22 12:49:44 +01:00
GridCartesian *Coarse5d = SpaceTimeGrid::makeFiveDimGrid(1,Coarse4d);
std::vector<int> seeds4({1,2,3,4});
std::vector<int> seeds5({5,6,7,8});
std::vector<int> cseeds({5,6,7,8});
GridParallelRNG RNG5(FGrid); RNG5.SeedFixedIntegers(seeds5);
GridParallelRNG RNG4(UGrid); RNG4.SeedFixedIntegers(seeds4);
GridParallelRNG CRNG(Coarse5d);CRNG.SeedFixedIntegers(cseeds);
Gamma g5(Gamma::Algebra::Gamma5);
LatticeFermion src(FGrid); gaussian(RNG5,src);// src=src+g5*src;
2018-01-27 23:46:02 +00:00
LatticeFermion result(FGrid); result=Zero();
LatticeFermion ref(FGrid); ref=Zero();
2015-06-22 12:49:44 +01:00
LatticeFermion tmp(FGrid);
LatticeFermion err(FGrid);
LatticeGaugeField Umu(UGrid);
LatticeGaugeField UmuDD(UGrid);
LatticeColourMatrix U(UGrid);
LatticeColourMatrix zz(UGrid);
2015-06-22 12:49:44 +01:00
2017-06-19 22:54:18 +01:00
FieldMetaData header;
std::string file("./ckpoint_lat.4000");
NerscIO::readConfiguration(Umu,header,file);
2015-06-22 12:49:44 +01:00
if ( params.domaindecompose ) {
Lattice<iScalar<vInteger> > coor(UGrid);
2018-01-27 23:46:02 +00:00
zz=Zero();
for(int mu=0;mu<Nd;mu++){
LatticeCoordinate(coor,mu);
U = PeekIndex<LorentzIndex>(Umu,mu);
U = where(mod(coor,params.domainsize)==(Integer)0,zz,U);
PokeIndex<LorentzIndex>(UmuDD,U,mu);
}
} else {
UmuDD = Umu;
}
2015-06-22 12:49:44 +01:00
// SU3::ColdConfiguration(RNG4,Umu);
// SU3::TepidConfiguration(RNG4,Umu);
// SU3::HotConfiguration(RNG4,Umu);
2018-01-27 23:46:02 +00:00
// Umu=Zero();
2015-06-22 12:49:44 +01:00
2016-03-16 09:31:04 +00:00
RealD mass=params.mq;
2015-06-22 12:49:44 +01:00
RealD M5=1.8;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "Building g5R5 hermitian DWF operator" <<std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
DomainWallFermionR Ddwf(Umu,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5);
DomainWallFermionR DdwfDD(UmuDD,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5);
2015-06-22 12:49:44 +01:00
typedef Aggregation<vSpinColourVector,vTComplex,nbasis> Subspace;
typedef CoarsenedMatrix<vSpinColourVector,vTComplex,nbasis> CoarseOperator;
typedef CoarseOperator::CoarseVector CoarseVector;
2019-12-28 15:32:35 +00:00
typedef CoarseOperator::siteVector siteVector;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "Calling Aggregation class to build subspace" <<std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
MdagMLinearOperator<DomainWallFermionR,LatticeFermion> HermDefOp(Ddwf);
2019-12-28 15:32:35 +00:00
2017-10-30 00:05:21 +00:00
Subspace Aggregates(Coarse5d,FGrid,0);
2019-12-28 15:32:35 +00:00
assert ( (nbasis & 0x1)==0);
int nb=nbasis/2;
std::cout<<GridLogMessage << " nbasis/2 = "<<nb<<std::endl;
2019-12-09 07:55:45 +00:00
// Aggregates.CreateSubspace(RNG5,HermDefOp,nb);
2017-10-30 00:05:21 +00:00
// Aggregates.CreateSubspaceLanczos(RNG5,HermDefOp,nb);
2019-12-28 15:32:35 +00:00
double f_first = 0.03;
double f_div = 1.2;
std::vector<double> f_lo(nb);
f_lo[0] = f_first;
for(int b=1;b<nb;b++) {
f_lo[b] = f_lo[b-1]/f_div;
}
std::vector<int> f_ord(nb,200);
f_ord[0]=500;
Aggregates.CreateSubspaceChebyshev(RNG5,HermDefOp,nb,60.0,f_lo,f_ord);
for(int n=0;n<nb;n++){
G5R5(Aggregates.subspace[n+nb],Aggregates.subspace[n]);
2019-12-28 15:32:35 +00:00
// std::cout<<GridLogMessage<<n<<" subspace "<<norm2(Aggregates.subspace[n+nb])<<" "<<norm2(Aggregates.subspace[n]) <<std::endl;
}
2019-12-28 15:32:35 +00:00
// for(int n=0;n<nbasis;n++){
// std::cout<<GridLogMessage << "vec["<<n<<"] = "<<norm2(Aggregates.subspace[n]) <<std::endl;
// }
// for(int i=0;i<nbasis;i++){
// result = Aggregates.subspace[i];
// Aggregates.subspace[i]=result+g5*result;
// }
2018-01-27 23:46:02 +00:00
result=Zero();
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "Building coarse representation of Indef operator" <<std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
Gamma5R5HermitianLinearOperator<DomainWallFermionR,LatticeFermion> HermIndefOp(Ddwf);
Gamma5R5HermitianLinearOperator<DomainWallFermionR,LatticeFermion> HermIndefOpDD(DdwfDD);
2019-12-09 07:55:45 +00:00
CoarsenedMatrix<vSpinColourVector,vTComplex,nbasis> LDOp(*Coarse5d,1); // Hermitian matrix
2015-06-22 12:49:44 +01:00
LDOp.CoarsenOperator(FGrid,HermIndefOp,Aggregates);
CoarseVector c_src (Coarse5d);
CoarseVector c_res (Coarse5d);
gaussian(CRNG,c_src);
2018-01-27 23:46:02 +00:00
c_res=Zero();
2015-06-22 12:49:44 +01:00
2019-12-28 15:32:35 +00:00
//////////////////////////////////////////////////
// Deflate the course space. Recursive multigrid?
//////////////////////////////////////////////////
2015-06-22 12:49:44 +01:00
2019-12-28 15:32:35 +00:00
typedef CoarsenedMatrix<vSpinColourVector,vTComplex,nbasis> Level1Op;
typedef CoarsenedMatrix<siteVector,iScalar<vTComplex>,nbasis> Level2Op;
auto cclatt = clatt;
for(int d=0;d<clatt.size();d++){
cclatt[d] = clatt[d]/block[d];
}
GridCartesian *CoarseCoarse4d = SpaceTimeGrid::makeFourDimGrid(cclatt, GridDefaultSimd(Nd,vComplex::Nsimd()),GridDefaultMpi());;
GridCartesian *CoarseCoarse5d = SpaceTimeGrid::makeFiveDimGrid(1,CoarseCoarse4d);
typedef Aggregation<siteVector,iScalar<vTComplex>,nbasis> CoarseSubspace;
CoarseSubspace CoarseAggregates(CoarseCoarse5d,Coarse5d,0);
double c_first = 0.2;
double c_div = 1.2;
std::vector<double> c_lo(nb);
c_lo[0] = c_first;
for(int b=1;b<nb;b++) {
c_lo[b] = c_lo[b-1]/c_div;
}
std::vector<int> c_ord(nb,200);
c_ord[0]=500;
2015-06-22 12:49:44 +01:00
2019-12-28 15:32:35 +00:00
#define RECURSIVE_MULTIGRID
#ifdef RECURSIVE_MULTIGRID
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
2019-12-28 15:32:35 +00:00
std::cout<<GridLogMessage << "Build deflation space in coarse operator "<< std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
2015-06-22 12:49:44 +01:00
2019-12-28 15:32:35 +00:00
MdagMLinearOperator<CoarseOperator,CoarseVector> PosdefLdop(LDOp);
// CoarseAggregates.CreateSubspaceChebyshev(CRNG,PosdefLdop,nbasis,14.0,c_lo,c_ord);
// CoarseAggregates.CreateSubspaceRandom(CRNG);
2019-12-28 15:32:35 +00:00
// Level2Op L2Op(*CoarseCoarse5d,1); // Hermitian matrix
// HermitianLinearOperator<Level1Op,CoarseVector> L1LinOp(LDOp);
// L2Op.CoarsenOperator(Coarse5d,L1LinOp,CoarseAggregates);
#endif
2019-12-09 07:55:45 +00:00
// std::cout<<GridLogMessage << "**************************************************"<< std::endl;
// std::cout<<GridLogMessage << "Unprec CG "<< std::endl;
// std::cout<<GridLogMessage << "**************************************************"<< std::endl;
// TrivialPrecon<LatticeFermion> simple;
2015-07-21 16:30:05 +01:00
// ConjugateGradient<LatticeFermion> fCG(1.0e-8,100000);
// fCG(HermDefOp,src,result);
2015-06-22 12:49:44 +01:00
2019-12-28 15:32:35 +00:00
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "Red Black Prec CG "<< std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
LatticeFermion src_o(FrbGrid);
LatticeFermion result_o(FrbGrid);
pickCheckerboard(Odd,src_o,src);
result_o=Zero();
SchurDiagMooeeOperator<DomainWallFermionR,LatticeFermion> HermOpEO(Ddwf);
ConjugateGradient<LatticeFermion> pCG(1.0e-8,10000);
// pCG(HermOpEO,src_o,result_o);
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << " Running coarse grid Lanczos "<< std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
MdagMLinearOperator<Level1Op,CoarseVector> IRLHermOp(LDOp);
Chebyshev<CoarseVector> IRLCheby(0.01,14,161);
FunctionHermOp<CoarseVector> IRLOpCheby(IRLCheby,IRLHermOp);
PlainHermOp<CoarseVector> IRLOp (IRLHermOp);
int Nstop=32;
int Nk=32;
int Nm=48;
ImplicitlyRestartedLanczos<CoarseVector> IRL(IRLOpCheby,IRLOp,Nstop,Nk,Nm,1.0e-4,20);
int Nconv;
std::vector<RealD> eval(Nm);
std::vector<CoarseVector> evec(Nm,Coarse5d);
IRL.calc(eval,evec,c_src,Nconv);
2019-12-09 07:55:45 +00:00
2015-06-22 12:49:44 +01:00
2019-12-28 15:32:35 +00:00
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "coarse grid CG "<< std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
// ConjugateGradient<CoarseVector> CG(3.0e-3,100000);
// CG(PosdefLdop,c_src,c_res);
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "coarse grid Deflated CG with "<< eval.size() << " evecs" << std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
2019-12-28 15:32:35 +00:00
c_res=Zero();
DeflatedGuesser<CoarseVector> DeflCoarseGuesser(evec,eval);
DeflCoarseGuesser(c_src,c_res);
// CG(PosdefLdop,c_src,c_res);
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
2019-12-28 15:32:35 +00:00
std::cout<<GridLogMessage <<" Applying Fine power method to find spectral range "<<std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
2019-12-28 15:32:35 +00:00
ZeroGuesser<CoarseVector> CoarseZeroGuesser;
2019-12-28 15:32:35 +00:00
MultiGridPreconditioner <vSpinColourVector,vTComplex,nbasis,DomainWallFermionR,
ZeroGuesser<CoarseVector> >
Precon (Aggregates, LDOp,
HermIndefOp,Ddwf,
HermIndefOp,Ddwf,
CoarseZeroGuesser,
params);
// Precon.PowerMethod(src);
/*
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage <<" Applying Coarse power method to find spectral range "<<std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
cparams = params;
cparams.hi = 20.0;
cparams.lo = 0.2;
cparams.order= 20;
MultiGridPreconditioner <siteVector,iScalar<vTComplex>,nbasis,Level1Op,ZeroGuesser<CoarseVector> >
CoarsePrecon (CoarseAggregates,
L2Op,
L1LinOp,LDOp,
L1LinOp,LDOp,
CoarseZeroGuesser,
cparams);
CoarsePrecon.PowerMethod(c_src);
*/
2019-12-28 15:32:35 +00:00
/*
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "Building a two level PGCR "<< std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
2016-03-16 09:31:04 +00:00
PrecGeneralisedConjugateResidual<LatticeFermion> PGCR(1.0e-8,100000,Precon,8,8);
std::cout<<GridLogMessage<<"checking norm src "<<norm2(src)<<std::endl;
2018-01-27 23:46:02 +00:00
result=Zero();
2016-03-16 09:31:04 +00:00
PGCR(HermIndefOp,src,result);
2019-12-28 15:32:35 +00:00
*/
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "Building a two level deflated PGCR "<< std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
MultiGridPreconditioner <vSpinColourVector,vTComplex,nbasis,DomainWallFermionR, DeflatedGuesser<CoarseVector> >
DeflatedPrecon (Aggregates, LDOp,
HermIndefOp,Ddwf,
HermIndefOp,Ddwf,
DeflCoarseGuesser,
params);
PrecGeneralisedConjugateResidual<LatticeFermion> deflPGCR(1.0e-8,100000,DeflatedPrecon,16,16);
std::cout<<GridLogMessage<<"checking norm src "<<norm2(src)<<std::endl;
result=Zero();
deflPGCR(HermIndefOp,src,result);
/*
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "Building deflation preconditioner "<< std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
PrecGeneralisedConjugateResidual<CoarseVector> CPGCR(1.0e-3,10000,CoarsePrecon,8,8);
std::cout<<GridLogMessage<<"checking norm src "<<norm2(c_src)<<std::endl;
c_res=Zero();
CPGCR(L1LinOp,c_src,c_res);
*/
2015-06-22 12:49:44 +01:00
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
std::cout<<GridLogMessage << "Done "<< std::endl;
std::cout<<GridLogMessage << "**************************************************"<< std::endl;
2015-06-22 12:49:44 +01:00
Grid_finalize();
}