mirror of
https://github.com/paboyle/Grid.git
synced 2025-04-10 06:00:45 +01:00
twist and boundary conditions for free propagator
This commit is contained in:
parent
84940fbdf0
commit
00963a7499
@ -43,7 +43,7 @@ namespace Grid {
|
|||||||
INHERIT_IMPL_TYPES(Impl);
|
INHERIT_IMPL_TYPES(Impl);
|
||||||
public:
|
public:
|
||||||
|
|
||||||
void FreePropagator(const FermionField &in,FermionField &out,RealD mass, std::vector<double> twist, bool fiveD) {
|
void FreePropagator(const FermionField &in,FermionField &out,RealD mass,std::vector<Complex> boundary, std::vector<double> twist, bool fiveD) {
|
||||||
FermionField in_k(in._grid);
|
FermionField in_k(in._grid);
|
||||||
FermionField prop_k(in._grid);
|
FermionField prop_k(in._grid);
|
||||||
|
|
||||||
@ -55,15 +55,19 @@ namespace Grid {
|
|||||||
FermionField in_buf(in._grid); in_buf = zero;
|
FermionField in_buf(in._grid); in_buf = zero;
|
||||||
Complex ci(0.0,1.0);
|
Complex ci(0.0,1.0);
|
||||||
assert(twist.size() == Nd);//check that twist is Nd
|
assert(twist.size() == Nd);//check that twist is Nd
|
||||||
|
assert(boundary.size() == Nd);//check that boundary conditions is Nd
|
||||||
int shift = 0;
|
int shift = 0;
|
||||||
if(fiveD) shift = 1;
|
if(fiveD) shift = 1;
|
||||||
for(unsigned int nu = 0; nu < Nd; nu++)
|
for(unsigned int nu = 0; nu < Nd; nu++)
|
||||||
{
|
{
|
||||||
// Shift coordinate lattice index by 1 to account for 5th dimension.
|
// Shift coordinate lattice index by 1 to account for 5th dimension.
|
||||||
LatticeCoordinate(coor, nu + shift);
|
LatticeCoordinate(coor, nu + shift);
|
||||||
ph = ph + twist[nu]*coor*((1./(in._grid->_fdimensions[nu+shift])));
|
double boundary_phase = ::acos(real(boundary[nu]));
|
||||||
|
ph = ph + boundary_phase*coor*((1./(in._grid->_fdimensions[nu+shift])));
|
||||||
|
//momenta for propagator shifted by twist+boundary
|
||||||
|
twist[nu] = twist[nu] + boundary_phase/((2.0*M_PI));
|
||||||
}
|
}
|
||||||
in_buf = exp((Real)(2.0*M_PI)*ci*ph*(-1.0))*in;
|
in_buf = exp(ci*ph*(-1.0))*in;
|
||||||
|
|
||||||
if(fiveD){//FFT only on temporal and spatial dimensions
|
if(fiveD){//FFT only on temporal and spatial dimensions
|
||||||
std::vector<int> mask(Nd+1,1); mask[0] = 0;
|
std::vector<int> mask(Nd+1,1); mask[0] = 0;
|
||||||
@ -76,25 +80,28 @@ namespace Grid {
|
|||||||
this->MomentumSpacePropagatorHt(prop_k,in_k,mass,twist);
|
this->MomentumSpacePropagatorHt(prop_k,in_k,mass,twist);
|
||||||
theFFT.FFT_all_dim(out,prop_k,FFT::backward);
|
theFFT.FFT_all_dim(out,prop_k,FFT::backward);
|
||||||
}
|
}
|
||||||
|
|
||||||
//phase for boundary condition
|
//phase for boundary condition
|
||||||
out = out * exp((Real)(2.0*M_PI)*ci*ph);
|
out = out * exp(ci*ph);
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass,std::vector<double> twist) {
|
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass,std::vector<Complex> boundary,std::vector<double> twist) {
|
||||||
bool fiveD = true; //5d propagator by default
|
bool fiveD = true; //5d propagator by default
|
||||||
FreePropagator(in,out,mass,twist,fiveD);
|
FreePropagator(in,out,mass,boundary,twist,fiveD);
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass, bool fiveD) {
|
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass, bool fiveD) {
|
||||||
std::vector<double> twist(Nd,0.0); //default: periodic boundarys in all directions
|
std::vector<double> twist(Nd,0.0); //default: periodic boundarys in all directions
|
||||||
FreePropagator(in,out,mass,twist,fiveD);
|
std::vector<Complex> boundary;
|
||||||
|
for(int i=0;i<Nd;i++) boundary.push_back(1);//default: periodic boundary conditions
|
||||||
|
FreePropagator(in,out,mass,boundary,twist,fiveD);
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass) {
|
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass) {
|
||||||
bool fiveD = true; //5d propagator by default
|
bool fiveD = true; //5d propagator by default
|
||||||
std::vector<double> twist(Nd,0.0); //default: periodic boundarys in all directions
|
std::vector<double> twist(Nd,0.0); //default: twist angle 0
|
||||||
FreePropagator(in,out,mass,twist,fiveD);
|
std::vector<Complex> boundary;
|
||||||
|
for(int i=0;i<Nd;i++) boundary.push_back(1); //default: periodic boundary conditions
|
||||||
|
FreePropagator(in,out,mass,boundary,twist,fiveD);
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void Instantiatable(void) {};
|
virtual void Instantiatable(void) {};
|
||||||
|
@ -96,7 +96,7 @@ namespace Grid {
|
|||||||
|
|
||||||
virtual void MomentumSpacePropagator(FermionField &out,const FermionField &in,RealD _m,std::vector<double> twist) { assert(0);};
|
virtual void MomentumSpacePropagator(FermionField &out,const FermionField &in,RealD _m,std::vector<double> twist) { assert(0);};
|
||||||
|
|
||||||
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass,std::vector<double> twist) {
|
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass,std::vector<Complex> boundary,std::vector<double> twist) {
|
||||||
FFT theFFT((GridCartesian *) in._grid);
|
FFT theFFT((GridCartesian *) in._grid);
|
||||||
|
|
||||||
FermionField in_k(in._grid);
|
FermionField in_k(in._grid);
|
||||||
@ -108,24 +108,31 @@ namespace Grid {
|
|||||||
FermionField in_buf(in._grid); in_buf = zero;
|
FermionField in_buf(in._grid); in_buf = zero;
|
||||||
Complex ci(0.0,1.0);
|
Complex ci(0.0,1.0);
|
||||||
assert(twist.size() == Nd);//check that twist is Nd
|
assert(twist.size() == Nd);//check that twist is Nd
|
||||||
|
assert(boundary.size() == Nd);//check that boundary conditions is Nd
|
||||||
for(unsigned int nu = 0; nu < Nd; nu++)
|
for(unsigned int nu = 0; nu < Nd; nu++)
|
||||||
{
|
{
|
||||||
LatticeCoordinate(coor, nu);
|
LatticeCoordinate(coor, nu);
|
||||||
ph = ph + twist[nu]*coor*((1./(in._grid->_fdimensions[nu])));
|
double boundary_phase = ::acos(real(boundary[nu]));
|
||||||
|
ph = ph + boundary_phase*coor*((1./(in._grid->_fdimensions[nu])));
|
||||||
|
//momenta for propagator shifted by twist+boundary
|
||||||
|
twist[nu] = twist[nu] + boundary_phase/((2.0*M_PI));
|
||||||
}
|
}
|
||||||
in_buf = exp((Real)(2.0*M_PI)*ci*ph*(-1.0))*in;
|
in_buf = exp(ci*ph*(-1.0))*in;
|
||||||
|
|
||||||
theFFT.FFT_all_dim(in_k,in_buf,FFT::forward);
|
theFFT.FFT_all_dim(in_k,in_buf,FFT::forward);
|
||||||
this->MomentumSpacePropagator(prop_k,in_k,mass,twist);
|
this->MomentumSpacePropagator(prop_k,in_k,mass,twist);
|
||||||
theFFT.FFT_all_dim(out,prop_k,FFT::backward);
|
theFFT.FFT_all_dim(out,prop_k,FFT::backward);
|
||||||
|
|
||||||
//phase for boundary condition
|
//phase for boundary condition
|
||||||
out = out * exp((Real)(2.0*M_PI)*ci*ph);
|
out = out * exp(ci*ph);
|
||||||
|
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass) {
|
virtual void FreePropagator(const FermionField &in,FermionField &out,RealD mass) {
|
||||||
|
std::vector<Complex> boundary;
|
||||||
|
for(int i=0;i<Nd;i++) boundary.push_back(1);//default: periodic boundary conditions
|
||||||
std::vector<double> twist(Nd,0.0); //default: periodic boundarys in all directions
|
std::vector<double> twist(Nd,0.0); //default: periodic boundarys in all directions
|
||||||
FreePropagator(in,out,mass,twist);
|
FreePropagator(in,out,mass,boundary,twist);
|
||||||
};
|
};
|
||||||
|
|
||||||
///////////////////////////////////////////////
|
///////////////////////////////////////////////
|
||||||
|
@ -72,6 +72,7 @@ public:
|
|||||||
std::string, action,
|
std::string, action,
|
||||||
std::string, emField,
|
std::string, emField,
|
||||||
double, mass,
|
double, mass,
|
||||||
|
std::string , boundary,
|
||||||
std::string, twist,
|
std::string, twist,
|
||||||
unsigned int, deltat);
|
unsigned int, deltat);
|
||||||
};
|
};
|
||||||
@ -163,8 +164,7 @@ void TEMLepton<FImpl>::execute(void)
|
|||||||
envGetTmp(FermionField, sol);
|
envGetTmp(FermionField, sol);
|
||||||
envGetTmp(FermionField, tmp);
|
envGetTmp(FermionField, tmp);
|
||||||
LOG(Message) << "Calculating a lepton Propagator with sequential Aslash insertion with lepton mass "
|
LOG(Message) << "Calculating a lepton Propagator with sequential Aslash insertion with lepton mass "
|
||||||
<< mass << " and twist ("
|
<< mass << " using the action '" << par().action
|
||||||
<< par().twist << ") using the action '" << par().action
|
|
||||||
<< "' for fixed source-sink separation of " << par().deltat << std::endl;
|
<< "' for fixed source-sink separation of " << par().deltat << std::endl;
|
||||||
|
|
||||||
envGetTmp(Lattice<iScalar<vInteger>>, tlat);
|
envGetTmp(Lattice<iScalar<vInteger>>, tlat);
|
||||||
@ -176,6 +176,11 @@ void TEMLepton<FImpl>::execute(void)
|
|||||||
{
|
{
|
||||||
HADRONS_ERROR(Size, "number of twist angles does not match number of dimensions");
|
HADRONS_ERROR(Size, "number of twist angles does not match number of dimensions");
|
||||||
}
|
}
|
||||||
|
std::vector<Complex> boundary = strToVec<Complex>(par().boundary);
|
||||||
|
if(boundary.size() != Nd)
|
||||||
|
{
|
||||||
|
HADRONS_ERROR(Size, "number of boundary conditions does not match number of dimensions");
|
||||||
|
}
|
||||||
|
|
||||||
auto &stoch_photon = envGet(EmField, par().emField);
|
auto &stoch_photon = envGet(EmField, par().emField);
|
||||||
unsigned int nt = env().getDim(Tp);
|
unsigned int nt = env().getDim(Tp);
|
||||||
@ -212,7 +217,7 @@ void TEMLepton<FImpl>::execute(void)
|
|||||||
mat.ImportPhysicalFermionSource(tmp, source);
|
mat.ImportPhysicalFermionSource(tmp, source);
|
||||||
}
|
}
|
||||||
sol = zero;
|
sol = zero;
|
||||||
mat.FreePropagator(source,sol,mass,twist);
|
mat.FreePropagator(source,sol,mass,boundary,twist);
|
||||||
if (Ls_ == 1)
|
if (Ls_ == 1)
|
||||||
{
|
{
|
||||||
FermToProp<FImpl>(freetmp, sol, s, 0);
|
FermToProp<FImpl>(freetmp, sol, s, 0);
|
||||||
@ -229,11 +234,7 @@ void TEMLepton<FImpl>::execute(void)
|
|||||||
|
|
||||||
//shift free propagator to different source positions
|
//shift free propagator to different source positions
|
||||||
proptmp = Cshift(freetmp,Tp, -tl);
|
proptmp = Cshift(freetmp,Tp, -tl);
|
||||||
//take anti-periodic boundary conditions into account, if used
|
proptmp = where( tlat < tl, boundary[Tp]*proptmp, proptmp);
|
||||||
if(twist[Tp]==0.5)
|
|
||||||
{
|
|
||||||
proptmp = where( tlat < tl, (-1.0)*proptmp, proptmp);
|
|
||||||
}
|
|
||||||
|
|
||||||
// free propagator for fixed source-sink separation
|
// free propagator for fixed source-sink separation
|
||||||
lep = where(tlat == (tl-par().deltat+nt)%nt, proptmp, lep);
|
lep = where(tlat == (tl-par().deltat+nt)%nt, proptmp, lep);
|
||||||
@ -264,7 +265,7 @@ void TEMLepton<FImpl>::execute(void)
|
|||||||
mat.ImportPhysicalFermionSource(tmp, source);
|
mat.ImportPhysicalFermionSource(tmp, source);
|
||||||
}
|
}
|
||||||
sol = zero;
|
sol = zero;
|
||||||
mat.FreePropagator(source,sol,mass,twist);
|
mat.FreePropagator(source,sol,mass,boundary,twist);
|
||||||
if (Ls_ == 1)
|
if (Ls_ == 1)
|
||||||
{
|
{
|
||||||
FermToProp<FImpl>(proptmp, sol, s, 0);
|
FermToProp<FImpl>(proptmp, sol, s, 0);
|
||||||
|
@ -49,6 +49,7 @@ public:
|
|||||||
std::string, source,
|
std::string, source,
|
||||||
std::string, action,
|
std::string, action,
|
||||||
double, mass,
|
double, mass,
|
||||||
|
std::string , boundary,
|
||||||
std::string, twist);
|
std::string, twist);
|
||||||
};
|
};
|
||||||
|
|
||||||
@ -172,7 +173,12 @@ void TFreeProp<FImpl>::execute(void)
|
|||||||
{
|
{
|
||||||
HADRONS_ERROR(Size, "number of twist angles does not match number of dimensions");
|
HADRONS_ERROR(Size, "number of twist angles does not match number of dimensions");
|
||||||
}
|
}
|
||||||
mat.FreePropagator(source,sol,mass,twist);
|
std::vector<Complex> boundary = strToVec<Complex>(par().boundary);
|
||||||
|
if(boundary.size() != Nd)
|
||||||
|
{
|
||||||
|
HADRONS_ERROR(Size, "number of boundary conditions does not match number of dimensions");
|
||||||
|
}
|
||||||
|
mat.FreePropagator(source,sol,mass,boundary,twist);
|
||||||
FermToProp<FImpl>(prop, sol, s, c);
|
FermToProp<FImpl>(prop, sol, s, c);
|
||||||
// create 4D propagators from 5D one if necessary
|
// create 4D propagators from 5D one if necessary
|
||||||
if (Ls_ > 1)
|
if (Ls_ > 1)
|
||||||
|
Loading…
x
Reference in New Issue
Block a user