cg_dot_kernel.cu 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010 Université de Bordeaux 1
  4. * Copyright (C) 2010, 2012 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 "cg.h"
  19. #define MAXNBLOCKS 128
  20. #define MAXTHREADSPERBLOCK 256
  21. static __global__ void dot_device(TYPE *vx, TYPE *vy, unsigned n, TYPE *dot_array)
  22. {
  23. __shared__ TYPE 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] = (TYPE)0.0;
  30. __syncthreads();
  31. int ind;
  32. for (ind = tid; ind < n; ind += nthreads)
  33. {
  34. TYPE x = vx[ind];
  35. TYPE y = vy[ind];
  36. scnt[threadIdx.x] += (x*y);
  37. }
  38. __syncthreads();
  39. /* Perform a reduction to compute the sum on each thread within that block */
  40. /* NB: We assume that the number of threads per block is a power of 2 ! */
  41. unsigned s;
  42. for (s = blockDim.x/2; s!=0; s>>=1)
  43. {
  44. if (threadIdx.x < s)
  45. scnt[threadIdx.x] += scnt[threadIdx.x + s];
  46. __syncthreads();
  47. }
  48. /* report the number of successful shots in the block */
  49. if (threadIdx.x == 0)
  50. dot_array[blockIdx.x] = scnt[0];
  51. __syncthreads();
  52. }
  53. static __global__ void gather_dot_device(TYPE *dot_array, TYPE *dot)
  54. {
  55. __shared__ TYPE accumulator[MAXNBLOCKS];
  56. unsigned i;
  57. /* Load the values from global mem */
  58. for (i = 0; i < blockDim.x; i++)
  59. accumulator[i] = dot_array[i];
  60. __syncthreads();
  61. /* Perform a reduction in shared memory */
  62. unsigned s;
  63. for (s = blockDim.x/2; s!=0; s>>=1)
  64. {
  65. if (threadIdx.x < s)
  66. accumulator[threadIdx.x] += accumulator[threadIdx.x + s];
  67. __syncthreads();
  68. }
  69. /* Save the result in global memory */
  70. if (threadIdx.x == 0)
  71. *dot = *dot + accumulator[0];
  72. }
  73. extern "C" void dot_host(TYPE *x, TYPE *y, unsigned nelems, TYPE *dot)
  74. {
  75. /* How many blocks do we use ? */
  76. unsigned nblocks = 128; // TODO
  77. STARPU_ASSERT(nblocks <= MAXNBLOCKS);
  78. TYPE *per_block_sum;
  79. cudaMalloc((void **)&per_block_sum, nblocks*sizeof(TYPE));
  80. STARPU_ASSERT((nelems % nblocks) == 0);
  81. /* How many threads per block ? At most 256, but no more threads than
  82. * there are entries to process per block. */
  83. unsigned nthread_per_block = STARPU_MIN(MAXTHREADSPERBLOCK, (nelems / nblocks));
  84. /* each entry of per_block_sum contains the number of successful shots
  85. * in the corresponding block. */
  86. dot_device<<<nblocks, nthread_per_block, 0, starpu_cuda_get_local_stream()>>>(x, y, nelems, per_block_sum);
  87. /* Note that we do not synchronize between kernel calls because there
  88. * is an implicit serialization */
  89. gather_dot_device<<<1, nblocks, 0, starpu_cuda_get_local_stream()>>>(per_block_sum, dot);
  90. cudaError_t cures;
  91. cures = cudaStreamSynchronize(starpu_cuda_get_local_stream());
  92. if (cures)
  93. STARPU_CUDA_REPORT_ERROR(cures);
  94. cudaFree(per_block_sum);
  95. }
  96. static __global__ void zero_vector_device(TYPE *x, unsigned nelems, unsigned nelems_per_thread)
  97. {
  98. unsigned i;
  99. unsigned first_i = blockDim.x * blockIdx.x + threadIdx.x;
  100. for (i = first_i; i < nelems; i += nelems_per_thread)
  101. x[i] = 0.0;
  102. }
  103. extern "C" void zero_vector(TYPE *x, unsigned nelems)
  104. {
  105. unsigned nblocks = STARPU_MIN(128, nelems);
  106. unsigned nthread_per_block = STARPU_MIN(MAXTHREADSPERBLOCK, (nelems / nblocks));
  107. unsigned nelems_per_thread = nelems / (nblocks * nthread_per_block);
  108. zero_vector_device<<<nblocks, nthread_per_block, 0, starpu_cuda_get_local_stream()>>>(x, nelems, nelems_per_thread);
  109. }