cholesky_implicit.c 9.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009-2017 Université de Bordeaux
  4. * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
  5. * Copyright (C) 2010, 2011, 2012, 2013, 2016, 2017 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. double start;
  38. double end;
  39. unsigned i,j,k;
  40. unsigned long n = starpu_matrix_get_nx(dataA);
  41. unsigned long nn = n/nblocks;
  42. unsigned unbound_prio = STARPU_MAX_PRIO == INT_MAX && STARPU_MIN_PRIO == INT_MIN;
  43. if (bound_p || bound_lp_p || bound_mps_p)
  44. starpu_bound_start(bound_deps_p, 0);
  45. starpu_fxt_start_profiling();
  46. start = starpu_timing_now();
  47. /* create all the DAG nodes */
  48. for (k = 0; k < nblocks; k++)
  49. {
  50. int ret;
  51. starpu_data_handle_t sdatakk = starpu_data_get_sub_data(dataA, 2, k, k);
  52. ret = starpu_task_insert(&cl11,
  53. STARPU_PRIORITY, noprio_p ? STARPU_DEFAULT_PRIO : unbound_prio ? 2*nblocks - 2*k : STARPU_MAX_PRIO,
  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. 0);
  59. if (ret == -ENODEV) return 77;
  60. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_insert");
  61. for (j = k+1; j<nblocks; j++)
  62. {
  63. starpu_data_handle_t sdatakj = starpu_data_get_sub_data(dataA, 2, k, j);
  64. ret = starpu_task_insert(&cl21,
  65. STARPU_PRIORITY, noprio_p ? STARPU_DEFAULT_PRIO : unbound_prio ? 2*nblocks - 2*k - j : (j == k+1)?STARPU_MAX_PRIO:STARPU_DEFAULT_PRIO,
  66. STARPU_R, sdatakk,
  67. STARPU_RW, sdatakj,
  68. STARPU_FLOPS, (double) FLOPS_STRSM(nn, nn),
  69. STARPU_TAG_ONLY, TAG21(k,j),
  70. 0);
  71. if (ret == -ENODEV) return 77;
  72. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_insert");
  73. }
  74. starpu_data_wont_use(sdatakk);
  75. for (j = k+1; j<nblocks; j++)
  76. {
  77. starpu_data_handle_t sdatakj = starpu_data_get_sub_data(dataA, 2, k, j);
  78. for (i = k+1; i<nblocks; i++)
  79. {
  80. if (i <= j)
  81. {
  82. starpu_data_handle_t sdataki = starpu_data_get_sub_data(dataA, 2, k, i);
  83. starpu_data_handle_t sdataij = starpu_data_get_sub_data(dataA, 2, i, j);
  84. ret = starpu_task_insert(&cl22,
  85. STARPU_PRIORITY, noprio_p ? STARPU_DEFAULT_PRIO : unbound_prio ? 2*nblocks - 2*k - j - i : ((i == k+1) && (j == k+1))?STARPU_MAX_PRIO:STARPU_DEFAULT_PRIO,
  86. STARPU_R, sdataki,
  87. STARPU_R, sdatakj,
  88. cl22.modes[2], sdataij,
  89. STARPU_FLOPS, (double) FLOPS_SGEMM(nn, nn, nn),
  90. STARPU_TAG_ONLY, TAG22(k,i,j),
  91. 0);
  92. if (ret == -ENODEV) return 77;
  93. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_insert");
  94. }
  95. }
  96. starpu_data_wont_use(sdatakj);
  97. }
  98. }
  99. starpu_task_wait_for_all();
  100. end = starpu_timing_now();
  101. starpu_fxt_stop_profiling();
  102. if (bound_p || bound_lp_p || bound_mps_p)
  103. starpu_bound_stop();
  104. double timing = end - start;
  105. double flop = FLOPS_SPOTRF(n);
  106. if(with_ctxs_p || with_noctxs_p || chole1_p || chole2_p)
  107. update_sched_ctx_timing_results((flop/timing/1000.0f), (timing/1000000.0f));
  108. else
  109. {
  110. PRINTF("# size\tms\tGFlops");
  111. if (bound_p)
  112. PRINTF("\tTms\tTGFlops");
  113. PRINTF("\n");
  114. PRINTF("%lu\t%.0f\t%.1f", n, timing/1000, (flop/timing/1000.0f));
  115. if (bound_lp_p)
  116. {
  117. FILE *f = fopen("cholesky.lp", "w");
  118. starpu_bound_print_lp(f);
  119. fclose(f);
  120. }
  121. if (bound_mps_p)
  122. {
  123. FILE *f = fopen("cholesky.mps", "w");
  124. starpu_bound_print_mps(f);
  125. fclose(f);
  126. }
  127. if (bound_p)
  128. {
  129. double res;
  130. starpu_bound_compute(&res, NULL, 0);
  131. PRINTF("\t%.0f\t%.1f", res, (flop/res/1000000.0f));
  132. }
  133. PRINTF("\n");
  134. }
  135. return 0;
  136. }
  137. static int cholesky(float *matA, unsigned size, unsigned ld, unsigned nblocks)
  138. {
  139. starpu_data_handle_t dataA;
  140. unsigned x, y;
  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. for (x = 0; x < nblocks; x++)
  156. for (y = 0; y < nblocks; y++)
  157. {
  158. starpu_data_handle_t data = starpu_data_get_sub_data(dataA, 2, x, y);
  159. starpu_data_set_coordinates(data, 2, x, y);
  160. }
  161. int ret = _cholesky(dataA, nblocks);
  162. starpu_data_unpartition(dataA, STARPU_MAIN_RAM);
  163. starpu_data_unregister(dataA);
  164. return ret;
  165. }
  166. static void execute_cholesky(unsigned size, unsigned nblocks)
  167. {
  168. float *mat = NULL;
  169. unsigned i,j;
  170. starpu_malloc_flags((void **)&mat, (size_t)size*size*sizeof(float), STARPU_MALLOC_PINNED|STARPU_MALLOC_SIMULATION_FOLDED);
  171. #ifndef STARPU_SIMGRID
  172. for (i = 0; i < size; i++)
  173. {
  174. for (j = 0; j < size; j++)
  175. {
  176. mat[j +i*size] = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);
  177. /* mat[j +i*size] = ((i == j)?1.0f*size:0.0f); */
  178. }
  179. }
  180. #endif
  181. /* #define PRINT_OUTPUT */
  182. #ifdef PRINT_OUTPUT
  183. FPRINTF(stdout, "Input :\n");
  184. for (j = 0; j < size; j++)
  185. {
  186. for (i = 0; i < size; i++)
  187. {
  188. if (i <= j)
  189. {
  190. FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);
  191. }
  192. else
  193. {
  194. FPRINTF(stdout, ".\t");
  195. }
  196. }
  197. FPRINTF(stdout, "\n");
  198. }
  199. #endif
  200. cholesky(mat, size, size, nblocks);
  201. #ifdef PRINT_OUTPUT
  202. FPRINTF(stdout, "Results :\n");
  203. for (j = 0; j < size; j++)
  204. {
  205. for (i = 0; i < size; i++)
  206. {
  207. if (i <= j)
  208. {
  209. FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);
  210. }
  211. else
  212. {
  213. FPRINTF(stdout, ".\t");
  214. mat[j+i*size] = 0.0f; /* debug */
  215. }
  216. }
  217. FPRINTF(stdout, "\n");
  218. }
  219. #endif
  220. if (check_p)
  221. {
  222. FPRINTF(stderr, "compute explicit LLt ...\n");
  223. for (j = 0; j < size; j++)
  224. {
  225. for (i = 0; i < size; i++)
  226. {
  227. if (i > j)
  228. {
  229. mat[j+i*size] = 0.0f; /* debug */
  230. }
  231. }
  232. }
  233. float *test_mat = malloc(size*size*sizeof(float));
  234. STARPU_ASSERT(test_mat);
  235. STARPU_SSYRK("L", "N", size, size, 1.0f,
  236. mat, size, 0.0f, test_mat, size);
  237. FPRINTF(stderr, "comparing results ...\n");
  238. #ifdef PRINT_OUTPUT
  239. for (j = 0; j < size; j++)
  240. {
  241. for (i = 0; i < size; i++)
  242. {
  243. if (i <= j)
  244. {
  245. FPRINTF(stdout, "%2.2f\t", test_mat[j +i*size]);
  246. }
  247. else
  248. {
  249. FPRINTF(stdout, ".\t");
  250. }
  251. }
  252. FPRINTF(stdout, "\n");
  253. }
  254. #endif
  255. for (j = 0; j < size; j++)
  256. {
  257. for (i = 0; i < size; i++)
  258. {
  259. if (i <= j)
  260. {
  261. float orig = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);
  262. float err = fabsf(test_mat[j +i*size] - orig) / orig;
  263. if (err > 0.00001)
  264. {
  265. FPRINTF(stderr, "Error[%u, %u] --> %2.6f != %2.6f (err %2.6f)\n", i, j, test_mat[j +i*size], orig, err);
  266. assert(0);
  267. }
  268. }
  269. }
  270. }
  271. free(test_mat);
  272. }
  273. starpu_free_flags(mat, (size_t)size*size*sizeof(float), STARPU_MALLOC_PINNED|STARPU_MALLOC_SIMULATION_FOLDED);
  274. }
  275. int main(int argc, char **argv)
  276. {
  277. /* create a simple definite positive symetric matrix example
  278. *
  279. * Hilbert matrix : h(i,j) = 1/(i+j+1)
  280. * */
  281. #ifdef STARPU_HAVE_MAGMA
  282. magma_init();
  283. #endif
  284. int ret;
  285. ret = starpu_init(NULL);
  286. if (ret == -ENODEV) return 77;
  287. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  288. //starpu_fxt_stop_profiling();
  289. init_sizes();
  290. parse_args(argc, argv);
  291. if(with_ctxs_p || with_noctxs_p || chole1_p || chole2_p)
  292. parse_args_ctx(argc, argv);
  293. #ifdef STARPU_USE_CUDA
  294. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,cuda_chol_task_11_cost);
  295. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,cuda_chol_task_21_cost);
  296. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,cuda_chol_task_22_cost);
  297. #else
  298. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,NULL);
  299. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,NULL);
  300. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,NULL);
  301. #endif
  302. starpu_cublas_init();
  303. if(with_ctxs_p)
  304. {
  305. construct_contexts(execute_cholesky);
  306. start_2benchs(execute_cholesky);
  307. }
  308. else if(with_noctxs_p)
  309. start_2benchs(execute_cholesky);
  310. else if(chole1_p)
  311. start_1stbench(execute_cholesky);
  312. else if(chole2_p)
  313. start_2ndbench(execute_cholesky);
  314. else
  315. execute_cholesky(size_p, nblocks_p);
  316. starpu_cublas_shutdown();
  317. starpu_shutdown();
  318. return 0;
  319. }