cudax_kernels.cu 4.8 KB

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