cholesky_implicit.c 7.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009, 2010, 2011 Université de Bordeaux 1
  4. * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
  5. * Copyright (C) 2010, 2011 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. /*
  20. * Create the codelets
  21. */
  22. static struct starpu_codelet cl11 =
  23. {
  24. .where = STARPU_CPU|STARPU_CUDA,
  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. #endif
  30. .nbuffers = 1,
  31. .model = &chol_model_11
  32. };
  33. static struct starpu_codelet cl21 =
  34. {
  35. .where = STARPU_CPU|STARPU_CUDA,
  36. .type = STARPU_SEQ,
  37. .cpu_funcs = {chol_cpu_codelet_update_u21, NULL},
  38. #ifdef STARPU_USE_CUDA
  39. .cuda_funcs = {chol_cublas_codelet_update_u21, NULL},
  40. #endif
  41. .nbuffers = 2,
  42. .model = &chol_model_21
  43. };
  44. static struct starpu_codelet cl22 =
  45. {
  46. .where = STARPU_CPU|STARPU_CUDA,
  47. .type = STARPU_SEQ,
  48. .max_parallelism = INT_MAX,
  49. .cpu_funcs = {chol_cpu_codelet_update_u22, NULL},
  50. #ifdef STARPU_USE_CUDA
  51. .cuda_funcs = {chol_cublas_codelet_update_u22, NULL},
  52. #endif
  53. .nbuffers = 3,
  54. .model = &chol_model_22
  55. };
  56. /*
  57. * code to bootstrap the factorization
  58. * and construct the DAG
  59. */
  60. static void callback_turn_spmd_on(void *arg __attribute__ ((unused)))
  61. {
  62. cl22.type = STARPU_SPMD;
  63. }
  64. static void _cholesky(starpu_data_handle_t dataA, unsigned nblocks)
  65. {
  66. struct timeval start;
  67. struct timeval end;
  68. unsigned i,j,k;
  69. int prio_level = noprio?STARPU_DEFAULT_PRIO:STARPU_MAX_PRIO;
  70. gettimeofday(&start, NULL);
  71. if (bound)
  72. starpu_bound_start(0, 0);
  73. /* create all the DAG nodes */
  74. for (k = 0; k < nblocks; k++)
  75. {
  76. starpu_data_handle_t sdatakk = starpu_data_get_sub_data(dataA, 2, k, k);
  77. starpu_insert_task(&cl11,
  78. STARPU_PRIORITY, prio_level,
  79. STARPU_RW, sdatakk,
  80. STARPU_CALLBACK, (k == 3*nblocks/4)?callback_turn_spmd_on:NULL,
  81. 0);
  82. for (j = k+1; j<nblocks; j++)
  83. {
  84. starpu_data_handle_t sdatakj = starpu_data_get_sub_data(dataA, 2, k, j);
  85. starpu_insert_task(&cl21,
  86. STARPU_PRIORITY, (j == k+1)?prio_level:STARPU_DEFAULT_PRIO,
  87. STARPU_R, sdatakk,
  88. STARPU_RW, sdatakj,
  89. 0);
  90. for (i = k+1; i<nblocks; i++)
  91. {
  92. if (i <= j)
  93. {
  94. starpu_data_handle_t sdataki = starpu_data_get_sub_data(dataA, 2, k, i);
  95. starpu_data_handle_t sdataij = starpu_data_get_sub_data(dataA, 2, i, j);
  96. starpu_insert_task(&cl22,
  97. STARPU_PRIORITY, ((i == k+1) && (j == k+1))?prio_level:STARPU_DEFAULT_PRIO,
  98. STARPU_R, sdataki,
  99. STARPU_R, sdatakj,
  100. STARPU_RW, sdataij,
  101. 0);
  102. }
  103. }
  104. }
  105. }
  106. starpu_task_wait_for_all();
  107. if (bound)
  108. starpu_bound_stop();
  109. starpu_data_unpartition(dataA, 0);
  110. gettimeofday(&end, NULL);
  111. double timing = (double)((end.tv_sec - start.tv_sec)*1000000 + (end.tv_usec - start.tv_usec));
  112. FPRINTF(stderr, "Computation took (in ms)\n");
  113. FPRINTF(stdout, "%2.2f\n", timing/1000);
  114. unsigned long n = starpu_matrix_get_nx(dataA);
  115. double flop = (1.0f*n*n*n)/3.0f;
  116. FPRINTF(stderr, "Synthetic GFlops : %2.2f\n", (flop/timing/1000.0f));
  117. if (bound) {
  118. double res;
  119. starpu_bound_compute(&res, NULL, 0);
  120. FPRINTF(stderr, "Theoretical GFlops: %2.2f\n", (flop/res/1000000.0f));
  121. }
  122. }
  123. static void 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, 0, (uintptr_t)matA, ld, size, size, sizeof(float));
  129. struct starpu_data_filter f = {
  130. .filter_func = starpu_vertical_block_filter_func,
  131. .nchildren = nblocks
  132. };
  133. struct starpu_data_filter f2 = {
  134. .filter_func = starpu_block_filter_func,
  135. .nchildren = nblocks
  136. };
  137. starpu_data_map_filters(dataA, 2, &f, &f2);
  138. _cholesky(dataA, nblocks);
  139. starpu_data_unregister(dataA);
  140. }
  141. int main(int argc, char **argv)
  142. {
  143. /* create a simple definite positive symetric matrix example
  144. *
  145. * Hilbert matrix : h(i,j) = 1/(i+j+1)
  146. * */
  147. parse_args(argc, argv);
  148. starpu_init(NULL);
  149. starpu_helper_cublas_init();
  150. float *mat;
  151. starpu_malloc((void **)&mat, (size_t)size*size*sizeof(float));
  152. unsigned i,j;
  153. for (i = 0; i < size; i++)
  154. {
  155. for (j = 0; j < size; j++)
  156. {
  157. mat[j +i*size] = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);
  158. /* mat[j +i*size] = ((i == j)?1.0f*size:0.0f); */
  159. }
  160. }
  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. FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);
  170. }
  171. else {
  172. FPRINTF(stdout, ".\t");
  173. }
  174. }
  175. FPRINTF(stdout, "\n");
  176. }
  177. #endif
  178. cholesky(mat, size, size, nblocks);
  179. #ifdef PRINT_OUTPUT
  180. FPRINTF(stdout, "Results :\n");
  181. for (j = 0; j < size; j++)
  182. {
  183. for (i = 0; i < size; i++)
  184. {
  185. if (i <= j) {
  186. FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);
  187. }
  188. else {
  189. FPRINTF(stdout, ".\t");
  190. mat[j+i*size] = 0.0f; /* debug */
  191. }
  192. }
  193. FPRINTF(stdout, "\n");
  194. }
  195. #endif
  196. if (check)
  197. {
  198. FPRINTF(stderr, "compute explicit LLt ...\n");
  199. for (j = 0; j < size; j++)
  200. {
  201. for (i = 0; i < size; i++)
  202. {
  203. if (i > j) {
  204. mat[j+i*size] = 0.0f; /* debug */
  205. }
  206. }
  207. }
  208. float *test_mat = malloc(size*size*sizeof(float));
  209. STARPU_ASSERT(test_mat);
  210. SSYRK("L", "N", size, size, 1.0f,
  211. mat, size, 0.0f, test_mat, size);
  212. FPRINTF(stderr, "comparing results ...\n");
  213. #ifdef PRINT_OUTPUT
  214. for (j = 0; j < size; j++)
  215. {
  216. for (i = 0; i < size; i++)
  217. {
  218. if (i <= j) {
  219. FPRINTF(stdout, "%2.2f\t", test_mat[j +i*size]);
  220. }
  221. else {
  222. FPRINTF(stdout, ".\t");
  223. }
  224. }
  225. FPRINTF(stdout, "\n");
  226. }
  227. #endif
  228. for (j = 0; j < size; j++)
  229. {
  230. for (i = 0; i < size; i++)
  231. {
  232. if (i <= j) {
  233. float orig = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);
  234. float err = abs(test_mat[j +i*size] - orig);
  235. if (err > 0.00001) {
  236. FPRINTF(stderr, "Error[%u, %u] --> %2.2f != %2.2f (err %2.2f)\n", i, j, test_mat[j +i*size], orig, err);
  237. assert(0);
  238. }
  239. }
  240. }
  241. }
  242. }
  243. starpu_helper_cublas_shutdown();
  244. starpu_shutdown();
  245. return 0;
  246. }