cudax_kernels.cu 4.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010-2012,2014,2015,2017 CNRS
  4. * Copyright (C) 2009,2010,2014,2016 Université de Bordeaux
  5. * Copyright (C) 2012 Inria
  6. *
  7. * StarPU is free software; you can redistribute it and/or modify
  8. * it under the terms of the GNU Lesser General Public License as published by
  9. * the Free Software Foundation; either version 2.1 of the License, or (at
  10. * your option) any later version.
  11. *
  12. * StarPU is distributed in the hope that it will be useful, but
  13. * WITHOUT ANY WARRANTY; without even the implied warranty of
  14. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  15. *
  16. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  17. */
  18. #define _externC extern "C"
  19. #include "cudax_kernels.h"
  20. /* Note: these assume that the sizes are powers of two */
  21. #define VARS_1d \
  22. unsigned start = threadIdx.x + blockIdx.x * blockDim.x; \
  23. unsigned numthreads = blockDim.x * gridDim.x;
  24. #define DISTRIB_1d(n, func,args) \
  25. unsigned threads_per_block = 128; \
  26. \
  27. if (n < threads_per_block) \
  28. { \
  29. dim3 dimGrid(n); \
  30. func <<<dimGrid, 1, 0, starpu_cuda_get_local_stream()>>> args; \
  31. } \
  32. else \
  33. { \
  34. dim3 dimGrid(n / threads_per_block); \
  35. dim3 dimBlock(threads_per_block); \
  36. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  37. } \
  38. cudaStreamSynchronize(starpu_cuda_get_local_stream()); \
  39. extern "C" __global__ void 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 STARPUFFT(cuda_twist1_1d_host)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned n1, unsigned n2)
  48. {
  49. DISTRIB_1d(n2, STARPUFFT(cuda_twist1_1d), (in, twisted1, i, n1, n2));
  50. }
  51. extern "C" __global__ void STARPUFFT(cuda_twiddle_1d)(_cuComplex * out, const _cuComplex * roots, unsigned n, unsigned i)
  52. {
  53. unsigned j;
  54. VARS_1d
  55. unsigned end = n;
  56. for (j = start; j < end; j += numthreads)
  57. out[j] = _cuCmul(out[j], roots[i*j]);
  58. return;
  59. }
  60. extern "C" void 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. { \
  74. if (m < threads_per_dim) \
  75. { \
  76. dim3 dimGrid(n, m); \
  77. func <<<dimGrid, 1, 0, starpu_cuda_get_local_stream()>>> args; \
  78. } \
  79. else \
  80. { \
  81. dim3 dimGrid(1, m / threads_per_dim); \
  82. dim3 dimBlock(n, threads_per_dim); \
  83. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  84. } \
  85. } \
  86. else \
  87. { \
  88. if (m < threads_per_dim) \
  89. { \
  90. dim3 dimGrid(n / threads_per_dim, 1); \
  91. dim3 dimBlock(threads_per_dim, m); \
  92. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  93. } \
  94. else \
  95. { \
  96. dim3 dimGrid(n / threads_per_dim, m / threads_per_dim); \
  97. dim3 dimBlock(threads_per_dim, threads_per_dim); \
  98. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  99. } \
  100. } \
  101. cudaStreamSynchronize(starpu_cuda_get_local_stream()); \
  102. extern "C" __global__ void STARPUFFT(cuda_twist1_2d)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned j, unsigned n1, unsigned n2, unsigned m1, unsigned m2)
  103. {
  104. unsigned k, l;
  105. VARS_2d
  106. unsigned endx = n2;
  107. unsigned endy = m2;
  108. unsigned m = m1*m2;
  109. for (k = startx; k < endx; k += numthreadsx)
  110. for (l = starty; l < endy; l += numthreadsy)
  111. twisted1[k*m2+l] = in[i*m+j+k*m*n1+l*m1];
  112. }
  113. extern "C" void STARPUFFT(cuda_twist1_2d_host)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned j, unsigned n1, unsigned n2, unsigned m1, unsigned m2)
  114. {
  115. DISTRIB_2d(n2, m2, STARPUFFT(cuda_twist1_2d), (in, twisted1, i, j, n1, n2, m1, m2));
  116. }
  117. extern "C" __global__ void STARPUFFT(cuda_twiddle_2d)(_cuComplex * out, const _cuComplex * roots0, const _cuComplex * roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
  118. {
  119. unsigned k, l;
  120. VARS_2d
  121. unsigned endx = n2;
  122. unsigned endy = m2;
  123. for (k = startx; k < endx ; k += numthreadsx)
  124. for (l = starty; l < endy ; l += numthreadsy)
  125. out[k*m2 + l] = _cuCmul(_cuCmul(out[k*m2 + l], roots0[i*k]), roots1[j*l]);
  126. return;
  127. }
  128. extern "C" void STARPUFFT(cuda_twiddle_2d_host)(_cuComplex *out, const _cuComplex *roots0, const _cuComplex *roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
  129. {
  130. DISTRIB_2d(n2, m2, STARPUFFT(cuda_twiddle_2d), (out, roots0, roots1, n2, m2, i, j));
  131. }