pi_kernel.cu 4.5 KB

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