/************************************************************************************* 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; /* 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" < block ({2,2,2,2}); const int nbasis= 8; 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); 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< SubspaceOp(Dw); Subspace Aggregates4D(Coarse4d,UGrid,0); Subspace Aggregates5D(Coarse5d,FGrid,0); assert ( (nbasis & 0x1)==0); std::cout< Level1Op; NonHermitianLinearOperator LinOpDwf(Ddwf); Level1Op LDOp (*Coarse5d,0); std::cout< CoarseMdagM(LDOp); BiCGSTAB CoarseBiCGSTAB(tol,MaxIt); ConjugateGradient CoarseCG(tol,MaxIt); c_res=Zero(); CoarseCG(CoarseMdagM,c_src,c_res); std::cout<