pi_kernel.cu 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149
  1. /*
  2. * StarPU
  3. * Copyright (C) Université Bordeaux 1, CNRS 2008-2010 (see AUTHORS file)
  4. *
  5. * This program 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. * This program 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 "SobolQRNG/sobol_gpu.h"
  17. #include "pi.h"
  18. #define MAXNBLOCKS 128
  19. #define MAXTHREADSPERBLOCK 256
  20. static __global__ void monte_carlo(TYPE *random_numbers_x, TYPE *random_numbers_y,
  21. unsigned n, unsigned *output_cnt)
  22. {
  23. __shared__ unsigned scnt[MAXTHREADSPERBLOCK];
  24. /* Do we have a successful shot ? */
  25. const int tid = threadIdx.x + blockIdx.x*blockDim.x;
  26. const int nthreads = gridDim.x * blockDim.x;
  27. /* Blank the shared mem buffer */
  28. if (threadIdx.x < MAXTHREADSPERBLOCK)
  29. scnt[threadIdx.x] = 0;
  30. __syncthreads();
  31. int ind;
  32. for (ind = tid; ind < n; ind += nthreads)
  33. {
  34. TYPE x = random_numbers_x[ind];
  35. TYPE y = random_numbers_y[ind];
  36. TYPE dist = (x*x + y*y);
  37. unsigned success = (dist <= 1.0f)?1:0;
  38. scnt[threadIdx.x] += success;
  39. }
  40. __syncthreads();
  41. /* Perform a reduction to compute the sum on each thread within that block */
  42. /* NB: We assume that the number of threads per block is a power of 2 ! */
  43. unsigned s;
  44. for (s = blockDim.x/2; s!=0; s>>=1)
  45. {
  46. if (threadIdx.x < s)
  47. scnt[threadIdx.x] += scnt[threadIdx.x + s];
  48. __syncthreads();
  49. }
  50. /* report the number of successful shots in the block */
  51. if (threadIdx.x == 0)
  52. output_cnt[blockIdx.x] = scnt[0];
  53. __syncthreads();
  54. }
  55. static __global__ void sum_per_block_cnt(unsigned *output_cnt, unsigned *cnt)
  56. {
  57. __shared__ unsigned accumulator[MAXNBLOCKS];
  58. unsigned i;
  59. /* Load the values from global mem */
  60. for (i = 0; i < blockDim.x; i++)
  61. accumulator[i] = output_cnt[i];
  62. __syncthreads();
  63. /* Perform a reduction in shared memory */
  64. unsigned s;
  65. for (s = blockDim.x/2; s!=0; s>>=1)
  66. {
  67. if (threadIdx.x < s)
  68. accumulator[threadIdx.x] += accumulator[threadIdx.x + s];
  69. __syncthreads();
  70. }
  71. /* Save the result in global memory */
  72. if (threadIdx.x == 0)
  73. *cnt = accumulator[0];
  74. }
  75. extern "C" void cuda_kernel(void *descr[], void *cl_arg)
  76. {
  77. cudaError_t cures;
  78. unsigned *directions = (unsigned *)STARPU_VECTOR_GET_PTR(descr[0]);
  79. unsigned nx = NSHOT_PER_TASK;
  80. /* Generate Random numbers */
  81. float *random_numbers;
  82. cudaMalloc((void **)&random_numbers, 2*nx*sizeof(float));
  83. STARPU_ASSERT(random_numbers);
  84. sobolGPU(2*nx/n_dimensions, n_dimensions, directions, random_numbers);
  85. cudaThreadSynchronize();
  86. TYPE *random_numbers_x = &random_numbers[0];
  87. TYPE *random_numbers_y = &random_numbers[nx];
  88. unsigned *cnt = (unsigned *)STARPU_VECTOR_GET_PTR(descr[1]);
  89. /* How many blocks do we use ? */
  90. unsigned nblocks = 128; // TODO
  91. STARPU_ASSERT(nblocks <= MAXNBLOCKS);
  92. unsigned *per_block_cnt;
  93. cudaMalloc((void **)&per_block_cnt, nblocks*sizeof(unsigned));
  94. STARPU_ASSERT((nx % nblocks) == 0);
  95. /* How many threads per block ? At most 256, but no more threads than
  96. * there are entries to process per block. */
  97. unsigned nthread_per_block = STARPU_MIN(MAXTHREADSPERBLOCK, (nx / nblocks));
  98. /* each entry of per_block_cnt contains the number of successful shots
  99. * in the corresponding block. */
  100. monte_carlo<<<nblocks, nthread_per_block>>>(random_numbers_x, random_numbers_y, nx, per_block_cnt);
  101. /* Note that we do not synchronize between kernel calls because there is an implicit serialization */
  102. /* compute the total number of successful shots by adding the elements
  103. * of the per_block_cnt array */
  104. sum_per_block_cnt<<<1, nblocks>>>(per_block_cnt, cnt);
  105. cures = cudaThreadSynchronize();
  106. if (cures)
  107. STARPU_CUDA_REPORT_ERROR(cures);
  108. cudaFree(per_block_cnt);
  109. cudaFree(random_numbers);
  110. }