dw_spmv.c 6.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301
  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. /*
  17. * Conjugate gradients for Sparse matrices
  18. */
  19. #include "dw_spmv.h"
  20. struct timeval start;
  21. struct timeval end;
  22. unsigned nblocks = 1;
  23. #ifdef STARPU_USE_CUDA
  24. extern void spmv_kernel_cpu_wrapper(uint32_t nnz, uint32_t nrow, float *nzval,
  25. uint32_t *colind, uint32_t *rowptr, uint32_t firstentry,
  26. float *vecin, uint32_t nx_in,
  27. float * vecout, uint32_t nx_out);
  28. void spmv_kernel_cuda(void *descr[], void *args)
  29. {
  30. uint32_t nnz = STARPU_GET_CSR_NNZ(descr[0]);
  31. uint32_t nrow = STARPU_GET_CSR_NROW(descr[0]);
  32. float *nzval = (float *)STARPU_GET_CSR_NZVAL(descr[0]);
  33. uint32_t *colind = STARPU_GET_CSR_COLIND(descr[0]);
  34. uint32_t *rowptr = STARPU_GET_CSR_ROWPTR(descr[0]);
  35. uint32_t firstentry = STARPU_GET_CSR_FIRSTENTRY(descr[0]);
  36. float *vecin = (float *)STARPU_GET_VECTOR_PTR(descr[1]);
  37. uint32_t nx_in = STARPU_GET_VECTOR_NX(descr[1]);
  38. float *vecout = (float *)STARPU_GET_VECTOR_PTR(descr[2]);
  39. uint32_t nx_out = STARPU_GET_VECTOR_NX(descr[2]);
  40. spmv_kernel_cpu_wrapper(nnz, nrow, nzval, colind, rowptr, firstentry, vecin, nx_in, vecout, nx_out);
  41. }
  42. #endif // STARPU_USE_CUDA
  43. uint32_t size = 4194304;
  44. starpu_data_handle sparse_matrix;
  45. starpu_data_handle vector_in, vector_out;
  46. float *sparse_matrix_nzval;
  47. uint32_t *sparse_matrix_colind;
  48. uint32_t *sparse_matrix_rowptr;
  49. float *vector_in_ptr;
  50. float *vector_out_ptr;
  51. unsigned usecpu = 0;
  52. void parse_args(int argc, char **argv)
  53. {
  54. int i;
  55. for (i = 1; i < argc; i++) {
  56. if (strcmp(argv[i], "-size") == 0) {
  57. char *argptr;
  58. size = strtol(argv[++i], &argptr, 10);
  59. }
  60. if (strcmp(argv[i], "-nblocks") == 0) {
  61. char *argptr;
  62. nblocks = strtol(argv[++i], &argptr, 10);
  63. }
  64. if (strcmp(argv[i], "-cpu") == 0) {
  65. usecpu = 1;
  66. }
  67. }
  68. }
  69. void cpu_spmv(void *descr[], __attribute__((unused)) void *arg)
  70. {
  71. float *nzval = (float *)STARPU_GET_CSR_NZVAL(descr[0]);
  72. uint32_t *colind = STARPU_GET_CSR_COLIND(descr[0]);
  73. uint32_t *rowptr = STARPU_GET_CSR_ROWPTR(descr[0]);
  74. float *vecin = (float *)STARPU_GET_VECTOR_PTR(descr[1]);
  75. float *vecout = (float *)STARPU_GET_VECTOR_PTR(descr[2]);
  76. uint32_t firstelem = STARPU_GET_CSR_FIRSTENTRY(descr[0]);
  77. uint32_t nnz;
  78. uint32_t nrow;
  79. nnz = STARPU_GET_CSR_NNZ(descr[0]);
  80. nrow = STARPU_GET_CSR_NROW(descr[0]);
  81. //STARPU_ASSERT(nrow == STARPU_GET_VECTOR_NX(descr[1]));
  82. STARPU_ASSERT(nrow == STARPU_GET_VECTOR_NX(descr[2]));
  83. unsigned row;
  84. for (row = 0; row < nrow; row++)
  85. {
  86. float tmp = 0.0f;
  87. unsigned index;
  88. unsigned firstindex = rowptr[row] - firstelem;
  89. unsigned lastindex = rowptr[row+1] - firstelem;
  90. for (index = firstindex; index < lastindex; index++)
  91. {
  92. unsigned col;
  93. col = colind[index];
  94. tmp += nzval[index]*vecin[col];
  95. }
  96. vecout[row] = tmp;
  97. }
  98. }
  99. void create_data(void)
  100. {
  101. /* we need a sparse symetric (definite positive ?) matrix and a "dense" vector */
  102. /* example of 3-band matrix */
  103. float *nzval;
  104. uint32_t nnz;
  105. uint32_t *colind;
  106. uint32_t *rowptr;
  107. nnz = 3*size-2;
  108. nzval = malloc(nnz*sizeof(float));
  109. colind = malloc(nnz*sizeof(uint32_t));
  110. rowptr = malloc((size+1)*sizeof(uint32_t));
  111. assert(nzval);
  112. assert(colind);
  113. assert(rowptr);
  114. /* fill the matrix */
  115. unsigned row;
  116. unsigned pos = 0;
  117. for (row = 0; row < size; row++)
  118. {
  119. rowptr[row] = pos;
  120. if (row > 0) {
  121. nzval[pos] = 1.0f;
  122. colind[pos] = row-1;
  123. pos++;
  124. }
  125. nzval[pos] = 5.0f;
  126. colind[pos] = row;
  127. pos++;
  128. if (row < size - 1) {
  129. nzval[pos] = 1.0f;
  130. colind[pos] = row+1;
  131. pos++;
  132. }
  133. }
  134. STARPU_ASSERT(pos == nnz);
  135. rowptr[size] = nnz;
  136. starpu_register_csr_data(&sparse_matrix, 0, nnz, size, (uintptr_t)nzval, colind, rowptr, 0, sizeof(float));
  137. sparse_matrix_nzval = nzval;
  138. sparse_matrix_colind = colind;
  139. sparse_matrix_rowptr = rowptr;
  140. /* initiate the 2 vectors */
  141. float *invec, *outvec;
  142. invec = malloc(size*sizeof(float));
  143. assert(invec);
  144. outvec = malloc(size*sizeof(float));
  145. assert(outvec);
  146. /* fill those */
  147. unsigned ind;
  148. for (ind = 0; ind < size; ind++)
  149. {
  150. invec[ind] = 2.0f;
  151. outvec[ind] = 0.0f;
  152. }
  153. starpu_register_vector_data(&vector_in, 0, (uintptr_t)invec, size, sizeof(float));
  154. starpu_register_vector_data(&vector_out, 0, (uintptr_t)outvec, size, sizeof(float));
  155. vector_in_ptr = invec;
  156. vector_out_ptr = outvec;
  157. }
  158. void call_spmv_codelet_filters(void)
  159. {
  160. starpu_codelet *cl = calloc(1, sizeof(starpu_codelet));
  161. /* partition the data along a block distribution */
  162. starpu_filter csr_f, vector_f;
  163. csr_f.filter_func = starpu_vertical_block_filter_func_csr;
  164. csr_f.filter_arg = nblocks;
  165. vector_f.filter_func = starpu_block_filter_func_vector;
  166. vector_f.filter_arg = nblocks;
  167. starpu_partition_data(sparse_matrix, &csr_f);
  168. starpu_partition_data(vector_out, &vector_f);
  169. cl->where = STARPU_CPU|STARPU_CUDA;
  170. cl->cpu_func = cpu_spmv;
  171. #ifdef STARPU_USE_CUDA
  172. cl->cuda_func = spmv_kernel_cuda;
  173. #endif
  174. cl->nbuffers = 3;
  175. cl->model = NULL;
  176. gettimeofday(&start, NULL);
  177. unsigned part;
  178. for (part = 0; part < nblocks; part++)
  179. {
  180. struct starpu_task *task = starpu_task_create();
  181. task->callback_func = NULL;
  182. task->cl = cl;
  183. task->cl_arg = NULL;
  184. task->buffers[0].handle = starpu_get_sub_data(sparse_matrix, 1, part);
  185. task->buffers[0].mode = STARPU_R;
  186. task->buffers[1].handle = vector_in;
  187. task->buffers[1].mode = STARPU_R;
  188. task->buffers[2].handle = starpu_get_sub_data(vector_out, 1, part);
  189. task->buffers[2].mode = STARPU_W;
  190. starpu_submit_task(task);
  191. }
  192. starpu_wait_all_tasks();
  193. gettimeofday(&end, NULL);
  194. starpu_unpartition_data(sparse_matrix, 0);
  195. starpu_unpartition_data(vector_out, 0);
  196. }
  197. void init_problem(void)
  198. {
  199. /* create the sparse input matrix */
  200. create_data();
  201. /* create a new codelet that will perform a SpMV on it */
  202. call_spmv_codelet_filters();
  203. }
  204. void print_results(void)
  205. {
  206. unsigned row;
  207. for (row = 0; row < STARPU_MIN(size, 16); row++)
  208. {
  209. printf("%2.2f\t%2.2f\n", vector_in_ptr[row], vector_out_ptr[row]);
  210. }
  211. }
  212. int main(__attribute__ ((unused)) int argc,
  213. __attribute__ ((unused)) char **argv)
  214. {
  215. parse_args(argc, argv);
  216. /* start the runtime */
  217. starpu_init(NULL);
  218. init_problem();
  219. print_results();
  220. double timing = (double)((end.tv_sec - start.tv_sec)*1000000 + (end.tv_usec - start.tv_usec));
  221. fprintf(stderr, "Computation took (in ms)\n");
  222. printf("%2.2f\n", timing/1000);
  223. return 0;
  224. }