OpenFPM  5.2.0
Project that contain the implementation of distributed structures
vector_gpu_unit_tests.cu
1 /*
2  * map_vector_cuda_funcs_tests.cu
3  *
4  * Created on: Aug 17, 2018
5  * Author: i-bird
6  */
7 
8 
9 #include "util/cuda_util.hpp"
10 #include "config.h"
11 #define BOOST_TEST_DYN_LINK
12 #include <boost/test/unit_test.hpp>
13 #include "Vector/map_vector.hpp"
14 
15 template<typename vector_vector_type, typename vector_out_type>
16 __global__ void vv_test_size(vector_vector_type vvt, vector_out_type out)
17 {
18  int p = threadIdx.x + blockIdx.x * blockDim.x;
19 
20  if (p >= vvt.size()) return;
21 
22  out.template get<0>(p) = vvt.template get<0>(p).size();
23 }
24 
25 template<typename vector_vector_type, typename vector_out_type>
26 __global__ void vv_test_pointer(vector_vector_type vvt, vector_out_type out)
27 {
28  int p = threadIdx.x + blockIdx.x * blockDim.x;
29 
30  if (p >= vvt.size()) return;
31 
32  out.template get<0>(p) = (size_t)vvt.template get<0>(p).template getPointer<0>();
33  out.template get<1>(p) = (size_t)vvt.template get<0>(p).template getPointer<1>();
34 }
35 
36 template<typename vector_vector_type, typename vector_out_type>
37 __global__ void vv_test_data_get(vector_vector_type vvt, vector_out_type out, int i_sz)
38 {
39  int p = threadIdx.x + blockIdx.x * blockDim.x;
40 
41  if (p >= out.size()) return;
42 
43  int id1 = p/i_sz;
44  int id2 = p%i_sz;
45 
46  out.template get<0>(p)[0] = (size_t)vvt.template get<0>(id1).template get<0>(id2)[0];
47  out.template get<0>(p)[1] = (size_t)vvt.template get<0>(id1).template get<0>(id2)[1];
48  out.template get<0>(p)[2] = (size_t)vvt.template get<0>(id1).template get<0>(id2)[2];
49 
50  out.template get<1>(p)[0] = (size_t)vvt.template get<0>(id1).template get<1>(id2)[0];
51  out.template get<1>(p)[1] = (size_t)vvt.template get<0>(id1).template get<1>(id2)[1];
52  out.template get<1>(p)[2] = (size_t)vvt.template get<0>(id1).template get<1>(id2)[2];
53 }
54 
55 BOOST_AUTO_TEST_SUITE( vector_cuda_tests )
56 
57 
58 
59 BOOST_AUTO_TEST_CASE ( test_vector_of_vector_gpu )
60 {
62 
64 
65  vb_int_proc.resize_no_device(5);
66 
68  ptr_dev.resize(5);
69 
70  for (size_t i = 0 ; i< vb_int_proc.size() ; i++)
71  {
72  vb_int_proc.template get<0>(i).resize(7);
73 
74  for (size_t j = 0 ; j < vb_int_proc.template get<0>(i).size() ; j++)
75  {
76  for (size_t k = 0 ; k < 3 ; k++)
77  {
78  vb_int_proc.template get<0>(i).template get<0>(j)[k] = i+j+k;
79  vb_int_proc.template get<0>(i).template get<1>(j)[k] = 100+i+j+k;
80  }
81  }
82 
83  vb_int_proc.template get<0>(i).template hostToDevice<0,1>();
84  ptr_dev.get(i).first = vb_int_proc.template get<0>(i).template getDeviceBuffer<0>();
85  ptr_dev.get(i).second = vb_int_proc.template get<0>(i).template getDeviceBuffer<1>();
86  }
87 
88  vb_int_proc.template hostToDevice<0>();
89 
91  out.resize(vb_int_proc.size());
92 
93  auto ite = vb_int_proc.getGPUIterator();
94 
95  CUDA_LAUNCH_DIM3((vv_test_size<decltype(vb_int_proc.toKernel()),decltype(out.toKernel())>),ite.wthr,ite.thr,vb_int_proc.toKernel(),out.toKernel());
96 
97  out.deviceToHost<0>();
98 
99  for (size_t i = 0 ; i < out.size() ; i++)
100  {
101  BOOST_REQUIRE_EQUAL(out.template get<0>(i),7);
102  }
103 
105  out_pointer.resize(vb_int_proc.size());
106 
107  CUDA_LAUNCH_DIM3((vv_test_pointer<decltype(vb_int_proc.toKernel()),decltype(out_pointer.toKernel())>),ite.wthr,ite.thr,vb_int_proc.toKernel(),out_pointer.toKernel());
108 
109  out_pointer.deviceToHost<0,1>();
110 
111  for (size_t i = 0 ; i < out_pointer.size() ; i++)
112  {
113  BOOST_REQUIRE_EQUAL((size_t)out_pointer.template get<0>(i),(size_t)ptr_dev.get(i).first);
114  BOOST_REQUIRE_EQUAL((size_t)out_pointer.template get<1>(i),(size_t)ptr_dev.get(i).second);
115  }
116 
118  out_data.resize(vb_int_proc.size()*7);
119 
120  auto ite2 = out_data.getGPUIterator();
121 
122  CUDA_LAUNCH_DIM3((vv_test_data_get<decltype(vb_int_proc.toKernel()),decltype(out_data.toKernel())>),ite2.wthr,ite2.thr,vb_int_proc.toKernel(),out_data.toKernel(),7);
123 
124  out_data.template deviceToHost<0,1>();
125 
126  size_t i_sz = 7;
127 
128  for (size_t p = 0 ; p < out_data.size() ; p++)
129  {
130  int id1 = p/i_sz;
131  int id2 = p%i_sz;
132 
133  BOOST_REQUIRE_EQUAL(out_data.template get<0>(p)[0],vb_int_proc.template get<0>(id1).template get<0>(id2)[0] );
134  BOOST_REQUIRE_EQUAL(out_data.template get<0>(p)[1],vb_int_proc.template get<0>(id1).template get<0>(id2)[1] );
135  BOOST_REQUIRE_EQUAL(out_data.template get<0>(p)[2],vb_int_proc.template get<0>(id1).template get<0>(id2)[2] );
136 
137  BOOST_REQUIRE_EQUAL(out_data.template get<1>(p)[0],vb_int_proc.template get<0>(id1).template get<1>(id2)[0] );
138  BOOST_REQUIRE_EQUAL(out_data.template get<1>(p)[1],vb_int_proc.template get<0>(id1).template get<1>(id2)[1] );
139  BOOST_REQUIRE_EQUAL(out_data.template get<1>(p)[2],vb_int_proc.template get<0>(id1).template get<1>(id2)[2] );
140  }
141 }
142 
143 
144 BOOST_AUTO_TEST_SUITE_END()
Implementation of 1-D std::vector like structure.
Definition: map_vector.hpp:204
size_t size()
Stub size.
Definition: map_vector.hpp:212
Transform the boost::fusion::vector into memory specification (memory_traits)
Definition: memory_conf.hpp:84