cg_dot_kernel.cu 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010 Université de Bordeaux 1
  4. * Copyright (C) 2010 Centre National de la Recherche Scientifique
  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. #include <starpu.h>
  18. #include <starpu_cuda.h>
  19. #include "cg.h"
  20. #define MAXNBLOCKS 128
  21. #define MAXTHREADSPERBLOCK 256
  22. static __global__ void dot_device(TYPE *vx, TYPE *vy, unsigned n, TYPE *dot_array)
  23. {
  24. __shared__ TYPE scnt[MAXTHREADSPERBLOCK];
  25. /* Do we have a successful shot ? */
  26. const int tid = threadIdx.x + blockIdx.x*blockDim.x;
  27. const int nthreads = gridDim.x * blockDim.x;
  28. /* Blank the shared mem buffer */
  29. if (threadIdx.x < MAXTHREADSPERBLOCK)
  30. scnt[threadIdx.x] = (TYPE)0.0;
  31. __syncthreads();
  32. int ind;
  33. for (ind = tid; ind < n; ind += nthreads)
  34. {
  35. TYPE x = vx[ind];
  36. TYPE y = vy[ind];
  37. scnt[threadIdx.x] += (x*y);
  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 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. dot_array[blockIdx.x] = scnt[0];
  52. __syncthreads();
  53. }
  54. static __global__ void gather_dot_device(TYPE *dot_array, TYPE *dot)
  55. {
  56. __shared__ TYPE accumulator[MAXNBLOCKS];
  57. unsigned i;
  58. /* Load the values from global mem */
  59. for (i = 0; i < blockDim.x; i++)
  60. accumulator[i] = dot_array[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. *dot = *dot + accumulator[0];
  73. }
  74. extern "C" void dot_host(TYPE *x, TYPE *y, unsigned nelems, TYPE *dot)
  75. {
  76. /* How many blocks do we use ? */
  77. unsigned nblocks = 128; // TODO
  78. STARPU_ASSERT(nblocks <= MAXNBLOCKS);
  79. TYPE *per_block_sum;
  80. cudaMalloc((void **)&per_block_sum, nblocks*sizeof(TYPE));
  81. STARPU_ASSERT((nelems % nblocks) == 0);
  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, (nelems / nblocks));
  85. /* each entry of per_block_sum contains the number of successful shots
  86. * in the corresponding block. */
  87. dot_device<<<nblocks, nthread_per_block, 0, starpu_cuda_get_local_stream()>>>(x, y, nelems, per_block_sum);
  88. /* Note that we do not synchronize between kernel calls because there
  89. * is an implicit serialization */
  90. gather_dot_device<<<1, nblocks, 0, starpu_cuda_get_local_stream()>>>(per_block_sum, dot);
  91. cudaError_t cures;
  92. cures = cudaStreamSynchronize(starpu_cuda_get_local_stream());
  93. if (cures)
  94. STARPU_CUDA_REPORT_ERROR(cures);
  95. cudaFree(per_block_sum);
  96. }
  97. static __global__ void zero_vector_device(TYPE *x, unsigned nelems, unsigned nelems_per_thread)
  98. {
  99. unsigned i;
  100. unsigned first_i = blockDim.x * blockIdx.x + threadIdx.x;
  101. for (i = first_i; i < nelems; i += nelems_per_thread)
  102. x[i] = 0.0;
  103. }
  104. extern "C" void zero_vector(TYPE *x, unsigned nelems)
  105. {
  106. unsigned nblocks = STARPU_MIN(128, nelems);
  107. unsigned nthread_per_block = STARPU_MIN(MAXTHREADSPERBLOCK, (nelems / nblocks));
  108. unsigned nelems_per_thread = nelems / (nblocks * nthread_per_block);
  109. zero_vector_device<<<nblocks, nthread_per_block, 0, starpu_cuda_get_local_stream()>>>(x, nelems, nelems_per_thread);
  110. }