pi_kernel.cu 4.3 KB

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