pi_redux_kernel.cu 3.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010, 2015 Université de Bordeaux
  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. /* This counts how many fall inside the circle quarter */
  17. #include <starpu.h>
  18. #define MAXNBLOCKS 128
  19. #define MAXTHREADSPERBLOCK 256
  20. static __global__ void monte_carlo(float *x, float *y, unsigned n, unsigned long *output_cnt)
  21. {
  22. __shared__ unsigned scnt[MAXTHREADSPERBLOCK];
  23. /* Do we have a successful shot ? */
  24. const int tid = threadIdx.x + blockIdx.x*blockDim.x;
  25. const int nthreads = gridDim.x * blockDim.x;
  26. /* Blank the shared mem buffer */
  27. if (threadIdx.x < MAXTHREADSPERBLOCK)
  28. scnt[threadIdx.x] = 0;
  29. __syncthreads();
  30. int ind;
  31. for (ind = tid; ind < n; ind += nthreads)
  32. {
  33. float xval = (2.0f * x[ind] - 1.0f);
  34. float yval = (2.0f * y[ind] - 1.0f);
  35. float dist = (xval*xval + yval*yval);
  36. unsigned long success = (dist <= 1.0f)?1:0;
  37. scnt[threadIdx.x] += success;
  38. }
  39. __syncthreads();
  40. /* Perform a reduction to compute the sum on each thread within that block */
  41. /* NB: We assume that the number of threads per block is a power of 2 ! */
  42. unsigned long s;
  43. for (s = blockDim.x/2; s!=0; s>>=1)
  44. {
  45. if (threadIdx.x < s)
  46. scnt[threadIdx.x] += scnt[threadIdx.x + s];
  47. __syncthreads();
  48. }
  49. /* report the number of successful shots in the block */
  50. if (threadIdx.x == 0)
  51. output_cnt[blockIdx.x] = scnt[0];
  52. __syncthreads();
  53. }
  54. static __global__ void sum_per_block_cnt(unsigned long *output_cnt, unsigned long *cnt)
  55. {
  56. __shared__ unsigned long accumulator[MAXNBLOCKS];
  57. unsigned i;
  58. /* Load the values from global mem */
  59. for (i = 0; i < blockDim.x; i++)
  60. accumulator[i] = output_cnt[i];
  61. __syncthreads();
  62. /* Perform a reduction in shared memory */
  63. unsigned s;
  64. for (s = blockDim.x/2; s!=0; s>>=1)
  65. {
  66. if (threadIdx.x < s)
  67. accumulator[threadIdx.x] += accumulator[threadIdx.x + s];
  68. __syncthreads();
  69. }
  70. /* Save the result in global memory */
  71. if (threadIdx.x == 0)
  72. *cnt = *cnt + accumulator[0];
  73. }
  74. extern "C" void pi_redux_cuda_kernel(float *x, float *y, unsigned n, unsigned long *shot_cnt)
  75. {
  76. cudaError_t cures;
  77. /* How many blocks do we use ? */
  78. unsigned nblocks = 128; // TODO
  79. STARPU_ASSERT(nblocks <= MAXNBLOCKS);
  80. STARPU_ASSERT((n % nblocks) == 0);
  81. unsigned long *per_block_cnt;
  82. cudaMalloc((void **)&per_block_cnt, nblocks*sizeof(unsigned long));
  83. /* How many threads per block ? At most 256, but no more threads than
  84. * there are entries to process per block. */
  85. unsigned nthread_per_block = STARPU_MIN(MAXTHREADSPERBLOCK, (n / nblocks));
  86. /* each entry of per_block_cnt contains the number of successful shots
  87. * in the corresponding block. */
  88. monte_carlo<<<nblocks, nthread_per_block, 0, starpu_cuda_get_local_stream()>>>(x, y, n, per_block_cnt);
  89. /* Note that we do not synchronize between kernel calls because there is an implicit serialization */
  90. /* compute the total number of successful shots by adding the elements
  91. * of the per_block_cnt array */
  92. sum_per_block_cnt<<<1, nblocks, 0, starpu_cuda_get_local_stream()>>>(per_block_cnt, shot_cnt);
  93. cures = cudaStreamSynchronize(starpu_cuda_get_local_stream());
  94. if (cures)
  95. STARPU_CUDA_REPORT_ERROR(cures);
  96. cudaFree(per_block_cnt);
  97. }