mirror of
https://github.com/paboyle/Grid.git
synced 2025-06-14 13:57:07 +01:00
Hadrons: Aslash field, tested
This commit is contained in:
@ -60,6 +60,14 @@ public:
|
||||
const FermionField *vj,
|
||||
int orthogdim);
|
||||
|
||||
template <typename TensorType> // output: rank 5 tensor, e.g. Eigen::Tensor<ComplexD, 5>
|
||||
static void AslashField(TensorType &mat,
|
||||
const FermionField *lhs_wi,
|
||||
const FermionField *rhs_vj,
|
||||
const std::vector<ComplexField> &emB0,
|
||||
const std::vector<ComplexField> &emB1,
|
||||
int orthogdim, double *t_kernel = nullptr, double *t_gsum = nullptr);
|
||||
|
||||
static void ContractWWVV(std::vector<PropagatorField> &WWVV,
|
||||
const Eigen::Tensor<ComplexD,3> &WW_sd,
|
||||
const FermionField *vs,
|
||||
@ -617,6 +625,189 @@ void A2Autils<FImpl>::PionFieldVV(Eigen::Tensor<ComplexD,3> &mat,
|
||||
PionFieldXX(mat,vi,vj,orthogdim,nog5);
|
||||
}
|
||||
|
||||
// "A-slash" field w_i(x)^dag * i * A_mu * gamma_mu * v_j(x)
|
||||
//
|
||||
// With:
|
||||
//
|
||||
// B_0 = A_0 + i A_1
|
||||
// B_1 = A_2 + i A_3
|
||||
//
|
||||
// then in spin space
|
||||
//
|
||||
// ( 0 0 -conj(B_1) -B_0 )
|
||||
// i * A_mu g_mu = ( 0 0 -conj(B_0) B_1 )
|
||||
// ( B_1 B_0 0 0 )
|
||||
// ( conj(B_0) -conj(B_1) 0 0 )
|
||||
template <class FImpl>
|
||||
template <typename TensorType>
|
||||
void A2Autils<FImpl>::AslashField(TensorType &mat,
|
||||
const FermionField *lhs_wi,
|
||||
const FermionField *rhs_vj,
|
||||
const std::vector<ComplexField> &emB0,
|
||||
const std::vector<ComplexField> &emB1,
|
||||
int orthogdim, double *t_kernel, double *t_gsum)
|
||||
{
|
||||
typedef typename FermionField::vector_object vobj;
|
||||
typedef typename vobj::scalar_object sobj;
|
||||
typedef typename vobj::scalar_type scalar_type;
|
||||
typedef typename vobj::vector_type vector_type;
|
||||
|
||||
typedef iSpinMatrix<vector_type> SpinMatrix_v;
|
||||
typedef iSpinMatrix<scalar_type> SpinMatrix_s;
|
||||
typedef iSinglet<vector_type> Singlet_v;
|
||||
typedef iSinglet<scalar_type> Singlet_s;
|
||||
|
||||
int Lblock = mat.dimension(3);
|
||||
int Rblock = mat.dimension(4);
|
||||
|
||||
GridBase *grid = lhs_wi[0]._grid;
|
||||
|
||||
const int Nd = grid->_ndimension;
|
||||
const int Nsimd = grid->Nsimd();
|
||||
|
||||
int Nt = grid->GlobalDimensions()[orthogdim];
|
||||
int Nem = emB0.size();
|
||||
assert(emB1.size() == Nem);
|
||||
|
||||
int fd=grid->_fdimensions[orthogdim];
|
||||
int ld=grid->_ldimensions[orthogdim];
|
||||
int rd=grid->_rdimensions[orthogdim];
|
||||
|
||||
// will locally sum vectors first
|
||||
// sum across these down to scalars
|
||||
// splitting the SIMD
|
||||
int MFrvol = rd*Lblock*Rblock*Nem;
|
||||
int MFlvol = ld*Lblock*Rblock*Nem;
|
||||
|
||||
Vector<vector_type> lvSum(MFrvol);
|
||||
parallel_for (int r = 0; r < MFrvol; r++)
|
||||
{
|
||||
lvSum[r] = zero;
|
||||
}
|
||||
|
||||
Vector<scalar_type> lsSum(MFlvol);
|
||||
parallel_for (int r = 0; r < MFlvol; r++)
|
||||
{
|
||||
lsSum[r] = scalar_type(0.0);
|
||||
}
|
||||
|
||||
int e1= grid->_slice_nblock[orthogdim];
|
||||
int e2= grid->_slice_block [orthogdim];
|
||||
int stride=grid->_slice_stride[orthogdim];
|
||||
|
||||
// Nested parallelism would be ok
|
||||
// Wasting cores here. Test case r
|
||||
if (t_kernel) *t_kernel = -usecond();
|
||||
parallel_for(int r=0;r<rd;r++)
|
||||
{
|
||||
int so=r*grid->_ostride[orthogdim]; // base offset for start of plane
|
||||
|
||||
for(int n=0;n<e1;n++)
|
||||
for(int b=0;b<e2;b++)
|
||||
{
|
||||
int ss= so+n*stride+b;
|
||||
|
||||
for(int i=0;i<Lblock;i++)
|
||||
{
|
||||
auto left = conjugate(lhs_wi[i]._odata[ss]);
|
||||
|
||||
for(int j=0;j<Rblock;j++)
|
||||
{
|
||||
SpinMatrix_v vv;
|
||||
auto right = rhs_vj[j]._odata[ss];
|
||||
|
||||
for(int s1=0;s1<Ns;s1++)
|
||||
for(int s2=0;s2<Ns;s2++)
|
||||
{
|
||||
vv()(s1,s2)() = left()(s2)(0) * right()(s1)(0)
|
||||
+ left()(s2)(1) * right()(s1)(1)
|
||||
+ left()(s2)(2) * right()(s1)(2);
|
||||
}
|
||||
|
||||
// After getting the sitewise product do the mom phase loop
|
||||
int base = Nem*i+Nem*Lblock*j+Nem*Lblock*Rblock*r;
|
||||
|
||||
for ( int m=0;m<Nem;m++)
|
||||
{
|
||||
int idx = m+base;
|
||||
auto b0 = emB0[m]._odata[ss];
|
||||
auto b1 = emB1[m]._odata[ss];
|
||||
auto cb0 = conjugate(b0);
|
||||
auto cb1 = conjugate(b1);
|
||||
|
||||
lvSum[idx] += - vv()(3,0)()*b0()()() - vv()(2,0)()*cb1()()()
|
||||
+ vv()(3,1)()*b1()()() - vv()(2,1)()*cb0()()()
|
||||
+ vv()(0,2)()*b1()()() + vv()(1,2)()*b0()()()
|
||||
+ vv()(0,3)()*cb0()()() - vv()(1,3)()*cb1()()();
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
// Sum across simd lanes in the plane, breaking out orthog dir.
|
||||
parallel_for(int rt=0;rt<rd;rt++)
|
||||
{
|
||||
std::vector<int> icoor(Nd);
|
||||
std::vector<scalar_type> extracted(Nsimd);
|
||||
|
||||
for(int i=0;i<Lblock;i++)
|
||||
for(int j=0;j<Rblock;j++)
|
||||
for(int m=0;m<Nem;m++)
|
||||
{
|
||||
|
||||
int ij_rdx = m+Nem*i+Nem*Lblock*j+Nem*Lblock*Rblock*rt;
|
||||
|
||||
extract<vector_type,scalar_type>(lvSum[ij_rdx],extracted);
|
||||
for(int idx=0;idx<Nsimd;idx++)
|
||||
{
|
||||
grid->iCoorFromIindex(icoor,idx);
|
||||
|
||||
int ldx = rt+icoor[orthogdim]*rd;
|
||||
int ij_ldx = m+Nem*i+Nem*Lblock*j+Nem*Lblock*Rblock*ldx;
|
||||
|
||||
lsSum[ij_ldx]=lsSum[ij_ldx]+extracted[idx];
|
||||
}
|
||||
}
|
||||
}
|
||||
if (t_kernel) *t_kernel += usecond();
|
||||
|
||||
// ld loop and local only??
|
||||
int pd = grid->_processors[orthogdim];
|
||||
int pc = grid->_processor_coor[orthogdim];
|
||||
parallel_for_nest2(int lt=0;lt<ld;lt++)
|
||||
{
|
||||
for(int pt=0;pt<pd;pt++)
|
||||
{
|
||||
int t = lt + pt*ld;
|
||||
if (pt == pc)
|
||||
{
|
||||
for(int i=0;i<Lblock;i++)
|
||||
for(int j=0;j<Rblock;j++)
|
||||
for(int m=0;m<Nem;m++)
|
||||
{
|
||||
int ij_dx = m+Nem*i + Nem*Lblock * j + Nem*Lblock * Rblock * lt;
|
||||
|
||||
mat(m,0,t,i,j) = lsSum[ij_dx];
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
const scalar_type zz(0.0);
|
||||
|
||||
for(int i=0;i<Lblock;i++)
|
||||
for(int j=0;j<Rblock;j++)
|
||||
for(int m=0;m<Nem;m++)
|
||||
{
|
||||
mat(m,0,t,i,j) = zz;
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
if (t_gsum) *t_gsum = -usecond();
|
||||
grid->GlobalSumVector(&mat(0,0,0,0,0),Nem*Nt*Lblock*Rblock);
|
||||
if (t_gsum) *t_gsum += usecond();
|
||||
}
|
||||
|
||||
////////////////////////////////////////////
|
||||
// Schematic thoughts about more generalised four quark insertion
|
||||
|
Reference in New Issue
Block a user