/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: Copyright (C) 2015-2016 Author: Peter Boyle Author: Guido Cossu Author: David Murphy Author: Chulwoo Jung 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 #ifdef GRID_DEFAULT_PRECISION_DOUBLE #define MIXED_PRECISION #endif // second level EOFA #undef EOFA_H #define USE_OBC #undef DO_IMPLICIT NAMESPACE_BEGIN(Grid); /* * Need a plan for gauge field update for mixed precision in HMC (2x speed up) * -- Store the single prec action operator. * -- Clone the gauge field from the operator function argument. * -- Build the mixed precision operator dynamically from the passed operator and single prec clone. */ template class MixedPrecisionConjugateGradientOperatorFunction : public OperatorFunction { public: typedef typename FermionOperatorD::FermionField FieldD; typedef typename FermionOperatorF::FermionField FieldF; using OperatorFunction::operator(); RealD Tolerance; RealD InnerTolerance; //Initial tolerance for inner CG. Defaults to Tolerance but can be changed Integer MaxInnerIterations; Integer MaxOuterIterations; GridBase* SinglePrecGrid4; //Grid for single-precision fields GridBase* SinglePrecGrid5; //Grid for single-precision fields RealD OuterLoopNormMult; //Stop the outer loop and move to a final double prec solve when the residual is OuterLoopNormMult * Tolerance FermionOperatorF &FermOpF; FermionOperatorD &FermOpD;; SchurOperatorF &LinOpF; SchurOperatorD &LinOpD; Integer TotalInnerIterations; //Number of inner CG iterations Integer TotalOuterIterations; //Number of restarts Integer TotalFinalStepIterations; //Number of CG iterations in final patch-up step MixedPrecisionConjugateGradientOperatorFunction(RealD tol, Integer maxinnerit, Integer maxouterit, GridBase* _sp_grid4, GridBase* _sp_grid5, FermionOperatorF &_FermOpF, FermionOperatorD &_FermOpD, SchurOperatorF &_LinOpF, SchurOperatorD &_LinOpD): LinOpF(_LinOpF), LinOpD(_LinOpD), FermOpF(_FermOpF), FermOpD(_FermOpD), Tolerance(tol), InnerTolerance(tol), MaxInnerIterations(maxinnerit), MaxOuterIterations(maxouterit), SinglePrecGrid4(_sp_grid4), SinglePrecGrid5(_sp_grid5), OuterLoopNormMult(100.) { /* Debugging instances of objects; references are stored std::cout << GridLogMessage << " Mixed precision CG wrapper LinOpF " < &LinOpU, const FieldD &src, FieldD &psi) { std::cout << GridLogMessage << " Mixed precision CG wrapper operator() "<(&LinOpU); // std::cout << GridLogMessage << " Mixed precision CG wrapper operator() FermOpU " <_Mat)<_Mat)==&(LinOpD._Mat)); //////////////////////////////////////////////////////////////////////////////////// // Must snarf a single precision copy of the gauge field in Linop_d argument //////////////////////////////////////////////////////////////////////////////////// typedef typename FermionOperatorF::GaugeField GaugeFieldF; typedef typename FermionOperatorF::GaugeLinkField GaugeLinkFieldF; typedef typename FermionOperatorD::GaugeField GaugeFieldD; typedef typename FermionOperatorD::GaugeLinkField GaugeLinkFieldD; GridBase * GridPtrF = SinglePrecGrid4; GridBase * GridPtrD = FermOpD.Umu.Grid(); GaugeFieldF U_f (GridPtrF); GaugeLinkFieldF Umu_f(GridPtrF); // std::cout << " Dim gauge field "<Nd()<Nd()<(FermOpD.Umu, mu); precisionChange(Umu_f,Umu_d); PokeIndex(FermOpF.Umu, Umu_f, mu); } pickCheckerboard(Even,FermOpF.UmuEven,FermOpF.Umu); pickCheckerboard(Odd ,FermOpF.UmuOdd ,FermOpF.Umu); //////////////////////////////////////////////////////////////////////////////////// // Make a mixed precision conjugate gradient //////////////////////////////////////////////////////////////////////////////////// MixedPrecisionConjugateGradient MPCG(Tolerance,MaxInnerIterations,MaxOuterIterations,SinglePrecGrid5,LinOpF,LinOpD); std::cout << GridLogMessage << "Calling mixed precision Conjugate Gradient" < HMCWrapper; typedef GenericHMCRunner HMCWrapper; HMCparams.MD.name =std::string("ImplicitMinimumNorm2"); #else // typedef GenericHMCRunner HMCWrapper; // typedef GenericHMCRunner HMCWrapper; typedef GenericHMCRunner HMCWrapper; HMCparams.MD.name =std::string("MinimumNorm2"); #endif std::cout << GridLogMessage<< HMCparams < PlaqObs; TheHMC.Resources.AddObservable(); ////////////////////////////////////////////// const int Ls = 12; Real beta = 5.983; std::cout << GridLogMessage << " beta "<< beta << std::endl; Real light_mass = 0.00049; Real strange_mass = 0.0158; Real charm_mass = 0.191; Real pv_mass = 1.0; RealD M5 = 1.4; RealD b = 2.0; RealD c = 1.0; // Copied from paper // std::vector hasenbusch({ 0.045 }); // Paper values from F1 incorrect run std::vector hasenbusch({ 0.0038, 0.0145, 0.045, 0.108 , 0.25, 0.51 }); // Paper values from F1 incorrect run std::vector hasenbusch2({ 0.4 }); // Paper values from F1 incorrect run // RealD eofa_mass=0.05 ; /////////////////////////////////////////////////////////////////////////////////////////////// //Bad choices with large dH. Equalising force L2 norm was not wise. /////////////////////////////////////////////////////////////////////////////////////////////// //std::vector hasenbusch({ 0.03, 0.2, 0.3, 0.5, 0.8 }); auto GridPtr = TheHMC.Resources.GetCartesian(); auto GridRBPtr = TheHMC.Resources.GetRBCartesian(); auto FGrid = SpaceTimeGrid::makeFiveDimGrid(Ls,GridPtr); auto FrbGrid = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls,GridPtr); Coordinate latt = GridDefaultLatt(); Coordinate mpi = GridDefaultMpi(); Coordinate simdF = GridDefaultSimd(Nd,vComplexF::Nsimd()); Coordinate simdD = GridDefaultSimd(Nd,vComplexD::Nsimd()); auto GridPtrF = SpaceTimeGrid::makeFourDimGrid(latt,simdF,mpi); //auto UGrid_f = SpaceTimeGrid::makeFourDimGrid(latt,simdF,mpi); auto GridRBPtrF = SpaceTimeGrid::makeFourDimRedBlackGrid(GridPtrF); auto FGridF = SpaceTimeGrid::makeFiveDimGrid(Ls,GridPtrF); auto FrbGridF = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls,GridPtrF); #ifndef USE_OBC // IwasakiGaugeActionR GaugeAction(beta); WilsonGaugeActionR GaugeAction(beta); #else std::vector boundaryG = {1,1,1,0}; WilsonGaugeActionR::ImplParams ParamsG(boundaryG); WilsonGaugeActionR GaugeAction(beta,ParamsG); #endif // temporarily need a gauge field LatticeGaugeField U(GridPtr); LatticeGaugeFieldF UF(GridPtrF); // These lines are unecessary if BC are all periodic #ifndef USE_OBC std::vector boundary = {1,1,1,-1}; #else std::vector boundary = {1,1,1,0}; #endif FermionAction::ImplParams Params(boundary); FermionActionF::ImplParams ParamsF(boundary); double ActionStoppingCondition = 1e-8; double DerivativeStoppingCondition = 1e-8; double MaxCGIterations = 100000; //////////////////////////////////// // Collect actions //////////////////////////////////// ActionLevel Level1(1); ActionLevel Level2(HMCparams.SW); //////////////////////////////////// // Strange action //////////////////////////////////// typedef SchurDiagMooeeOperator LinearOperatorF; typedef SchurDiagMooeeOperator LinearOperatorD; typedef SchurDiagMooeeOperator LinearOperatorEOFAF; typedef SchurDiagMooeeOperator LinearOperatorEOFAD; typedef MixedPrecisionConjugateGradientOperatorFunction MxPCG; typedef MixedPrecisionConjugateGradientOperatorFunction MxPCG_EOFA; // DJM: setup for EOFA ratio (Mobius) OneFlavourRationalParams OFRp; OFRp.lo = 0.99; // How do I know this on F1? OFRp.hi = 20; OFRp.MaxIter = 100000; OFRp.tolerance= 1.0e-12; OFRp.degree = 12; OFRp.precision= 50; MobiusEOFAFermionD Strange_Op_L (U , *FGrid , *FrbGrid , *GridPtr , *GridRBPtr , strange_mass, strange_mass, charm_mass, 0.0, -1, M5, b, c); MobiusEOFAFermionF Strange_Op_LF(UF, *FGridF, *FrbGridF, *GridPtrF, *GridRBPtrF, strange_mass, strange_mass, charm_mass, 0.0, -1, M5, b, c); MobiusEOFAFermionD Strange_Op_R (U , *FGrid , *FrbGrid , *GridPtr , *GridRBPtr , charm_mass, strange_mass, charm_mass, -1.0, 1, M5, b, c); MobiusEOFAFermionF Strange_Op_RF(UF, *FGridF, *FrbGridF, *GridPtrF, *GridRBPtrF, charm_mass, strange_mass, charm_mass, -1.0, 1, M5, b, c); #ifdef EOFA_H MobiusEOFAFermionD Strange2_Op_L (U , *FGrid , *FrbGrid , *GridPtr , *GridRBPtr , eofa_mass, eofa_mass, charm_mass , 0.0, -1, M5, b, c); MobiusEOFAFermionF Strange2_Op_LF(UF, *FGridF, *FrbGridF, *GridPtrF, *GridRBPtrF, eofa_mass, eofa_mass, charm_mass , 0.0, -1, M5, b, c); MobiusEOFAFermionD Strange2_Op_R (U , *FGrid , *FrbGrid , *GridPtr , *GridRBPtr , charm_mass , eofa_mass, charm_mass , -1.0, 1, M5, b, c); MobiusEOFAFermionF Strange2_Op_RF(UF, *FGridF, *FrbGridF, *GridPtrF, *GridRBPtrF, charm_mass , eofa_mass, charm_mass , -1.0, 1, M5, b, c); #endif ConjugateGradient ActionCG(ActionStoppingCondition,MaxCGIterations); ConjugateGradient DerivativeCG(DerivativeStoppingCondition,MaxCGIterations); #ifdef MIXED_PRECISION const int MX_inner = 5000; // Mixed precision EOFA LinearOperatorEOFAD Strange_LinOp_L (Strange_Op_L); LinearOperatorEOFAD Strange_LinOp_R (Strange_Op_R); LinearOperatorEOFAF Strange_LinOp_LF(Strange_Op_LF); LinearOperatorEOFAF Strange_LinOp_RF(Strange_Op_RF); #ifdef EOFA_H // Mixed precision EOFA LinearOperatorEOFAD Strange2_LinOp_L (Strange2_Op_L); LinearOperatorEOFAD Strange2_LinOp_R (Strange2_Op_R); LinearOperatorEOFAF Strange2_LinOp_LF(Strange2_Op_LF); LinearOperatorEOFAF Strange2_LinOp_RF(Strange2_Op_RF); #endif MxPCG_EOFA ActionCGL(ActionStoppingCondition, MX_inner, MaxCGIterations, GridPtrF, FrbGridF, Strange_Op_LF,Strange_Op_L, Strange_LinOp_LF,Strange_LinOp_L); #ifdef EOFA_H MxPCG_EOFA ActionCGL2(ActionStoppingCondition, MX_inner, MaxCGIterations, GridPtrF, FrbGridF, Strange2_Op_LF,Strange2_Op_L, Strange2_LinOp_LF,Strange2_LinOp_L); #endif MxPCG_EOFA DerivativeCGL(DerivativeStoppingCondition, MX_inner, MaxCGIterations, GridPtrF, FrbGridF, Strange_Op_LF,Strange_Op_L, Strange_LinOp_LF,Strange_LinOp_L); #ifdef EOFA_H MxPCG_EOFA DerivativeCGL2(DerivativeStoppingCondition, MX_inner, MaxCGIterations, GridPtrF, FrbGridF, Strange2_Op_LF,Strange2_Op_L, Strange2_LinOp_LF,Strange2_LinOp_L); #endif MxPCG_EOFA ActionCGR(ActionStoppingCondition, MX_inner, MaxCGIterations, GridPtrF, FrbGridF, Strange_Op_RF,Strange_Op_R, Strange_LinOp_RF,Strange_LinOp_R); #ifdef EOFA_H MxPCG_EOFA ActionCGR2(ActionStoppingCondition, MX_inner, MaxCGIterations, GridPtrF, FrbGridF, Strange2_Op_RF,Strange2_Op_R, Strange2_LinOp_RF,Strange2_LinOp_R); #endif MxPCG_EOFA DerivativeCGR(DerivativeStoppingCondition, MX_inner, MaxCGIterations, GridPtrF, FrbGridF, Strange_Op_RF,Strange_Op_R, Strange_LinOp_RF,Strange_LinOp_R); #ifdef EOFA_H MxPCG_EOFA DerivativeCGR2(DerivativeStoppingCondition, MX_inner, MaxCGIterations, GridPtrF, FrbGridF, Strange2_Op_RF,Strange2_Op_R, Strange2_LinOp_RF,Strange2_LinOp_R); #endif ExactOneFlavourRatioPseudoFermionAction EOFA(Strange_Op_L, Strange_Op_R, ActionCG, ActionCGL, ActionCGR, DerivativeCGL, DerivativeCGR, OFRp, true); #ifdef EOFA_H ExactOneFlavourRatioPseudoFermionAction EOFA2(Strange2_Op_L, Strange2_Op_R, ActionCG, ActionCGL2, ActionCGR2, DerivativeCGL2, DerivativeCGR2, OFRp, true); #endif Level1.push_back(&EOFA); #ifdef EOFA_H Level1.push_back(&EOFA2); #endif #else ExactOneFlavourRatioPseudoFermionAction EOFA(Strange_Op_L, Strange_Op_R, ActionCG, ActionCG, ActionCG, ActionCG, ActionCG, // DerivativeCG, DerivativeCG, OFRp, true); Level1.push_back(&EOFA); #endif //////////////////////////////////// // up down action //////////////////////////////////// std::vector light_den; std::vector light_num; int n_hasenbusch = hasenbusch.size(); light_den.push_back(light_mass); for(int h=0;h Numerators; std::vector Denominators; std::vector *> Quotients; std::vector ActionMPCG; std::vector MPCG; std::vector DenominatorsF; std::vector LinOpD; std::vector LinOpF; for(int h=0;h(*Numerators[h],*Denominators[h],*MPCG[h],*ActionMPCG[h],ActionCG)); #else //////////////////////////////////////////////////////////////////////////// // Standard CG for 2f force //////////////////////////////////////////////////////////////////////////// Quotients.push_back (new TwoFlavourEvenOddRatioPseudoFermionAction(*Numerators[h],*Denominators[h],DerivativeCG,ActionCG)); #endif } for(int h=0;h S; #ifndef DO_IMPLICIT TrivialMetric Mtr; #else LaplacianRatParams gpar(2),mpar(2); gpar.offset = 1.; gpar.a0[0] = 500.; gpar.a1[0] = 0.; gpar.b0[0] = 0.25; gpar.b1[0] = 1.; gpar.a0[1] = -500.; gpar.a1[1] = 0.; gpar.b0[1] = 0.36; gpar.b1[1] = 1.2; gpar.b2=1.; mpar.offset = 1.; mpar.a0[0] = -0.850891906532; mpar.a1[0] = -1.54707654538; mpar. b0[0] = 2.85557166137; mpar. b1[0] = 5.74194794773; mpar.a0[1] = -13.5120056831218384729709214298; mpar.a1[1] = 1.54707654538396877086370295729; mpar.b0[1] = 19.2921090880640520026645390317; mpar.b1[1] = -3.54194794773029020262811172870; mpar.b2=1.; for(int i=0;i<2;i++){ gpar.a1[i] *=16.; gpar.b1[i] *=16.; mpar.a1[i] *=16.; mpar.b1[i] *=16.; } gpar.b2 *= 16.*16.; mpar.b2 *= 16.*16.; ConjugateGradient CG(1.0e-8,10000); LaplacianParams LapPar(0.0001, 1.0, 10000, 1e-8, 12, 64); std::cout << GridLogMessage << "LaplacianRat " << std::endl; gpar.tolerance=HMCparams.MD.RMHMCCGTol; mpar.tolerance=HMCparams.MD.RMHMCCGTol; std::cout << GridLogMessage << "gpar offset= " << gpar.offset <