cholesky_tile_tag.c 7.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009-2017 Université de Bordeaux
  4. * Copyright (C) 2012,2013 Inria
  5. * Copyright (C) 2010-2013,2015-2017 CNRS
  6. * Copyright (C) 2013 Thibaut Lambert
  7. *
  8. * StarPU is free software; you can redistribute it and/or modify
  9. * it under the terms of the GNU Lesser General Public License as published by
  10. * the Free Software Foundation; either version 2.1 of the License, or (at
  11. * your option) any later version.
  12. *
  13. * StarPU is distributed in the hope that it will be useful, but
  14. * WITHOUT ANY WARRANTY; without even the implied warranty of
  15. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  16. *
  17. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  18. */
  19. /*
  20. * This version of the Cholesky factorization uses explicit dependency
  21. * declaration through dependency tags.
  22. * It also directly registers matrix tiles instead of using partitioning.
  23. */
  24. #include "cholesky.h"
  25. #if defined(STARPU_USE_CUDA) && defined(STARPU_HAVE_MAGMA)
  26. #include "magma.h"
  27. #endif
  28. /* A [ y ] [ x ] */
  29. float *A[NMAXBLOCKS][NMAXBLOCKS];
  30. starpu_data_handle_t A_state[NMAXBLOCKS][NMAXBLOCKS];
  31. /*
  32. * Some useful functions
  33. */
  34. static struct starpu_task *create_task(starpu_tag_t id)
  35. {
  36. struct starpu_task *task = starpu_task_create();
  37. task->cl_arg = NULL;
  38. task->use_tag = 1;
  39. task->tag_id = id;
  40. return task;
  41. }
  42. /*
  43. * Create the codelets
  44. */
  45. static struct starpu_task * create_task_11(unsigned k, unsigned nblocks)
  46. {
  47. (void)nblocks;
  48. /* FPRINTF(stdout, "task 11 k = %d TAG = %llx\n", k, (TAG11(k))); */
  49. struct starpu_task *task = create_task(TAG11(k));
  50. task->cl = &cl11;
  51. /* which sub-data is manipulated ? */
  52. task->handles[0] = A_state[k][k];
  53. /* this is an important task */
  54. task->priority = STARPU_MAX_PRIO;
  55. /* enforce dependencies ... */
  56. if (k > 0)
  57. {
  58. starpu_tag_declare_deps(TAG11(k), 1, TAG22(k-1, k, k));
  59. }
  60. int n = starpu_matrix_get_nx(task->handles[0]);
  61. task->flops = FLOPS_SPOTRF(n);
  62. return task;
  63. }
  64. static int create_task_21(unsigned k, unsigned j)
  65. {
  66. int ret;
  67. struct starpu_task *task = create_task(TAG21(k, j));
  68. task->cl = &cl21;
  69. /* which sub-data is manipulated ? */
  70. task->handles[0] = A_state[k][k];
  71. task->handles[1] = A_state[j][k];
  72. if (j == k+1)
  73. {
  74. task->priority = STARPU_MAX_PRIO;
  75. }
  76. /* enforce dependencies ... */
  77. if (k > 0)
  78. {
  79. starpu_tag_declare_deps(TAG21(k, j), 2, TAG11(k), TAG22(k-1, k, j));
  80. }
  81. else
  82. {
  83. starpu_tag_declare_deps(TAG21(k, j), 1, TAG11(k));
  84. }
  85. int n = starpu_matrix_get_nx(task->handles[0]);
  86. task->flops = FLOPS_STRSM(n, n);
  87. ret = starpu_task_submit(task);
  88. if (ret != -ENODEV) STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_submit");
  89. return ret;
  90. }
  91. static int create_task_22(unsigned k, unsigned i, unsigned j)
  92. {
  93. int ret;
  94. /* FPRINTF(stdout, "task 22 k,i,j = %d,%d,%d TAG = %llx\n", k,i,j, TAG22(k,i,j)); */
  95. struct starpu_task *task = create_task(TAG22(k, i, j));
  96. task->cl = &cl22;
  97. /* which sub-data is manipulated ? */
  98. task->handles[0] = A_state[i][k];
  99. task->handles[1] = A_state[j][k];
  100. task->handles[2] = A_state[j][i];
  101. if ( (i == k + 1) && (j == k +1) )
  102. {
  103. task->priority = STARPU_MAX_PRIO;
  104. }
  105. /* enforce dependencies ... */
  106. if (k > 0)
  107. {
  108. starpu_tag_declare_deps(TAG22(k, i, j), 3, TAG22(k-1, i, j), TAG21(k, i), TAG21(k, j));
  109. }
  110. else
  111. {
  112. starpu_tag_declare_deps(TAG22(k, i, j), 2, TAG21(k, i), TAG21(k, j));
  113. }
  114. int n = starpu_matrix_get_nx(task->handles[0]);
  115. task->flops = FLOPS_SGEMM(n, n, n);
  116. ret = starpu_task_submit(task);
  117. if (ret != -ENODEV) STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_submit");
  118. return ret;
  119. }
  120. /*
  121. * code to bootstrap the factorization
  122. * and construct the DAG
  123. */
  124. static int cholesky_no_stride(void)
  125. {
  126. int ret;
  127. double start;
  128. double end;
  129. struct starpu_task *entry_task = NULL;
  130. /* create all the DAG nodes */
  131. unsigned i,j,k;
  132. for (k = 0; k < nblocks_p; k++)
  133. {
  134. starpu_iteration_push(k);
  135. struct starpu_task *task = create_task_11(k, nblocks_p);
  136. /* we defer the launch of the first task */
  137. if (k == 0)
  138. {
  139. entry_task = task;
  140. }
  141. else
  142. {
  143. ret = starpu_task_submit(task);
  144. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_submit");
  145. }
  146. for (j = k+1; j<nblocks_p; j++)
  147. {
  148. ret = create_task_21(k, j);
  149. if (ret == -ENODEV) return 77;
  150. for (i = k+1; i<nblocks_p; i++)
  151. {
  152. if (i <= j)
  153. {
  154. ret = create_task_22(k, i, j);
  155. if (ret == -ENODEV) return 77;
  156. }
  157. }
  158. }
  159. starpu_iteration_pop();
  160. }
  161. /* schedule the codelet */
  162. start = starpu_timing_now();
  163. ret = starpu_task_submit(entry_task);
  164. if (ret == -ENODEV) return 77;
  165. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_submit");
  166. /* stall the application until the end of computations */
  167. starpu_tag_wait(TAG11(nblocks_p-1));
  168. end = starpu_timing_now();
  169. double timing = end - start;
  170. double flop = (1.0f*size_p*size_p*size_p)/3.0f;
  171. PRINTF("# size\tms\tGFlops\n");
  172. PRINTF("%u\t%.0f\t%.1f\n", size_p, timing/1000, (flop/timing/1000.0f));
  173. return 0;
  174. }
  175. int main(int argc, char **argv)
  176. {
  177. unsigned x, y;
  178. int ret;
  179. #ifdef STARPU_HAVE_MAGMA
  180. magma_init();
  181. #endif
  182. ret = starpu_init(NULL);
  183. if (ret == -ENODEV)
  184. return 77;
  185. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  186. init_sizes();
  187. parse_args(argc, argv);
  188. assert(nblocks_p <= NMAXBLOCKS);
  189. FPRINTF(stderr, "BLOCK SIZE = %u\n", size_p / nblocks_p);
  190. #ifdef STARPU_USE_CUDA
  191. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,cuda_chol_task_11_cost);
  192. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,cuda_chol_task_21_cost);
  193. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,cuda_chol_task_22_cost);
  194. #else
  195. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,NULL);
  196. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,NULL);
  197. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,NULL);
  198. #endif
  199. /* Disable sequential consistency */
  200. starpu_data_set_default_sequential_consistency_flag(0);
  201. starpu_cublas_init();
  202. for (y = 0; y < nblocks_p; y++)
  203. for (x = 0; x < nblocks_p; x++)
  204. {
  205. if (x <= y)
  206. {
  207. starpu_malloc_flags((void **)&A[y][x], BLOCKSIZE*BLOCKSIZE*sizeof(float), STARPU_MALLOC_PINNED|STARPU_MALLOC_SIMULATION_FOLDED);
  208. assert(A[y][x]);
  209. }
  210. }
  211. #ifndef STARPU_SIMGRID
  212. /* create a simple definite positive symetric matrix example
  213. *
  214. * Hilbert matrix : h(i,j) = 1/(i+j+1) ( + n In to make is stable )
  215. * */
  216. for (y = 0; y < nblocks_p; y++)
  217. for (x = 0; x < nblocks_p; x++)
  218. if (x <= y)
  219. {
  220. unsigned i, j;
  221. for (i = 0; i < BLOCKSIZE; i++)
  222. for (j = 0; j < BLOCKSIZE; j++)
  223. {
  224. A[y][x][i*BLOCKSIZE + j] =
  225. (float)(1.0f/((float) (1.0+(x*BLOCKSIZE+i)+(y*BLOCKSIZE+j))));
  226. /* make it a little more numerically stable ... ;) */
  227. if ((x == y) && (i == j))
  228. A[y][x][i*BLOCKSIZE + j] += (float)(2*size_p);
  229. }
  230. }
  231. #endif
  232. for (y = 0; y < nblocks_p; y++)
  233. for (x = 0; x < nblocks_p; x++)
  234. {
  235. if (x <= y)
  236. {
  237. starpu_matrix_data_register(&A_state[y][x], STARPU_MAIN_RAM, (uintptr_t)A[y][x],
  238. BLOCKSIZE, BLOCKSIZE, BLOCKSIZE, sizeof(float));
  239. starpu_data_set_coordinates(A_state[y][x], 2, x, y);
  240. }
  241. }
  242. ret = cholesky_no_stride();
  243. for (y = 0; y < nblocks_p; y++)
  244. for (x = 0; x < nblocks_p; x++)
  245. {
  246. if (x <= y)
  247. {
  248. starpu_data_unregister(A_state[y][x]);
  249. starpu_free_flags(A[y][x], BLOCKSIZE*BLOCKSIZE*sizeof(float), STARPU_MALLOC_PINNED|STARPU_MALLOC_SIMULATION_FOLDED);
  250. }
  251. }
  252. starpu_cublas_shutdown();
  253. starpu_shutdown();
  254. return ret;
  255. }