mirror of
				https://github.com/paboyle/Grid.git
				synced 2025-11-03 21:44:33 +00:00 
			
		
		
		
	Fixes to fermion force terms after sign of gamma_mu (0...3) change.
Thought I had already committed these.
Believe I have got the Gparity fermion force working.
* tests/Test_gpdwf_force.cc     -- correctly predicts dS for two flavour pseudofermion
                                   based on a small dt update of U field.
* tests/Test_hmc_EODWFRatio_Gparity.cc -- ran 1 trajectory on 8^4 with dH=0.21.
Need to accumulate a full plaquette log to believe fully which will take some hours of run time.
			
			
This commit is contained in:
		@@ -1,5 +1,5 @@
 | 
			
		||||
 | 
			
		||||
bin_PROGRAMS = Test_GaugeAction Test_cayley_cg Test_cayley_coarsen_support Test_cayley_even_odd Test_cayley_ldop_cr Test_cf_coarsen_support Test_cf_cr_unprec Test_cheby Test_contfrac_cg Test_contfrac_even_odd Test_contfrac_force Test_cshift Test_cshift_red_black Test_dwf_cg_prec Test_dwf_cg_schur Test_dwf_cg_unprec Test_dwf_cr_unprec Test_dwf_even_odd Test_dwf_force Test_dwf_fpgcr Test_dwf_hdcr Test_dwf_lanczos Test_gamma Test_hmc_EODWFRatio Test_hmc_EOWilsonFermionGauge Test_hmc_EOWilsonRatio Test_hmc_WilsonFermionGauge Test_hmc_WilsonGauge Test_hmc_WilsonRatio Test_lie_generators Test_main Test_multishift_sqrt Test_nersc_io Test_partfrac_force Test_quenched_update Test_remez Test_rhmc_EOWilson1p1 Test_rhmc_EOWilsonRatio Test_rhmc_Wilson1p1 Test_rhmc_WilsonRatio Test_rng Test_rng_fixed Test_serialisation Test_simd Test_stencil Test_synthetic_lanczos Test_wilson_cg_prec Test_wilson_cg_schur Test_wilson_cg_unprec Test_wilson_cr_unprec Test_wilson_even_odd Test_wilson_force Test_wilson_force_phiMdagMphi Test_wilson_force_phiMphi 
 | 
			
		||||
bin_PROGRAMS = Test_GaugeAction Test_cayley_cg Test_cayley_coarsen_support Test_cayley_even_odd Test_cayley_ldop_cr Test_cf_coarsen_support Test_cf_cr_unprec Test_cheby Test_contfrac_cg Test_contfrac_even_odd Test_contfrac_force Test_cshift Test_cshift_red_black Test_dwf_cg_prec Test_dwf_cg_schur Test_dwf_cg_unprec Test_dwf_cr_unprec Test_dwf_even_odd Test_dwf_force Test_dwf_fpgcr Test_dwf_gpforce Test_dwf_hdcr Test_dwf_lanczos Test_gamma Test_gparity Test_gpdwf_force Test_gpwilson_even_odd Test_hmc_EODWFRatio Test_hmc_EODWFRatio_Gparity Test_hmc_EOWilsonFermionGauge Test_hmc_EOWilsonRatio Test_hmc_WilsonFermionGauge Test_hmc_WilsonGauge Test_hmc_WilsonRatio Test_lie_generators Test_main Test_multishift_sqrt Test_nersc_io Test_partfrac_force Test_quenched_update Test_remez Test_rhmc_EOWilson1p1 Test_rhmc_EOWilsonRatio Test_rhmc_Wilson1p1 Test_rhmc_WilsonRatio Test_rng Test_rng_fixed Test_serialisation Test_simd Test_stencil Test_synthetic_lanczos Test_wilson_cg_prec Test_wilson_cg_schur Test_wilson_cg_unprec Test_wilson_cr_unprec Test_wilson_even_odd Test_wilson_force Test_wilson_force_phiMdagMphi Test_wilson_force_phiMphi 
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
Test_GaugeAction_SOURCES=Test_GaugeAction.cc
 | 
			
		||||
@@ -82,6 +82,10 @@ Test_dwf_fpgcr_SOURCES=Test_dwf_fpgcr.cc
 | 
			
		||||
Test_dwf_fpgcr_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
Test_dwf_gpforce_SOURCES=Test_dwf_gpforce.cc
 | 
			
		||||
Test_dwf_gpforce_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
Test_dwf_hdcr_SOURCES=Test_dwf_hdcr.cc
 | 
			
		||||
Test_dwf_hdcr_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
@@ -98,6 +102,10 @@ Test_gparity_SOURCES=Test_gparity.cc
 | 
			
		||||
Test_gparity_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
Test_gpdwf_force_SOURCES=Test_gpdwf_force.cc
 | 
			
		||||
Test_gpdwf_force_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
Test_gpwilson_even_odd_SOURCES=Test_gpwilson_even_odd.cc
 | 
			
		||||
Test_gpwilson_even_odd_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
@@ -106,6 +114,10 @@ Test_hmc_EODWFRatio_SOURCES=Test_hmc_EODWFRatio.cc
 | 
			
		||||
Test_hmc_EODWFRatio_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
Test_hmc_EODWFRatio_Gparity_SOURCES=Test_hmc_EODWFRatio_Gparity.cc
 | 
			
		||||
Test_hmc_EODWFRatio_Gparity_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
Test_hmc_EOWilsonFermionGauge_SOURCES=Test_hmc_EOWilsonFermionGauge.cc
 | 
			
		||||
Test_hmc_EOWilsonFermionGauge_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
@@ -225,10 +237,7 @@ Test_wilson_force_phiMdagMphi_LDADD=-lGrid
 | 
			
		||||
Test_wilson_force_phiMphi_SOURCES=Test_wilson_force_phiMphi.cc
 | 
			
		||||
Test_wilson_force_phiMphi_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
Test_zmm_SOURCES=Test_zmm.cc
 | 
			
		||||
Test_zmm_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
Test_RectPlaq_SOURCES=Test_RectPlaq.cc
 | 
			
		||||
Test_RectPlaq_LDADD=-lGrid
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 
 | 
			
		||||
@@ -8,5 +8,8 @@ endif
 | 
			
		||||
AM_CXXFLAGS = -I$(top_srcdir)/lib
 | 
			
		||||
AM_LDFLAGS = -L$(top_builddir)/lib
 | 
			
		||||
 | 
			
		||||
if BUILD_ZMM
 | 
			
		||||
  bin_PROGRAMS=Test_zmm
 | 
			
		||||
endif
 | 
			
		||||
 | 
			
		||||
include Make.inc
 | 
			
		||||
 
 | 
			
		||||
							
								
								
									
										191
									
								
								tests/Test_dwf_gpforce.cc
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										191
									
								
								tests/Test_dwf_gpforce.cc
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,191 @@
 | 
			
		||||
#include <Grid.h>
 | 
			
		||||
 | 
			
		||||
using namespace std;
 | 
			
		||||
using namespace Grid;
 | 
			
		||||
using namespace Grid::QCD;
 | 
			
		||||
 | 
			
		||||
#define parallel_for PARALLEL_FOR_LOOP for
 | 
			
		||||
 | 
			
		||||
int main (int argc, char ** argv)
 | 
			
		||||
{
 | 
			
		||||
  Grid_init(&argc,&argv);
 | 
			
		||||
 | 
			
		||||
  std::vector<int> latt_size   = GridDefaultLatt();
 | 
			
		||||
  std::vector<int> simd_layout = GridDefaultSimd(Nd,vComplex::Nsimd());
 | 
			
		||||
  std::vector<int> mpi_layout  = GridDefaultMpi();
 | 
			
		||||
 | 
			
		||||
  const int Ls=8;
 | 
			
		||||
  GridCartesian         * UGrid   = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd,vComplex::Nsimd()),GridDefaultMpi());
 | 
			
		||||
  GridRedBlackCartesian * UrbGrid = SpaceTimeGrid::makeFourDimRedBlackGrid(UGrid);
 | 
			
		||||
  GridCartesian         * FGrid   = SpaceTimeGrid::makeFiveDimGrid(Ls,UGrid);
 | 
			
		||||
  GridRedBlackCartesian * FrbGrid = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls,UGrid);
 | 
			
		||||
 | 
			
		||||
  typedef typename GparityDomainWallFermionR::FermionField FermionField;
 | 
			
		||||
 | 
			
		||||
  int threads = GridThread::GetThreads();
 | 
			
		||||
  std::cout<<GridLogMessage << "Grid is setup to use "<<threads<<" threads"<<std::endl;
 | 
			
		||||
 | 
			
		||||
  std::vector<int> seeds({1,2,3,4});
 | 
			
		||||
 | 
			
		||||
  GridParallelRNG          RNG5(FGrid);  RNG5.SeedRandomDevice();
 | 
			
		||||
  GridParallelRNG          RNG4(UGrid);  RNG4.SeedRandomDevice();
 | 
			
		||||
  
 | 
			
		||||
  FermionField phi        (FGrid); gaussian(RNG5,phi);
 | 
			
		||||
  FermionField Mphi       (FGrid); 
 | 
			
		||||
  FermionField MphiPrime  (FGrid); 
 | 
			
		||||
 | 
			
		||||
  LatticeGaugeField U(UGrid);
 | 
			
		||||
 | 
			
		||||
  SU3::HotConfiguration(RNG4,U);
 | 
			
		||||
  //  SU3::ColdConfiguration(pRNG,U);
 | 
			
		||||
  
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  // Unmodified matrix element
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  RealD mass=0.2; //kills the diagonal term
 | 
			
		||||
  RealD M5=1.8;
 | 
			
		||||
  //  const int nu = 3;
 | 
			
		||||
  //  std::vector<int> twists(Nd,0); // twists[nu] = 1;
 | 
			
		||||
  //  GparityDomainWallFermionR::ImplParams params;  params.twists = twists;
 | 
			
		||||
  //  GparityDomainWallFermionR Ddwf(U,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5,params);
 | 
			
		||||
 | 
			
		||||
  //  DomainWallFermionR Dw     (U,     Grid,RBGrid,mass,M5);
 | 
			
		||||
 | 
			
		||||
  const int nu = 3;
 | 
			
		||||
  std::vector<int> twists(Nd,0);
 | 
			
		||||
  twists[nu] = 1;
 | 
			
		||||
  GparityDomainWallFermionR::ImplParams params;
 | 
			
		||||
  params.twists = twists;
 | 
			
		||||
 | 
			
		||||
  GparityDomainWallFermionR Dw(U,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5,params);
 | 
			
		||||
 | 
			
		||||
  Dw.M   (phi,Mphi);
 | 
			
		||||
 | 
			
		||||
  ComplexD S    = innerProduct(Mphi,Mphi); // pdag MdagM p
 | 
			
		||||
 | 
			
		||||
  // get the deriv of phidag MdagM phi with respect to "U"
 | 
			
		||||
  LatticeGaugeField UdSdU(UGrid);
 | 
			
		||||
  LatticeGaugeField tmp(UGrid);
 | 
			
		||||
 | 
			
		||||
  Dw.MDeriv(tmp , Mphi,  phi,DaggerNo );  UdSdU=tmp;
 | 
			
		||||
  Dw.MDeriv(tmp , phi,  Mphi,DaggerYes ); UdSdU=(UdSdU+tmp);
 | 
			
		||||
  
 | 
			
		||||
  FermionField Ftmp      (FGrid);
 | 
			
		||||
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  // Modify the gauge field a little 
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  RealD dt = 0.0001;
 | 
			
		||||
  RealD Hmom = 0.0;
 | 
			
		||||
  RealD Hmomprime = 0.0;
 | 
			
		||||
  RealD Hmompp    = 0.0;
 | 
			
		||||
  LatticeColourMatrix mommu(UGrid); 
 | 
			
		||||
  LatticeColourMatrix forcemu(UGrid); 
 | 
			
		||||
  LatticeGaugeField mom(UGrid); 
 | 
			
		||||
  LatticeGaugeField Uprime(UGrid); 
 | 
			
		||||
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
 | 
			
		||||
    SU3::GaussianLieAlgebraMatrix(RNG4, mommu); // Traceless antihermitian momentum; gaussian in lie alg
 | 
			
		||||
 | 
			
		||||
    Hmom -= real(sum(trace(mommu*mommu)));
 | 
			
		||||
 | 
			
		||||
    PokeIndex<LorentzIndex>(mom,mommu,mu);
 | 
			
		||||
 | 
			
		||||
    // fourth order exponential approx
 | 
			
		||||
    parallel_for(auto i=mom.begin();i<mom.end();i++){
 | 
			
		||||
      Uprime[i](mu) =
 | 
			
		||||
	  U[i](mu)
 | 
			
		||||
	+ mom[i](mu)*U[i](mu)*dt 
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt/2.0)
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt*dt/6.0)
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt*dt*dt/24.0)
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt*dt*dt*dt/120.0)
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt*dt*dt*dt*dt/720.0)
 | 
			
		||||
	;
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage <<"Initial mom hamiltonian is "<< Hmom <<std::endl;
 | 
			
		||||
  Dw.ImportGauge(Uprime);
 | 
			
		||||
  Dw.M          (phi,MphiPrime);
 | 
			
		||||
 | 
			
		||||
  ComplexD Sprime    = innerProduct(MphiPrime   ,MphiPrime);
 | 
			
		||||
 | 
			
		||||
  //////////////////////////////////////////////
 | 
			
		||||
  // Use derivative to estimate dS
 | 
			
		||||
  //////////////////////////////////////////////
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
    std::cout << "" <<std::endl;
 | 
			
		||||
    mommu   = PeekIndex<LorentzIndex>(mom,mu);
 | 
			
		||||
    std::cout << GridLogMessage<< " Mommu  " << norm2(mommu)<<std::endl;
 | 
			
		||||
    mommu   = mommu+adj(mommu);
 | 
			
		||||
    std::cout << GridLogMessage<< " Mommu + Mommudag " << norm2(mommu)<<std::endl;
 | 
			
		||||
    mommu   = PeekIndex<LorentzIndex>(UdSdU,mu);
 | 
			
		||||
    std::cout << GridLogMessage<< " dsdumu  " << norm2(mommu)<<std::endl;
 | 
			
		||||
    mommu   = mommu+adj(mommu);
 | 
			
		||||
    std::cout << GridLogMessage<< " dsdumu + dag  " << norm2(mommu)<<std::endl;
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  LatticeComplex dS(UGrid); dS = zero;
 | 
			
		||||
  LatticeComplex dSmom(UGrid); dSmom = zero;
 | 
			
		||||
  LatticeComplex dSmom2(UGrid); dSmom2 = zero;
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
    mommu   = PeekIndex<LorentzIndex>(UdSdU,mu);
 | 
			
		||||
    mommu=Ta(mommu)*2.0;
 | 
			
		||||
    PokeIndex<LorentzIndex>(UdSdU,mommu,mu);
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
    mommu   = PeekIndex<LorentzIndex>(mom,mu);
 | 
			
		||||
    std::cout << GridLogMessage<< " Mommu  " << norm2(mommu)<<std::endl;
 | 
			
		||||
    mommu   = mommu+adj(mommu);
 | 
			
		||||
    std::cout << GridLogMessage<< " Mommu + Mommudag " << norm2(mommu)<<std::endl;
 | 
			
		||||
    mommu   = PeekIndex<LorentzIndex>(UdSdU,mu);
 | 
			
		||||
    std::cout << GridLogMessage<< " dsdumu  " << norm2(mommu)<<std::endl;
 | 
			
		||||
    mommu   = mommu+adj(mommu);
 | 
			
		||||
    std::cout << GridLogMessage<< " dsdumu + dag  " << norm2(mommu)<<std::endl;
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
    forcemu = PeekIndex<LorentzIndex>(UdSdU,mu);
 | 
			
		||||
    mommu   = PeekIndex<LorentzIndex>(mom,mu);
 | 
			
		||||
 | 
			
		||||
    // Update PF action density
 | 
			
		||||
    dS = dS+trace(mommu*forcemu)*dt;
 | 
			
		||||
 | 
			
		||||
    dSmom  = dSmom  - trace(mommu*forcemu) * dt;
 | 
			
		||||
    dSmom2 = dSmom2 - trace(forcemu*forcemu) *(0.25* dt*dt);
 | 
			
		||||
 | 
			
		||||
    // Update mom action density
 | 
			
		||||
    mommu = mommu + forcemu*(dt*0.5);
 | 
			
		||||
 | 
			
		||||
    Hmomprime -= real(sum(trace(mommu*mommu)));
 | 
			
		||||
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  Complex dSpred    = sum(dS);
 | 
			
		||||
  Complex dSm       = sum(dSmom);
 | 
			
		||||
  Complex dSm2      = sum(dSmom2);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage <<"Initial mom hamiltonian is "<< Hmom <<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage <<"Final   mom hamiltonian is "<< Hmomprime <<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage <<"Delta   mom hamiltonian is "<< Hmomprime-Hmom <<std::endl;
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage << " S      "<<S<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << " Sprime "<<Sprime<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << "dS      "<<Sprime-S<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << "predict dS    "<< dSpred <<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage <<"dSm "<< dSm<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage <<"dSm2"<< dSm2<<std::endl;
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage << "Total dS    "<< Hmomprime - Hmom + Sprime - S <<std::endl;
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  std::cout<< GridLogMessage << "Done" <<std::endl;
 | 
			
		||||
  Grid_finalize();
 | 
			
		||||
}
 | 
			
		||||
@@ -13,8 +13,6 @@ int main (int argc, char ** argv)
 | 
			
		||||
 | 
			
		||||
  Grid_init(&argc,&argv);
 | 
			
		||||
 | 
			
		||||
  std::vector<int> twists(Nd,0);
 | 
			
		||||
  twists[nu] = 1;
 | 
			
		||||
 | 
			
		||||
  const int Ls=4;
 | 
			
		||||
  const int L =4;
 | 
			
		||||
@@ -99,6 +97,9 @@ int main (int argc, char ** argv)
 | 
			
		||||
  ConjugateGradient<LatticeFermion> CG(1.0e-8,10000);
 | 
			
		||||
  CG(HermOpEO,src_o_1f,result_o_1f);
 | 
			
		||||
  
 | 
			
		||||
  //  const int nu = 3;
 | 
			
		||||
  std::vector<int> twists(Nd,0);
 | 
			
		||||
  twists[nu] = 1;
 | 
			
		||||
  GparityDomainWallFermionR::ImplParams params;
 | 
			
		||||
  params.twists = twists;
 | 
			
		||||
  GparityDomainWallFermionR GPDdwf(Umu_2f,*FGrid_2f,*FrbGrid_2f,*UGrid_2f,*UrbGrid_2f,mass,M5,params);
 | 
			
		||||
 
 | 
			
		||||
							
								
								
									
										179
									
								
								tests/Test_gpdwf_force.cc
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										179
									
								
								tests/Test_gpdwf_force.cc
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,179 @@
 | 
			
		||||
#include <Grid.h>
 | 
			
		||||
 | 
			
		||||
using namespace std;
 | 
			
		||||
using namespace Grid;
 | 
			
		||||
using namespace Grid::QCD;
 | 
			
		||||
 | 
			
		||||
#define parallel_for PARALLEL_FOR_LOOP for
 | 
			
		||||
 | 
			
		||||
int main (int argc, char ** argv)
 | 
			
		||||
{
 | 
			
		||||
  Grid_init(&argc,&argv);
 | 
			
		||||
 | 
			
		||||
  std::vector<int> latt_size   = GridDefaultLatt();
 | 
			
		||||
  std::vector<int> simd_layout = GridDefaultSimd(Nd,vComplex::Nsimd());
 | 
			
		||||
  std::vector<int> mpi_layout  = GridDefaultMpi();
 | 
			
		||||
 | 
			
		||||
  const int Ls=8;
 | 
			
		||||
 | 
			
		||||
  GridCartesian         * UGrid   = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd,vComplex::Nsimd()),GridDefaultMpi());
 | 
			
		||||
  GridRedBlackCartesian * UrbGrid = SpaceTimeGrid::makeFourDimRedBlackGrid(UGrid);
 | 
			
		||||
  GridCartesian         * FGrid   = SpaceTimeGrid::makeFiveDimGrid(Ls,UGrid);
 | 
			
		||||
  GridRedBlackCartesian * FrbGrid = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls,UGrid);
 | 
			
		||||
 | 
			
		||||
  std::vector<int> seeds4({1,2,3,4});
 | 
			
		||||
  std::vector<int> seeds5({5,6,7,8});
 | 
			
		||||
  GridParallelRNG          RNG5(FGrid);  RNG5.SeedFixedIntegers(seeds5);
 | 
			
		||||
  GridParallelRNG          RNG4(UGrid);  RNG4.SeedFixedIntegers(seeds4);
 | 
			
		||||
 | 
			
		||||
  int threads = GridThread::GetThreads();
 | 
			
		||||
  std::cout<<GridLogMessage << "Grid is setup to use "<<threads<<" threads"<<std::endl;
 | 
			
		||||
 | 
			
		||||
  typedef typename GparityDomainWallFermionR::FermionField FermionField;
 | 
			
		||||
  FermionField phi        (FGrid); gaussian(RNG5,phi);
 | 
			
		||||
  FermionField Mphi       (FGrid); 
 | 
			
		||||
  FermionField MphiPrime  (FGrid); 
 | 
			
		||||
 | 
			
		||||
  LatticeGaugeField U(UGrid);
 | 
			
		||||
 | 
			
		||||
  SU3::HotConfiguration(RNG4,U);
 | 
			
		||||
  
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  // Unmodified matrix element
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  RealD mass=0.01; 
 | 
			
		||||
  RealD M5=1.8; 
 | 
			
		||||
 | 
			
		||||
  const int nu = 3;
 | 
			
		||||
  std::vector<int> twists(Nd,0);  twists[nu] = 1;
 | 
			
		||||
  GparityDomainWallFermionR::ImplParams params;  params.twists = twists;
 | 
			
		||||
  GparityDomainWallFermionR Ddwf(U,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5,params);
 | 
			
		||||
  Ddwf.M   (phi,Mphi);
 | 
			
		||||
 | 
			
		||||
  ComplexD S    = innerProduct(Mphi,Mphi); // pdag MdagM p
 | 
			
		||||
 | 
			
		||||
  // get the deriv of phidag MdagM phi with respect to "U"
 | 
			
		||||
  LatticeGaugeField UdSdU(UGrid);
 | 
			
		||||
  LatticeGaugeField tmp(UGrid);
 | 
			
		||||
 | 
			
		||||
  Ddwf.MDeriv(tmp , Mphi,  phi,DaggerNo );  UdSdU=tmp;
 | 
			
		||||
  Ddwf.MDeriv(tmp , phi,  Mphi,DaggerYes ); UdSdU=(UdSdU+tmp);  
 | 
			
		||||
  
 | 
			
		||||
  FermionField Ftmp      (FGrid);
 | 
			
		||||
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  // Modify the gauge field a little 
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  RealD dt = 0.001;
 | 
			
		||||
 | 
			
		||||
  LatticeColourMatrix mommu(UGrid); 
 | 
			
		||||
  LatticeColourMatrix forcemu(UGrid); 
 | 
			
		||||
  LatticeGaugeField mom(UGrid); 
 | 
			
		||||
  LatticeGaugeField Uprime(UGrid); 
 | 
			
		||||
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
 | 
			
		||||
    SU3::GaussianLieAlgebraMatrix(RNG4, mommu); // Traceless antihermitian momentum; gaussian in lie alg
 | 
			
		||||
 | 
			
		||||
    PokeIndex<LorentzIndex>(mom,mommu,mu);
 | 
			
		||||
 | 
			
		||||
    // fourth order exponential approx
 | 
			
		||||
    parallel_for(auto i=mom.begin();i<mom.end();i++){
 | 
			
		||||
      Uprime[i](mu) =
 | 
			
		||||
	  U[i](mu)
 | 
			
		||||
	+ mom[i](mu)*U[i](mu)*dt 
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt/2.0)
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt*dt/6.0)
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt*dt*dt/24.0)
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt*dt*dt*dt/120.0)
 | 
			
		||||
	+ mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *mom[i](mu) *U[i](mu)*(dt*dt*dt*dt*dt*dt/720.0)
 | 
			
		||||
	;
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
  }
 | 
			
		||||
  
 | 
			
		||||
  Ddwf.ImportGauge(Uprime);
 | 
			
		||||
  Ddwf.M          (phi,MphiPrime);
 | 
			
		||||
 | 
			
		||||
  ComplexD Sprime    = innerProduct(MphiPrime   ,MphiPrime);
 | 
			
		||||
 | 
			
		||||
  //////////////////////////////////////////////
 | 
			
		||||
  // Use derivative to estimate dS
 | 
			
		||||
  //////////////////////////////////////////////
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  LatticeComplex dS(UGrid); dS = zero;
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
    mommu   = PeekIndex<LorentzIndex>(UdSdU,mu);
 | 
			
		||||
    mommu=Ta(mommu)*2.0;
 | 
			
		||||
    PokeIndex<LorentzIndex>(UdSdU,mommu,mu);
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
    forcemu = PeekIndex<LorentzIndex>(UdSdU,mu);
 | 
			
		||||
    mommu   = PeekIndex<LorentzIndex>(mom,mu);
 | 
			
		||||
 | 
			
		||||
    // Update PF action density
 | 
			
		||||
    dS = dS+trace(mommu*forcemu)*dt;
 | 
			
		||||
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  Complex dSpred    = sum(dS);
 | 
			
		||||
 | 
			
		||||
  // From TwoFlavourPseudoFermion:
 | 
			
		||||
  //////////////////////////////////////////////////////
 | 
			
		||||
  // dS/du = - phi^dag  (Mdag M)^-1 [ Mdag dM + dMdag M ]  (Mdag M)^-1 phi
 | 
			
		||||
  //       = - phi^dag M^-1 dM (MdagM)^-1 phi -  phi^dag (MdagM)^-1 dMdag dM (Mdag)^-1 phi 
 | 
			
		||||
  //
 | 
			
		||||
  //       = - Ydag dM X  - Xdag dMdag Y
 | 
			
		||||
  //
 | 
			
		||||
  //////////////////////////////////////////////////////
 | 
			
		||||
  // Our conventions really make this UdSdU; We do not differentiate wrt Udag here.
 | 
			
		||||
  // So must take dSdU - adj(dSdU) and left multiply by mom to get dS/dt.
 | 
			
		||||
  //
 | 
			
		||||
  //  When we have Gparity -- U and Uconj enter.
 | 
			
		||||
  //
 | 
			
		||||
  // dU/dt  = p U
 | 
			
		||||
  // dUc/dt = p* Uc     // Is p real, traceless, etc..
 | 
			
		||||
  //
 | 
			
		||||
  // dS/dt = dUdt dSdU = p U dSdU
 | 
			
		||||
  //
 | 
			
		||||
  // Gparity --- deriv is pc Uc dSdUc + p U dSdU 
 | 
			
		||||
  //
 | 
			
		||||
  //	Pmu = zero;
 | 
			
		||||
  //	for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
  //	  SU<Ncol>::GaussianLieAlgebraMatrix(pRNG, Pmu);
 | 
			
		||||
  //	  PokeIndex<LorentzIndex>(P, Pmu, mu);
 | 
			
		||||
  //	}
 | 
			
		||||
  //
 | 
			
		||||
  //
 | 
			
		||||
  //    GridBase *grid = out._grid;
 | 
			
		||||
  //    LatticeReal ca (grid);
 | 
			
		||||
  //    LatticeMatrix  la (grid);
 | 
			
		||||
  //    Complex ci(0.0,scale);
 | 
			
		||||
  //    Matrix ta;
 | 
			
		||||
  //    out=zero;
 | 
			
		||||
  //    for(int a=0;a<generators();a++){
 | 
			
		||||
  //      gaussian(pRNG,ca); 
 | 
			
		||||
  //      generator(a,ta);
 | 
			
		||||
  //      la=toComplex(ca)*ci*ta;  // i t_a Lambda_a c_a // c_a is gaussian
 | 
			
		||||
  //      out += la; 
 | 
			
		||||
  //    }
 | 
			
		||||
  //  p = sum_a i gauss_a t_a 
 | 
			
		||||
  //
 | 
			
		||||
  // dU = p U dt
 | 
			
		||||
  //
 | 
			
		||||
  // dUc = p^c Uc dt = -i gauss_a t_a^c Uc 
 | 
			
		||||
  // 
 | 
			
		||||
  //
 | 
			
		||||
  // For Gparity the dS /dt from Uc links 
 | 
			
		||||
  // 
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage << " S      "<<S<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << " Sprime "<<Sprime<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << "dS      "<<Sprime-S<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << "predict dS    "<< dSpred <<std::endl;
 | 
			
		||||
 | 
			
		||||
  std::cout<< GridLogMessage << "Done" <<std::endl;
 | 
			
		||||
  Grid_finalize();
 | 
			
		||||
}
 | 
			
		||||
@@ -52,7 +52,10 @@ int main (int argc, char ** argv)
 | 
			
		||||
 | 
			
		||||
  RealD mass=0.1;
 | 
			
		||||
 | 
			
		||||
  GparityWilsonFermionR Dw(Umu,Grid,RBGrid,mass);
 | 
			
		||||
  GparityWilsonFermionR::ImplParams params;
 | 
			
		||||
  std::vector<int> twists(Nd,0);  twists[1] = 1;
 | 
			
		||||
  params.twists = twists;
 | 
			
		||||
  GparityWilsonFermionR Dw(Umu,Grid,RBGrid,mass,params);
 | 
			
		||||
 | 
			
		||||
  FermionField src_e   (&RBGrid);
 | 
			
		||||
  FermionField src_o   (&RBGrid);
 | 
			
		||||
 
 | 
			
		||||
							
								
								
									
										70
									
								
								tests/Test_hmc_EODWFRatio_Gparity.cc
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										70
									
								
								tests/Test_hmc_EODWFRatio_Gparity.cc
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,70 @@
 | 
			
		||||
#include "Grid.h"
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
using namespace std;
 | 
			
		||||
using namespace Grid;
 | 
			
		||||
using namespace Grid::QCD;
 | 
			
		||||
 | 
			
		||||
int main (int argc, char ** argv)
 | 
			
		||||
{
 | 
			
		||||
  Grid_init(&argc,&argv);
 | 
			
		||||
 | 
			
		||||
  const int Ls = 8;
 | 
			
		||||
 | 
			
		||||
  GridCartesian         * UGrid   = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd,vComplex::Nsimd()),GridDefaultMpi());
 | 
			
		||||
  GridRedBlackCartesian * UrbGrid = SpaceTimeGrid::makeFourDimRedBlackGrid(UGrid);
 | 
			
		||||
 | 
			
		||||
  GridCartesian         * FGrid   = SpaceTimeGrid::makeFiveDimGrid(Ls,UGrid);
 | 
			
		||||
  GridRedBlackCartesian * FrbGrid = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls,UGrid);
 | 
			
		||||
 | 
			
		||||
  GridSerialRNG    sRNG;
 | 
			
		||||
  GridParallelRNG  pRNG(UGrid);
 | 
			
		||||
  sRNG.SeedRandomDevice();
 | 
			
		||||
  pRNG.SeedRandomDevice();
 | 
			
		||||
 | 
			
		||||
  LatticeLorentzColourMatrix     U(UGrid);
 | 
			
		||||
 | 
			
		||||
  SU3::HotConfiguration(pRNG, U);
 | 
			
		||||
 | 
			
		||||
  // simplify template declaration? Strip the lorentz from the second template
 | 
			
		||||
  WilsonGaugeActionR Waction(5.6);
 | 
			
		||||
 | 
			
		||||
  Real mass=0.04;
 | 
			
		||||
  Real pv  =1.0;
 | 
			
		||||
  RealD M5=1.5;
 | 
			
		||||
 | 
			
		||||
  typedef typename GparityDomainWallFermionR::FermionField FermionField;
 | 
			
		||||
 | 
			
		||||
  const int nu = 3;
 | 
			
		||||
  std::vector<int> twists(Nd,0);
 | 
			
		||||
  twists[nu] = 1;
 | 
			
		||||
  GparityDomainWallFermionR::ImplParams params;
 | 
			
		||||
  params.twists = twists;
 | 
			
		||||
 | 
			
		||||
  GparityDomainWallFermionR DenOp(U,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5,params);
 | 
			
		||||
  GparityDomainWallFermionR NumOp(U,*FGrid,*FrbGrid,*UGrid,*UrbGrid,pv,M5,params);
 | 
			
		||||
  
 | 
			
		||||
  ConjugateGradient<FermionField>  CG(1.0e-8,10000);
 | 
			
		||||
 | 
			
		||||
  TwoFlavourEvenOddRatioPseudoFermionAction<GparityWilsonImplR> Nf2(NumOp, DenOp,CG,CG);
 | 
			
		||||
  
 | 
			
		||||
  //Collect actions
 | 
			
		||||
  ActionLevel<LatticeGaugeField> Level1;
 | 
			
		||||
  Level1.push_back(&Nf2);
 | 
			
		||||
  Level1.push_back(&Waction);
 | 
			
		||||
  ActionSet<LatticeGaugeField> FullSet;
 | 
			
		||||
  FullSet.push_back(Level1);
 | 
			
		||||
 | 
			
		||||
  // Create integrator
 | 
			
		||||
  typedef MinimumNorm2<LatticeGaugeField>  IntegratorType;// change here to change the algorithm
 | 
			
		||||
  IntegratorParameters MDpar(20);
 | 
			
		||||
  IntegratorType MDynamics(UGrid,MDpar, FullSet);
 | 
			
		||||
 | 
			
		||||
  // Create HMC
 | 
			
		||||
  HMCparameters HMCpar;
 | 
			
		||||
  HybridMonteCarlo<LatticeGaugeField,IntegratorType>  HMC(HMCpar, MDynamics,sRNG,pRNG);
 | 
			
		||||
 | 
			
		||||
  // Run it
 | 
			
		||||
  HMC.evolve(U);
 | 
			
		||||
 | 
			
		||||
}
 | 
			
		||||
		Reference in New Issue
	
	Block a user