cholesky_implicit.c 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2011-2013,2015 Inria
  4. * Copyright (C) 2009-2017,2020 Université de Bordeaux
  5. * Copyright (C) 2010-2013,2015-2017 CNRS
  6. * Copyright (C) 2013 Thibaut Lambert
  7. * Copyright (C) 2010 Mehdi Juhoor
  8. *
  9. * StarPU is free software; you can redistribute it and/or modify
  10. * it under the terms of the GNU Lesser General Public License as published by
  11. * the Free Software Foundation; either version 2.1 of the License, or (at
  12. * your option) any later version.
  13. *
  14. * StarPU is distributed in the hope that it will be useful, but
  15. * WITHOUT ANY WARRANTY; without even the implied warranty of
  16. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  17. *
  18. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  19. */
  20. /*
  21. * This version of the Cholesky factorization uses implicit dependency computation.
  22. * The whole algorithm thus appears clearly in the task submission loop in _cholesky().
  23. */
  24. /* Note: this is using fortran ordering, i.e. column-major ordering, i.e.
  25. * elements with consecutive row number are consecutive in memory */
  26. #include "cholesky.h"
  27. #include "../sched_ctx_utils/sched_ctx_utils.h"
  28. #if defined(STARPU_USE_CUDA) && defined(STARPU_HAVE_MAGMA)
  29. #include "magma.h"
  30. #endif
  31. /*
  32. * code to bootstrap the factorization
  33. * and construct the DAG
  34. */
  35. static void callback_turn_spmd_on(void *arg)
  36. {
  37. (void)arg;
  38. cl22.type = STARPU_SPMD;
  39. }
  40. static int _cholesky(starpu_data_handle_t dataA, unsigned nblocks)
  41. {
  42. double start;
  43. double end;
  44. unsigned k,m,n;
  45. unsigned long nx = starpu_matrix_get_nx(dataA);
  46. unsigned long nn = nx/nblocks;
  47. unsigned unbound_prio = STARPU_MAX_PRIO == INT_MAX && STARPU_MIN_PRIO == INT_MIN;
  48. if (bound_p || bound_lp_p || bound_mps_p)
  49. starpu_bound_start(bound_deps_p, 0);
  50. starpu_fxt_start_profiling();
  51. start = starpu_timing_now();
  52. /* create all the DAG nodes */
  53. for (k = 0; k < nblocks; k++)
  54. {
  55. int ret;
  56. starpu_iteration_push(k);
  57. starpu_data_handle_t sdatakk = starpu_data_get_sub_data(dataA, 2, k, k);
  58. ret = starpu_task_insert(&cl11,
  59. STARPU_PRIORITY, noprio_p ? STARPU_DEFAULT_PRIO : unbound_prio ? (int)(2*nblocks - 2*k) : STARPU_MAX_PRIO,
  60. STARPU_RW, sdatakk,
  61. STARPU_CALLBACK, (k == 3*nblocks/4)?callback_turn_spmd_on:NULL,
  62. STARPU_FLOPS, (double) FLOPS_SPOTRF(nn),
  63. STARPU_TAG_ONLY, TAG11(k),
  64. 0);
  65. if (ret == -ENODEV) return 77;
  66. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_insert");
  67. for (m = k+1; m<nblocks; m++)
  68. {
  69. starpu_data_handle_t sdatamk = starpu_data_get_sub_data(dataA, 2, m, k);
  70. ret = starpu_task_insert(&cl21,
  71. STARPU_PRIORITY, noprio_p ? STARPU_DEFAULT_PRIO : unbound_prio ? (int)(2*nblocks - 2*k - m) : (m == k+1)?STARPU_MAX_PRIO:STARPU_DEFAULT_PRIO,
  72. STARPU_R, sdatakk,
  73. STARPU_RW, sdatamk,
  74. STARPU_FLOPS, (double) FLOPS_STRSM(nn, nn),
  75. STARPU_TAG_ONLY, TAG21(m,k),
  76. 0);
  77. if (ret == -ENODEV) return 77;
  78. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_insert");
  79. }
  80. starpu_data_wont_use(sdatakk);
  81. for (m = k+1; m<nblocks; m++)
  82. {
  83. starpu_data_handle_t sdatamk = starpu_data_get_sub_data(dataA, 2, m, k);
  84. for (n = k+1; n<nblocks; n++)
  85. {
  86. if (n <= m)
  87. {
  88. starpu_data_handle_t sdatank = starpu_data_get_sub_data(dataA, 2, n, k);
  89. starpu_data_handle_t sdatamn = starpu_data_get_sub_data(dataA, 2, m, n);
  90. ret = starpu_task_insert(&cl22,
  91. STARPU_PRIORITY, noprio_p ? STARPU_DEFAULT_PRIO : unbound_prio ? (int)(2*nblocks - 2*k - m - n) : ((n == k+1) && (m == k+1))?STARPU_MAX_PRIO:STARPU_DEFAULT_PRIO,
  92. STARPU_R, sdatamk,
  93. STARPU_R, sdatank,
  94. cl22.modes[2], sdatamn,
  95. STARPU_FLOPS, (double) FLOPS_SGEMM(nn, nn, nn),
  96. STARPU_TAG_ONLY, TAG22(k,m,n),
  97. 0);
  98. if (ret == -ENODEV) return 77;
  99. STARPU_CHECK_RETURN_VALUE(ret, "starpu_task_insert");
  100. }
  101. }
  102. starpu_data_wont_use(sdatamk);
  103. }
  104. starpu_iteration_pop();
  105. }
  106. starpu_task_wait_for_all();
  107. end = starpu_timing_now();
  108. starpu_fxt_stop_profiling();
  109. if (bound_p || bound_lp_p || bound_mps_p)
  110. starpu_bound_stop();
  111. double timing = end - start;
  112. double flop = FLOPS_SPOTRF(nx);
  113. if(with_ctxs_p || with_noctxs_p || chole1_p || chole2_p)
  114. update_sched_ctx_timing_results((flop/timing/1000.0f), (timing/1000000.0f));
  115. else
  116. {
  117. PRINTF("# size\tms\tGFlops");
  118. if (bound_p)
  119. PRINTF("\tTms\tTGFlops");
  120. PRINTF("\n");
  121. PRINTF("%lu\t%.0f\t%.1f", nx, timing/1000, (flop/timing/1000.0f));
  122. if (bound_lp_p)
  123. {
  124. FILE *f = fopen("cholesky.lp", "w");
  125. starpu_bound_print_lp(f);
  126. fclose(f);
  127. }
  128. if (bound_mps_p)
  129. {
  130. FILE *f = fopen("cholesky.mps", "w");
  131. starpu_bound_print_mps(f);
  132. fclose(f);
  133. }
  134. if (bound_p)
  135. {
  136. double res;
  137. starpu_bound_compute(&res, NULL, 0);
  138. PRINTF("\t%.0f\t%.1f", res, (flop/res/1000000.0f));
  139. }
  140. PRINTF("\n");
  141. }
  142. return 0;
  143. }
  144. static int cholesky(float *matA, unsigned size, unsigned ld, unsigned nblocks)
  145. {
  146. starpu_data_handle_t dataA;
  147. unsigned m, n;
  148. /* monitor and partition the A matrix into blocks :
  149. * one block is now determined by 2 unsigned (m,n) */
  150. starpu_matrix_data_register(&dataA, STARPU_MAIN_RAM, (uintptr_t)matA, ld, size, size, sizeof(float));
  151. /* Split into blocks of complete rows first */
  152. struct starpu_data_filter f =
  153. {
  154. .filter_func = starpu_matrix_filter_block,
  155. .nchildren = nblocks
  156. };
  157. /* Then split rows into tiles */
  158. struct starpu_data_filter f2 =
  159. {
  160. /* Note: here "vertical" is for row-major, we are here using column-major. */
  161. .filter_func = starpu_matrix_filter_vertical_block,
  162. .nchildren = nblocks
  163. };
  164. starpu_data_map_filters(dataA, 2, &f, &f2);
  165. for (m = 0; m < nblocks; m++)
  166. for (n = 0; n < nblocks; n++)
  167. {
  168. starpu_data_handle_t data = starpu_data_get_sub_data(dataA, 2, m, n);
  169. starpu_data_set_coordinates(data, 2, m, n);
  170. }
  171. int ret = _cholesky(dataA, nblocks);
  172. starpu_data_unpartition(dataA, STARPU_MAIN_RAM);
  173. starpu_data_unregister(dataA);
  174. return ret;
  175. }
  176. static void execute_cholesky(unsigned size, unsigned nblocks)
  177. {
  178. float *mat = NULL;
  179. #ifndef STARPU_SIMGRID
  180. unsigned m,n;
  181. starpu_malloc_flags((void **)&mat, (size_t)size*size*sizeof(float), STARPU_MALLOC_PINNED|STARPU_MALLOC_SIMULATION_FOLDED);
  182. for (n = 0; n < size; n++)
  183. {
  184. for (m = 0; m < size; m++)
  185. {
  186. mat[m +n*size] = (1.0f/(1.0f+m+n)) + ((m == n)?1.0f*size:0.0f);
  187. /* mat[m +n*size] = ((m == n)?1.0f*size:0.0f); */
  188. }
  189. }
  190. /* #define PRINT_OUTPUT */
  191. #ifdef PRINT_OUTPUT
  192. FPRINTF(stdout, "Input :\n");
  193. for (m = 0; m < size; m++)
  194. {
  195. for (n = 0; n < size; n++)
  196. {
  197. if (n <= m)
  198. {
  199. FPRINTF(stdout, "%2.2f\t", mat[m +n*size]);
  200. }
  201. else
  202. {
  203. FPRINTF(stdout, ".\t");
  204. }
  205. }
  206. FPRINTF(stdout, "\n");
  207. }
  208. #endif
  209. #endif
  210. cholesky(mat, size, size, nblocks);
  211. #ifndef STARPU_SIMGRID
  212. #ifdef PRINT_OUTPUT
  213. FPRINTF(stdout, "Results :\n");
  214. for (m = 0; m < size; m++)
  215. {
  216. for (n = 0; n < size; n++)
  217. {
  218. if (n <= m)
  219. {
  220. FPRINTF(stdout, "%2.2f\t", mat[m +n*size]);
  221. }
  222. else
  223. {
  224. FPRINTF(stdout, ".\t");
  225. }
  226. }
  227. FPRINTF(stdout, "\n");
  228. }
  229. #endif
  230. if (check_p)
  231. {
  232. FPRINTF(stderr, "compute explicit LLt ...\n");
  233. for (m = 0; m < size; m++)
  234. {
  235. for (n = 0; n < size; n++)
  236. {
  237. if (n > m)
  238. {
  239. mat[m+n*size] = 0.0f; /* debug */
  240. }
  241. }
  242. }
  243. float *test_mat = malloc(size*size*sizeof(float));
  244. STARPU_ASSERT(test_mat);
  245. STARPU_SSYRK("L", "N", size, size, 1.0f,
  246. mat, size, 0.0f, test_mat, size);
  247. FPRINTF(stderr, "comparing results ...\n");
  248. #ifdef PRINT_OUTPUT
  249. for (m = 0; m < size; m++)
  250. {
  251. for (n = 0; n < size; n++)
  252. {
  253. if (n <= m)
  254. {
  255. FPRINTF(stdout, "%2.2f\t", test_mat[m +n*size]);
  256. }
  257. else
  258. {
  259. FPRINTF(stdout, ".\t");
  260. }
  261. }
  262. FPRINTF(stdout, "\n");
  263. }
  264. #endif
  265. for (m = 0; m < size; m++)
  266. {
  267. for (n = 0; n < size; n++)
  268. {
  269. if (n <= m)
  270. {
  271. float orig = (1.0f/(1.0f+m+n)) + ((m == n)?1.0f*size:0.0f);
  272. float err = fabsf(test_mat[m +n*size] - orig) / orig;
  273. if (err > 0.0001)
  274. {
  275. FPRINTF(stderr, "Error[%u, %u] --> %2.6f != %2.6f (err %2.6f)\n", m, n, test_mat[m +n*size], orig, err);
  276. assert(0);
  277. }
  278. }
  279. }
  280. }
  281. free(test_mat);
  282. }
  283. starpu_free_flags(mat, (size_t)size*size*sizeof(float), STARPU_MALLOC_PINNED|STARPU_MALLOC_SIMULATION_FOLDED);
  284. #endif
  285. }
  286. int main(int argc, char **argv)
  287. {
  288. /* create a simple definite positive symetric matrix example
  289. *
  290. * Hilbert matrix : h(i,j) = 1/(i+j+1)
  291. * */
  292. #ifdef STARPU_HAVE_MAGMA
  293. magma_init();
  294. #endif
  295. int ret;
  296. ret = starpu_init(NULL);
  297. if (ret == -ENODEV) return 77;
  298. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  299. //starpu_fxt_stop_profiling();
  300. init_sizes();
  301. parse_args(argc, argv);
  302. if(with_ctxs_p || with_noctxs_p || chole1_p || chole2_p)
  303. parse_args_ctx(argc, argv);
  304. #ifdef STARPU_USE_CUDA
  305. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,cuda_chol_task_11_cost);
  306. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,cuda_chol_task_21_cost);
  307. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,cuda_chol_task_22_cost);
  308. #else
  309. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,NULL);
  310. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,NULL);
  311. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,NULL);
  312. #endif
  313. starpu_cublas_init();
  314. if(with_ctxs_p)
  315. {
  316. construct_contexts();
  317. start_2benchs(execute_cholesky);
  318. }
  319. else if(with_noctxs_p)
  320. start_2benchs(execute_cholesky);
  321. else if(chole1_p)
  322. start_1stbench(execute_cholesky);
  323. else if(chole2_p)
  324. start_2ndbench(execute_cholesky);
  325. else
  326. execute_cholesky(size_p, nblocks_p);
  327. starpu_cublas_shutdown();
  328. starpu_shutdown();
  329. return 0;
  330. }