cholesky_implicit.c 8.6 KB

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