/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./tests/Test_dwf_hdcr.cc Copyright (C) 2015 Author: Peter Boyle 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 */ #include #include #include #include using namespace std; using namespace Grid; // TODO // // Coarse Grid axpby_ssp_pminus // Inherit from spProj5pm // Coarse Grid axpby_ssp_pplus template class CayleyBase : public SparseMatrixBase { public: int Ls; // protected: RealD mass; RealD M5; // Save arguments to SetCoefficientsInternal Vector _gamma; RealD _zolo_hi; RealD _b; RealD _c; // Cayley form Moebius (tanh and zolotarev) Vector omega; Vector bs; // S dependent coeffs Vector cs; Vector as; // For preconditioning Cayley form Vector bee; Vector cee; Vector aee; Vector beo; Vector ceo; Vector aeo; // LDU factorisation of the eeoo matrix Vector lee; Vector leem; Vector uee; Vector ueem; Vector dee; public: CayleyBase(RealD _M5, RealD _mass, int _Ls, RealD b_, RealD c_) : M5(_M5), mass(_mass), Ls(_Ls), _b(b_), _c(c_) { RealD eps = 1.0; Approx::zolotarev_data *zdata = Approx::higham(eps,this->Ls);// eps is ignored for higham this->SetCoefficientsTanh(zdata,1.0,0.0); Approx::zolotarev_free(zdata); } ///////////////////////////////////////////////////////// // Replicates functionality // Use a common base class approach ///////////////////////////////////////////////////////// // Tanh void SetCoefficientsTanh(Approx::zolotarev_data *zdata,RealD b,RealD c) { Vector gamma(this->Ls); for(int s=0;sLs;s++) gamma[s] = zdata->gamma[s]; SetCoefficientsInternal(1.0,gamma,b,c); } //Zolo void SetCoefficientsZolotarev(RealD zolo_hi,Approx::zolotarev_data *zdata,RealD b,RealD c) { Vector gamma(this->Ls); for(int s=0;sLs;s++) gamma[s] = zdata->gamma[s]; SetCoefficientsInternal(zolo_hi,gamma,b,c); } //Zolo void SetCoefficientsInternal(RealD zolo_hi,Vector & gamma,RealD b,RealD c) { int Ls=this->Ls; /////////////////////////////////////////////////////////// // The Cayley coeffs (unprec) /////////////////////////////////////////////////////////// assert(gamma.size()==Ls); omega.resize(Ls); bs.resize(Ls); cs.resize(Ls); as.resize(Ls); double bpc = b+c; double bmc = b-c; _b = b; _c = c; _gamma = gamma; // Save the parameters so we can change mass later. _zolo_hi= zolo_hi; for(int i=0; i < Ls; i++){ as[i] = 1.0; omega[i] = _gamma[i]*_zolo_hi; //NB reciprocal relative to Chroma NEF code assert(omega[i]!=Coeff_t(0.0)); bs[i] = 0.5*(bpc/omega[i] + bmc); cs[i] = 0.5*(bpc/omega[i] - bmc); } //////////////////////////////////////////////////////// // Constants for the preconditioned matrix Cayley form //////////////////////////////////////////////////////// bee.resize(Ls); cee.resize(Ls); beo.resize(Ls); ceo.resize(Ls); for(int i=0;iM5) +1.0); assert(bee[i]!=Coeff_t(0.0)); cee[i]=as[i]*(1.0-cs[i]*(4.0-this->M5)); beo[i]=as[i]*bs[i]; ceo[i]=-as[i]*cs[i]; } aee.resize(Ls); aeo.resize(Ls); for(int i=0;i &out){}; virtual void DW (const Field &psi, Field &chi)=0; virtual void DWDag (const Field &psi, Field &chi)=0; void M (const Field &psi, Field &chi) { Field Din(psi.Grid()); Meooe5D(psi,Din); DW(Din,chi); axpby(chi,1.0,1.0,chi,psi); M5D(psi,chi); } void Mdag (const Field &psi, Field &chi) { Field Din(psi.Grid()); DWDag(psi,Din); MeooeDag5D(Din,chi); M5Ddag(psi,chi); axpby (chi,1.0,1.0,chi,psi); } ///////////////////////////////// // P and Pdag - might be needed ///////////////////////////////// void P(const Field &psi, Field &chi) { int Ls= this->Ls; chi=Zero(); for(int s=0;sLs; chi=Zero(); for(int s=0;sLs; Vector diag (Ls,1.0); Vector upper(Ls,-1.0); upper[Ls-1]=mass; Vector lower(Ls,-1.0); lower[0] =mass; M5D(psi,chi,chi,lower,diag,upper); } void M5Ddag (const Field &psi, Field &chi) { int Ls=this->Ls; Vector diag(Ls,1.0); Vector upper(Ls,-1.0); Vector lower(Ls,-1.0); upper[Ls-1]=-mass*upper[Ls-1]; lower[0] =-mass*lower[0]; M5Ddag(psi,chi,chi,lower,diag,upper); } void Meooe5D (const Field &psi, Field &Din) { int Ls=this->Ls; Vector diag = bs; Vector upper= cs; Vector lower= cs; upper[Ls-1]=-mass*upper[Ls-1]; lower[0] =-mass*lower[0]; M5D(psi,psi,Din,lower,diag,upper); } void MeooeDag5D (const Field &psi, Field &Din) { int Ls=this->Ls; Vector diag =bs; Vector upper=cs; Vector lower=cs; for (int s=0;s &lower, Vector &diag, Vector &upper) { chi_i.Checkerboard()=psi_i.Checkerboard(); GridBase *grid=psi_i.Grid(); autoView(psi , psi_i,AcceleratorRead); autoView(phi , phi_i,AcceleratorRead); autoView(chi , chi_i,AcceleratorWrite); assert(phi.Checkerboard() == psi.Checkerboard()); auto pdiag = &diag[0]; auto pupper = &upper[0]; auto plower = &lower[0]; int Ls =this->Ls; // 10 = 3 complex mult + 2 complex add // Flops = 10.0*(Nc*Ns) *Ls*vol (/2 for red black counting) uint64_t nloop = grid->oSites()/Ls; const int Nsimd = Field::vector_type::Nsimd(); accelerator_for(sss,nloop,Nsimd,{ uint64_t ss= sss*Ls; typedef decltype(coalescedRead(psi[0])) spinor; spinor tmp1, tmp2; for(int s=0;s &lower, Vector &diag, Vector &upper) { chi_i.Checkerboard()=psi_i.Checkerboard(); GridBase *grid=psi_i.Grid(); autoView(psi , psi_i,AcceleratorRead); autoView(phi , phi_i,AcceleratorRead); autoView(chi , chi_i,AcceleratorWrite); assert(phi.Checkerboard() == psi.Checkerboard()); auto pdiag = &diag[0]; auto pupper = &upper[0]; auto plower = &lower[0]; int Ls=this->Ls; uint64_t nloop = grid->oSites()/Ls; const int Nsimd = Field::vector_type::Nsimd(); accelerator_for(sss,nloop,Nsimd,{ uint64_t ss=sss*Ls; typedef decltype(coalescedRead(psi[0])) spinor; spinor tmp1,tmp2; for(int s=0;s class CoarseCayleyFermion : public CayleyBase< Lattice > , ComplexD > { public: typedef iVector siteVector; typedef Lattice CoarseComplexField; typedef Lattice CoarseVector; typedef Lattice > CoarseMatrix; typedef iMatrix Cobj; typedef Lattice< CComplex > CoarseScalar; // used for inner products on fine field typedef Lattice FineField; // Similar to the CoarseOperator but add 5D support. Geometry geom; GridBase *Coarse5D; GridBase *Coarse4D; CartesianStencil Stencil; CoarsenedMatrix &Dw; Vector NNsv; int *NNs; GridBase * Grid(void) { return Coarse5D; }; // this is all the linalg routines need to know CoarseCayleyFermion(GridCartesian &CoarseGrid4, GridCartesian &CoarseGrid5, CoarsenedMatrix &_Dw, RealD M5, RealD mass, int Ls, RealD b, RealD c) : CayleyBase(M5,mass,Ls,b,c), Coarse4D(&CoarseGrid4), Coarse5D(&CoarseGrid5), Dw(_Dw), geom(CoarseGrid5._ndimension), Stencil( &CoarseGrid5,geom.npoint,Even,geom.directions,geom.displacements,0), NNsv(Ls) { int surface=4; NNs =&NNsv[0]; for(int s=0;s=(nbasis/2)) && (b< ((nbasis/2)+NNs[s]) ) ) ) }; public: void Project( CoarseVector &C ) { const int Nsimd = CComplex::Nsimd(); autoView(Cv,C, AcceleratorWrite); int nb2=nbasis/2; int Ls = this->Ls; for(int s=0;soSites(), Nsimd, { int sF= sU*Ls+s; auto tmp = coalescedRead(Cv[sF]); for(int b=NN;b compressor; Stencil.HaloExchange(in,compressor); typedef LatticeView Aview; const int Nsimd = CComplex::Nsimd(); // Ls loop for2D int Ls=this->Ls; Vector AcceleratorViewContainer; for(int p=0;poSites(), b, nbasis, Nsimd, { typedef decltype(coalescedRead(in_v[0])) calcVector; typedef decltype(coalescedRead(in_v[0](0))) calcComplex; int sU = sF/Ls; int s = sF%Ls; calcComplex res = Zero(); if ( IS_INCLUDED(s,b) ) { calcVector nbr; int ptype; for(int point=0;point_is_local) { nbr = coalescedReadPermute(in_v[SE->_offset],ptype,SE->_permute); } else { nbr = coalescedRead(CBp[SE->_offset]); } acceleratorSynchronise(); for(int bb=0;bb Aggregates; void PromoteFromSubspace(Aggregates &_Aggregates,CoarseVector &C,FineField &F) { auto FineGrid4 = _Aggregates.FineGrid; FineField F4(FineGrid4); CoarseVector C4(Coarse4D); for(int s=0;sLs;s++){ ExtractSlice(C4,C,s,0); _Aggregates.PromoteFromSubspace(C4,F4); InsertSlice(F4,F,s,0); } } void ProjectToSubspace(Aggregates &_Aggregates,CoarseVector &C,FineField &F) { auto FineGrid4 = _Aggregates.FineGrid; FineField F4(FineGrid4); CoarseVector C4(Coarse4D); for(int s=0;sLs;s++){ ExtractSlice(F4,F,s,0); _Aggregates.ProjectToSubspace (C4,F4); InsertSlice(C4,C,s,0); } Project(C); } template void Test(Aggregates &_Aggregates,GridBase *FineGrid, Ddwf &_Ddwf) { typedef Lattice FineField; CoarseVector Cin(Coarse5D); CoarseVector Cout(Coarse5D); CoarseVector CFout(Coarse5D); FineField Fin(FineGrid); FineField Fout(FineGrid); std::vector seeds({1,2,3,4,5}); GridParallelRNG RNG(Coarse5D); RNG.SeedFixedIntegers(seeds); gaussian(RNG,Cin); PromoteFromSubspace(_Aggregates,Cin,Fin); ProjectToSubspace(_Aggregates,Cin,Fin); std::cout << "************ "<M(Cin,Cout); this->Project(Cout); std::cout << " Cout "<Mdag(Cin,Cout); this->Project(Cout); std::cout << " Cout "< class SolverWrapper : public LinearFunction { private: LinearOperatorBase & _Matrix; OperatorFunction & _Solver; LinearFunction & _Guess; public: ///////////////////////////////////////////////////// // Wrap the usual normal equations trick ///////////////////////////////////////////////////// SolverWrapper(LinearOperatorBase &Matrix, OperatorFunction &Solver, LinearFunction &Guess) : _Matrix(Matrix), _Solver(Solver), _Guess(Guess) {}; void operator() (const Field &in, Field &out){ _Guess(in,out); _Solver(_Matrix,in,out); // Mdag M out = Mdag in } }; // Must use a non-hermitian solver template class PVdagMLinearOperator : public LinearOperatorBase { Matrix &_Mat; Matrix &_PV; public: PVdagMLinearOperator(Matrix &Mat,Matrix &PV): _Mat(Mat),_PV(PV){}; void OpDiag (const Field &in, Field &out) { assert(0); } void OpDir (const Field &in, Field &out,int dir,int disp) { assert(0); } void OpDirAll (const Field &in, std::vector &out){ assert(0); }; void Op (const Field &in, Field &out){ Field tmp(in.Grid()); _Mat.M(in,tmp); _PV.Mdag(tmp,out); } void AdjOp (const Field &in, Field &out){ Field tmp(in.Grid()); _PV.M(tmp,out); _Mat.Mdag(in,tmp); } void HermOpAndNorm(const Field &in, Field &out,RealD &n1,RealD &n2){ assert(0); } void HermOp(const Field &in, Field &out){ assert(0); } }; RealD InverseApproximation(RealD x){ return 1.0/x; } template class ChebyshevSmoother : public LinearFunction { public: typedef LinearOperatorBase FineOperator; Matrix & _SmootherMatrix; FineOperator & _SmootherOperator; Chebyshev Cheby; ChebyshevSmoother(RealD _lo,RealD _hi,int _ord, FineOperator &SmootherOperator,Matrix &SmootherMatrix) : _SmootherOperator(SmootherOperator), _SmootherMatrix(SmootherMatrix), Cheby(_lo,_hi,_ord,InverseApproximation) {}; void operator() (const Field &in, Field &out) { Field tmp(in.Grid()); MdagMLinearOperator MdagMOp(_SmootherMatrix); _SmootherOperator.AdjOp(in,tmp); Cheby(MdagMOp,tmp,out); } }; template class MGPreconditioner : public LinearFunction< Lattice > { public: typedef Aggregation Aggregates; typedef CoarseCayleyFermion CoarseOperator; typedef typename Aggregation::CoarseVector CoarseVector; typedef typename Aggregation::CoarseMatrix CoarseMatrix; typedef typename Aggregation::FineField FineField; typedef LinearOperatorBase FineOperator; typedef LinearFunction FineSmoother; Aggregates & _Aggregates; FineOperator & _FineOperator; FineSmoother & _PreSmoother; FineSmoother & _PostSmoother; CoarseOperator & _CoarseOperator; CoarseSolver & _CoarseSolve; int level; void Level(int lv) {level = lv; }; MGPreconditioner(Aggregates &Agg, FineOperator &Fine, FineSmoother &PreSmoother, FineSmoother &PostSmoother, CoarseOperator &CoarseOperator_, CoarseSolver &CoarseSolve_) : _Aggregates(Agg), _FineOperator(Fine), _PreSmoother(PreSmoother), _PostSmoother(PostSmoother), _CoarseOperator(CoarseOperator_), _CoarseSolve(CoarseSolve_), level(1) { } virtual void operator()(const FineField &in, FineField & out) { auto CoarseGrid = _CoarseOperator.Grid(); CoarseVector Csrc(CoarseGrid); CoarseVector Csol(CoarseGrid); FineField vec1(in.Grid()); FineField vec2(in.Grid()); double t; // Fine Smoother t=-usecond(); _PreSmoother(in,out); t+=usecond(); std::cout< block ({2,2,2,2}); const int nbasis= 24; auto clatt = GridDefaultLatt(); for(int d=0;d seeds({1,2,3,4}); GridParallelRNG RNG5(FGrid); RNG5.SeedFixedIntegers(seeds); GridParallelRNG RNG4(UGrid); RNG4.SeedFixedIntegers(seeds); GridParallelRNG CRNG(Coarse5d);CRNG.SeedFixedIntegers(seeds); LatticeGaugeField Umu(UGrid); #if 0 SU3::TepidConfiguration(RNG4,Umu); #else FieldMetaData header; std::string file("./ckpoint_lat.4000"); NerscIO::readConfiguration(Umu,header,file); #endif std::cout< Subspace; typedef CoarsenedMatrix CoarseOperator; typedef CoarseOperator::CoarseVector CoarseVector; typedef CoarseOperator::siteVector siteVector; std::cout< MdagM_Dw(Dw_null); std::cout< WilsonCG(1.0e-10,40000); LatticeFermion w_src(UGrid); w_src=1.0; LatticeFermion w_res(UGrid); WilsonCG(MdagM_Dw,w_src,w_res); exit(0); */ std::cout< Level1Op4; typedef CoarseCayleyFermion Level1Op5; Level1Op4 c_Dw (*Coarse4d,0); std::cout< LinOpDw(Dw); std::cout< CoarseBiCGSTAB(tol,MaxIt); ConjugateGradient CoarseCG(tol,MaxIt); BiCGSTAB FineBiCGSTAB(tol,MaxIt); ConjugateGradient FineCG(tol,MaxIt); NonHermitianLinearOperator FineM(Ddwf); MdagMLinearOperator FineMdagM(Ddwf); // M^\dag M PVdagMLinearOperator FinePVdagM(Ddwf,Dpv);// M_{pv}^\dag M SchurDiagMooeeOperator FineDiagMooee(Ddwf); // M_ee - Meo Moo^-1 Moe SchurDiagOneOperator FineDiagOne(Ddwf); // 1 - M_ee^{-1} Meo Moo^{-1} Moe e MdagMLinearOperator CoarseMdagM(c_Dwf); PVdagMLinearOperator CoarsePVdagM(c_Dwf,c_Dpv); std::cout< FineDiagMooee(Ddwf); // M_ee - Meo Moo^-1 Moe f_res_e=Zero(); FineCG(FineDiagMooee,f_src_e,f_res_e); */ std::cout< PM; PM(MdagM_Dw,w_src); std::cout< cPM; cPM(CoarseMdagM,c_src); std::cout< CoarseZeroGuesser; NormalEquations CoarseCGNE (c_Dwf,CoarseCG,CoarseZeroGuesser); ConjugateGradient CoarseMgridCG(0.001,1000); ChebyshevSmoother FineSmoother(0.5,60.0,14,FineM,Ddwf); typedef MGPreconditioner > TwoLevelHDCR; TwoLevelHDCR TwoLevelPrecon(Aggregates4D, FineM, FineSmoother, FineSmoother, c_Dwf, CoarseCGNE); TwoLevelPrecon.Level(1); PrecGeneralisedConjugateResidualNonHermitian l1PGCR(1.0e-8,100,FineM,TwoLevelPrecon,16,16); l1PGCR.Level(1); f_res=Zero(); CoarseCG.Tolerance=0.02; l1PGCR(f_src,f_res); Grid_finalize(); exit(0); std::cout< CoarseMgridBiCGSTAB(0.01,1000,false); BiCGSTAB FineMgridPreBiCGSTAB(0.0,24,false); BiCGSTAB FineMgridPostBiCGSTAB(0.0,24,false); ZeroGuesser FineZeroGuesser; SolverWrapper FineBiCGPreSmoother( FinePVdagM, FineMgridPreBiCGSTAB, FineZeroGuesser); SolverWrapper FineBiCGPostSmoother( FinePVdagM, FineMgridPostBiCGSTAB, FineZeroGuesser); SolverWrapper CoarsePVdagMSolver(CoarsePVdagM,CoarseMgridBiCGSTAB,CoarseZeroGuesser); typedef MGPreconditioner > TwoLevelMG; TwoLevelMG _TwoLevelMG(Aggregates4D, FinePVdagM, FineBiCGPreSmoother, FineBiCGPostSmoother, c_Dwf, CoarsePVdagMSolver); _TwoLevelMG.Level(1); PrecGeneralisedConjugateResidualNonHermitian pvPGCR(1.0e-8,100,FinePVdagM,_TwoLevelMG,16,16); pvPGCR.Level(1); f_res=Zero(); pvPGCR(f_src,f_res); std::cout<