cudax_kernels.cu 4.8 KB

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