spmv_kernels.c 3.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009, 2010, 2011 Université de Bordeaux 1
  4. * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
  5. * Copyright (C) 2010, 2011 Centre National de la Recherche Scientifique
  6. *
  7. * StarPU is free software; you can redistribute it and/or modify
  8. * it under the terms of the GNU Lesser General Public License as published by
  9. * the Free Software Foundation; either version 2.1 of the License, or (at
  10. * your option) any later version.
  11. *
  12. * StarPU is distributed in the hope that it will be useful, but
  13. * WITHOUT ANY WARRANTY; without even the implied warranty of
  14. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  15. *
  16. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  17. */
  18. #include "spmv.h"
  19. #ifdef STARPU_USE_OPENCL
  20. struct starpu_opencl_program opencl_codelet;
  21. void spmv_kernel_opencl(void *descr[], void *args)
  22. {
  23. cl_kernel kernel;
  24. cl_command_queue queue;
  25. cl_event event;
  26. int id, devid, err, n;
  27. uint32_t nnz = STARPU_CSR_GET_NNZ(descr[0]);
  28. uint32_t nrow = STARPU_CSR_GET_NROW(descr[0]);
  29. cl_mem nzval = (cl_mem)STARPU_CSR_GET_NZVAL(descr[0]);
  30. uint32_t *colind = STARPU_CSR_GET_COLIND(descr[0]);
  31. uint32_t *rowptr = STARPU_CSR_GET_ROWPTR(descr[0]);
  32. uint32_t firstentry = STARPU_CSR_GET_FIRSTENTRY(descr[0]);
  33. cl_mem vecin = (cl_mem)STARPU_VECTOR_GET_PTR(descr[1]);
  34. uint32_t nx_in = STARPU_VECTOR_GET_NX(descr[1]);
  35. cl_mem vecout = (cl_mem)STARPU_VECTOR_GET_PTR(descr[2]);
  36. uint32_t nx_out = STARPU_VECTOR_GET_NX(descr[2]);
  37. id = starpu_worker_get_id();
  38. devid = starpu_worker_get_devid(id);
  39. err = starpu_opencl_load_kernel(&kernel, &queue, &opencl_codelet, "spvm", devid);
  40. if (err != CL_SUCCESS) STARPU_OPENCL_REPORT_ERROR(err);
  41. err = 0;
  42. n=0;
  43. err = clSetKernelArg(kernel, n++, sizeof(nnz), &nnz);
  44. err = clSetKernelArg(kernel, n++, sizeof(nrow), &nrow);
  45. err = clSetKernelArg(kernel, n++, sizeof(nzval), &nzval);
  46. err = clSetKernelArg(kernel, n++, sizeof(colind), &colind);
  47. err = clSetKernelArg(kernel, n++, sizeof(rowptr), &rowptr);
  48. err = clSetKernelArg(kernel, n++, sizeof(firstentry), &firstentry);
  49. err = clSetKernelArg(kernel, n++, sizeof(vecin), &vecin);
  50. err = clSetKernelArg(kernel, n++, sizeof(nx_in), &nx_in);
  51. err = clSetKernelArg(kernel, n++, sizeof(vecout), &vecout);
  52. err = clSetKernelArg(kernel, n++, sizeof(nx_out), &nx_out);
  53. if (err) STARPU_OPENCL_REPORT_ERROR(err);
  54. {
  55. size_t global=1024;
  56. err = clEnqueueNDRangeKernel(queue, kernel, 1, NULL, &global, NULL, 0, NULL, &event);
  57. if (err != CL_SUCCESS) STARPU_OPENCL_REPORT_ERROR(err);
  58. }
  59. clFinish(queue);
  60. starpu_opencl_collect_stats(event);
  61. clReleaseEvent(event);
  62. starpu_opencl_release_kernel(kernel);
  63. }
  64. void compile_spmv_opencl_kernel(void)
  65. {
  66. int ret;
  67. ret = starpu_opencl_load_opencl_from_file("examples/spmv/spmv_opencl.cl", &opencl_codelet, NULL);
  68. if (ret)
  69. {
  70. FPRINTF(stderr, "Failed to compile OpenCL codelet\n");
  71. exit(ret);
  72. }
  73. }
  74. #endif
  75. void spmv_kernel_cpu(void *descr[], __attribute__((unused)) void *arg)
  76. {
  77. float *nzval = (float *)STARPU_CSR_GET_NZVAL(descr[0]);
  78. uint32_t *colind = STARPU_CSR_GET_COLIND(descr[0]);
  79. uint32_t *rowptr = STARPU_CSR_GET_ROWPTR(descr[0]);
  80. float *vecin = (float *)STARPU_VECTOR_GET_PTR(descr[1]);
  81. float *vecout = (float *)STARPU_VECTOR_GET_PTR(descr[2]);
  82. uint32_t firstelem = STARPU_CSR_GET_FIRSTENTRY(descr[0]);
  83. uint32_t nrow;
  84. nrow = STARPU_CSR_GET_NROW(descr[0]);
  85. STARPU_ASSERT(nrow == STARPU_VECTOR_GET_NX(descr[2]));
  86. unsigned row;
  87. for (row = 0; row < nrow; row++)
  88. {
  89. float tmp = 0.0f;
  90. unsigned index;
  91. unsigned firstindex = rowptr[row] - firstelem;
  92. unsigned lastindex = rowptr[row+1] - firstelem;
  93. for (index = firstindex; index < lastindex; index++)
  94. {
  95. unsigned col;
  96. col = colind[index];
  97. tmp += nzval[index]*vecin[col];
  98. }
  99. vecout[row] = tmp;
  100. }
  101. }