cholesky_implicit.c 9.7 KB

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