cholesky_implicit.c 9.7 KB

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