cudax_kernels.cu 5.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009-2021 Université de Bordeaux, CNRS (LaBRI UMR 5800), Inria
  4. *
  5. * StarPU is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU Lesser General Public License as published by
  7. * the Free Software Foundation; either version 2.1 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * StarPU is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  13. *
  14. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  15. */
  16. #define _externC extern "C"
  17. #include "cudax_kernels.h"
  18. /* Note: these assume that the sizes are powers of two */
  19. #define VARS_1d \
  20. unsigned start = threadIdx.x + blockIdx.x * blockDim.x; \
  21. unsigned numthreads = blockDim.x * gridDim.x;
  22. #define DISTRIB_1d(n, func,args) \
  23. unsigned threads_per_block = 128; \
  24. \
  25. if (n < threads_per_block) \
  26. { \
  27. dim3 dimGrid(n); \
  28. func <<<dimGrid, 1, 0, starpu_cuda_get_local_stream()>>> args; \
  29. cudaError_t status = cudaGetLastError(); \
  30. if (status != cudaSuccess) STARPU_CUDA_REPORT_ERROR(status); \
  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. cudaError_t status = cudaGetLastError(); \
  38. if (status != cudaSuccess) STARPU_CUDA_REPORT_ERROR(status); \
  39. } \
  40. cudaStreamSynchronize(starpu_cuda_get_local_stream()); \
  41. extern "C" __global__ void STARPUFFT(cuda_twist1_1d)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned n1, unsigned n2)
  42. {
  43. unsigned j;
  44. VARS_1d
  45. unsigned end = n2;
  46. for (j = start; j < end; j += numthreads)
  47. twisted1[j] = in[i+j*n1];
  48. }
  49. extern "C" void STARPUFFT(cuda_twist1_1d_host)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned n1, unsigned n2)
  50. {
  51. DISTRIB_1d(n2, STARPUFFT(cuda_twist1_1d), (in, twisted1, i, n1, n2));
  52. }
  53. extern "C" __global__ void 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 STARPUFFT(cuda_twiddle_1d_host)(_cuComplex *out, const _cuComplex *roots, unsigned n, unsigned i)
  63. {
  64. DISTRIB_1d(n, STARPUFFT(cuda_twiddle_1d), (out, roots, n, i));
  65. }
  66. #define VARS_2d \
  67. unsigned startx = threadIdx.x + blockIdx.x * blockDim.x; \
  68. unsigned starty = threadIdx.y + blockIdx.y * blockDim.y; \
  69. unsigned numthreadsx = blockDim.x * gridDim.x; \
  70. unsigned numthreadsy = blockDim.y * gridDim.y;
  71. /* FIXME: introduce threads_per_dim_n / m instead */
  72. #define DISTRIB_2d(n, m, func, args) \
  73. unsigned threads_per_dim = 16; \
  74. if (n < threads_per_dim) \
  75. { \
  76. if (m < threads_per_dim) \
  77. { \
  78. dim3 dimGrid(n, m); \
  79. func <<<dimGrid, 1, 0, starpu_cuda_get_local_stream()>>> args; \
  80. cudaError_t status = cudaGetLastError(); \
  81. if (status != cudaSuccess) STARPU_CUDA_REPORT_ERROR(status); \
  82. } \
  83. else \
  84. { \
  85. dim3 dimGrid(1, m / threads_per_dim); \
  86. dim3 dimBlock(n, threads_per_dim); \
  87. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  88. cudaError_t status = cudaGetLastError(); \
  89. if (status != cudaSuccess) STARPU_CUDA_REPORT_ERROR(status); \
  90. } \
  91. } \
  92. else \
  93. { \
  94. if (m < threads_per_dim) \
  95. { \
  96. dim3 dimGrid(n / threads_per_dim, 1); \
  97. dim3 dimBlock(threads_per_dim, m); \
  98. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  99. cudaError_t status = cudaGetLastError(); \
  100. if (status != cudaSuccess) STARPU_CUDA_REPORT_ERROR(status); \
  101. } \
  102. else \
  103. { \
  104. dim3 dimGrid(n / threads_per_dim, m / threads_per_dim); \
  105. dim3 dimBlock(threads_per_dim, threads_per_dim); \
  106. func <<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>> args; \
  107. cudaError_t status = cudaGetLastError(); \
  108. if (status != cudaSuccess) STARPU_CUDA_REPORT_ERROR(status); \
  109. } \
  110. } \
  111. cudaStreamSynchronize(starpu_cuda_get_local_stream()); \
  112. 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)
  113. {
  114. unsigned k, l;
  115. VARS_2d
  116. unsigned endx = n2;
  117. unsigned endy = m2;
  118. unsigned m = m1*m2;
  119. for (k = startx; k < endx; k += numthreadsx)
  120. for (l = starty; l < endy; l += numthreadsy)
  121. twisted1[k*m2+l] = in[i*m+j+k*m*n1+l*m1];
  122. }
  123. 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)
  124. {
  125. DISTRIB_2d(n2, m2, STARPUFFT(cuda_twist1_2d), (in, twisted1, i, j, n1, n2, m1, m2));
  126. }
  127. extern "C" __global__ void STARPUFFT(cuda_twiddle_2d)(_cuComplex * out, const _cuComplex * roots0, const _cuComplex * roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
  128. {
  129. unsigned k, l;
  130. VARS_2d
  131. unsigned endx = n2;
  132. unsigned endy = m2;
  133. for (k = startx; k < endx ; k += numthreadsx)
  134. for (l = starty; l < endy ; l += numthreadsy)
  135. out[k*m2 + l] = _cuCmul(_cuCmul(out[k*m2 + l], roots0[i*k]), roots1[j*l]);
  136. return;
  137. }
  138. extern "C" void STARPUFFT(cuda_twiddle_2d_host)(_cuComplex *out, const _cuComplex *roots0, const _cuComplex *roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
  139. {
  140. DISTRIB_2d(n2, m2, STARPUFFT(cuda_twiddle_2d), (out, roots0, roots1, n2, m2, i, j));
  141. }