spmv_cuda.cu 3.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2008-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. /* CUDA kernel for SPMV */
  17. #include <starpu.h>
  18. #define MIN(a,b) ((a)<(b)?(a):(b))
  19. extern "C" __global__ void spmv_kernel(uint32_t nnz, uint32_t nrow, float *nzval, uint32_t *colind, uint32_t *rowptr,
  20. uint32_t firstentry, uint32_t elemsize,
  21. float *vecin, uint32_t nx_in, uint32_t elemsize1, float * vecout, uint32_t nx_out, uint32_t elemsize2)
  22. {
  23. /* only one dimension is used here */
  24. unsigned nthreads = gridDim.x*blockDim.x;
  25. unsigned threadid = threadIdx.x + blockIdx.x*blockDim.x;
  26. unsigned rowstart = threadid * ((nrow + (nthreads - 1))/nthreads);
  27. unsigned rowend = MIN(nrow, (threadid+1) * ((nrow + (nthreads - 1))/nthreads));
  28. unsigned row;
  29. for (row = rowstart; row < rowend; row++)
  30. {
  31. float tmp = 0.0f;
  32. unsigned index;
  33. unsigned firstindex = rowptr[row] - firstentry;
  34. unsigned lastindex = rowptr[row+1] - firstentry;
  35. for (index = firstindex; index < lastindex; index++)
  36. {
  37. tmp += nzval[index]*vecin[colind[index]];
  38. }
  39. vecout[row] = tmp;
  40. }
  41. }
  42. extern "C" __global__ void spmv_kernel_3(uint32_t nnz, uint32_t nrow, float *nzval, uint32_t *colind, uint32_t *rowptr,
  43. uint32_t firstentry,
  44. float *vecin, uint32_t nx_in, float * vecout, uint32_t nx_out)
  45. {
  46. /* only one dimension is used here */
  47. unsigned block_rowstart = blockIdx.x*( (nrow + gridDim.x - 1)/gridDim.x );
  48. unsigned block_rowend = MIN((blockIdx.x+1)*( (nrow + gridDim.x - 1)/gridDim.x ), nrow);
  49. unsigned row;
  50. for (row = block_rowstart + threadIdx.x; row < block_rowend; row+=blockDim.x)
  51. {
  52. float tmp = 0.0f;
  53. unsigned index;
  54. unsigned firstindex = rowptr[row] - firstentry;
  55. unsigned lastindex = rowptr[row+1] - firstentry;
  56. for (index = firstindex; index < lastindex; index++)
  57. {
  58. tmp += nzval[index]*vecin[colind[index]];
  59. }
  60. vecout[row] = tmp;
  61. }
  62. }
  63. extern "C" void spmv_kernel_cuda(void *descr[], void *args)
  64. {
  65. uint32_t nnz = STARPU_CSR_GET_NNZ(descr[0]);
  66. uint32_t nrow = STARPU_CSR_GET_NROW(descr[0]);
  67. float *nzval = (float *)STARPU_CSR_GET_NZVAL(descr[0]);
  68. uint32_t *colind = STARPU_CSR_GET_COLIND(descr[0]);
  69. uint32_t *rowptr = STARPU_CSR_GET_ROWPTR(descr[0]);
  70. uint32_t firstentry = STARPU_CSR_GET_FIRSTENTRY(descr[0]);
  71. float *vecin = (float *)STARPU_VECTOR_GET_PTR(descr[1]);
  72. uint32_t nx_in = STARPU_VECTOR_GET_NX(descr[1]);
  73. float *vecout = (float *)STARPU_VECTOR_GET_PTR(descr[2]);
  74. uint32_t nx_out = STARPU_VECTOR_GET_NX(descr[2]);
  75. dim3 dimBlock(8, 1);
  76. dim3 dimGrid(512, 1);
  77. spmv_kernel_3<<<dimGrid, dimBlock, 0, starpu_cuda_get_local_stream()>>>
  78. (nnz, nrow, nzval, colind, rowptr, firstentry, vecin, nx_in, vecout, nx_out);
  79. cudaError_t status = cudaGetLastError();
  80. if (status != cudaSuccess) STARPU_CUDA_REPORT_ERROR(status);
  81. }