cholesky_implicit.c 8.7 KB

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