cholesky_implicit.c 8.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009-2016 Université de Bordeaux
  4. * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
  5. * Copyright (C) 2010, 2011, 2012, 2013, 2016 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. int ret;
  38. double start;
  39. double end;
  40. unsigned i,j,k;
  41. unsigned long n = starpu_matrix_get_nx(dataA);
  42. unsigned long nn = n/nblocks;
  43. int prio_level = noprio?STARPU_DEFAULT_PRIO:STARPU_MAX_PRIO;
  44. if (bound || bound_lp || bound_mps)
  45. starpu_bound_start(bound_deps, 0);
  46. starpu_fxt_start_profiling();
  47. start = starpu_timing_now();
  48. /* create all the DAG nodes */
  49. for (k = 0; k < nblocks; k++)
  50. {
  51. starpu_data_handle_t sdatakk = starpu_data_get_sub_data(dataA, 2, k, k);
  52. ret = starpu_task_insert(&cl11,
  53. STARPU_PRIORITY, prio_level,
  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, (j == k+1)?prio_level: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, ((i == k+1) && (j == k+1))?prio_level: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 || bound_lp || bound_mps)
  103. starpu_bound_stop();
  104. double timing = end - start;
  105. double flop = FLOPS_SPOTRF(n);
  106. if(with_ctxs || with_noctxs || chole1 || chole2)
  107. update_sched_ctx_timing_results((flop/timing/1000.0f), (timing/1000000.0f));
  108. else
  109. {
  110. PRINTF("# size\tms\tGFlops");
  111. if (bound)
  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)
  116. {
  117. FILE *f = fopen("cholesky.lp", "w");
  118. starpu_bound_print_lp(f);
  119. }
  120. if (bound_mps)
  121. {
  122. FILE *f = fopen("cholesky.mps", "w");
  123. starpu_bound_print_mps(f);
  124. }
  125. if (bound)
  126. {
  127. double res;
  128. starpu_bound_compute(&res, NULL, 0);
  129. PRINTF("\t%.0f\t%.1f", res, (flop/res/1000000.0f));
  130. }
  131. PRINTF("\n");
  132. }
  133. return 0;
  134. }
  135. static int cholesky(float *matA, unsigned size, unsigned ld, unsigned nblocks)
  136. {
  137. starpu_data_handle_t dataA;
  138. /* monitor and partition the A matrix into blocks :
  139. * one block is now determined by 2 unsigned (i,j) */
  140. starpu_matrix_data_register(&dataA, STARPU_MAIN_RAM, (uintptr_t)matA, ld, size, size, sizeof(float));
  141. struct starpu_data_filter f =
  142. {
  143. .filter_func = starpu_matrix_filter_vertical_block,
  144. .nchildren = nblocks
  145. };
  146. struct starpu_data_filter f2 =
  147. {
  148. .filter_func = starpu_matrix_filter_block,
  149. .nchildren = nblocks
  150. };
  151. starpu_data_map_filters(dataA, 2, &f, &f2);
  152. int ret = _cholesky(dataA, nblocks);
  153. starpu_data_unpartition(dataA, STARPU_MAIN_RAM);
  154. starpu_data_unregister(dataA);
  155. return ret;
  156. }
  157. static void execute_cholesky(unsigned size, unsigned nblocks)
  158. {
  159. float *mat = NULL;
  160. unsigned i,j;
  161. starpu_malloc_flags((void **)&mat, (size_t)size*size*sizeof(float), STARPU_MALLOC_PINNED|STARPU_MALLOC_SIMULATION_FOLDED);
  162. #ifndef STARPU_SIMGRID
  163. for (i = 0; i < size; i++)
  164. {
  165. for (j = 0; j < size; j++)
  166. {
  167. mat[j +i*size] = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);
  168. /* mat[j +i*size] = ((i == j)?1.0f*size:0.0f); */
  169. }
  170. }
  171. #endif
  172. /* #define PRINT_OUTPUT */
  173. #ifdef PRINT_OUTPUT
  174. FPRINTF(stdout, "Input :\n");
  175. for (j = 0; j < size; j++)
  176. {
  177. for (i = 0; i < size; i++)
  178. {
  179. if (i <= j)
  180. {
  181. FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);
  182. }
  183. else
  184. {
  185. FPRINTF(stdout, ".\t");
  186. }
  187. }
  188. FPRINTF(stdout, "\n");
  189. }
  190. #endif
  191. cholesky(mat, size, size, nblocks);
  192. #ifdef PRINT_OUTPUT
  193. FPRINTF(stdout, "Results :\n");
  194. for (j = 0; j < size; j++)
  195. {
  196. for (i = 0; i < size; i++)
  197. {
  198. if (i <= j)
  199. {
  200. FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);
  201. }
  202. else
  203. {
  204. FPRINTF(stdout, ".\t");
  205. mat[j+i*size] = 0.0f; /* debug */
  206. }
  207. }
  208. FPRINTF(stdout, "\n");
  209. }
  210. #endif
  211. if (check)
  212. {
  213. FPRINTF(stderr, "compute explicit LLt ...\n");
  214. for (j = 0; j < size; j++)
  215. {
  216. for (i = 0; i < size; i++)
  217. {
  218. if (i > j)
  219. {
  220. mat[j+i*size] = 0.0f; /* debug */
  221. }
  222. }
  223. }
  224. float *test_mat = malloc(size*size*sizeof(float));
  225. STARPU_ASSERT(test_mat);
  226. STARPU_SSYRK("L", "N", size, size, 1.0f,
  227. mat, size, 0.0f, test_mat, size);
  228. FPRINTF(stderr, "comparing results ...\n");
  229. #ifdef PRINT_OUTPUT
  230. for (j = 0; j < size; j++)
  231. {
  232. for (i = 0; i < size; i++)
  233. {
  234. if (i <= j)
  235. {
  236. FPRINTF(stdout, "%2.2f\t", test_mat[j +i*size]);
  237. }
  238. else
  239. {
  240. FPRINTF(stdout, ".\t");
  241. }
  242. }
  243. FPRINTF(stdout, "\n");
  244. }
  245. #endif
  246. for (j = 0; j < size; j++)
  247. {
  248. for (i = 0; i < size; i++)
  249. {
  250. if (i <= j)
  251. {
  252. float orig = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);
  253. float err = fabsf(test_mat[j +i*size] - orig);
  254. if (err > 0.00001)
  255. {
  256. FPRINTF(stderr, "Error[%u, %u] --> %2.2f != %2.2f (err %2.2f)\n", i, j, test_mat[j +i*size], orig, err);
  257. assert(0);
  258. }
  259. }
  260. }
  261. }
  262. free(test_mat);
  263. }
  264. starpu_free_flags(mat, (size_t)size*size*sizeof(float), STARPU_MALLOC_PINNED|STARPU_MALLOC_SIMULATION_FOLDED);
  265. }
  266. int main(int argc, char **argv)
  267. {
  268. /* create a simple definite positive symetric matrix example
  269. *
  270. * Hilbert matrix : h(i,j) = 1/(i+j+1)
  271. * */
  272. parse_args(argc, argv);
  273. if(with_ctxs || with_noctxs || chole1 || chole2)
  274. parse_args_ctx(argc, argv);
  275. #ifdef STARPU_HAVE_MAGMA
  276. magma_init();
  277. #endif
  278. int ret;
  279. ret = starpu_init(NULL);
  280. //starpu_fxt_stop_profiling();
  281. if (ret == -ENODEV) return 77;
  282. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  283. #ifdef STARPU_USE_CUDA
  284. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,cuda_chol_task_11_cost);
  285. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,cuda_chol_task_21_cost);
  286. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,cuda_chol_task_22_cost);
  287. #else
  288. initialize_chol_model(&chol_model_11,"chol_model_11",cpu_chol_task_11_cost,NULL);
  289. initialize_chol_model(&chol_model_21,"chol_model_21",cpu_chol_task_21_cost,NULL);
  290. initialize_chol_model(&chol_model_22,"chol_model_22",cpu_chol_task_22_cost,NULL);
  291. #endif
  292. starpu_cublas_init();
  293. if(with_ctxs)
  294. {
  295. construct_contexts(execute_cholesky);
  296. start_2benchs(execute_cholesky);
  297. }
  298. else if(with_noctxs)
  299. start_2benchs(execute_cholesky);
  300. else if(chole1)
  301. start_1stbench(execute_cholesky);
  302. else if(chole2)
  303. start_2ndbench(execute_cholesky);
  304. else
  305. execute_cholesky(size, nblocks);
  306. starpu_cublas_shutdown();
  307. starpu_shutdown();
  308. return 0;
  309. }