cudax_kernels.cu 4.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144
  1. /*
  2. * StarPU
  3. * Copyright (C) Université Bordeaux 1, CNRS 2009 (see AUTHORS file)
  4. *
  5. * This program 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. * This program 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. dim3 dimGrid(n); \
  27. func <<<dimGrid, 1>>> args; \
  28. } else { \
  29. dim3 dimGrid(n / threads_per_block); \
  30. dim3 dimBlock(threads_per_block); \
  31. func <<<dimGrid, dimBlock>>> args; \
  32. } \
  33. cudaThreadSynchronize(); \
  34. extern "C" __global__ void
  35. STARPUFFT(cuda_twist1_1d)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned n1, unsigned n2)
  36. {
  37. unsigned j;
  38. VARS_1d
  39. unsigned end = n2;
  40. for (j = start; j < end; j += numthreads)
  41. twisted1[j] = in[i+j*n1];
  42. }
  43. extern "C" void
  44. STARPUFFT(cuda_twist1_1d_host)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned n1, unsigned n2)
  45. {
  46. DISTRIB_1d(n2, STARPUFFT(cuda_twist1_1d), (in, twisted1, i, n1, n2));
  47. }
  48. extern "C" __global__ void
  49. 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
  59. 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. if (m < threads_per_dim) { \
  73. dim3 dimGrid(n, m); \
  74. func <<<dimGrid, 1>>> args; \
  75. } else { \
  76. dim3 dimGrid(1, m / threads_per_dim); \
  77. dim3 dimBlock(n, threads_per_dim); \
  78. func <<<dimGrid, dimBlock>>> args; \
  79. } \
  80. } else { \
  81. if (m < threads_per_dim) { \
  82. dim3 dimGrid(n / threads_per_dim, 1); \
  83. dim3 dimBlock(threads_per_dim, m); \
  84. func <<<dimGrid, dimBlock>>> args; \
  85. } else { \
  86. dim3 dimGrid(n / threads_per_dim, m / threads_per_dim); \
  87. dim3 dimBlock(threads_per_dim, threads_per_dim); \
  88. func <<<dimGrid, dimBlock>>> args; \
  89. } \
  90. } \
  91. cudaThreadSynchronize(); \
  92. extern "C" __global__ void
  93. STARPUFFT(cuda_twist1_2d)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned j, unsigned n1, unsigned n2, unsigned m1, unsigned m2)
  94. {
  95. unsigned k, l;
  96. VARS_2d
  97. unsigned endx = n2;
  98. unsigned endy = m2;
  99. unsigned m = m1*m2;
  100. for (k = startx; k < endx; k += numthreadsx)
  101. for (l = starty; l < endy; l += numthreadsy)
  102. twisted1[k*m2+l] = in[i*m+j+k*m*n1+l*m1];
  103. }
  104. extern "C" void
  105. STARPUFFT(cuda_twist1_2d_host)(const _cuComplex *in, _cuComplex *twisted1, unsigned i, unsigned j, unsigned n1, unsigned n2, unsigned m1, unsigned m2)
  106. {
  107. DISTRIB_2d(n2, m2, STARPUFFT(cuda_twist1_2d), (in, twisted1, i, j, n1, n2, m1, m2));
  108. }
  109. extern "C" __global__ void
  110. STARPUFFT(cuda_twiddle_2d)(_cuComplex * out, const _cuComplex * roots0, const _cuComplex * roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
  111. {
  112. unsigned k, l;
  113. VARS_2d
  114. unsigned endx = n2;
  115. unsigned endy = m2;
  116. for (k = startx; k < endx ; k += numthreadsx)
  117. for (l = starty; l < endy ; l += numthreadsy)
  118. out[k*m2 + l] = _cuCmul(_cuCmul(out[k*m2 + l], roots0[i*k]), roots1[j*l]);
  119. return;
  120. }
  121. extern "C" void
  122. STARPUFFT(cuda_twiddle_2d_host)(_cuComplex *out, const _cuComplex *roots0, const _cuComplex *roots1, unsigned n2, unsigned m2, unsigned i, unsigned j)
  123. {
  124. DISTRIB_2d(n2, m2, STARPUFFT(cuda_twiddle_2d), (out, roots0, roots1, n2, m2, i, j));
  125. }