cholesky_implicit.c 8.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009-2015 Université de Bordeaux
  4. * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
  5. * Copyright (C) 2010, 2011, 2012, 2013 CNRS
  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. /*
  19. * This version of the Cholesky factorization uses implicit dependency computation.
  20. * The whole algorithm thus appears clearly in the task submission loop in _cholesky().
  21. */
  22. #include "cholesky.h"
  23. #include "../sched_ctx_utils/sched_ctx_utils.h"
  24. #if defined(STARPU_USE_CUDA) && defined(STARPU_HAVE_MAGMA)
  25. #include "magma.h"
  26. #endif
  27. /*
  28. * code to bootstrap the factorization
  29. * and construct the DAG
  30. */
  31. static void callback_turn_spmd_on(void *arg STARPU_ATTRIBUTE_UNUSED)
  32. {
  33. cl22.type = STARPU_SPMD;
  34. }
  35. static int _cholesky(starpu_data_handle_t dataA, unsigned nblocks)
  36. {
  37. int ret;
  38. double start;
  39. double end;
  40. unsigned i,j,k;
  41. unsigned long n = starpu_matrix_get_nx(dataA);
  42. unsigned long nn = n/nblocks;
  43. int prio_level = noprio?STARPU_DEFAULT_PRIO:STARPU_MAX_PRIO;
  44. if (bound || bound_lp || bound_mps)
  45. starpu_bound_start(bound_deps, 0);
  46. starpu_fxt_start_profiling();
  47. start = starpu_timing_now();
  48. /* create all the DAG nodes */
  49. for (k = 0; k < nblocks; k++)
  50. {
  51. starpu_data_handle_t sdatakk = starpu_data_get_sub_data(dataA, 2, k, k);
  52. ret = starpu_task_insert(&cl11,
  53. STARPU_PRIORITY, prio_level,
  54. STARPU_RW, sdatakk,
  55. STARPU_CALLBACK, (k == 3*nblocks/4)?callback_turn_spmd_on:NULL,
  56. STARPU_FLOPS, (double) FLOPS_SPOTRF(nn),
  57. STARPU_TAG_ONLY, TAG11(k),
  58. STARPU_NAME, "cl11",
  59. 0);
  60. if (ret == -ENODEV) return 77;
  61. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_insert");
  62. for (j = k+1; j<nblocks; j++)
  63. {
  64. starpu_data_handle_t sdatakj = starpu_data_get_sub_data(dataA, 2, k, j);
  65. ret = starpu_task_insert(&cl21,
  66. STARPU_PRIORITY, (j == k+1)?prio_level:STARPU_DEFAULT_PRIO,
  67. STARPU_R, sdatakk,
  68. STARPU_RW, sdatakj,
  69. STARPU_FLOPS, (double) FLOPS_STRSM(nn, nn),
  70. STARPU_TAG_ONLY, TAG21(k,j),
  71. STARPU_NAME, "cl21",
  72. 0);
  73. if (ret == -ENODEV) return 77;
  74. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_insert");
  75. }
  76. starpu_data_wont_use(sdatakk);
  77. for (j = k+1; j<nblocks; j++)
  78. {
  79. starpu_data_handle_t sdatakj = starpu_data_get_sub_data(dataA, 2, k, j);
  80. for (i = k+1; i<nblocks; i++)
  81. {
  82. if (i <= j)
  83. {
  84. starpu_data_handle_t sdataki = starpu_data_get_sub_data(dataA, 2, k, i);
  85. starpu_data_handle_t sdataij = starpu_data_get_sub_data(dataA, 2, i, j);
  86. ret = starpu_task_insert(&cl22,
  87. STARPU_PRIORITY, ((i == k+1) && (j == k+1))?prio_level:STARPU_DEFAULT_PRIO,
  88. STARPU_R, sdataki,
  89. STARPU_R, sdatakj,
  90. cl22.modes[2], sdataij,
  91. STARPU_FLOPS, (double) FLOPS_SGEMM(nn, nn, nn),
  92. STARPU_TAG_ONLY, TAG22(k,i,j),
  93. STARPU_NAME, ((i == k+1) && (j == k+1))?"cl22_p":"cl22",
  94. 0);
  95. if (ret == -ENODEV) return 77;
  96. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_insert");
  97. }
  98. }
  99. starpu_data_wont_use(sdatakj);
  100. }
  101. }
  102. starpu_task_wait_for_all();
  103. end = starpu_timing_now();
  104. starpu_fxt_stop_profiling();
  105. if (bound || bound_lp || bound_mps)
  106. starpu_bound_stop();
  107. double timing = end - start;
  108. double flop = FLOPS_SPOTRF(n);
  109. if(with_ctxs || with_noctxs || chole1 || chole2)
  110. update_sched_ctx_timing_results((flop/timing/1000.0f), (timing/1000000.0f));
  111. else
  112. {
  113. PRINTF("# size\tms\tGFlops");
  114. if (bound)
  115. PRINTF("\tTms\tTGFlops");
  116. PRINTF("\n");
  117. PRINTF("%lu\t%.0f\t%.1f", n, timing/1000, (flop/timing/1000.0f));
  118. if (bound_lp)
  119. {
  120. FILE *f = fopen("cholesky.lp", "w");
  121. starpu_bound_print_lp(f);
  122. }
  123. if (bound_mps)
  124. {
  125. FILE *f = fopen("cholesky.mps", "w");
  126. starpu_bound_print_mps(f);
  127. }
  128. if (bound)
  129. {
  130. double res;
  131. starpu_bound_compute(&res, NULL, 0);
  132. PRINTF("\t%.0f\t%.1f", res, (flop/res/1000000.0f));
  133. }
  134. PRINTF("\n");
  135. }
  136. return 0;
  137. }
  138. static int cholesky(float *matA, unsigned size, unsigned ld, unsigned nblocks)
  139. {
  140. starpu_data_handle_t dataA;
  141. /* monitor and partition the A matrix into blocks :
  142. * one block is now determined by 2 unsigned (i,j) */
  143. starpu_matrix_data_register(&dataA, STARPU_MAIN_RAM, (uintptr_t)matA, ld, size, size, sizeof(float));
  144. struct starpu_data_filter f =
  145. {
  146. .filter_func = starpu_matrix_filter_vertical_block,
  147. .nchildren = nblocks
  148. };
  149. struct starpu_data_filter f2 =
  150. {
  151. .filter_func = starpu_matrix_filter_block,
  152. .nchildren = nblocks
  153. };
  154. starpu_data_map_filters(dataA, 2, &f, &f2);
  155. int ret = _cholesky(dataA, nblocks);
  156. starpu_data_unpartition(dataA, STARPU_MAIN_RAM);
  157. starpu_data_unregister(dataA);
  158. return ret;
  159. }
  160. static void execute_cholesky(unsigned size, unsigned nblocks)
  161. {
  162. int ret;
  163. float *mat = NULL;
  164. unsigned i,j;
  165. #ifndef STARPU_SIMGRID
  166. starpu_malloc((void **)&mat, (size_t)size*size*sizeof(float));
  167. for (i = 0; i < size; i++)
  168. {
  169. for (j = 0; j < size; j++)
  170. {
  171. mat[j +i*size] = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);
  172. /* mat[j +i*size] = ((i == j)?1.0f*size:0.0f); */
  173. }
  174. }
  175. #endif
  176. /* #define PRINT_OUTPUT */
  177. #ifdef PRINT_OUTPUT
  178. FPRINTF(stdout, "Input :\n");
  179. for (j = 0; j < size; j++)
  180. {
  181. for (i = 0; i < size; i++)
  182. {
  183. if (i <= j)
  184. {
  185. FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);
  186. }
  187. else
  188. {
  189. FPRINTF(stdout, ".\t");
  190. }
  191. }
  192. FPRINTF(stdout, "\n");
  193. }
  194. #endif
  195. ret = cholesky(mat, size, size, nblocks);
  196. #ifdef PRINT_OUTPUT
  197. FPRINTF(stdout, "Results :\n");
  198. for (j = 0; j < size; j++)
  199. {
  200. for (i = 0; i < size; i++)
  201. {
  202. if (i <= j)
  203. {
  204. FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);
  205. }
  206. else
  207. {
  208. FPRINTF(stdout, ".\t");
  209. mat[j+i*size] = 0.0f; /* debug */
  210. }
  211. }
  212. FPRINTF(stdout, "\n");
  213. }
  214. #endif
  215. if (check)
  216. {
  217. FPRINTF(stderr, "compute explicit LLt ...\n");
  218. for (j = 0; j < size; j++)
  219. {
  220. for (i = 0; i < size; i++)
  221. {
  222. if (i > j)
  223. {
  224. mat[j+i*size] = 0.0f; /* debug */
  225. }
  226. }
  227. }
  228. float *test_mat = malloc(size*size*sizeof(float));
  229. STARPU_ASSERT(test_mat);
  230. STARPU_SSYRK("L", "N", size, size, 1.0f,
  231. mat, size, 0.0f, test_mat, size);
  232. FPRINTF(stderr, "comparing results ...\n");
  233. #ifdef PRINT_OUTPUT
  234. for (j = 0; j < size; j++)
  235. {
  236. for (i = 0; i < size; i++)
  237. {
  238. if (i <= j)
  239. {
  240. FPRINTF(stdout, "%2.2f\t", test_mat[j +i*size]);
  241. }
  242. else
  243. {
  244. FPRINTF(stdout, ".\t");
  245. }
  246. }
  247. FPRINTF(stdout, "\n");
  248. }
  249. #endif
  250. for (j = 0; j < size; j++)
  251. {
  252. for (i = 0; i < size; i++)
  253. {
  254. if (i <= j)
  255. {
  256. float orig = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);
  257. float err = abs(test_mat[j +i*size] - orig);
  258. if (err > 0.00001)
  259. {
  260. FPRINTF(stderr, "Error[%u, %u] --> %2.2f != %2.2f (err %2.2f)\n", i, j, test_mat[j +i*size], orig, err);
  261. assert(0);
  262. }
  263. }
  264. }
  265. }
  266. free(test_mat);
  267. }
  268. starpu_free(mat);
  269. }
  270. int main(int argc, char **argv)
  271. {
  272. /* create a simple definite positive symetric matrix example
  273. *
  274. * Hilbert matrix : h(i,j) = 1/(i+j+1)
  275. * */
  276. parse_args(argc, argv);
  277. if(with_ctxs || with_noctxs || chole1 || chole2)
  278. parse_args_ctx(argc, argv);
  279. #ifdef STARPU_HAVE_MAGMA
  280. magma_init();
  281. #endif
  282. int ret;
  283. ret = starpu_init(NULL);
  284. starpu_fxt_stop_profiling();
  285. if (ret == -ENODEV)
  286. return 77;
  287. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  288. #ifdef STARPU_USE_CUDA
  289. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,cuda_chol_task_11_cost);
  290. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,cuda_chol_task_21_cost);
  291. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,cuda_chol_task_22_cost);
  292. #else
  293. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,NULL);
  294. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,NULL);
  295. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,NULL);
  296. #endif
  297. starpu_cublas_init();
  298. if(with_ctxs)
  299. {
  300. construct_contexts(execute_cholesky);
  301. start_2benchs(execute_cholesky);
  302. }
  303. else if(with_noctxs)
  304. start_2benchs(execute_cholesky);
  305. else if(chole1)
  306. start_1stbench(execute_cholesky);
  307. else if(chole2)
  308. start_2ndbench(execute_cholesky);
  309. else
  310. execute_cholesky(size, nblocks);
  311. starpu_cublas_shutdown();
  312. starpu_shutdown();
  313. return ret;
  314. }