cholesky_implicit.c 8.4 KB

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