StatMech
Loading...
Searching...
No Matches
ManyBodySpaceBase.cuh File Reference

Classes

class  ObjectOnGPU< ManyBodySpaceBase< Derived > >
 

Functions

template<class Derived >
void copyTransEqClass (ObjectOnGPU< ManyBodySpaceBase< Derived > > &obj, ManyBodySpaceBase< Derived > const &src)
 
template<typename Derived >
__global__ void setTransEqDim_kernel (ManyBodySpaceBase< Derived > *obj, int transEqDim)
 
template<typename Derived >
__global__ void copyTransEqClass_kernel (ManyBodySpaceBase< Derived > *obj, int *srcTransEqClassRep, int *srcTransPeriod)
 

Function Documentation

◆ copyTransEqClass()

template<class Derived >
void copyTransEqClass ( ObjectOnGPU< ManyBodySpaceBase< Derived > > &  obj,
ManyBodySpaceBase< Derived > const &  src 
)
89 {
90 debug_print(__PRETTY_FUNCTION__ << "\n\tsrc.transEqDim() = " << src.transEqDim() << "\n");
91 if(src.transEqDim() == -1) return;
92
93 int* srcTransEqClassRep = nullptr;
94 int* srcTransPeriod = nullptr;
95 cuCHECK(cudaMalloc(&srcTransEqClassRep, src.transEqDim() * sizeof(int)));
96 cuCHECK(cudaMalloc(&srcTransPeriod, src.transEqDim() * sizeof(int)));
97 cuCHECK(cudaMemcpy(srcTransEqClassRep, src.transEqClassRep().data(),
98 src.transEqDim() * sizeof(int), cudaMemcpyHostToDevice));
99 cuCHECK(cudaMemcpy(srcTransPeriod, src.transPeriod().data(), src.transEqDim() * sizeof(int),
100 cudaMemcpyHostToDevice));
101 setTransEqDim_kernel<<<1, 1>>>(obj.ptr(), src.transEqDim());
102
103 cudaDeviceProp deviceProp;
104 cudaGetDeviceProperties(&deviceProp, 0); // 0-th device
105 constexpr int warpSize = 32;
106 int const nWarps = (src.transEqDim() % warpSize == 0 ? src.transEqDim() / warpSize
107 : src.transEqDim() / warpSize + 1);
108 int const nWarpsPerSM = (nWarps % deviceProp.multiProcessorCount == 0
109 ? nWarps / deviceProp.multiProcessorCount
110 : nWarps / deviceProp.multiProcessorCount + 1);
111 int const nThreads = min(nWarpsPerSM * warpSize, deviceProp.maxThreadsPerBlock);
112 int const nBlocks = (src.transEqDim() % nThreads == 0 ? src.transEqDim() / nThreads
113 : src.transEqDim() / nThreads + 1);
114 // std::cout << "\tdeviceProp.multiProcessorCount = " << deviceProp.multiProcessorCount
115 // << ", deviceProp.sharedMemPerBlock = " << deviceProp.sharedMemPerBlock << std::endl;
116 // std::cout << "\tnBlocks = " << nBlocks << ", nThreads = " << nThreads << std::endl;
117 copyTransEqClass_kernel<<<nBlocks, nThreads>>>(obj.ptr(), srcTransEqClassRep, srcTransPeriod);
118 cuCHECK(cudaGetLastError());
119 cuCHECK(cudaFree(srcTransEqClassRep));
120 cuCHECK(cudaFree(srcTransPeriod));
121}
__host__ __device__ Vector_t const & transPeriod() const
Definition HilbertSpace.hpp:317
__host__ __device__ int transEqDim() const
Definition HilbertSpace.hpp:309
__host__ __device__ Vector_t const & transEqClassRep() const
Definition HilbertSpace.hpp:311
Object_t * ptr() const
Definition ObjectOnGPU.cuh:144
cuCHECK(cudaFuncGetAttributes(&attr, MatrixElementsInSector))
debug_print("# Determining GPU configuration.")

◆ copyTransEqClass_kernel()

template<typename Derived >
__global__ void copyTransEqClass_kernel ( ManyBodySpaceBase< Derived > *  obj,
int *  srcTransEqClassRep,
int *  srcTransPeriod 
)
77 {
78 int const idx = blockIdx.x * blockDim.x + threadIdx.x;
79 if(idx >= obj->transEqDim()) return;
80 obj->m_transEqClassRep(idx) = srcTransEqClassRep[idx];
81 obj->m_transPeriod(idx) = srcTransPeriod[idx];
82 for(auto trans = 0; trans != srcTransPeriod[idx]; ++trans) {
83 auto translated = obj->translate(srcTransEqClassRep[idx], trans);
84 obj->m_stateToTransEqClass(translated) = idx;
85 }
86}
__host__ __device__ int translate(int state, int trans) const
Translate the input state to the left by one.
Definition HilbertSpace.hpp:300
Vector_t m_transPeriod
Definition HilbertSpace.hpp:145
Vector_t m_stateToTransEqClass
Definition HilbertSpace.hpp:146
Vector_t m_transEqClassRep
Definition HilbertSpace.hpp:144

◆ setTransEqDim_kernel()

template<typename Derived >
__global__ void setTransEqDim_kernel ( ManyBodySpaceBase< Derived > *  obj,
int  transEqDim 
)
68 {
69 debug_printf("%s\n\ttransEqDim=%d\n", __PRETTY_FUNCTION__, transEqDim);
70 obj->m_transEqDim = transEqDim;
71 obj->m_transEqClassRep.resize(transEqDim);
72 obj->m_transPeriod.resize(transEqDim);
73 obj->m_stateToTransEqClass.resize(obj->dim());
74}
__host__ __device__ int dim() const
Definition HilbertSpace.hpp:34
int m_transEqDim
Definition HilbertSpace.hpp:143