/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./tests/Test_wilson_force.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 using namespace std; using namespace Grid; using namespace Grid::QCD; int main (int argc, char ** argv) { Grid_init(&argc,&argv); std::vector latt_size = GridDefaultLatt(); std::vector simd_layout = GridDefaultSimd(Nd,vComplex::Nsimd()); std::vector mpi_layout = GridDefaultMpi(); GridCartesian Grid(latt_size,simd_layout,mpi_layout); GridRedBlackCartesian RBGrid(latt_size,simd_layout,mpi_layout); int threads = GridThread::GetThreads(); std::cout< seeds({1,2,3,4}); GridParallelRNG pRNG(&Grid); pRNG.SeedFixedIntegers(std::vector({45,12,81,9})); LatticeFermion phi (&Grid); gaussian(pRNG,phi); LatticeFermion Mphi (&Grid); LatticeFermion MphiPrime (&Grid); LatticeGaugeField U(&Grid); //SU2::HotConfiguration(pRNG,U); SU3::ColdConfiguration(pRNG,U); //////////////////////////////////// // Unmodified matrix element //////////////////////////////////// RealD mass=-4.0; //kills the diagonal term WilsonFermionR Dw (U, Grid,RBGrid,mass); 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(&Grid); LatticeGaugeField tmp(&Grid); Dw.MDeriv(tmp , Mphi, phi,DaggerNo ); UdSdU=tmp; Dw.MDeriv(tmp , phi, Mphi,DaggerYes ); UdSdU=(UdSdU+tmp); // Take the trace UdSdU = Ta(UdSdU); LatticeFermion Ftmp (&Grid); //////////////////////////////////// // Modify the gauge field a little //////////////////////////////////// RealD dt = 0.0001; RealD Hmom = 0.0; RealD Hmomprime = 0.0; RealD Hmompp = 0.0; LatticeColourMatrix mommu(&Grid); LatticeColourMatrix forcemu(&Grid); LatticeGaugeField mom(&Grid); LatticeGaugeField Uprime(&Grid); for(int mu=0;mu(mom,mommu,mu); // fourth order exponential approx parallel_for(auto i=mom.begin();i(mom,mu); std::cout << GridLogMessage<< " Mommu " << norm2(mommu)<(UdSdU,mu); std::cout << GridLogMessage<< " dsdumu " << norm2(mommu)<(UdSdU,mu); mommu=Ta(mommu)*2.0; PokeIndex(UdSdU,mommu,mu); } for(int mu=0;mu(mom,mu); std::cout << GridLogMessage<< " Mommu " << norm2(mommu)<(UdSdU,mu); std::cout << GridLogMessage<< " dsdumu " << norm2(mommu)<(UdSdU,mu); mommu = PeekIndex(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))); } ComplexD dSpred = sum(dS); ComplexD dSm = sum(dSmom); ComplexD dSm2 = sum(dSmom2); std::cout << GridLogMessage <<"Initial mom hamiltonian is "<< Hmom <