spmv_cuda.cu 3.1 KB

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