OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
trash_bin.hpp
1 /*
2  * trash_bin.hpp
3  *
4  * Created on: Nov 24, 2014
5  * Author: i-bird
6  */
7 
8 #ifndef TRASH_BIN_HPP_
9 #define TRASH_BIN_HPP_
10 
11 
20 template<typename T>
21 class E_p
22 {
23 public:
25  typedef typename boost::fusion::result_of::push_front<typename T::type, size_t>::type E_a;
26 
28  typedef E_a type;
29 
30  E_a data;
31 };
32 
41 template<typename T>
42 class V_p
43 {
44 public:
46  typedef typename boost::fusion::result_of::push_front<typename T::type, size_t>::type V_a;
47 
49  typedef V_a type;
50 
51  V_a data;
52 };
53 
54 
55 #endif /* TRASH_BIN_HPP_ */
56 
57 
66 /* template <unsigned int p>inline typename type_cpu_prop<p,T>::type & getBoostVector(grid_key_dx<dim> & v1)
67 {
68 #ifdef MEMLEAK_CHECK
69  check_valid(&boost::fusion::at_c<p>(data.mem_r->operator[](g1.LinId(v1))));
70 #endif
71  return boost::fusion::at_c<p>(data.mem_r->operator[](g1.LinId(v1)));
72 }*/
73 
74 
75 template<unsigned int prp>
76 struct red_max
77 {
78  template<typename aggr_vect, typename encap_type>
79  static void red(aggr_vect & aggr, encap_type & ec)
80  {
81  auto op1 = ec.template get<prp>();
82  auto op2 = aggr.template get<prp>();
83 
84  aggr.template get<prp>() = (op2 < op1)?op1:op2;
85  }
86 
87  template<typename red_type, typename aggr_vect, typename BlockReduce>
88  static red_type red_final(red_type * rt, aggr_vect & aggr)
89  {
90  return BlockReduce(rt).Max(aggr.template get<prp>());
91  }
92 };
93 
94 template<unsigned int prp>
95 struct red_min
96 {
97  template<typename aggr_vect, typename encap_type>
98  static void red(aggr_vect & aggr, encap_type & ec)
99  {
100  auto op1 = ec.template get<prp>();
101  auto op2 = aggr.template get<prp>();
102 
103  aggr.template get<prp>() = (op2 > op1)?op1:op2;
104  }
105 
106 
107  template<typename red_type, typename aggr_vect, typename BlockReduce>
108  static red_type red_final(red_type * rt, aggr_vect & aggr)
109  {
110  return BlockReduce(rt).Min(aggr.template get<prp>());
111  }
112 };
113 
114 template<unsigned int prp>
115 struct red_sum
116 {
117  template<typename aggr_vect, typename encap_type>
118  static void red(aggr_vect & aggr, encap_type & ec)
119  {
120  aggr.template get<prp>() += ec.template get<prp>();
121  }
122 
123 
124  template<typename red_type, typename aggr_vect, typename BlockReduce>
125  static red_type red_final(red_type * rt, aggr_vect & aggr)
126  {
127  return BlockReduce(rt).Sum(aggr.template get<prp>());
128  }
129 };
130 
141 template<typename aggr_vect, typename reduction_vectors>
142 struct reduce_op
143 {
145  aggr_vect & red;
147  reduction_vectors & input;
148 
155  __device__ __host__ inline reduce_op(aggr_vect & red, reduction_vectors & input)
156  :red(red),input(input)
157  {
158  };
159 
161  template<typename T>
162  __device__ __host__ inline void operator()(T& t) const
163  {
164  boost::mpl::at<reduction_vectors,boost::mpl::int_<T::value>>::red(red,input);
165  }
166 };
167 
178 template<typename aggr_vect, typename red_type, typename reduction_vectors>
180 {
182  aggr_vect & red;
183 
185  red_type * red_space;
186 
187  typedef typename boost::mpl::size<reduction_vectors>::type nred;
188 
189  red_type red_final[nred::value];
190 
197  __device__ __host__ inline reduce_op_final(aggr_vect & red, red_type & red_space)
199  {
200  };
201 
203  template<typename T>
204  __device__ __host__ inline void operator()(T& t)
205  {
206  // fill the block
207  red_space[threadIdx.x] = red.template get<T::value>();
208 
209  red_final[T::value] = boost::mpl::at<reduction_vectors,boost::mpl::int_<T::value>>::red_final(red_space);
210  }
211 };
212 
213 
224 template<typename red_type, typename encap_type , typename reduction_vectors>
226 {
227  typedef typename boost::mpl::size<reduction_vectors>::type nred;
228 
229  red_type (& red_final)[nred::value];
230 
231  encap_type & destination;
232 
239  __device__ __host__ inline store_reduce_op_final(encap_type & destination, red_type (& red_final)[nred::value])
240  :red_final(red_final),destination(destination)
241  {
242  };
243 
245  template<typename T>
246  __device__ __host__ inline void operator()(T& t)
247  {
248  typedef typename boost::mpl::at<reduction_vectors,boost::mpl::int_<T::value>>::prop prp;
249 
250  destination.template get<prp::value>() = red_final[T::value];
251  }
252 };
253 
254 template<typename vector_index_type, typename vector_index_type2, typename vector_data_type, typename reduction_operation_vector,
255  typename red_type, unsigned int blockSize>
256 __global__ void reduce_data(vector_index_type v_offsets, vector_index_type2 v_ord, vector_data_type v_data, vector_data_type v_data_red)
257 {
258  int offset_start = v_offsets.template get<0>(blockIdx.x);
259  int offset_end = v_offsets.template get<0>(blockIdx.x+1);
260  int n_block_reduce = (offset_end - offset_start) / blockDim.x;
261 
262  // Specialize BlockReduce for a 1D block
263  typedef cub::BlockReduce<red_type, blockSize> BlockReduceT;
264 
265  __shared__ typename BlockReduceT::TempStorage temp_storage;
266 
267  typename vector_data_type::value_type red_v;
268 
269  typedef boost::mpl::size<reduction_operation_vector> n_reductions;
270 
271  int i = 0;
272  for ( ; i < n_block_reduce ; i++)
273  {
274  reduce_op<typename vector_data_type::value_type,reduction_operation_vector> ro(red_v,v_data.get_o(offset_start+i*blockDim.x+threadIdx.x));
275 
276  boost::mpl::for_each_ref<boost::mpl::range_c<int,0,n_reductions::value>>(ro);
277  }
278 
279  // within block reduction
281  boost::mpl::for_each_ref<boost::mpl::range_c<int,0,n_reductions::value>>(rof);
282 
283  if (threadIdx.x == 0)
284  {
285 
286  }
287 }
aggr_vect & red
reduction aggregate
Definition: trash_bin.hpp:182
boost::fusion::result_of::push_front< typename T::type, size_t >::type E_a
insert a size_t property to the boost::vector
Definition: trash_bin.hpp:25
The BlockReduce class provides collective methods for computing a parallel reduction of items partiti...
boost::fusion::result_of::push_front< typename T::type, size_t >::type V_a
insert a size_t property to the boost::vector
Definition: trash_bin.hpp:46
V_a type
type for the vertex
Definition: trash_bin.hpp:49
__device__ __host__ reduce_op(aggr_vect &red, reduction_vectors &input)
constructor
Definition: trash_bin.hpp:155
__device__ __host__ void operator()(T &t) const
It call the copy function for each property.
Definition: trash_bin.hpp:162
Get the reference of the selected element.
Definition: trash_bin.hpp:76
reduction_vectors & input
encapsulated input object
Definition: trash_bin.hpp:147
Vertex class that encapsulate an object T.
Definition: trash_bin.hpp:42
E_a type
type for the edge
Definition: trash_bin.hpp:28
this class is a functor for "for_each" algorithm
Definition: trash_bin.hpp:179
this class is a functor for "for_each" algorithm
Definition: trash_bin.hpp:225
__device__ __host__ reduce_op_final(aggr_vect &red, red_type &red_space)
constructor
Definition: trash_bin.hpp:197
__device__ __host__ void operator()(T &t)
It call the copy function for each property.
Definition: trash_bin.hpp:246
__device__ __host__ void operator()(T &t)
It call the copy function for each property.
Definition: trash_bin.hpp:204
aggr_vect & red
reduction aggregate
Definition: trash_bin.hpp:145
this class is a functor for "for_each" algorithm
Definition: trash_bin.hpp:142
__device__ __host__ store_reduce_op_final(encap_type &destination, red_type(&red_final)[nred::value])
constructor
Definition: trash_bin.hpp:239
temporal buffer for reductions
Edge class that encapsulate an object T.
Definition: trash_bin.hpp:21
red_type * red_space
block to reduce
Definition: trash_bin.hpp:185