OpenFPM  5.2.0
Project that contain the implementation of distributed structures
cudify_sequential.hpp
1 #ifndef CUDIFY_SEQUENCIAL_HPP_
2 #define CUDIFY_SEQUENCIAL_HPP_
3 
4 #define CUDA_ON_BACKEND CUDA_BACKEND_SEQUENTIAL
5 
6 #include "config.h"
7 
8 constexpr int default_kernel_wg_threads_ = 1024;
9 
10 #include "util/cudify/cudify_hardware_cpu.hpp"
11 
12 #ifdef HAVE_BOOST_CONTEXT
13 
14 #include <boost/bind/bind.hpp>
15 #include <type_traits>
16 #ifdef HAVE_BOOST_CONTEXT
17 #include <boost/context/continuation.hpp>
18 #endif
19 #include <vector>
20 #include <limits>
21 #include <string.h>
22 
23 
24 #ifndef CUDIFY_BOOST_CONTEXT_STACK_SIZE
25 #define CUDIFY_BOOST_CONTEXT_STACK_SIZE 8192
26 #endif
27 
28 extern std::vector<void *>mem_stack;
29 
30 extern thread_local dim3 threadIdx;
31 extern thread_local dim3 blockIdx;
32 
33 extern dim3 blockDim;
34 extern dim3 gridDim;
35 
36 extern std::vector<void *> mem_stack;
37 extern std::vector<boost::context::detail::fcontext_t> contexts;
38 extern thread_local void * par_glob;
39 extern thread_local boost::context::detail::fcontext_t main_ctx;
40 
41 static void __syncthreads()
42 {
43  boost::context::detail::jump_fcontext(main_ctx,par_glob);
44 };
45 
46 
47 
48 extern int thread_local vct_atomic_add;
49 extern int thread_local vct_atomic_rem;
50 
51 
52 namespace cub
53 {
54  template<typename T, unsigned int dim>
55  class BlockScan
56  {
57  public:
58  typedef std::array<T,dim> TempStorage;
59 
60  private:
61  TempStorage & tmp;
62 
63  public:
64 
65 
66 
67  BlockScan(TempStorage & tmp)
68  :tmp(tmp)
69  {};
70 
71  void ExclusiveSum(T & in, T & out)
72  {
73  tmp[threadIdx.x] = in;
74 
75  __syncthreads();
76 
77  if (threadIdx.x == 0)
78  {
79  T prec = tmp[0];
80  tmp[0] = 0;
81  for (int i = 1 ; i < dim ; i++)
82  {
83  auto next = tmp[i-1] + prec;
84  prec = tmp[i];
85  tmp[i] = next;
86  }
87  }
88 
89  __syncthreads();
90 
91  out = tmp[threadIdx.x];
92  return;
93  }
94  };
95 }
96 
97 
98 template<typename T, typename T2>
99 static T atomicAdd(T * address, T2 val)
100 {
101  T old = *address;
102  *address += val;
103  return old;
104 };
105 
106 
107 namespace gpu
108 {
109  template<typename type_t>
110  struct less_t : public std::binary_function<type_t, type_t, bool> {
111  bool operator()(type_t a, type_t b) const {
112  return a < b;
113  }
114 
115  template<typename type2_t, typename type3_t>
116  bool operator()(type2_t a, type3_t b) const {
117  return a < b;
118  }
119  };
120 /* template<typename type_t>
121  struct less_equal_t : public std::binary_function<type_t, type_t, bool> {
122  bool operator()(type_t a, type_t b) const {
123  return a <= b;
124  }
125  };*/
126  template<typename type_t>
127  struct greater_t : public std::binary_function<type_t, type_t, bool> {
128  bool operator()(type_t a, type_t b) const {
129  return a > b;
130  }
131 
132  template<typename type2_t, typename type3_t>
133  bool operator()(type2_t a, type3_t b) const {
134  return a > b;
135  }
136  };
137 /* template<typename type_t>
138  struct greater_equal_t : public std::binary_function<type_t, type_t, bool> {
139  bool operator()(type_t a, type_t b) const {
140  return a >= b;
141  }
142  };
143  template<typename type_t>
144  struct equal_to_t : public std::binary_function<type_t, type_t, bool> {
145  bool operator()(type_t a, type_t b) const {
146  return a == b;
147  }
148  };
149  template<typename type_t>
150  struct not_equal_to_t : public std::binary_function<type_t, type_t, bool> {
151  bool operator()(type_t a, type_t b) const {
152  return a != b;
153  }
154  };*/
155 
157  // Device-side arithmetic operators.
158 
159  template<typename type_t>
160  struct plus_t : public std::binary_function<type_t, type_t, type_t> {
161  type_t operator()(type_t a, type_t b) const {
162  return a + b;
163  }
164 
165  type_t reduceInitValue() const {
166  return 0;
167  }
168  };
169 
170 /* template<typename type_t>
171  struct minus_t : public std::binary_function<type_t, type_t, type_t> {
172  type_t operator()(type_t a, type_t b) const {
173  return a - b;
174  }
175  };
176 
177  template<typename type_t>
178  struct multiplies_t : public std::binary_function<type_t, type_t, type_t> {
179  type_t operator()(type_t a, type_t b) const {
180  return a * b;
181  }
182  };*/
183 
184  template<typename type_t>
185  struct maximum_t : public std::binary_function<type_t, type_t, type_t> {
186  type_t operator()(type_t a, type_t b) const {
187  return std::max(a, b);
188  }
189 
190  type_t reduceInitValue() const {
191  return std::numeric_limits<T>::min();
192  }
193  };
194 
195  template<typename type_t>
196  struct minimum_t : public std::binary_function<type_t, type_t, type_t> {
197  type_t operator()(type_t a, type_t b) const {
198  return std::min(a, b);
199  }
200 
201  type_t reduceInitValue() const {
202  return std::numeric_limits<T>::max();
203  }
204  };
205 }
206 
207 
208 namespace gpu
209 {
210  template<typename input_it,
211  typename segments_it, typename output_it, typename op_t, typename type_t, typename context_t>
212  void segreduce(input_it input, int count, segments_it segments,
213  int num_segments, output_it output, op_t op, type_t init,
214  context_t& context)
215  {
216  int i = 0;
217  for ( ; i < num_segments - 1; i++)
218  {
219  int j = segments[i];
220  output[i] = input[j];
221  ++j;
222  for ( ; j < segments[i+1] ; j++)
223  {
224  output[i] = op(output[i],input[j]);
225  }
226  }
227 
228  // Last segment
229  int j = segments[i];
230  output[i] = input[j];
231  ++j;
232  for ( ; j < count ; j++)
233  {
234  output[i] = op(output[i],input[j]);
235  }
236  }
237 
238  // Key-value merge.
239  template<typename a_keys_it, typename a_vals_it,
240  typename b_keys_it, typename b_vals_it,
241  typename c_keys_it, typename c_vals_it,
242  typename comp_t, typename context_t>
243  void merge(a_keys_it a_keys, a_vals_it a_vals, int a_count,
244  b_keys_it b_keys, b_vals_it b_vals, int b_count,
245  c_keys_it c_keys, c_vals_it c_vals, comp_t comp, context_t& context)
246  {
247  int a_it = 0;
248  int b_it = 0;
249  int c_it = 0;
250 
251  while (a_it < a_count || b_it < b_count)
252  {
253  if (a_it < a_count)
254  {
255  if (b_it < b_count)
256  {
257  if (comp(b_keys[b_it],a_keys[a_it]))
258  {
259  c_keys[c_it] = b_keys[b_it];
260  c_vals[c_it] = b_vals[b_it];
261  c_it++;
262  b_it++;
263  }
264  else
265  {
266  c_keys[c_it] = a_keys[a_it];
267  c_vals[c_it] = a_vals[a_it];
268  c_it++;
269  a_it++;
270  }
271  }
272  else
273  {
274  c_keys[c_it] = a_keys[a_it];
275  c_vals[c_it] = a_vals[a_it];
276  c_it++;
277  a_it++;
278  }
279  }
280  else
281  {
282  c_keys[c_it] = b_keys[b_it];
283  c_vals[c_it] = b_vals[b_it];
284  c_it++;
285  b_it++;
286  }
287  }
288  }
289 }
290 
291 static void init_wrappers()
292 {}
293 
294 template<typename lambda_f>
295 struct Fun_enc
296 {
297  lambda_f Fn;
298 
299  Fun_enc(lambda_f Fn)
300  :Fn(Fn)
301  {}
302 
303  void run()
304  {
305  Fn();
306  }
307 };
308 
309 template<typename lambda_f>
310 struct Fun_enc_bt
311 {
312  lambda_f Fn;
313  dim3 & blockIdx;
314  dim3 & threadIdx;
315 
316  Fun_enc_bt(lambda_f Fn,dim3 & blockIdx,dim3 & threadIdx)
317  :Fn(Fn),blockIdx(blockIdx),threadIdx(threadIdx)
318  {}
319 
320  void run()
321  {
322  Fn(blockIdx,threadIdx);
323  }
324 };
325 
326 template<typename Fun_enc_type>
327 void launch_kernel(boost::context::detail::transfer_t par)
328 {
329  main_ctx = par.fctx;
330  par_glob = par.data;
331  Fun_enc_type * ptr = (Fun_enc_type *)par.data;
332 
333  ptr->run();
334 
335  boost::context::detail::jump_fcontext(par.fctx,0);
336 }
337 
338 template<typename lambda_f, typename ite_type>
339 static void exe_kernel(lambda_f f, const ite_type & ite)
340 {
341  if (ite.nthrs() == 0 || ite.nblocks() == 0) {return;}
342 
343  if (mem_stack.size() < ite.nthrs())
344  {
345  int old_size = mem_stack.size();
346  mem_stack.resize(ite.nthrs());
347 
348  for (int i = old_size ; i < mem_stack.size() ; i++)
349  {
350  mem_stack[i] = new char [CUDIFY_BOOST_CONTEXT_STACK_SIZE];
351  }
352  }
353 
354  // Resize contexts
355  contexts.resize(mem_stack.size());
356 
357  Fun_enc<lambda_f> fe(f);
358 
359  for (int i = 0 ; i < ite.wthr.z ; i++)
360  {
361  blockIdx.z = i;
362  for (int j = 0 ; j < ite.wthr.y ; j++)
363  {
364  blockIdx.y = j;
365  for (int k = 0 ; k < ite.wthr.x ; k++)
366  {
367  blockIdx.x = k;
368  int nc = 0;
369  for (int it = 0 ; it < ite.thr.z ; it++)
370  {
371  for (int jt = 0 ; jt < ite.thr.y ; jt++)
372  {
373  for (int kt = 0 ; kt < ite.thr.x ; kt++)
374  {
375  contexts[nc] = boost::context::detail::make_fcontext((char *)mem_stack[nc]+CUDIFY_BOOST_CONTEXT_STACK_SIZE-16,CUDIFY_BOOST_CONTEXT_STACK_SIZE,launch_kernel<Fun_enc<lambda_f>>);
376  nc++;
377  }
378  }
379  }
380 
381  bool work_to_do = true;
382  while(work_to_do)
383  {
384  nc = 0;
385  // Work threads
386  for (int it = 0 ; it < ite.thr.z ; it++)
387  {
388  threadIdx.z = it;
389  for (int jt = 0 ; jt < ite.thr.y ; jt++)
390  {
391  threadIdx.y = jt;
392  for (int kt = 0 ; kt < ite.thr.x ; kt++)
393  {
394  threadIdx.x = kt;
395  auto t = boost::context::detail::jump_fcontext(contexts[nc],&fe);
396  contexts[nc] = t.fctx;
397  work_to_do &= (t.data != 0);
398  nc++;
399  }
400  }
401  }
402  }
403  }
404  }
405  }
406 }
407 
408 template<typename lambda_f, typename ite_type>
409 static void exe_kernel_lambda(lambda_f f, ite_type & ite)
410 {
411  if (ite.nthrs() == 0 || ite.nblocks() == 0) {return;}
412 
413  if (mem_stack.size() < ite.nthrs())
414  {
415  int old_size = mem_stack.size();
416  mem_stack.resize(ite.nthrs());
417 
418  for (int i = old_size ; i < mem_stack.size() ; i++)
419  {
420  mem_stack[i] = new char [CUDIFY_BOOST_CONTEXT_STACK_SIZE];
421  }
422  }
423 
424  // Resize contexts
425  contexts.resize(mem_stack.size());
426 
427  bool is_sync_free = true;
428 
429  bool first_block = true;
430 
431  for (int i = 0 ; i < ite.wthr.z ; i++)
432  {
433  for (int j = 0 ; j < ite.wthr.y ; j++)
434  {
435  for (int k = 0 ; k < ite.wthr.x ; k++)
436  {
437  dim3 blockIdx;
438  dim3 threadIdx;
439  Fun_enc_bt<lambda_f> fe(f,blockIdx,threadIdx);
440  if (first_block == true || is_sync_free == false)
441  {
442  blockIdx.z = i;
443  blockIdx.y = j;
444  blockIdx.x = k;
445  int nc = 0;
446  for (int it = 0 ; it < ite.thr.z ; it++)
447  {
448  for (int jt = 0 ; jt < ite.thr.y ; jt++)
449  {
450  for (int kt = 0 ; kt < ite.thr.x ; kt++)
451  {
452  contexts[nc] = boost::context::detail::make_fcontext((char *)mem_stack[nc]+CUDIFY_BOOST_CONTEXT_STACK_SIZE-16,CUDIFY_BOOST_CONTEXT_STACK_SIZE,launch_kernel<Fun_enc_bt<lambda_f>>);
453 
454  nc++;
455  }
456  }
457  }
458 
459  bool work_to_do = true;
460  while(work_to_do)
461  {
462  nc = 0;
463  // Work threads
464  for (int it = 0 ; it < ite.thr.z ; it++)
465  {
466  threadIdx.z = it;
467  for (int jt = 0 ; jt < ite.thr.y ; jt++)
468  {
469  threadIdx.y = jt;
470  for (int kt = 0 ; kt < ite.thr.x ; kt++)
471  {
472  threadIdx.x = kt;
473  auto t = boost::context::detail::jump_fcontext(contexts[nc],&fe);
474  contexts[nc] = t.fctx;
475 
476  work_to_do &= (t.data != 0);
477  is_sync_free &= !(work_to_do);
478  nc++;
479  }
480  }
481  }
482  }
483  }
484  else
485  {
486  blockIdx.z = i;
487  blockIdx.y = j;
488  blockIdx.x = k;
489  int fb = 0;
490  // Work threads
491  for (int it = 0 ; it < ite.thr.z ; it++)
492  {
493  threadIdx.z = it;
494  for (int jt = 0 ; jt < ite.thr.y ; jt++)
495  {
496  threadIdx.y = jt;
497  for (int kt = 0 ; kt < ite.thr.x ; kt++)
498  {
499  threadIdx.x = kt;
500  f(blockIdx,threadIdx);
501  }
502  }
503  }
504  }
505 
506  first_block = false;
507  }
508  }
509  }
510 }
511 
512 template<typename lambda_f, typename ite_type>
513 static void exe_kernel_lambda_tls(lambda_f f, ite_type & ite)
514 {
515  if (ite.nthrs() == 0 || ite.nblocks() == 0) {return;}
516 
517  if (mem_stack.size() < ite.nthrs())
518  {
519  int old_size = mem_stack.size();
520  mem_stack.resize(ite.nthrs());
521 
522  for (int i = old_size ; i < mem_stack.size() ; i++)
523  {
524  mem_stack[i] = new char [CUDIFY_BOOST_CONTEXT_STACK_SIZE];
525  }
526  }
527 
528  // Resize contexts
529  contexts.resize(mem_stack.size());
530 
531  bool is_sync_free = true;
532 
533  bool first_block = true;
534 
535  for (int i = 0 ; i < ite.wthr.z ; i++)
536  {
537  for (int j = 0 ; j < ite.wthr.y ; j++)
538  {
539  for (int k = 0 ; k < ite.wthr.x ; k++)
540  {
541  Fun_enc<lambda_f> fe(f);
542  if (first_block == true || is_sync_free == false)
543  {
544  blockIdx.z = i;
545  blockIdx.y = j;
546  blockIdx.x = k;
547  int nc = 0;
548  for (int it = 0 ; it < ite.thr.z ; it++)
549  {
550  for (int jt = 0 ; jt < ite.thr.y ; jt++)
551  {
552  for (int kt = 0 ; kt < ite.thr.x ; kt++)
553  {
554  contexts[nc] = boost::context::detail::make_fcontext((char *)mem_stack[nc]+CUDIFY_BOOST_CONTEXT_STACK_SIZE-16,CUDIFY_BOOST_CONTEXT_STACK_SIZE,launch_kernel<Fun_enc<lambda_f>>);
555 
556  nc++;
557  }
558  }
559  }
560 
561  bool work_to_do = true;
562  while(work_to_do)
563  {
564  nc = 0;
565  // Work threads
566  for (int it = 0 ; it < ite.thr.z ; it++)
567  {
568  threadIdx.z = it;
569  for (int jt = 0 ; jt < ite.thr.y ; jt++)
570  {
571  threadIdx.y = jt;
572  for (int kt = 0 ; kt < ite.thr.x ; kt++)
573  {
574  threadIdx.x = kt;
575  auto t = boost::context::detail::jump_fcontext(contexts[nc],&fe);
576  contexts[nc] = t.fctx;
577 
578  work_to_do &= (t.data != 0);
579  is_sync_free &= !(work_to_do);
580  nc++;
581  }
582  }
583  }
584  }
585  }
586  else
587  {
588  blockIdx.z = i;
589  blockIdx.y = j;
590  blockIdx.x = k;
591  int fb = 0;
592  // Work threads
593  for (int it = 0 ; it < ite.thr.z ; it++)
594  {
595  threadIdx.z = it;
596  for (int jt = 0 ; jt < ite.thr.y ; jt++)
597  {
598  threadIdx.y = jt;
599  for (int kt = 0 ; kt < ite.thr.x ; kt++)
600  {
601  threadIdx.x = kt;
602  f();
603  }
604  }
605  }
606  }
607 
608  first_block = false;
609  }
610  }
611  }
612 }
613 
614 template<typename lambda_f, typename ite_type>
615 static void exe_kernel_no_sync(lambda_f f, ite_type & ite)
616 {
617  for (int i = 0 ; i < ite.wthr.z ; i++)
618  {
619  blockIdx.z = i;
620  for (int j = 0 ; j < ite.wthr.y ; j++)
621  {
622  blockIdx.y = j;
623  for (int k = 0 ; k < ite.wthr.x ; k++)
624  {
625  blockIdx.x = k;
626  int fb = 0;
627  // Work threads
628  for (int it = 0 ; it < ite.wthr.z ; it++)
629  {
630  threadIdx.z = it;
631  for (int jt = 0 ; jt < ite.wthr.y ; jt++)
632  {
633  threadIdx.y = jt;
634  for (int kt = 0 ; kt < ite.wthr.x ; kt++)
635  {
636  threadIdx.x = kt;
637  f();
638  }
639  }
640  }
641  }
642  }
643  }
644 }
645 
646 #ifdef PRINT_CUDA_LAUNCHES
647 
648 #define CUDA_LAUNCH(cuda_call,ite, ...)\
649  \
650  gridDim.x = ite.wthr.x;\
651  gridDim.y = ite.wthr.y;\
652  gridDim.z = ite.wthr.z;\
653  \
654  blockDim.x = ite.thr.x;\
655  blockDim.y = ite.thr.y;\
656  blockDim.z = ite.thr.z;\
657  \
658  CHECK_SE_CLASS1_PRE\
659  \
660  std::cout << "Launching: " << #cuda_call << std::endl;\
661  \
662  exe_kernel(\
663  [&](boost::context::fiber && main) -> void {\
664  \
665  \
666 \
667  cuda_call(__VA_ARGS__);\
668  },ite);\
669  CHECK_SE_CLASS1_POST(#cuda_call,__VA_ARGS__)\
670  }
671 
672 
673 #define CUDA_LAUNCH_DIM3(cuda_call,wthr_,thr_, ...)\
674  {\
675  dim3 wthr__(wthr_);\
676  dim3 thr__(thr_);\
677  \
678  ite_gpu<1> itg;\
679  itg.wthr = wthr;\
680  itg.thr = thr;\
681  \
682  gridDim.x = wthr__.x;\
683  gridDim.y = wthr__.y;\
684  gridDim.z = wthr__.z;\
685  \
686  blockDim.x = thr__.x;\
687  blockDim.y = thr__.y;\
688  blockDim.z = thr__.z;\
689  \
690  CHECK_SE_CLASS1_PRE\
691  std::cout << "Launching: " << #cuda_call << std::endl;\
692  \
693  exe_kernel(\
694  [&] (boost::context::fiber && main) -> void {\
695  \
696  \
697 \
698  cuda_call(__VA_ARGS__);\
699  \
700  \
701  });\
702  CHECK_SE_CLASS1_POST(#cuda_call,__VA_ARGS__)\
703  }
704 
705 #define CUDA_CHECK()
706 
707 #else
708 
709 #define CUDA_LAUNCH(cuda_call,ite, ...) \
710  {\
711  gridDim.x = ite.wthr.x;\
712  gridDim.y = ite.wthr.y;\
713  gridDim.z = ite.wthr.z;\
714  \
715  blockDim.x = ite.thr.x;\
716  blockDim.y = ite.thr.y;\
717  blockDim.z = ite.thr.z;\
718  \
719  CHECK_SE_CLASS1_PRE\
720  \
721  exe_kernel([&]() -> void {\
722  \
723  \
724  cuda_call(__VA_ARGS__);\
725  \
726  },ite);\
727  \
728  CHECK_SE_CLASS1_POST(#cuda_call,__VA_ARGS__)\
729  }
730 
731 #define CUDA_LAUNCH_LAMBDA(ite,lambda_f) \
732  {\
733  gridDim.x = ite.wthr.x;\
734  gridDim.y = ite.wthr.y;\
735  gridDim.z = ite.wthr.z;\
736  \
737  blockDim.x = ite.thr.x;\
738  blockDim.y = ite.thr.y;\
739  blockDim.z = ite.thr.z;\
740  \
741  CHECK_SE_CLASS1_PRE\
742  \
743  exe_kernel_lambda(lambda_f,ite);\
744  \
745  CHECK_SE_CLASS1_POST("lambda",0)\
746  }
747 
748 #define CUDA_LAUNCH_LAMBDA_TLS(ite,lambda_f) \
749  {\
750  gridDim.x = ite.wthr.x;\
751  gridDim.y = ite.wthr.y;\
752  gridDim.z = ite.wthr.z;\
753  \
754  blockDim.x = ite.thr.x;\
755  blockDim.y = ite.thr.y;\
756  blockDim.z = ite.thr.z;\
757  \
758  CHECK_SE_CLASS1_PRE\
759  \
760  exe_kernel_lambda_tls(lambda_f,ite);\
761  \
762  CHECK_SE_CLASS1_POST("lambda",0)\
763  }
764 
765 #define CUDA_LAUNCH_DIM3(cuda_call,wthr_,thr_, ...)\
766  {\
767  dim3 wthr__(wthr_);\
768  dim3 thr__(thr_);\
769  \
770  ite_gpu<1> itg;\
771  itg.wthr = wthr_;\
772  itg.thr = thr_;\
773  \
774  gridDim.x = wthr__.x;\
775  gridDim.y = wthr__.y;\
776  gridDim.z = wthr__.z;\
777  \
778  blockDim.x = thr__.x;\
779  blockDim.y = thr__.y;\
780  blockDim.z = thr__.z;\
781  \
782  CHECK_SE_CLASS1_PRE\
783  \
784  exe_kernel([&]() -> void {\
785  \
786  cuda_call(__VA_ARGS__);\
787  \
788  },itg);\
789  \
790  CHECK_SE_CLASS1_POST(#cuda_call,__VA_ARGS__)\
791  }
792 
793 #define CUDA_LAUNCH_LAMBDA_DIM3_TLS(wthr_,thr_,lambda_f) \
794  {\
795  dim3 wthr__(wthr_);\
796  dim3 thr__(thr_);\
797  \
798  ite_gpu<1> itg;\
799  itg.wthr = wthr_;\
800  itg.thr = thr_;\
801  gridDim.x = itg.wthr.x;\
802  gridDim.y = itg.wthr.y;\
803  gridDim.z = itg.wthr.z;\
804  \
805  blockDim.x = itg.thr.x;\
806  blockDim.y = itg.thr.y;\
807  blockDim.z = itg.thr.z;\
808  \
809  CHECK_SE_CLASS1_PRE\
810  \
811  exe_kernel_lambda_tls(lambda_f,itg);\
812  \
813  }
814 
815 #define CUDA_LAUNCH_DIM3_DEBUG_SE1(cuda_call,wthr_,thr_, ...)\
816  {\
817  dim3 wthr__(wthr_);\
818  dim3 thr__(thr_);\
819  \
820  ite_gpu<1> itg;\
821  itg.wthr = wthr_;\
822  itg.thr = thr_;\
823  \
824  gridDim.x = wthr__.x;\
825  gridDim.y = wthr__.y;\
826  gridDim.z = wthr__.z;\
827  \
828  blockDim.x = thr__.x;\
829  blockDim.y = thr__.y;\
830  blockDim.z = thr__.z;\
831  \
832  CHECK_SE_CLASS1_PRE\
833  \
834  exe_kernel([&]() -> void {\
835  \
836  cuda_call(__VA_ARGS__);\
837  \
838  },itg);\
839  \
840  }
841 
842 #define CUDA_CHECK()
843 
844 #endif
845 
846 #endif
847 
848 #else
849 
850 constexpr int default_kernel_wg_threads_ = 1024;
851 
852 #endif /* CUDIFY_SEQUENCIAL_HPP_ */
__device__ __forceinline__ void ExclusiveSum(T input, T &output)
Computes an exclusive block-wide prefix scan using addition (+) as the scan operator....
Definition: block_scan.cuh:333
__device__ __forceinline__ BlockScan()
Collective constructor using a private static allocation of shared memory as temporary storage.
Definition: block_scan.cuh:271
Optional outer namespace(s)
OutputIteratorT OffsetT ReductionOpT OuputT init
< [in] The initial value of the reduction