cholesky_implicit.c 8.7 KB

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