cholesky_implicit.c 8.6 KB

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