cholesky_implicit.c 8.6 KB

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