mirror of
https://github.com/paboyle/Grid.git
synced 2025-06-16 06:47:06 +01:00
Merge pull request #453 from dbollweg/feature/sliceSum_gpu
Feature/slice sum gpu
This commit is contained in:
@ -31,6 +31,7 @@ Author: Christoph Lehner <christoph@lhnr.de>
|
||||
#if defined(GRID_SYCL)
|
||||
#include <Grid/lattice/Lattice_reduction_sycl.h>
|
||||
#endif
|
||||
#include <Grid/lattice/Lattice_slicesum_core.h>
|
||||
|
||||
NAMESPACE_BEGIN(Grid);
|
||||
|
||||
@ -448,19 +449,10 @@ template<class vobj> inline void sliceSum(const Lattice<vobj> &Data,std::vector<
|
||||
int e1= grid->_slice_nblock[orthogdim];
|
||||
int e2= grid->_slice_block [orthogdim];
|
||||
int stride=grid->_slice_stride[orthogdim];
|
||||
|
||||
// sum over reduced dimension planes, breaking out orthog dir
|
||||
// Parallel over orthog direction
|
||||
autoView( Data_v, Data, CpuRead);
|
||||
thread_for( r,rd, {
|
||||
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;
|
||||
lvSum[r]=lvSum[r]+Data_v[ss];
|
||||
}
|
||||
}
|
||||
});
|
||||
int ostride=grid->_ostride[orthogdim];
|
||||
|
||||
//Reduce Data down to lvSum
|
||||
sliceSumReduction(Data,lvSum,rd, e1,e2,stride,ostride,Nsimd);
|
||||
|
||||
// Sum across simd lanes in the plane, breaking out orthog dir.
|
||||
Coordinate icoor(Nd);
|
||||
@ -504,6 +496,7 @@ sliceSum(const Lattice<vobj> &Data,int orthogdim)
|
||||
return result;
|
||||
}
|
||||
|
||||
|
||||
template<class vobj>
|
||||
static void sliceInnerProductVector( std::vector<ComplexD> & result, const Lattice<vobj> &lhs,const Lattice<vobj> &rhs,int orthogdim)
|
||||
{
|
||||
|
213
Grid/lattice/Lattice_slicesum_core.h
Normal file
213
Grid/lattice/Lattice_slicesum_core.h
Normal file
@ -0,0 +1,213 @@
|
||||
#pragma once
|
||||
#include <type_traits>
|
||||
#if defined(GRID_CUDA)
|
||||
|
||||
#include <cub/cub.cuh>
|
||||
#define gpucub cub
|
||||
#define gpuError_t cudaError_t
|
||||
#define gpuSuccess cudaSuccess
|
||||
|
||||
#elif defined(GRID_HIP)
|
||||
|
||||
#include <hipcub/hipcub.hpp>
|
||||
#define gpucub hipcub
|
||||
#define gpuError_t hipError_t
|
||||
#define gpuSuccess hipSuccess
|
||||
|
||||
#endif
|
||||
|
||||
|
||||
NAMESPACE_BEGIN(Grid);
|
||||
|
||||
|
||||
#if defined(GRID_CUDA) || defined(GRID_HIP)
|
||||
template<class vobj> inline void sliceSumReduction_cub_small(const vobj *Data, Vector<vobj> &lvSum, const int rd, const int e1, const int e2, const int stride, const int ostride, const int Nsimd) {
|
||||
size_t subvol_size = e1*e2;
|
||||
commVector<vobj> reduction_buffer(rd*subvol_size);
|
||||
auto rb_p = &reduction_buffer[0];
|
||||
vobj zero_init;
|
||||
zeroit(zero_init);
|
||||
|
||||
|
||||
void *temp_storage_array = NULL;
|
||||
size_t temp_storage_bytes = 0;
|
||||
vobj *d_out;
|
||||
int* d_offsets;
|
||||
|
||||
std::vector<int> offsets(rd+1,0);
|
||||
|
||||
for (int i = 0; i < offsets.size(); i++) {
|
||||
offsets[i] = i*subvol_size;
|
||||
}
|
||||
|
||||
//Allocate memory for output and offset arrays on device
|
||||
d_out = static_cast<vobj*>(acceleratorAllocDevice(rd*sizeof(vobj)));
|
||||
|
||||
d_offsets = static_cast<int*>(acceleratorAllocDevice((rd+1)*sizeof(int)));
|
||||
|
||||
//copy offsets to device
|
||||
acceleratorCopyToDeviceAsync(&offsets[0],d_offsets,sizeof(int)*(rd+1),computeStream);
|
||||
|
||||
|
||||
gpuError_t gpuErr = gpucub::DeviceSegmentedReduce::Reduce(temp_storage_array, temp_storage_bytes, rb_p,d_out, rd, d_offsets, d_offsets+1, ::gpucub::Sum(), zero_init, computeStream);
|
||||
if (gpuErr!=gpuSuccess) {
|
||||
std::cout << GridLogError << "Lattice_slicesum_gpu.h: Encountered error during gpucub::DeviceSegmentedReduce::Reduce (setup)! Error: " << gpuErr <<std::endl;
|
||||
exit(EXIT_FAILURE);
|
||||
}
|
||||
|
||||
//allocate memory for temp_storage_array
|
||||
temp_storage_array = acceleratorAllocDevice(temp_storage_bytes);
|
||||
|
||||
//prepare buffer for reduction
|
||||
//use non-blocking accelerator_for to avoid syncs (ok because we submit to same computeStream)
|
||||
//use 2d accelerator_for to avoid launch latencies found when serially looping over rd
|
||||
accelerator_for2dNB( s,subvol_size, r,rd, Nsimd,{
|
||||
|
||||
int n = s / e2;
|
||||
int b = s % e2;
|
||||
int so=r*ostride; // base offset for start of plane
|
||||
int ss= so+n*stride+b;
|
||||
|
||||
coalescedWrite(rb_p[r*subvol_size+s], coalescedRead(Data[ss]));
|
||||
|
||||
});
|
||||
|
||||
//issue segmented reductions in computeStream
|
||||
gpuErr = gpucub::DeviceSegmentedReduce::Reduce(temp_storage_array, temp_storage_bytes, rb_p, d_out, rd, d_offsets, d_offsets+1,::gpucub::Sum(), zero_init, computeStream);
|
||||
if (gpuErr!=gpuSuccess) {
|
||||
std::cout << GridLogError << "Lattice_slicesum_gpu.h: Encountered error during gpucub::DeviceSegmentedReduce::Reduce! Error: " << gpuErr <<std::endl;
|
||||
exit(EXIT_FAILURE);
|
||||
}
|
||||
|
||||
acceleratorCopyFromDeviceAsync(d_out,&lvSum[0],rd*sizeof(vobj),computeStream);
|
||||
|
||||
//sync after copy
|
||||
accelerator_barrier();
|
||||
|
||||
acceleratorFreeDevice(temp_storage_array);
|
||||
acceleratorFreeDevice(d_out);
|
||||
acceleratorFreeDevice(d_offsets);
|
||||
|
||||
|
||||
}
|
||||
|
||||
template<class vobj> inline void sliceSumReduction_cub_large(const vobj *Data, Vector<vobj> &lvSum, const int rd, const int e1, const int e2, const int stride, const int ostride, const int Nsimd) {
|
||||
typedef typename vobj::vector_type vector;
|
||||
const int words = sizeof(vobj)/sizeof(vector);
|
||||
const int osites = rd*e1*e2;
|
||||
commVector<vector>buffer(osites);
|
||||
vector *dat = (vector *)Data;
|
||||
vector *buf = &buffer[0];
|
||||
Vector<vector> lvSum_small(rd);
|
||||
vector *lvSum_ptr = (vector *)&lvSum[0];
|
||||
|
||||
for (int w = 0; w < words; w++) {
|
||||
accelerator_for(ss,osites,1,{
|
||||
buf[ss] = dat[ss*words+w];
|
||||
});
|
||||
|
||||
sliceSumReduction_cub_small(buf,lvSum_small,rd,e1,e2,stride, ostride,Nsimd);
|
||||
|
||||
for (int r = 0; r < rd; r++) {
|
||||
lvSum_ptr[w+words*r]=lvSum_small[r];
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
|
||||
}
|
||||
|
||||
template<class vobj> inline void sliceSumReduction_cub(const Lattice<vobj> &Data, Vector<vobj> &lvSum, const int rd, const int e1, const int e2, const int stride, const int ostride, const int Nsimd)
|
||||
{
|
||||
autoView(Data_v, Data, AcceleratorRead); //hipcub/cub cannot deal with large vobjs so we split into small/large case.
|
||||
if constexpr (sizeof(vobj) <= 256) {
|
||||
sliceSumReduction_cub_small(&Data_v[0], lvSum, rd, e1, e2, stride, ostride, Nsimd);
|
||||
}
|
||||
else {
|
||||
sliceSumReduction_cub_large(&Data_v[0], lvSum, rd, e1, e2, stride, ostride, Nsimd);
|
||||
}
|
||||
}
|
||||
#endif
|
||||
|
||||
|
||||
#if defined(GRID_SYCL)
|
||||
template<class vobj> inline void sliceSumReduction_sycl(const Lattice<vobj> &Data, Vector <vobj> &lvSum, const int &rd, const int &e1, const int &e2, const int &stride, const int &ostride, const int &Nsimd)
|
||||
{
|
||||
typedef typename vobj::scalar_object sobj;
|
||||
size_t subvol_size = e1*e2;
|
||||
|
||||
vobj *mysum = (vobj *) malloc_shared(sizeof(vobj),*theGridAccelerator);
|
||||
vobj vobj_zero;
|
||||
zeroit(vobj_zero);
|
||||
|
||||
commVector<vobj> reduction_buffer(rd*subvol_size);
|
||||
|
||||
auto rb_p = &reduction_buffer[0];
|
||||
|
||||
autoView(Data_v, Data, AcceleratorRead);
|
||||
|
||||
//prepare reduction buffer
|
||||
accelerator_for2d( s,subvol_size, r,rd, (size_t)Nsimd,{
|
||||
|
||||
int n = s / e2;
|
||||
int b = s % e2;
|
||||
int so=r*ostride; // base offset for start of plane
|
||||
int ss= so+n*stride+b;
|
||||
|
||||
coalescedWrite(rb_p[r*subvol_size+s], coalescedRead(Data_v[ss]));
|
||||
|
||||
});
|
||||
|
||||
for (int r = 0; r < rd; r++) {
|
||||
mysum[0] = vobj_zero; //dirty hack: cannot pass vobj_zero as identity to sycl::reduction as its not device_copyable
|
||||
theGridAccelerator->submit([&](cl::sycl::handler &cgh) {
|
||||
auto Reduction = cl::sycl::reduction(mysum,std::plus<>());
|
||||
cgh.parallel_for(cl::sycl::range<1>{subvol_size},
|
||||
Reduction,
|
||||
[=](cl::sycl::id<1> item, auto &sum) {
|
||||
auto s = item[0];
|
||||
sum += rb_p[r*subvol_size+s];
|
||||
});
|
||||
});
|
||||
theGridAccelerator->wait();
|
||||
lvSum[r] = mysum[0];
|
||||
}
|
||||
|
||||
free(mysum,*theGridAccelerator);
|
||||
}
|
||||
#endif
|
||||
|
||||
template<class vobj> inline void sliceSumReduction_cpu(const Lattice<vobj> &Data, Vector<vobj> &lvSum, const int &rd, const int &e1, const int &e2, const int &stride, const int &ostride, const int &Nsimd)
|
||||
{
|
||||
// sum over reduced dimension planes, breaking out orthog dir
|
||||
// Parallel over orthog direction
|
||||
autoView( Data_v, Data, CpuRead);
|
||||
thread_for( r,rd, {
|
||||
int so=r*ostride; // 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;
|
||||
lvSum[r]=lvSum[r]+Data_v[ss];
|
||||
}
|
||||
}
|
||||
});
|
||||
}
|
||||
|
||||
template<class vobj> inline void sliceSumReduction(const Lattice<vobj> &Data, Vector<vobj> &lvSum, const int &rd, const int &e1, const int &e2, const int &stride, const int &ostride, const int &Nsimd)
|
||||
{
|
||||
#if defined(GRID_CUDA) || defined(GRID_HIP)
|
||||
|
||||
sliceSumReduction_cub(Data, lvSum, rd, e1, e2, stride, ostride, Nsimd);
|
||||
|
||||
#elif defined(GRID_SYCL)
|
||||
|
||||
sliceSumReduction_sycl(Data, lvSum, rd, e1, e2, stride, ostride, Nsimd);
|
||||
|
||||
#else
|
||||
sliceSumReduction_cpu(Data, lvSum, rd, e1, e2, stride, ostride, Nsimd);
|
||||
|
||||
#endif
|
||||
}
|
||||
|
||||
|
||||
NAMESPACE_END(Grid);
|
@ -34,7 +34,7 @@ class GridTracer {
|
||||
};
|
||||
inline void tracePush(const char *name) { roctxRangePushA(name); }
|
||||
inline void tracePop(const char *name) { roctxRangePop(); }
|
||||
inline int traceStart(const char *name) { roctxRangeStart(name); }
|
||||
inline int traceStart(const char *name) { return roctxRangeStart(name); }
|
||||
inline void traceStop(int ID) { roctxRangeStop(ID); }
|
||||
#endif
|
||||
|
||||
|
@ -225,6 +225,8 @@ inline void acceleratorFreeShared(void *ptr){ cudaFree(ptr);};
|
||||
inline void acceleratorFreeDevice(void *ptr){ cudaFree(ptr);};
|
||||
inline void acceleratorCopyToDevice(void *from,void *to,size_t bytes) { cudaMemcpy(to,from,bytes, cudaMemcpyHostToDevice);}
|
||||
inline void acceleratorCopyFromDevice(void *from,void *to,size_t bytes){ cudaMemcpy(to,from,bytes, cudaMemcpyDeviceToHost);}
|
||||
inline void acceleratorCopyToDeviceAsync(void *from, void *to, size_t bytes, cudaStream_t stream = copyStream) { cudaMemcpyAsync(to,from,bytes, cudaMemcpyHostToDevice, stream);}
|
||||
inline void acceleratorCopyFromDeviceAsync(void *from, void *to, size_t bytes, cudaStream_t stream = copyStream) { cudaMemcpyAsync(to,from,bytes, cudaMemcpyDeviceToHost, stream);}
|
||||
inline void acceleratorMemSet(void *base,int value,size_t bytes) { cudaMemset(base,value,bytes);}
|
||||
inline void acceleratorCopyDeviceToDeviceAsynch(void *from,void *to,size_t bytes) // Asynch
|
||||
{
|
||||
@ -443,6 +445,8 @@ inline void acceleratorFreeShared(void *ptr){ auto r=hipFree(ptr);};
|
||||
inline void acceleratorFreeDevice(void *ptr){ auto r=hipFree(ptr);};
|
||||
inline void acceleratorCopyToDevice(void *from,void *to,size_t bytes) { auto r=hipMemcpy(to,from,bytes, hipMemcpyHostToDevice);}
|
||||
inline void acceleratorCopyFromDevice(void *from,void *to,size_t bytes){ auto r=hipMemcpy(to,from,bytes, hipMemcpyDeviceToHost);}
|
||||
inline void acceleratorCopyToDeviceAsync(void *from, void *to, size_t bytes, hipStream_t stream = copyStream) { auto r = hipMemcpyAsync(to,from,bytes, hipMemcpyHostToDevice, stream);}
|
||||
inline void acceleratorCopyFromDeviceAsync(void *from, void *to, size_t bytes, hipStream_t stream = copyStream) { auto r = hipMemcpyAsync(to,from,bytes, hipMemcpyDeviceToHost, stream);}
|
||||
//inline void acceleratorCopyDeviceToDeviceAsynch(void *from,void *to,size_t bytes) { hipMemcpy(to,from,bytes, hipMemcpyDeviceToDevice);}
|
||||
//inline void acceleratorCopySynchronise(void) { }
|
||||
inline void acceleratorMemSet(void *base,int value,size_t bytes) { auto r=hipMemset(base,value,bytes);}
|
||||
|
Reference in New Issue
Block a user