cholesky_implicit.c 9.6 KB

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