pi_redux_kernel.cu 3.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010,2011,2014,2015 Université de Bordeaux
  4. * Copyright (C) 2011,2012,2015,2017 CNRS
  5. * Copyright (C) 2012 Inria
  6. *
  7. * StarPU is free software; you can redistribute it and/or modify
  8. * it under the terms of the GNU Lesser General Public License as published by
  9. * the Free Software Foundation; either version 2.1 of the License, or (at
  10. * your option) any later version.
  11. *
  12. * StarPU is distributed in the hope that it will be useful, but
  13. * WITHOUT ANY WARRANTY; without even the implied warranty of
  14. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  15. *
  16. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  17. */
  18. /* This counts how many fall inside the circle quarter */
  19. #include <starpu.h>
  20. #define MAXNBLOCKS 128
  21. #define MAXTHREADSPERBLOCK 256
  22. static __global__ void monte_carlo(float *x, float *y, unsigned n, unsigned long *output_cnt)
  23. {
  24. __shared__ unsigned scnt[MAXTHREADSPERBLOCK];
  25. /* Do we have a successful shot ? */
  26. const int tid = threadIdx.x + blockIdx.x*blockDim.x;
  27. const int nthreads = gridDim.x * blockDim.x;
  28. /* Blank the shared mem buffer */
  29. if (threadIdx.x < MAXTHREADSPERBLOCK)
  30. scnt[threadIdx.x] = 0;
  31. __syncthreads();
  32. int ind;
  33. for (ind = tid; ind < n; ind += nthreads)
  34. {
  35. float xval = (2.0f * x[ind] - 1.0f);
  36. float yval = (2.0f * y[ind] - 1.0f);
  37. float dist = (xval*xval + yval*yval);
  38. unsigned long success = (dist <= 1.0f)?1:0;
  39. scnt[threadIdx.x] += success;
  40. }
  41. __syncthreads();
  42. /* Perform a reduction to compute the sum on each thread within that block */
  43. /* NB: We assume that the number of threads per block is a power of 2 ! */
  44. unsigned long s;
  45. for (s = blockDim.x/2; s!=0; s>>=1)
  46. {
  47. if (threadIdx.x < s)
  48. scnt[threadIdx.x] += scnt[threadIdx.x + s];
  49. __syncthreads();
  50. }
  51. /* report the number of successful shots in the block */
  52. if (threadIdx.x == 0)
  53. output_cnt[blockIdx.x] = scnt[0];
  54. __syncthreads();
  55. }
  56. static __global__ void sum_per_block_cnt(unsigned long *output_cnt, unsigned long *cnt)
  57. {
  58. __shared__ unsigned long accumulator[MAXNBLOCKS];
  59. unsigned i;
  60. /* Load the values from global mem */
  61. for (i = 0; i < blockDim.x; i++)
  62. accumulator[i] = output_cnt[i];
  63. __syncthreads();
  64. /* Perform a reduction in shared memory */
  65. unsigned s;
  66. for (s = blockDim.x/2; s!=0; s>>=1)
  67. {
  68. if (threadIdx.x < s)
  69. accumulator[threadIdx.x] += accumulator[threadIdx.x + s];
  70. __syncthreads();
  71. }
  72. /* Save the result in global memory */
  73. if (threadIdx.x == 0)
  74. *cnt = *cnt + accumulator[0];
  75. }
  76. extern "C" void pi_redux_cuda_kernel(float *x, float *y, unsigned n, unsigned long *shot_cnt)
  77. {
  78. cudaError_t cures;
  79. /* How many blocks do we use ? */
  80. unsigned nblocks = 128; // TODO
  81. STARPU_ASSERT(nblocks <= MAXNBLOCKS);
  82. STARPU_ASSERT((n % nblocks) == 0);
  83. unsigned long *per_block_cnt;
  84. cudaMalloc((void **)&per_block_cnt, nblocks*sizeof(unsigned long));
  85. /* How many threads per block ? At most 256, but no more threads than
  86. * there are entries to process per block. */
  87. unsigned nthread_per_block = STARPU_MIN(MAXTHREADSPERBLOCK, (n / nblocks));
  88. /* each entry of per_block_cnt contains the number of successful shots
  89. * in the corresponding block. */
  90. monte_carlo<<<nblocks, nthread_per_block, 0, starpu_cuda_get_local_stream()>>>(x, y, n, per_block_cnt);
  91. /* Note that we do not synchronize between kernel calls because there is an implicit serialization */
  92. /* compute the total number of successful shots by adding the elements
  93. * of the per_block_cnt array */
  94. sum_per_block_cnt<<<1, nblocks, 0, starpu_cuda_get_local_stream()>>>(per_block_cnt, shot_cnt);
  95. cures = cudaStreamSynchronize(starpu_cuda_get_local_stream());
  96. if (cures)
  97. STARPU_CUDA_REPORT_ERROR(cures);
  98. cudaFree(per_block_cnt);
  99. }