/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./lib/lattice/Lattice_reduction.h Copyright (C) 2015 Author: Azusa Yamaguchi Author: Peter Boyle Author: paboyle 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 */ #ifndef GRID_LATTICE_REDUCTION_H #define GRID_LATTICE_REDUCTION_H namespace Grid { #ifdef GRID_WARN_SUBOPTIMAL #warning "Optimisation alert all these reduction loops are NOT threaded " #endif //////////////////////////////////////////////////////////////////////////////////////////////////// // Deterministic Reduction operations //////////////////////////////////////////////////////////////////////////////////////////////////// template inline RealD norm2(const Lattice &arg) { ComplexD nrm = innerProduct(arg, arg); return std::real(nrm); } template inline ComplexD innerProduct(const Lattice &left, const Lattice &right) { typedef typename vobj::scalar_type scalar_type; typedef typename vobj::vector_type vector_type; scalar_type nrm; GridBase *grid = left._grid; std::vector > sumarray( grid->SumArraySize()); for (int i = 0; i < grid->SumArraySize(); i++) { sumarray[i] = zero; } PARALLEL_FOR_LOOP for (int thr = 0; thr < grid->SumArraySize(); thr++) { int nwork, mywork, myoff; GridThread::GetWork(left._grid->oSites(), thr, mywork, myoff); decltype(innerProduct(left._odata[0], right._odata[0])) vnrm = zero; // private to thread; sub summation for (int ss = myoff; ss < mywork + myoff; ss++) { vnrm = vnrm + innerProduct(left._odata[ss], right._odata[ss]); } sumarray[thr] = TensorRemove(vnrm); } vector_type vvnrm; vvnrm = zero; // sum across threads for (int i = 0; i < grid->SumArraySize(); i++) { vvnrm = vvnrm + sumarray[i]; } nrm = Reduce(vvnrm); // sum across simd right._grid->GlobalSum(nrm); return nrm; } template inline auto sum(const LatticeUnaryExpression &expr) -> typename decltype( expr.first.func(eval(0, std::get<0>(expr.second))))::scalar_object { return sum(closure(expr)); } template inline auto sum(const LatticeBinaryExpression & expr) ->typename decltype(expr.first.func(eval(0,std::get<0>(expr.second)),eval(0,std::get<1>(expr.second))))::scalar_object { return sum(closure(expr)); } template inline auto sum(const LatticeTrinaryExpression & expr) ->typename decltype(expr.first.func(eval(0,std::get<0>(expr.second)), eval(0,std::get<1>(expr.second)), eval(0,std::get<2>(expr.second)) ))::scalar_object { return sum(closure(expr)); } template inline typename vobj::scalar_object sum(const Lattice &arg){ GridBase *grid=arg._grid; int Nsimd = grid->Nsimd(); std::vector > sumarray(grid->SumArraySize()); for(int i=0;iSumArraySize();i++){ sumarray[i]=zero; } PARALLEL_FOR_LOOP for(int thr=0;thrSumArraySize();thr++){ int nwork, mywork, myoff; GridThread::GetWork(grid->oSites(),thr,mywork,myoff); vobj vvsum=zero; for(int ss=myoff;ssSumArraySize();i++){ vsum = vsum+sumarray[i]; } typedef typename vobj::scalar_object sobj; sobj ssum=zero; std::vector buf(Nsimd); extract(vsum,buf); for(int i=0;iGlobalSum(ssum); return ssum; } template inline void sliceSum(const Lattice &Data,std::vector &result,int orthogdim) { typedef typename vobj::scalar_object sobj; GridBase *grid = Data._grid; assert(grid!=NULL); // FIXME // std::cout<SumArraySize()<<" threads "<_ndimension; const int Nsimd = grid->Nsimd(); assert(orthogdim >= 0); assert(orthogdim < Nd); int fd=grid->_fdimensions[orthogdim]; int ld=grid->_ldimensions[orthogdim]; int rd=grid->_rdimensions[orthogdim]; std::vector > lvSum(rd); // will locally sum vectors first std::vector lsSum(ld,zero); // sum across these down to scalars std::vector extracted(Nsimd); // splitting the SIMD result.resize(fd); // And then global sum to return the same vector to every node for IO to file for(int r=0;r coor(Nd); // sum over reduced dimension planes, breaking out orthog dir for(int ss=0;ssoSites();ss++){ Lexicographic::CoorFromIndex(coor,ss,grid->_rdimensions); int r = coor[orthogdim]; lvSum[r]=lvSum[r]+Data._odata[ss]; } // Sum across simd lanes in the plane, breaking out orthog dir. std::vector icoor(Nd); for(int rt=0;rtiCoorFromIindex(icoor,idx); int ldx =rt+icoor[orthogdim]*rd; lsSum[ldx]=lsSum[ldx]+extracted[idx]; } } // sum over nodes. sobj gsum; for(int t=0;t_processor_coor[orthogdim] ) { gsum=lsSum[lt]; } else { gsum=zero; } grid->GlobalSum(gsum); result[t]=gsum; } } } #endif