StatMech
Loading...
Searching...
No Matches
mBodyOperatorSpace_GPUtest.cu File Reference

Functions

template<typename Scalar >
__global__ void create_mBodyOperatorSpace_kernel (mBodyOperatorSpace< Scalar > **ptrptr, int m, int sysSize, int dimLoc)
 
int main (int argc, char **argv)
 

Function Documentation

◆ create_mBodyOperatorSpace_kernel()

template<typename Scalar >
__global__ void create_mBodyOperatorSpace_kernel ( mBodyOperatorSpace< Scalar > **  ptrptr,
int  m,
int  sysSize,
int  dimLoc 
)
11 {
12 int const idx = blockIdx.x * blockDim.x + threadIdx.x;
13 if(idx > 1) return;
14 *ptrptr = new mBodyOperatorSpace<Scalar>(m, std::move(ManyBodySpinSpace(sysSize, dimLoc)));
15 debug_printf("%s:\n\t m=%d, sysSize=%d,dimLoc=%d\n\t Created a ManyBodyOperatorSpace object at "
16 "ptr=%p\n ",
17 __PRETTY_FUNCTION__, m, sysSize, dimLoc, *ptrptr);
18}
Definition HilbertSpace.hpp:423
Definition OperatorSpace.hpp:430

◆ main()

int main ( int  argc,
char **  argv 
)
20 {
21 if(argc != 4) {
22 std::cerr << "Usage: 1.(This) 2.(m) 3.(sysSize) 4.(dimLoc)\n"
23 << " argc=" << argc << std::endl;
24 std::exit(EXIT_FAILURE);
25 }
26 int const m = std::atoi(argv[1]);
27 int const sysSize = std::atoi(argv[2]);
28 int const dimLoc = std::atoi(argv[3]);
29
31 cuCHECK(cudaMalloc((void***)&ptrTo_dPtr, sizeof(mBodyOperatorSpace< Complex_t<double>>*)));
32 create_mBodyOperatorSpace_kernel<<<1, 1>>>(ptrTo_dPtr, m, sysSize, dimLoc);
33 cuCHECK(cudaGetLastError());
34 cuCHECK(cudaDeviceSynchronize());
35
37 cuCHECK(cudaMemcpy(&dPtr, ptrTo_dPtr, sizeof(mBodyOperatorSpace< Complex_t<double>>*),
38 cudaMemcpyDeviceToHost));
39 std::cout << "dPtr=" << dPtr << std::endl;
40 cuCHECK(cudaFree(ptrTo_dPtr));
41
42 testManyBodySpaceBase_kernel<<<1, 1>>>(dPtr);
43 cuCHECK(cudaGetLastError());
44 cuCHECK(cudaDeviceSynchronize());
45 std::cout << "Passed \"testManyBodySpaceBase_kernel\"\n" << std::endl;
46
47 testOperatorSpaceBase_kernel<<<1, 1>>>(dPtr);
48 cuCHECK(cudaGetLastError());
49 cuCHECK(cudaDeviceSynchronize());
50 std::cout << "Passed \"testOperatorSpaceBase_kernel\"\n" << std::endl;
51
52 return EXIT_SUCCESS;
53}
Definition mytypes.hpp:147
cuCHECK(cudaFuncGetAttributes(&attr, MatrixElementsInSector))