cudax_kernels.cu 4.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009, 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. #define _externC extern "C"
  18. #include "cudax_kernels.h"
  19. /* Note: these assume that the sizes are powers of two */
  20. #define VARS_1d \
  21. unsigned start = threadIdx.x + blockIdx.x * blockDim.x; \
  22. unsigned numthreads = blockDim.x * gridDim.x;
  23. #define DISTRIB_1d(n, func,args) \
  24. unsigned threads_per_block = 128; \
  25. \
  26. if (n < threads_per_block) { \
  27. dim3 dimGrid(n); \
  28. func <<<dimGrid, 1, 0, starpu_cuda_get_local_stream()>>> args; \
  29. } else { \
  30. dim3 dimGrid(n / threads_per_block); \
  31. dim3 dimBlock(threads_per_block); \
  32. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  33. } \
  34. cudaStreamSynchronize(starpu_cuda_get_local_stream()); \
  35. extern "C" __global__ void
  36. STARPUFFT(cuda_twist1_1d)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned n1, unsigned n2)
  37. {
  38. unsigned j;
  39. VARS_1d
  40. unsigned end = n2;
  41. for (j = start; j < end; j += numthreads)
  42. twisted1[j] = in[i+j*n1];
  43. }
  44. extern "C" void
  45. STARPUFFT(cuda_twist1_1d_host)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned n1, unsigned n2)
  46. {
  47. DISTRIB_1d(n2, STARPUFFT(cuda_twist1_1d), (in, twisted1, i, n1, n2));
  48. }
  49. extern "C" __global__ void
  50. STARPUFFT(cuda_twiddle_1d)(_cuComplex * out, const _cuComplex * roots, unsigned n, unsigned i)
  51. {
  52. unsigned j;
  53. VARS_1d
  54. unsigned end = n;
  55. for (j = start; j < end; j += numthreads)
  56. out[j] = _cuCmul(out[j], roots[i*j]);
  57. return;
  58. }
  59. extern "C" void
  60. STARPUFFT(cuda_twiddle_1d_host)(_cuComplex *out, const _cuComplex *roots, unsigned n, unsigned i)
  61. {
  62. DISTRIB_1d(n, STARPUFFT(cuda_twiddle_1d), (out, roots, n, i));
  63. }
  64. #define VARS_2d \
  65. unsigned startx = threadIdx.x + blockIdx.x * blockDim.x; \
  66. unsigned starty = threadIdx.y + blockIdx.y * blockDim.y; \
  67. unsigned numthreadsx = blockDim.x * gridDim.x; \
  68. unsigned numthreadsy = blockDim.y * gridDim.y;
  69. /* FIXME: introduce threads_per_dim_n / m instead */
  70. #define DISTRIB_2d(n, m, func, args) \
  71. unsigned threads_per_dim = 16; \
  72. if (n < threads_per_dim) { \
  73. if (m < threads_per_dim) { \
  74. dim3 dimGrid(n, m); \
  75. func <<<dimGrid, 1, 0, starpu_cuda_get_local_stream()>>> args; \
  76. } else { \
  77. dim3 dimGrid(1, m / threads_per_dim); \
  78. dim3 dimBlock(n, threads_per_dim); \
  79. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  80. } \
  81. } else { \
  82. if (m < threads_per_dim) { \
  83. dim3 dimGrid(n / threads_per_dim, 1); \
  84. dim3 dimBlock(threads_per_dim, m); \
  85. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  86. } else { \
  87. dim3 dimGrid(n / threads_per_dim, m / threads_per_dim); \
  88. dim3 dimBlock(threads_per_dim, threads_per_dim); \
  89. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  90. } \
  91. } \
  92. cudaStreamSynchronize(starpu_cuda_get_local_stream()); \
  93. extern "C" __global__ void
  94. STARPUFFT(cuda_twist1_2d)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned j, unsigned n1, unsigned n2, unsigned m1, unsigned m2)
  95. {
  96. unsigned k, l;
  97. VARS_2d
  98. unsigned endx = n2;
  99. unsigned endy = m2;
  100. unsigned m = m1*m2;
  101. for (k = startx; k < endx; k += numthreadsx)
  102. for (l = starty; l < endy; l += numthreadsy)
  103. twisted1[k*m2+l] = in[i*m+j+k*m*n1+l*m1];
  104. }
  105. extern "C" void
  106. STARPUFFT(cuda_twist1_2d_host)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned j, unsigned n1, unsigned n2, unsigned m1, unsigned m2)
  107. {
  108. DISTRIB_2d(n2, m2, STARPUFFT(cuda_twist1_2d), (in, twisted1, i, j, n1, n2, m1, m2));
  109. }
  110. extern "C" __global__ void
  111. STARPUFFT(cuda_twiddle_2d)(_cuComplex * out, const _cuComplex * roots0, const _cuComplex * roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
  112. {
  113. unsigned k, l;
  114. VARS_2d
  115. unsigned endx = n2;
  116. unsigned endy = m2;
  117. for (k = startx; k < endx ; k += numthreadsx)
  118. for (l = starty; l < endy ; l += numthreadsy)
  119. out[k*m2 + l] = _cuCmul(_cuCmul(out[k*m2 + l], roots0[i*k]), roots1[j*l]);
  120. return;
  121. }
  122. extern "C" void
  123. STARPUFFT(cuda_twiddle_2d_host)(_cuComplex *out, const _cuComplex *roots0, const _cuComplex *roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
  124. {
  125. DISTRIB_2d(n2, m2, STARPUFFT(cuda_twiddle_2d), (out, roots0, roots1, n2, m2, i, j));
  126. }