plu_example.c 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474
  1. /*
  2. * StarPU
  3. * Copyright (C) INRIA 2008-2010 (see AUTHORS file)
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU Lesser General Public License as published by
  7. * the Free Software Foundation; either version 2.1 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  13. *
  14. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  15. */
  16. #include <stdlib.h>
  17. #include <stdio.h>
  18. #include <string.h>
  19. #include <time.h>
  20. #include <math.h>
  21. #include <starpu.h>
  22. #include "pxlu.h"
  23. //#include "pxlu_kernels.h"
  24. static unsigned long size = 16384;
  25. static unsigned nblocks = 16;
  26. static unsigned check = 0;
  27. static unsigned p = 1;
  28. static unsigned q = 1;
  29. static unsigned display = 0;
  30. static starpu_data_handle *dataA_handles;
  31. static TYPE **dataA;
  32. /* In order to implement the distributed LU decomposition, we allocate
  33. * temporary buffers */
  34. #ifdef SINGLE_TMP11
  35. static starpu_data_handle tmp_11_block_handle;
  36. static TYPE *tmp_11_block;
  37. #else
  38. static starpu_data_handle *tmp_11_block_handles;
  39. static TYPE **tmp_11_block;
  40. #endif
  41. #ifdef SINGLE_TMP1221
  42. static starpu_data_handle *tmp_12_block_handles;
  43. static TYPE **tmp_12_block;
  44. static starpu_data_handle *tmp_21_block_handles;
  45. static TYPE **tmp_21_block;
  46. #else
  47. static starpu_data_handle *(tmp_12_block_handles[2]);
  48. static TYPE **(tmp_12_block[2]);
  49. static starpu_data_handle *(tmp_21_block_handles[2]);
  50. static TYPE **(tmp_21_block[2]);
  51. #endif
  52. static void parse_args(int argc, char **argv)
  53. {
  54. int i;
  55. for (i = 1; i < argc; i++) {
  56. if (strcmp(argv[i], "-size") == 0) {
  57. char *argptr;
  58. size = strtol(argv[++i], &argptr, 10);
  59. }
  60. if (strcmp(argv[i], "-nblocks") == 0) {
  61. char *argptr;
  62. nblocks = strtol(argv[++i], &argptr, 10);
  63. }
  64. if (strcmp(argv[i], "-check") == 0) {
  65. check = 1;
  66. }
  67. if (strcmp(argv[i], "-display") == 0) {
  68. display = 1;
  69. }
  70. if (strcmp(argv[i], "-p") == 0) {
  71. char *argptr;
  72. p = strtol(argv[++i], &argptr, 10);
  73. }
  74. if (strcmp(argv[i], "-q") == 0) {
  75. char *argptr;
  76. q = strtol(argv[++i], &argptr, 10);
  77. }
  78. }
  79. }
  80. unsigned STARPU_PLU(display_flag)(void)
  81. {
  82. return display;
  83. }
  84. static void fill_block_with_random(TYPE *blockptr, unsigned size, unsigned nblocks)
  85. {
  86. const unsigned block_size = (size/nblocks);
  87. unsigned i, j;
  88. for (i = 0; i < block_size; i++)
  89. for (j = 0; j < block_size; j++)
  90. {
  91. blockptr[j+i*block_size] = (TYPE)drand48();
  92. }
  93. }
  94. #ifdef SINGLE_TMP11
  95. starpu_data_handle STARPU_PLU(get_tmp_11_block_handle)(void)
  96. {
  97. return tmp_11_block_handle;
  98. }
  99. #else
  100. starpu_data_handle STARPU_PLU(get_tmp_11_block_handle)(unsigned k)
  101. {
  102. return tmp_11_block_handles[k];
  103. }
  104. #endif
  105. #ifdef SINGLE_TMP1221
  106. starpu_data_handle STARPU_PLU(get_tmp_12_block_handle)(unsigned j)
  107. {
  108. return tmp_12_block_handles[j];
  109. }
  110. starpu_data_handle STARPU_PLU(get_tmp_21_block_handle)(unsigned i)
  111. {
  112. return tmp_21_block_handles[i];
  113. }
  114. #else
  115. starpu_data_handle STARPU_PLU(get_tmp_12_block_handle)(unsigned j, unsigned k)
  116. {
  117. return tmp_12_block_handles[k%2][j];
  118. }
  119. starpu_data_handle STARPU_PLU(get_tmp_21_block_handle)(unsigned i, unsigned k)
  120. {
  121. return tmp_21_block_handles[k%2][i];
  122. }
  123. #endif
  124. static void init_matrix(int rank)
  125. {
  126. /* Allocate a grid of data handles, not all of them have to be allocated later on */
  127. dataA_handles = calloc(nblocks*nblocks, sizeof(starpu_data_handle));
  128. dataA = calloc(nblocks*nblocks, sizeof(TYPE *));
  129. size_t blocksize = (size_t)(size/nblocks)*(size/nblocks)*sizeof(TYPE);
  130. /* Allocate all the blocks that belong to this mpi node */
  131. unsigned long i,j;
  132. for (j = 0; j < nblocks; j++)
  133. {
  134. for (i = 0; i < nblocks; i++)
  135. {
  136. TYPE **blockptr = &dataA[j+i*nblocks];
  137. // starpu_data_handle *handleptr = &dataA_handles[j+nblocks*i];
  138. starpu_data_handle *handleptr = &dataA_handles[j+nblocks*i];
  139. if (get_block_rank(i, j) == rank)
  140. {
  141. /* This blocks should be treated by the current MPI process */
  142. /* Allocate and fill it */
  143. starpu_malloc_pinned_if_possible((void **)blockptr, blocksize);
  144. //fprintf(stderr, "Rank %d : fill block (i = %d, j = %d)\n", rank, i, j);
  145. fill_block_with_random(*blockptr, size, nblocks);
  146. //fprintf(stderr, "Rank %d : fill block (i = %d, j = %d)\n", rank, i, j);
  147. if (i == j)
  148. {
  149. unsigned tmp;
  150. for (tmp = 0; tmp < size/nblocks; tmp++)
  151. {
  152. (*blockptr)[tmp*((size/nblocks)+1)] += (TYPE)10*nblocks;
  153. }
  154. }
  155. /* Register it to StarPU */
  156. starpu_register_blas_data(handleptr, 0,
  157. (uintptr_t)*blockptr, size/nblocks,
  158. size/nblocks, size/nblocks, sizeof(TYPE));
  159. }
  160. else {
  161. *blockptr = STARPU_POISON_PTR;
  162. *handleptr = STARPU_POISON_PTR;
  163. }
  164. }
  165. }
  166. /* Allocate the temporary buffers required for the distributed algorithm */
  167. unsigned k;
  168. /* tmp buffer 11 */
  169. #ifdef SINGLE_TMP11
  170. starpu_malloc_pinned_if_possible((void **)&tmp_11_block, blocksize);
  171. starpu_register_blas_data(&tmp_11_block_handle, 0, (uintptr_t)tmp_11_block,
  172. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  173. #else
  174. tmp_11_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  175. tmp_11_block = calloc(nblocks, sizeof(TYPE *));
  176. for (k = 0; k < nblocks; k++)
  177. {
  178. starpu_malloc_pinned_if_possible((void **)&tmp_11_block[k], blocksize);
  179. STARPU_ASSERT(tmp_11_block[k]);
  180. starpu_register_blas_data(&tmp_11_block_handles[k], 0,
  181. (uintptr_t)tmp_11_block[k],
  182. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  183. }
  184. #endif
  185. /* tmp buffers 12 and 21 */
  186. #ifdef SINGLE_TMP1221
  187. tmp_12_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  188. tmp_21_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  189. tmp_12_block = calloc(nblocks, sizeof(TYPE *));
  190. tmp_21_block = calloc(nblocks, sizeof(TYPE *));
  191. #else
  192. for (i = 0; i < 2; i++) {
  193. tmp_12_block_handles[i] = calloc(nblocks, sizeof(starpu_data_handle));
  194. tmp_21_block_handles[i] = calloc(nblocks, sizeof(starpu_data_handle));
  195. tmp_12_block[i] = calloc(nblocks, sizeof(TYPE *));
  196. tmp_21_block[i] = calloc(nblocks, sizeof(TYPE *));
  197. }
  198. #endif
  199. for (k = 0; k < nblocks; k++)
  200. {
  201. #ifdef SINGLE_TMP1221
  202. starpu_malloc_pinned_if_possible((void **)&tmp_12_block[k], blocksize);
  203. STARPU_ASSERT(tmp_12_block[k]);
  204. starpu_register_blas_data(&tmp_12_block_handles[k], 0,
  205. (uintptr_t)tmp_12_block[k],
  206. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  207. starpu_malloc_pinned_if_possible((void **)&tmp_21_block[k], blocksize);
  208. STARPU_ASSERT(tmp_21_block[k]);
  209. starpu_register_blas_data(&tmp_21_block_handles[k], 0,
  210. (uintptr_t)tmp_21_block[k],
  211. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  212. #else
  213. for (i = 0; i < 2; i++) {
  214. starpu_malloc_pinned_if_possible((void **)&tmp_12_block[i][k], blocksize);
  215. STARPU_ASSERT(tmp_12_block[i][k]);
  216. starpu_register_blas_data(&tmp_12_block_handles[i][k], 0,
  217. (uintptr_t)tmp_12_block[i][k],
  218. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  219. starpu_malloc_pinned_if_possible((void **)&tmp_21_block[i][k], blocksize);
  220. STARPU_ASSERT(tmp_21_block[i][k]);
  221. starpu_register_blas_data(&tmp_21_block_handles[i][k], 0,
  222. (uintptr_t)tmp_21_block[i][k],
  223. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  224. }
  225. #endif
  226. }
  227. //display_all_blocks(nblocks, size/nblocks);
  228. }
  229. TYPE *STARPU_PLU(get_block)(unsigned i, unsigned j)
  230. {
  231. return dataA[j+i*nblocks];
  232. }
  233. int get_block_rank(unsigned i, unsigned j)
  234. {
  235. /* Take a 2D block cyclic distribution */
  236. /* NB: p (resp. q) is for "direction" i (resp. j) */
  237. return (j % q) * p + (i % p);
  238. }
  239. starpu_data_handle STARPU_PLU(get_block_handle)(unsigned i, unsigned j)
  240. {
  241. return dataA_handles[j+i*nblocks];
  242. }
  243. static void display_grid(int rank, unsigned nblocks)
  244. {
  245. if (!display)
  246. return;
  247. //if (rank == 0)
  248. {
  249. fprintf(stderr, "2D grid layout (Rank %d): \n", rank);
  250. unsigned i, j;
  251. for (j = 0; j < nblocks; j++)
  252. {
  253. for (i = 0; i < nblocks; i++)
  254. {
  255. TYPE *blockptr = STARPU_PLU(get_block)(i, j);
  256. starpu_data_handle handle = STARPU_PLU(get_block_handle)(i, j);
  257. fprintf(stderr, "%d (data %p handle %p)", get_block_rank(i, j), blockptr, handle);
  258. }
  259. fprintf(stderr, "\n");
  260. }
  261. }
  262. }
  263. int main(int argc, char **argv)
  264. {
  265. int rank;
  266. int world_size;
  267. int thread_support;
  268. /*
  269. * Initialization
  270. */
  271. if (MPI_Init_thread(&argc, &argv, MPI_THREAD_SERIALIZED, &thread_support) != MPI_SUCCESS) {
  272. fprintf(stderr,"MPI_Init_thread failed\n");
  273. exit(1);
  274. }
  275. if (thread_support == MPI_THREAD_FUNNELED)
  276. fprintf(stderr,"Warning: MPI only has funneled thread support, not serialized, hoping this will work\n");
  277. if (thread_support < MPI_THREAD_FUNNELED)
  278. fprintf(stderr,"Warning: MPI does not have thread support!\n");
  279. MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  280. MPI_Comm_size(MPI_COMM_WORLD, &world_size);
  281. srand48((long int)time(NULL));
  282. parse_args(argc, argv);
  283. STARPU_ASSERT(p*q == world_size);
  284. starpu_init(NULL);
  285. starpu_mpi_initialize();
  286. starpu_helper_init_cublas();
  287. int barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  288. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  289. /*
  290. * Problem Init
  291. */
  292. init_matrix(rank);
  293. display_grid(rank, nblocks);
  294. TYPE *a_r;
  295. // STARPU_PLU(display_data_content)(a_r, size);
  296. TYPE *x, *y;
  297. if (check)
  298. {
  299. x = calloc(size, sizeof(TYPE));
  300. STARPU_ASSERT(x);
  301. y = calloc(size, sizeof(TYPE));
  302. STARPU_ASSERT(y);
  303. if (rank == 0)
  304. {
  305. unsigned ind;
  306. for (ind = 0; ind < size; ind++)
  307. x[ind] = (TYPE)drand48();
  308. }
  309. a_r = STARPU_PLU(reconstruct_matrix)(size, nblocks);
  310. if (rank == 0)
  311. STARPU_PLU(display_data_content)(a_r, size);
  312. // STARPU_PLU(compute_ax)(size, x, y, nblocks, rank);
  313. }
  314. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  315. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  316. double timing = STARPU_PLU(plu_main)(nblocks, rank, world_size);
  317. /*
  318. * Report performance
  319. */
  320. int reduce_ret;
  321. double min_timing = timing;
  322. double max_timing = timing;
  323. double sum_timing = timing;
  324. reduce_ret = MPI_Reduce(&timing, &min_timing, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD);
  325. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  326. reduce_ret = MPI_Reduce(&timing, &max_timing, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
  327. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  328. reduce_ret = MPI_Reduce(&timing, &sum_timing, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
  329. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  330. if (rank == 0)
  331. {
  332. fprintf(stderr, "Computation took: %lf ms\n", max_timing/1000);
  333. fprintf(stderr, "\tMIN : %lf ms\n", min_timing/1000);
  334. fprintf(stderr, "\tMAX : %lf ms\n", max_timing/1000);
  335. fprintf(stderr, "\tAVG : %lf ms\n", sum_timing/(world_size*1000));
  336. unsigned n = size;
  337. double flop = (2.0f*n*n*n)/3.0f;
  338. fprintf(stderr, "Synthetic GFlops : %2.2f\n", (flop/max_timing/1000.0f));
  339. }
  340. /*
  341. * Test Result Correctness
  342. */
  343. TYPE *y2;
  344. if (check)
  345. {
  346. /*
  347. * Compute || A - LU ||
  348. */
  349. STARPU_PLU(compute_lu_matrix)(size, nblocks, a_r);
  350. #if 0
  351. /*
  352. * Compute || Ax - LUx ||
  353. */
  354. unsigned ind;
  355. y2 = calloc(size, sizeof(TYPE));
  356. STARPU_ASSERT(y);
  357. if (rank == 0)
  358. {
  359. for (ind = 0; ind < size; ind++)
  360. {
  361. y2[ind] = (TYPE)0.0;
  362. }
  363. }
  364. STARPU_PLU(compute_lux)(size, x, y2, nblocks, rank);
  365. /* Compute y2 = y2 - y */
  366. CPU_AXPY(size, -1.0, y, 1, y2, 1);
  367. TYPE err = CPU_ASUM(size, y2, 1);
  368. int max = CPU_IAMAX(size, y2, 1);
  369. fprintf(stderr, "(A - LU)X Avg error : %e\n", err/(size*size));
  370. fprintf(stderr, "(A - LU)X Max error : %e\n", y2[max]);
  371. #endif
  372. }
  373. /*
  374. * Termination
  375. */
  376. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  377. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  378. starpu_helper_shutdown_cublas();
  379. starpu_mpi_shutdown();
  380. starpu_shutdown();
  381. MPI_Finalize();
  382. return 0;
  383. }