/************************************************************************************* 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 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 */ RealD InverseApproximation(RealD x){ return 1.0/x; } template class SolverWrapper : public LinearFunction { private: CheckerBoardedSparseMatrixBase & _Matrix; SchurRedBlackBase & _Solver; public: ///////////////////////////////////////////////////// // Wrap the usual normal equations trick ///////////////////////////////////////////////////// SolverWrapper(CheckerBoardedSparseMatrixBase &Matrix, SchurRedBlackBase &Solver) : _Matrix(Matrix), _Solver(Solver) {}; void operator() (const Field &in, Field &out){ _Solver(_Matrix,in,out); // Mdag M out = Mdag in } }; 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); } }; template 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; Matrix & _FineMatrix; FineOperator & _FineOperator; Guesser & _Guess; FineSmoother & _Smoother; CoarseSolver & _CoarseSolve; int level; void Level(int lv) {level = lv; }; #define GridLogLevel std::cout << GridLogMessage < 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"); //std::string file("./ckpoint_lat.1000"); 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; LatticeFermion A(FGrid); LatticeFermion B(FGrid); // Aggregates.CreateSubspaceChebyshev(RNG5,HermDefOp,nb,60.0,0.002,1000,800,100,0.0); // Aggregates.CreateSubspaceChebyshev(RNG5,HermDefOp,nb,60.0,0.02,1000,800,100,0.0); Aggregates.CreateSubspaceChebyshev(RNG5,HermDefOp,nb,60.0,0.05,500,200,150,0.0);// // Aggregates.CreateSubspaceChebyshev(RNG5,HermDefOp,nb,60.0,0.01,1000,100,100,0.0); // Slightly faster for(int n=0;n Level1Op; typedef CoarsenedMatrix,nbasisc> Level2Op; Gamma5R5HermitianLinearOperator HermIndefOp(Ddwf); GridRedBlackCartesian * Coarse4dRB = SpaceTimeGrid::makeFourDimRedBlackGrid(Coarse4d); std::cout << " Making 5D coarse RB grid " <,nbasisc> CoarseSubspace; // CoarseSubspace CoarseAggregates(CoarseCoarse5d,Coarse5d,0); std::cout< PosdefLdop(LDOp); typedef Level2Op::CoarseVector CoarseCoarseVector; CoarseVector c_src(Coarse5d); c_src=1.0; std::cout< , SolverWrapper > TwoLevelMG; typedef MultiGridPreconditioner,nbasisc,Level1Op, DeflatedGuesser, NormalEquations > CoarseMG; typedef MultiGridPreconditioner, LinearFunction > ThreeLevelMG; ChebyshevSmoother FineSmoother(0.5,60.0,12,HermIndefOp,Ddwf); std::cout< CoarseZeroGuesser; ConjugateGradient CoarseCG(0.005,1000); // SchurDiagMooeeOperator CoarseMpcDagMpc(LDOp); SchurRedBlackDiagMooeeSolve CoarseRBCG(CoarseCG); SolverWrapper CoarseSolver(LDOp,CoarseRBCG); // NormalEquations CoarseCGNE(LDOp,CoarseCG,CoarseZeroGuesser); TwoLevelMG TwoLevelPrecon(Aggregates, LDOp, HermIndefOp,Ddwf, FineSmoother, CoarseZeroGuesser, CoarseSolver); TwoLevelPrecon.Level(1); PrecGeneralisedConjugateResidual l1PGCR(1.0e-8,20,HermIndefOp,TwoLevelPrecon,16,16); l1PGCR.Level(1); l1PGCR(src,result); std::cout< pCG(1.0e-8,60000); result=Zero(); // pCG(HermDefOp,src,result); std::cout< HermOpEO(Ddwf); // pCG(HermOpEO,src_o,result_o); std::cout< PM; PM(HermDefOp,src); std::cout< cPM; cPM(PosdefLdop,c_src); // std::cout< ccPM; ccPM(IRLHermOpL2,cc_src); std::cout<