mirror of
https://github.com/paboyle/Grid.git
synced 2024-11-10 15:55:37 +00:00
a6e1ea216d
not even SU(3) for now) gauge field. Convergence history is correctly indepdendent of decomposition on 1,2,4,8,16 mpi tasks. Found a couple of simd bugs which required fixed and enhanced the Grid_simd.cc test suite. Implemented the Mdag, M, MdagM, Meooe Mooee schur type stuff in the wilson dop.
188 lines
6.2 KiB
C++
188 lines
6.2 KiB
C++
#ifndef GRID_LATTICE_ARITH_H
|
|
#define GRID_LATTICE_ARITH_H
|
|
|
|
namespace Grid {
|
|
|
|
|
|
//////////////////////////////////////////////////////////////////////////////////////////////////////
|
|
// avoid copy back routines for mult, mac, sub, add
|
|
//////////////////////////////////////////////////////////////////////////////////////////////////////
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void mult(Lattice<obj1> &ret,const Lattice<obj2> &lhs,const Lattice<obj3> &rhs){
|
|
conformable(lhs,rhs);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<lhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
mult(&tmp,&lhs._odata[ss],&rhs._odata[ss]);
|
|
vstream(ret._odata[ss],tmp);
|
|
// mult(&ret._odata[ss],&lhs._odata[ss],&rhs._odata[ss]);
|
|
}
|
|
}
|
|
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void mac(Lattice<obj1> &ret,const Lattice<obj2> &lhs,const Lattice<obj3> &rhs){
|
|
conformable(lhs,rhs);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<lhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
mac(&tmp,&lhs._odata[ss],&rhs._odata[ss]);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void sub(Lattice<obj1> &ret,const Lattice<obj2> &lhs,const Lattice<obj3> &rhs){
|
|
conformable(lhs,rhs);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<lhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
sub(&tmp,&lhs._odata[ss],&rhs._odata[ss]);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void add(Lattice<obj1> &ret,const Lattice<obj2> &lhs,const Lattice<obj3> &rhs){
|
|
conformable(lhs,rhs);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<lhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
add(&tmp,&lhs._odata[ss],&rhs._odata[ss]);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
|
|
//////////////////////////////////////////////////////////////////////////////////////////////////////
|
|
// avoid copy back routines for mult, mac, sub, add
|
|
//////////////////////////////////////////////////////////////////////////////////////////////////////
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void mult(Lattice<obj1> &ret,const Lattice<obj2> &lhs,const obj3 &rhs){
|
|
conformable(lhs,ret);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<lhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
mult(&tmp,&lhs._odata[ss],&rhs);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void mac(Lattice<obj1> &ret,const Lattice<obj2> &lhs,const obj3 &rhs){
|
|
conformable(lhs,ret);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<lhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
mac(&tmp,&lhs._odata[ss],&rhs);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void sub(Lattice<obj1> &ret,const Lattice<obj2> &lhs,const obj3 &rhs){
|
|
conformable(lhs,ret);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<lhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
sub(&tmp,&lhs._odata[ss],&rhs);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void add(Lattice<obj1> &ret,const Lattice<obj2> &lhs,const obj3 &rhs){
|
|
conformable(lhs,ret);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<lhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
add(&tmp,&lhs._odata[ss],&rhs);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
|
|
//////////////////////////////////////////////////////////////////////////////////////////////////////
|
|
// avoid copy back routines for mult, mac, sub, add
|
|
//////////////////////////////////////////////////////////////////////////////////////////////////////
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void mult(Lattice<obj1> &ret,const obj2 &lhs,const Lattice<obj3> &rhs){
|
|
conformable(ret,rhs);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<rhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
mult(&tmp,&lhs,&rhs._odata[ss]);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void mac(Lattice<obj1> &ret,const obj2 &lhs,const Lattice<obj3> &rhs){
|
|
conformable(ret,rhs);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<rhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
mac(&tmp,&lhs,&rhs._odata[ss]);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void sub(Lattice<obj1> &ret,const obj2 &lhs,const Lattice<obj3> &rhs){
|
|
conformable(ret,rhs);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<rhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
sub(&tmp,&lhs,&rhs._odata[ss]);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
template<class obj1,class obj2,class obj3> strong_inline
|
|
void add(Lattice<obj1> &ret,const obj2 &lhs,const Lattice<obj3> &rhs){
|
|
conformable(ret,rhs);
|
|
PARALLEL_FOR_LOOP
|
|
for(int ss=0;ss<rhs._grid->oSites();ss++){
|
|
obj1 tmp;
|
|
add(&tmp,&lhs,&rhs._odata[ss]);
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
|
|
template<class sobj,class vobj> strong_inline
|
|
void axpy(Lattice<vobj> &ret,sobj a,const Lattice<vobj> &x,const Lattice<vobj> &y){
|
|
conformable(x,y);
|
|
#pragma omp parallel for
|
|
for(int ss=0;ss<x._grid->oSites();ss++){
|
|
vobj tmp = a*x._odata[ss]+y._odata[ss];
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
template<class sobj,class vobj> strong_inline
|
|
void axpby(Lattice<vobj> &ret,sobj a,sobj b,const Lattice<vobj> &x,const Lattice<vobj> &y){
|
|
conformable(x,y);
|
|
#pragma omp parallel for
|
|
for(int ss=0;ss<x._grid->oSites();ss++){
|
|
vobj tmp = a*x._odata[ss]+b*y._odata[ss];
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
}
|
|
|
|
template<class sobj,class vobj> strong_inline
|
|
RealD axpy_norm(Lattice<vobj> &ret,sobj a,const Lattice<vobj> &x,const Lattice<vobj> &y){
|
|
conformable(x,y);
|
|
#pragma omp parallel for
|
|
for(int ss=0;ss<x._grid->oSites();ss++){
|
|
vobj tmp = a*x._odata[ss]+y._odata[ss];
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
return norm2(ret);
|
|
}
|
|
template<class sobj,class vobj> strong_inline
|
|
RealD axpby_norm(Lattice<vobj> &ret,sobj a,sobj b,const Lattice<vobj> &x,const Lattice<vobj> &y){
|
|
conformable(x,y);
|
|
#pragma omp parallel for
|
|
for(int ss=0;ss<x._grid->oSites();ss++){
|
|
vobj tmp = a*x._odata[ss]+b*y._odata[ss];
|
|
vstream(ret._odata[ss],tmp);
|
|
}
|
|
return norm2(ret); // FIXME implement parallel norm in ss loop
|
|
}
|
|
|
|
}
|
|
#endif
|