/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./tests/Test_dwf_hdcr.cc Copyright (C) 2015 Author: Antonin Portelli Author: Peter Boyle Author: paboyle 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; /* Params * Grid: * block1(4) * block2(4) * * Subspace * * Fine : Subspace(nbasis,hi,lo,order,first,step) -- 32, 60,0.02,500,100,100 * * Coarse: Subspace(nbasis,hi,lo,order,first,step) -- 32, 18,0.02,500,100,100 * Smoother: * * Fine: Cheby(hi, lo, order) -- 60,0.5,10 * * Coarse: Cheby(hi, lo, order) -- 12,0.1,4 * Lanczos: * CoarseCoarse IRL( Nk, Nm, Nstop, poly(lo,hi,order)) 24,36,24,0.002,4.0,61 */ template 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 MirsSmoother : public LinearFunction { public: typedef LinearOperatorBase FineOperator; Matrix & SmootherMatrix; FineOperator & SmootherOperator; RealD tol; RealD shift; int maxit; MirsSmoother(RealD _shift,RealD _tol,int _maxit,FineOperator &_SmootherOperator,Matrix &_SmootherMatrix) : shift(_shift),tol(_tol),maxit(_maxit), SmootherOperator(_SmootherOperator), SmootherMatrix(_SmootherMatrix) {}; void operator() (const Field &in, Field &out) { ZeroGuesser Guess; ConjugateGradient CG(tol,maxit,false); Field src(in.Grid()); ShiftedMdagMLinearOperator,Field> MdagMOp(SmootherMatrix,shift); SmootherOperator.AdjOp(in,src); Guess(src,out); CG(MdagMOp,src,out); } }; #define GridLogLevel std::cout << GridLogMessage < class HDCRPreconditioner : public LinearFunction< Lattice > { public: typedef Aggregation Aggregates; typedef CoarsenedMatrix 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 & _Smoother; CoarseSolver & _CoarseSolve; int level; void Level(int lv) {level = lv; }; HDCRPreconditioner(Aggregates &Agg, FineOperator &Fine, FineSmoother &Smoother, CoarseSolver &CoarseSolve_) : _Aggregates(Agg), _FineOperator(Fine), _Smoother(Smoother), _CoarseSolve(CoarseSolve_), level(1) { } virtual void operator()(const FineField &in, FineField & out) { auto CoarseGrid = _Aggregates.CoarseGrid; CoarseVector Csrc(CoarseGrid); CoarseVector Csol(CoarseGrid); FineField vec1(in.Grid()); FineField vec2(in.Grid()); double t; // Fine Smoother t=-usecond(); _Smoother(in,out); t+=usecond(); GridLogLevel << "Smoother took "<< t/1000.0<< "ms" < class MultiGridPreconditioner : public LinearFunction< Lattice > { public: typedef Aggregation Aggregates; typedef CoarsenedMatrix CoarseOperator; typedef typename Aggregation::CoarseVector CoarseVector; typedef typename Aggregation::CoarseMatrix CoarseMatrix; typedef typename Aggregation::FineField FineField; typedef LinearOperatorBase FineOperator; typedef LinearFunction FineSmoother; Aggregates & _Aggregates; CoarseOperator & _CoarseOperator; FineOperator & _FineOperator; Guesser & _Guess; FineSmoother & _Smoother; CoarseSolver & _CoarseSolve; int level; void Level(int lv) {level = lv; }; MultiGridPreconditioner(Aggregates &Agg, CoarseOperator &Coarse, FineOperator &Fine, FineSmoother &Smoother, Guesser &Guess_, CoarseSolver &CoarseSolve_) : _Aggregates(Agg), _CoarseOperator(Coarse), _FineOperator(Fine), _Smoother(Smoother), _Guess(Guess_), _CoarseSolve(CoarseSolve_), level(1) { } virtual void operator()(const FineField &in, FineField & out) { CoarseVector Csrc(_CoarseOperator.Grid()); CoarseVector Csol(_CoarseOperator.Grid()); FineField vec1(in.Grid()); FineField vec2(in.Grid()); double t; // Fine Smoother t=-usecond(); _Smoother(in,out); t+=usecond(); GridLogLevel << "Smoother took "<< t/1000.0<< "ms" < block ({2,2,2,2}); std::vector blockc ({2,2,2,2}); const int nbasis= 32; const int nbasisc= 32; auto clatt = GridDefaultLatt(); for(int d=0;d seeds4({1,2,3,4}); std::vector seeds5({5,6,7,8}); std::vector cseeds({5,6,7,8}); GridParallelRNG RNG5(FGrid); RNG5.SeedFixedIntegers(seeds5); GridParallelRNG RNG4(UGrid); RNG4.SeedFixedIntegers(seeds4); GridParallelRNG CRNG(Coarse5d);CRNG.SeedFixedIntegers(cseeds); LatticeFermion src(FGrid); gaussian(RNG5,src);// src=src+g5*src; LatticeFermion result(FGrid); LatticeGaugeField Umu(UGrid); FieldMetaData header; std::string file("./ckpoint_lat.4000"); NerscIO::readConfiguration(Umu,header,file); std::cout< Subspace; typedef CoarsenedMatrix CoarseOperator; typedef CoarseOperator::CoarseVector CoarseVector; typedef CoarseOperator::siteVector siteVector; std::cout< HermDefOp(Ddwf); Subspace Aggregates(Coarse5d,FGrid,0); assert ( (nbasis & 0x1)==0); { int nb=nbasis/2; Aggregates.CreateSubspaceChebyshev(RNG5,HermDefOp,nb,60.0,0.02,500,100,100,0.0); for(int n=0;n Level1Op; typedef CoarsenedMatrix,nbasisc> Level2Op; Gamma5R5HermitianLinearOperator HermIndefOp(Ddwf); Gamma5R5HermitianLinearOperator HermIndefOpPV(Dpv); std::cout< CoarseBiCGSTAB(tol,MaxIt); ConjugateGradient CoarseCG(tol,MaxIt); // GeneralisedMinimalResidual CoarseGMRES(tol,MaxIt,20); BiCGSTAB FineBiCGSTAB(tol,MaxIt); ConjugateGradient FineCG(tol,MaxIt); // GeneralisedMinimalResidual FineGMRES(tol,MaxIt,20); 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(LDOp); PVdagMLinearOperator CoarsePVdagM(LDOp,LDOpPV); std::cout< IRLCheby(0.03,12.0,71); // 1 iter FunctionHermOp IRLOpCheby(IRLCheby,CoarseMdagM); PlainHermOp IRLOp (CoarseMdagM); int Nk=64; int Nm=128; int Nstop=Nk; ImplicitlyRestartedLanczos IRL(IRLOpCheby,IRLOp,Nstop,Nk,Nm,1.0e-3,20); int Nconv; std::vector eval(Nm); std::vector evec(Nm,Coarse5d); IRL.calc(eval,evec,c_src,Nconv); std::cout< DeflCoarseGuesser(evec,eval); NormalEquations DeflCoarseCGNE (LDOp,CoarseCG,DeflCoarseGuesser); c_res=Zero(); DeflCoarseCGNE(c_src,c_res); std::cout< CoarseMgridCG(0.001,1000); ChebyshevSmoother FineSmoother(0.5,60.0,10,HermIndefOp,Ddwf); typedef HDCRPreconditioner > TwoLevelHDCR; TwoLevelHDCR TwoLevelPrecon(Aggregates, HermIndefOp, FineSmoother, DeflCoarseCGNE); TwoLevelPrecon.Level(1); // PrecGeneralisedConjugateResidual l1PGCR(1.0e-8,100,HermIndefOp,TwoLevelPrecon,16,16); PrecGeneralisedConjugateResidualNonHermitian l1PGCR(1.0e-8,100,HermIndefOp,TwoLevelPrecon,16,16); l1PGCR.Level(1); f_res=Zero(); CoarseCG.Tolerance=0.02; l1PGCR(f_src,f_res); std::cout< CoarseMgridBiCGSTAB(0.01,1000); BiCGSTAB FineMgridBiCGSTAB(0.0,24); ZeroGuesser CoarseZeroGuesser; ZeroGuesser FineZeroGuesser; SolverWrapper FineBiCGSmoother( FinePVdagM, FineMgridBiCGSTAB, FineZeroGuesser); SolverWrapper CoarsePVdagMSolver(CoarsePVdagM,CoarseMgridBiCGSTAB,CoarseZeroGuesser); typedef HDCRPreconditioner > TwoLevelMG; TwoLevelMG _TwoLevelMG(Aggregates, FinePVdagM, FineBiCGSmoother, 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<