cudax_kernels.cu 4.8 KB

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