xgemm.c 6.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009, 2010 Université de Bordeaux 1
  4. * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
  5. * Copyright (C) 2010 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 "dw_mult.h"
  19. TYPE *A, *B, *C;
  20. starpu_data_handle A_handle, B_handle, C_handle;
  21. /*
  22. * This program computes C = A * B
  23. *
  24. * A of size (z,y)
  25. * B of size (x,z)
  26. * C of size (x,y)
  27. |---------------|
  28. z | B |
  29. |---------------|
  30. z x
  31. |----| |---------------|
  32. | | | |
  33. | | | |
  34. | A | y | C |
  35. | | | |
  36. | | | |
  37. |----| |---------------|
  38. */
  39. static void check_output(void)
  40. {
  41. /* check results */
  42. /* compute C = C - AB */
  43. CPU_GEMM("N", "N", ydim, xdim, zdim, (TYPE)-1.0f, A, ydim, B, zdim, (TYPE)1.0f, C, ydim);
  44. /* make sure C = 0 */
  45. TYPE err;
  46. err = CPU_ASUM(xdim*ydim, C, 1);
  47. if (err < xdim*ydim*0.001) {
  48. fprintf(stderr, "Results are OK\n");
  49. }
  50. else {
  51. int max;
  52. max = CPU_IAMAX(xdim*ydim, C, 1);
  53. fprintf(stderr, "There were errors ... err = %f\n", err);
  54. fprintf(stderr, "Max error : %e\n", C[max]);
  55. }
  56. }
  57. void callback_func(void *arg)
  58. {
  59. /* do some accounting */
  60. int id = starpu_worker_get_id();
  61. flop_per_worker[id] += BLAS3_FLOP(conf.m, conf.n, conf.k);
  62. ls_per_worker[id] += BLAS3_LS(conf.m, conf.n, conf.k);
  63. }
  64. static void init_problem_data(void)
  65. {
  66. unsigned i,j;
  67. #ifdef STARPU_USE_CUDA
  68. if (pin) {
  69. starpu_data_malloc_pinned_if_possible((void **)&A, zdim*ydim*sizeof(TYPE));
  70. starpu_data_malloc_pinned_if_possible((void **)&B, xdim*zdim*sizeof(TYPE));
  71. starpu_data_malloc_pinned_if_possible((void **)&C, xdim*ydim*sizeof(TYPE));
  72. } else
  73. #endif
  74. {
  75. #ifdef STARPU_HAVE_POSIX_MEMALIGN
  76. posix_memalign((void **)&A, 4096, zdim*ydim*sizeof(TYPE));
  77. posix_memalign((void **)&B, 4096, xdim*zdim*sizeof(TYPE));
  78. posix_memalign((void **)&C, 4096, xdim*ydim*sizeof(TYPE));
  79. #else
  80. A = malloc(zdim*ydim*sizeof(TYPE));
  81. B = malloc(xdim*zdim*sizeof(TYPE));
  82. C = malloc(xdim*ydim*sizeof(TYPE));
  83. #endif
  84. }
  85. /* fill the A and B matrices */
  86. if (norandom) {
  87. for (j=0; j < ydim; j++) {
  88. for (i=0; i < zdim; i++) {
  89. A[j+i*ydim] = (TYPE)(i);
  90. }
  91. }
  92. for (j=0; j < zdim; j++) {
  93. for (i=0; i < xdim; i++) {
  94. B[j+i*zdim] = (TYPE)(j);
  95. }
  96. }
  97. }
  98. else {
  99. #ifdef NORANDOM
  100. srand(2008);
  101. STARPU_ABORT();
  102. #endif
  103. for (j=0; j < ydim; j++) {
  104. for (i=0; i < zdim; i++) {
  105. A[j+i*ydim] = (TYPE)(starpu_drand48());
  106. }
  107. }
  108. for (j=0; j < zdim; j++) {
  109. for (i=0; i < xdim; i++) {
  110. B[j+i*zdim] = (TYPE)(starpu_drand48());
  111. }
  112. }
  113. }
  114. for (j=0; j < ydim; j++) {
  115. for (i=0; i < xdim; i++) {
  116. C[j+i*ydim] = (TYPE)(0);
  117. }
  118. }
  119. display_memory_consumption();
  120. }
  121. static void partition_mult_data(void)
  122. {
  123. starpu_matrix_data_register(&A_handle, 0, (uintptr_t)A,
  124. ydim, ydim, zdim, sizeof(TYPE));
  125. starpu_matrix_data_register(&B_handle, 0, (uintptr_t)B,
  126. zdim, zdim, xdim, sizeof(TYPE));
  127. starpu_matrix_data_register(&C_handle, 0, (uintptr_t)C,
  128. ydim, ydim, xdim, sizeof(TYPE));
  129. starpu_data_set_wt_mask(C_handle, 1<<0);
  130. conf.k = zdim;
  131. conf.m = ydim/nslicesy;
  132. conf.n = xdim/nslicesx;
  133. struct starpu_data_filter f;
  134. f.filter_func = starpu_vertical_block_filter_func;
  135. f.nchildren = nslicesx;
  136. f.get_nchildren = NULL;
  137. f.get_child_ops = NULL;
  138. struct starpu_data_filter f2;
  139. f2.filter_func = starpu_block_filter_func;
  140. f2.nchildren = nslicesy;
  141. f2.get_nchildren = NULL;
  142. f2.get_child_ops = NULL;
  143. starpu_data_partition(B_handle, &f);
  144. starpu_data_partition(A_handle, &f2);
  145. starpu_data_map_filters(C_handle, 2, &f, &f2);
  146. }
  147. static void unpartition_mult_data(void)
  148. {
  149. fprintf(stderr, "unpartition !!\n");
  150. starpu_data_unpartition(C_handle, 0);
  151. starpu_data_unregister(C_handle);
  152. }
  153. static starpu_codelet cl = {
  154. .where = STARPU_CPU|STARPU_CUDA
  155. #ifdef SPU_FUNC_SGEMM
  156. |STARPU_GORDON
  157. #endif
  158. ,
  159. .cpu_func = STARPU_GEMM(cpu_mult),
  160. #ifdef STARPU_USE_CUDA
  161. .cuda_func = STARPU_GEMM(cublas_mult),
  162. #endif
  163. #ifdef STARPU_USE_GORDON
  164. #ifdef SPU_FUNC_SGEMM
  165. .gordon_func = SPU_FUNC_SGEMM,
  166. #else
  167. #warning SPU_FUNC_SGEMM is not available
  168. #endif
  169. #endif
  170. .nbuffers = 3
  171. };
  172. static struct starpu_task *construct_task(unsigned x, unsigned y, unsigned z, unsigned iter)
  173. {
  174. /* A B[task] = C[task] */
  175. struct starpu_task *task = starpu_task_create();
  176. task->cl = &cl;
  177. /* we have a callback to do some accounting */
  178. task->callback_func = callback_func;
  179. task->callback_arg = NULL;
  180. task->buffers[0].handle = starpu_data_get_sub_data(A_handle, 1, y);
  181. task->buffers[0].mode = STARPU_R;
  182. task->buffers[1].handle = starpu_data_get_sub_data(B_handle, 1, x);
  183. task->buffers[1].mode = STARPU_R;
  184. task->buffers[2].handle = starpu_data_get_sub_data(C_handle, 2, x, y);
  185. task->buffers[2].mode = STARPU_RW;
  186. task->cl_arg = &conf;
  187. task->cl_arg_size = sizeof(struct block_conf);
  188. return task;
  189. }
  190. static void submit_new_iter(unsigned x, unsigned y, unsigned iter)
  191. {
  192. unsigned z;
  193. z = 0;
  194. {
  195. struct starpu_task *task;
  196. task = construct_task(x, y, z, iter);
  197. starpu_task_submit(task);
  198. }
  199. }
  200. static void launch_codelets(void)
  201. {
  202. #ifdef STARPU_USE_FXT
  203. _starpu_fxt_register_thread(0);
  204. #endif
  205. /* partition the work into slices */
  206. unsigned taskx, tasky;
  207. srand(time(NULL));
  208. /* should we use a single performance model for all archs and use an
  209. * acceleration factor ? */
  210. if (use_common_model) {
  211. cl.model = &STARPU_GEMM(model_common);
  212. }
  213. else {
  214. cl.model = &STARPU_GEMM(model);
  215. }
  216. for (taskx = 0; taskx < nslicesx; taskx++)
  217. {
  218. for (tasky = 0; tasky < nslicesy; tasky++)
  219. {
  220. submit_new_iter(taskx, tasky, 0);
  221. }
  222. }
  223. }
  224. int main(__attribute__ ((unused)) int argc,
  225. __attribute__ ((unused)) char **argv)
  226. {
  227. parse_args(argc, argv);
  228. /* start the runtime */
  229. starpu_init(NULL);
  230. starpu_helper_cublas_init();
  231. init_problem_data();
  232. gettimeofday(&start, NULL);
  233. partition_mult_data();
  234. launch_codelets();
  235. starpu_task_wait_for_all();
  236. gettimeofday(&end, NULL);
  237. double timing = (double)((end.tv_sec - start.tv_sec)*1000000 + (end.tv_usec - start.tv_usec));
  238. display_stats(timing);
  239. unpartition_mult_data();
  240. if (check)
  241. check_output();
  242. starpu_helper_cublas_shutdown();
  243. starpu_shutdown();
  244. return 0;
  245. }