/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: Copyright (C) 2015-2016 Author: Peter Boyle Author: Guido Cossu Author: David Murphy 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 namespace Grid{ namespace QCD{ /* * 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; 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 //Option to speed up *inner single precision* solves using a LinearFunction that produces a guess LinearFunction *guesser; 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.), guesser(NULL) { std::cout << GridLogMessage << " Mixed precision CG wrapper LinOpF " < &g){ guesser = &g; } void operator()(LinearOperatorBase &LinOpU, const FieldD &src, FieldD &psi) { SchurOperatorD * SchurOpU = static_cast(&LinOpU); std::cout << GridLogMessage << " Mixed precision CG wrapper FermOpD " <_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); //////////////////////////////////////////////////////////////////////////////////// // Could test to make sure that LinOpF and LinOpD agree to single prec? //////////////////////////////////////////////////////////////////////////////////// GridBase *Fgrid = psi._grid; FieldD tmp2(Fgrid); FieldD tmp1(Fgrid); LinOpU.Op(src,tmp1); LinOpD.Op(src,tmp2); std::cout << " Double gauge field "<< norm2(FermOpD.Umu)< MPCG(Tolerance,MaxInnerIterations,MaxOuterIterations,SinglePrecGrid5,LinOpF,LinOpD); std::cout << GridLogMessage << "Calling mixed precision Conjugate Gradient" < HMCWrapper; // MD.name = std::string("Leap Frog"); typedef GenericHMCRunner HMCWrapper; MD.name = std::string("Force Gradient"); // typedef GenericHMCRunner HMCWrapper; // MD.name = std::string("MinimumNorm2"); MD.MDsteps = 6; MD.trajL = 1.0; HMCparameters HMCparams; HMCparams.StartTrajectory = 530; HMCparams.Trajectories = 1000; HMCparams.NoMetropolisUntil= 0; // "[HotStart, ColdStart, TepidStart, CheckpointStart]\n"; // HMCparams.StartingType =std::string("ColdStart"); HMCparams.StartingType =std::string("CheckpointStart"); HMCparams.MD = MD; HMCWrapper TheHMC(HMCparams); // Grid from the command line arguments --grid and --mpi TheHMC.Resources.AddFourDimGrid("gauge"); // use default simd lanes decomposition CheckpointerParameters CPparams; CPparams.config_prefix = "ckpoint_EODWF_lat"; CPparams.rng_prefix = "ckpoint_EODWF_rng"; CPparams.saveInterval = 10; CPparams.format = "IEEE64BIG"; TheHMC.Resources.LoadNerscCheckpointer(CPparams); RNGModuleParameters RNGpar; RNGpar.serial_seeds = "1 2 3 4 5"; RNGpar.parallel_seeds = "6 7 8 9 10"; TheHMC.Resources.SetRNGSeeds(RNGpar); // Construct observables // here there is too much indirection typedef PlaquetteMod PlaqObs; TheHMC.Resources.AddObservable(); ////////////////////////////////////////////// const int Ls = 16; Real beta = 2.13; Real light_mass = 0.01; Real strange_mass = 0.04; Real pv_mass = 1.0; RealD M5 = 1.8; RealD b = 1.0; RealD c = 0.0; std::vector hasenbusch({ 0.1, 0.3, 0.6 }); auto GridPtr = TheHMC.Resources.GetCartesian(); auto GridRBPtr = TheHMC.Resources.GetRBCartesian(); auto FGrid = SpaceTimeGrid::makeFiveDimGrid(Ls,GridPtr); auto FrbGrid = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls,GridPtr); std::vector latt = GridDefaultLatt(); std::vector mpi = GridDefaultMpi(); std::vector simdF = GridDefaultSimd(Nd,vComplexF::Nsimd()); std::vector simdD = GridDefaultSimd(Nd,vComplexD::Nsimd()); auto GridPtrF = SpaceTimeGrid::makeFourDimGrid(latt,simdF,mpi); auto GridRBPtrF = SpaceTimeGrid::makeFourDimRedBlackGrid(GridPtrF); auto FGridF = SpaceTimeGrid::makeFiveDimGrid(Ls,GridPtrF); auto FrbGridF = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls,GridPtrF); IwasakiGaugeActionR GaugeAction(beta); // temporarily need a gauge field LatticeGaugeField U(GridPtr); LatticeGaugeFieldF UF(GridPtrF); // These lines are unecessary if BC are all periodic std::vector boundary = {1,1,1,-1}; FermionAction::ImplParams Params(boundary); FermionActionF::ImplParams ParamsF(boundary); double ActionStoppingCondition = 1e-10; double DerivativeStoppingCondition = 1e-6; double MaxCGIterations = 30000; ConjugateGradient ActionCG(ActionStoppingCondition,MaxCGIterations); ConjugateGradient DerivativeCG(DerivativeStoppingCondition,MaxCGIterations); //////////////////////////////////// // Collect actions //////////////////////////////////// ActionLevel Level1(1); ActionLevel Level2(8); //////////////////////////////////// // Strange action //////////////////////////////////// // DJM: setup for EOFA ratio (Mobius) OneFlavourRationalParams OFRp; OFRp.lo = 0.1; OFRp.hi = 25.0; OFRp.MaxIter = 10000; OFRp.tolerance= 1.0e-9; OFRp.degree = 14; OFRp.precision= 50; MobiusEOFAFermionR Strange_Op_L(U, *FGrid, *FrbGrid, *GridPtr, *GridRBPtr, strange_mass, strange_mass, pv_mass, 0.0, -1, M5, b, c); MobiusEOFAFermionR Strange_Op_R(U, *FGrid, *FrbGrid, *GridPtr, *GridRBPtr, pv_mass, strange_mass, pv_mass, -1.0, 1, M5, b, c); ExactOneFlavourRatioPseudoFermionAction EOFA(Strange_Op_L, Strange_Op_R, ActionCG, DerivativeCG, OFRp, true); Level1.push_back(&EOFA); //////////////////////////////////// // 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 LinearOperatorF; typedef SchurDiagMooeeOperator LinearOperatorD; std::vector Numerators; std::vector Denominators; std::vector LinOpD; std::vector DenominatorsF; std::vector LinOpF; typedef MixedPrecisionConjugateGradientOperatorFunction MxPCG; std::vector MPCG; std::vector *> Quotients; for(int h=0;h(*Numerators[h],*Denominators[h],DerivativeCG,ActionCG)); #else //////////////////////////////////////////////////////////////////////////// // Mixed precision CG for 2f force //////////////////////////////////////////////////////////////////////////// DenominatorsF.push_back(new FermionActionF(UF,*FGridF,*FrbGridF,*GridPtrF,*GridRBPtrF,light_den[h],M5,b,c, ParamsF)); LinOpF.push_back(new LinearOperatorF(*DenominatorsF[h])); LinOpD.push_back(new LinearOperatorD(*Denominators[h])); MPCG.push_back(new MxPCG(DerivativeStoppingCondition, 200, MaxCGIterations, GridPtrF, FrbGridF, *DenominatorsF[h],*Denominators[h], *LinOpF[h], *LinOpD[h]) ); Quotients.push_back (new TwoFlavourEvenOddRatioPseudoFermionAction(*Numerators[h],*Denominators[h],*MPCG[h],ActionCG)); #endif } for(int h=0;h