pi_redux_kernel.cu 3.6 KB

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