plu_example.c 9.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407
  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 starpu_data_handle *dataA_handles;
  30. static TYPE **dataA;
  31. /* In order to implement the distributed LU decomposition, we allocate
  32. * temporary buffers */
  33. static starpu_data_handle tmp_11_block_handle;
  34. static TYPE *tmp_11_block;
  35. static starpu_data_handle *tmp_12_block_handles;
  36. static TYPE **tmp_12_block;
  37. static starpu_data_handle *tmp_21_block_handles;
  38. static TYPE **tmp_21_block;
  39. TYPE *STARPU_PLU(reconstruct_matrix)(unsigned size, unsigned nblocks);
  40. static void display_block_content(unsigned bi, unsigned bj, unsigned blocksize);
  41. static void display_all_blocks(unsigned nblocks, unsigned blocksize)
  42. {
  43. fprintf(stderr, "DISPLAY ALL\n");
  44. unsigned bi, bj;
  45. for (bj = 0; bj < nblocks; bj++)
  46. for (bi = 0; bi < nblocks; bi++)
  47. display_block_content(bi, bj, blocksize);
  48. fprintf(stderr, "*****************\n");
  49. }
  50. static void parse_args(int argc, char **argv)
  51. {
  52. int i;
  53. for (i = 1; i < argc; i++) {
  54. if (strcmp(argv[i], "-size") == 0) {
  55. char *argptr;
  56. size = strtol(argv[++i], &argptr, 10);
  57. }
  58. if (strcmp(argv[i], "-nblocks") == 0) {
  59. char *argptr;
  60. nblocks = strtol(argv[++i], &argptr, 10);
  61. }
  62. if (strcmp(argv[i], "-check") == 0) {
  63. check = 1;
  64. }
  65. if (strcmp(argv[i], "-p") == 0) {
  66. char *argptr;
  67. p = strtol(argv[++i], &argptr, 10);
  68. }
  69. if (strcmp(argv[i], "-q") == 0) {
  70. char *argptr;
  71. q = strtol(argv[++i], &argptr, 10);
  72. }
  73. }
  74. }
  75. static void fill_block_with_random(TYPE *blockptr, unsigned size, unsigned nblocks)
  76. {
  77. const unsigned block_size = (size/nblocks);
  78. unsigned i, j;
  79. for (j = 0; j < block_size; j++)
  80. for (i = 0; i < block_size; i++)
  81. {
  82. // blockptr[i+j*block_size] = (TYPE)drand48();
  83. blockptr[i+j*block_size] = (i == j)?2.0:(TYPE)j;
  84. }
  85. }
  86. starpu_data_handle STARPU_PLU(get_block_handle)(unsigned j, unsigned i)
  87. {
  88. return dataA_handles[i+j*nblocks];
  89. }
  90. starpu_data_handle STARPU_PLU(get_tmp_11_block_handle)(void)
  91. {
  92. return tmp_11_block_handle;
  93. }
  94. starpu_data_handle STARPU_PLU(get_tmp_12_block_handle)(unsigned j)
  95. {
  96. return tmp_12_block_handles[j];
  97. }
  98. starpu_data_handle STARPU_PLU(get_tmp_21_block_handle)(unsigned i)
  99. {
  100. return tmp_21_block_handles[i];
  101. }
  102. TYPE *STARPU_PLU(get_block)(unsigned j, unsigned i)
  103. {
  104. return dataA[i+j*nblocks];
  105. }
  106. static void init_matrix(int rank)
  107. {
  108. /* Allocate a grid of data handles, not all of them have to be allocated later on */
  109. dataA_handles = calloc(nblocks*nblocks, sizeof(starpu_data_handle));
  110. dataA = calloc(nblocks*nblocks, sizeof(TYPE *));
  111. size_t blocksize = (size_t)(size/nblocks)*(size/nblocks)*sizeof(TYPE);
  112. /* Allocate all the blocks that belong to this mpi node */
  113. unsigned long i,j;
  114. for (j = 0; j < nblocks; j++)
  115. {
  116. for (i = 0; i < nblocks; i++)
  117. {
  118. if (get_block_rank(i, j) == rank)
  119. {
  120. /* This blocks should be treated by the current MPI process */
  121. /* Allocate and fill it */
  122. starpu_malloc_pinned_if_possible((void **)&dataA[i+j*nblocks], blocksize);
  123. fill_block_with_random(STARPU_PLU(get_block)(j, i), size, nblocks);
  124. if (i == j)
  125. {
  126. TYPE *b = STARPU_PLU(get_block)(j, i);
  127. unsigned tmp;
  128. for (tmp = 0; tmp < size/nblocks; tmp++)
  129. {
  130. b[tmp*((size/nblocks)+1)] += (TYPE)10*nblocks;
  131. }
  132. }
  133. /* Register it to StarPU */
  134. starpu_register_blas_data(&dataA_handles[i+nblocks*j], 0,
  135. (uintptr_t)dataA[i+nblocks*j], size/nblocks,
  136. size/nblocks, size/nblocks, sizeof(TYPE));
  137. }
  138. else {
  139. dataA[i+j*nblocks] = STARPU_POISON_PTR;
  140. dataA_handles[i+j*nblocks] = STARPU_POISON_PTR;
  141. }
  142. }
  143. }
  144. /* Allocate the temporary buffers required for the distributed algorithm */
  145. /* tmp buffer 11 */
  146. starpu_malloc_pinned_if_possible((void **)&tmp_11_block, blocksize);
  147. starpu_register_blas_data(&tmp_11_block_handle, 0, (uintptr_t)tmp_11_block,
  148. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  149. /* tmp buffers 12 and 21 */
  150. tmp_12_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  151. tmp_21_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  152. tmp_12_block = calloc(nblocks, sizeof(TYPE *));
  153. tmp_21_block = calloc(nblocks, sizeof(TYPE *));
  154. unsigned k;
  155. for (k = 0; k < nblocks; k++)
  156. {
  157. starpu_malloc_pinned_if_possible((void **)&tmp_12_block[k], blocksize);
  158. STARPU_ASSERT(tmp_12_block[k]);
  159. starpu_register_blas_data(&tmp_12_block_handles[k], 0,
  160. (uintptr_t)tmp_12_block[k],
  161. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  162. starpu_malloc_pinned_if_possible((void **)&tmp_21_block[k], blocksize);
  163. STARPU_ASSERT(tmp_21_block[k]);
  164. starpu_register_blas_data(&tmp_21_block_handles[k], 0,
  165. (uintptr_t)tmp_21_block[k],
  166. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  167. }
  168. display_all_blocks(nblocks, size/nblocks);
  169. }
  170. int get_block_rank(unsigned i, unsigned j)
  171. {
  172. /* Take a 2D block cyclic distribution */
  173. /* NB: p (resp. q) is for "direction" i (resp. j) */
  174. return (j % q) * p + (i % p);
  175. }
  176. static void display_grid(int rank, unsigned nblocks)
  177. {
  178. if (rank == 0)
  179. {
  180. fprintf(stderr, "2D grid layout: \n");
  181. unsigned i, j;
  182. for (j = 0; j < nblocks; j++)
  183. {
  184. for (i = 0; i < nblocks; i++)
  185. {
  186. fprintf(stderr, "%d ", get_block_rank(i, j));
  187. }
  188. fprintf(stderr, "\n");
  189. }
  190. }
  191. }
  192. static void display_block_content(unsigned bi, unsigned bj, unsigned blocksize)
  193. {
  194. TYPE *data = STARPU_PLU(get_block)(bj, bi);
  195. fprintf(stderr, "BLOCK i = %d j = %d\n", bi, bj);
  196. unsigned i, j;
  197. for (j = 0; j < blocksize; j++)
  198. {
  199. for (i = 0; i < blocksize; i++)
  200. {
  201. fprintf(stderr, "%f ", data[j+i*blocksize]);
  202. }
  203. fprintf(stderr, "\n");
  204. }
  205. fprintf(stderr, "****\n");
  206. }
  207. int main(int argc, char **argv)
  208. {
  209. int rank;
  210. int world_size;
  211. /*
  212. * Initialization
  213. */
  214. MPI_Init(&argc, &argv);
  215. MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  216. MPI_Comm_size(MPI_COMM_WORLD, &world_size);
  217. srand48((long int)time(NULL));
  218. parse_args(argc, argv);
  219. STARPU_ASSERT(p*q == world_size);
  220. //display_grid(rank, nblocks);
  221. starpu_init(NULL);
  222. starpu_mpi_initialize();
  223. starpu_helper_init_cublas();
  224. int barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  225. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  226. /*
  227. * Problem Init
  228. */
  229. init_matrix(rank);
  230. TYPE *a_r = STARPU_PLU(reconstruct_matrix)(size, nblocks);
  231. STARPU_PLU(display_data_content)(a_r, size);
  232. TYPE *x, *y;
  233. if (check)
  234. {
  235. unsigned ind;
  236. x = calloc(size, sizeof(TYPE));
  237. STARPU_ASSERT(x);
  238. y = calloc(size, sizeof(TYPE));
  239. STARPU_ASSERT(y);
  240. if (rank == 0)
  241. {
  242. fprintf(stderr, "Compute AX = B\n");
  243. for (ind = 0; ind < size; ind++)
  244. {
  245. x[ind] = (TYPE)ind;
  246. // x[ind] = (TYPE)drand48();
  247. y[ind] = (TYPE)0.0;
  248. }
  249. }
  250. STARPU_PLU(compute_ax)(size, x, y, nblocks, rank);
  251. if (rank == 0)
  252. for (ind = 0; ind < STARPU_MIN(10, size); ind++)
  253. {
  254. fprintf(stderr, "y[%d] = %f\n", ind, (float)y[ind]);
  255. }
  256. }
  257. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  258. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  259. double timing = STARPU_PLU(plu_main)(nblocks, rank, world_size);
  260. /*
  261. * Report performance
  262. */
  263. int reduce_ret;
  264. double min_timing = timing;
  265. double max_timing = timing;
  266. double sum_timing = timing;
  267. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  268. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  269. reduce_ret = MPI_Reduce(&timing, &min_timing, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD);
  270. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  271. reduce_ret = MPI_Reduce(&timing, &max_timing, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
  272. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  273. reduce_ret = MPI_Reduce(&timing, &sum_timing, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
  274. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  275. if (rank == 0)
  276. {
  277. fprintf(stderr, "Computation took: %lf ms\n", max_timing/1000);
  278. fprintf(stderr, "\tMIN : %lf ms\n", min_timing/1000);
  279. fprintf(stderr, "\tMAX : %lf ms\n", max_timing/1000);
  280. fprintf(stderr, "\tAVG : %lf ms\n", sum_timing/(world_size*1000));
  281. unsigned n = size;
  282. double flop = (2.0f*n*n*n)/3.0f;
  283. fprintf(stderr, "Synthetic GFlops : %2.2f\n", (flop/max_timing/1000.0f));
  284. }
  285. /*
  286. * Test Result Correctness
  287. */
  288. STARPU_PLU(compute_lu_matrix)(size, nblocks);
  289. TYPE *y2;
  290. if (check)
  291. {
  292. unsigned ind;
  293. y2 = calloc(size, sizeof(TYPE));
  294. STARPU_ASSERT(y);
  295. if (rank == 0)
  296. {
  297. fprintf(stderr, "Compute LUX = B2\n");
  298. for (ind = 0; ind < size; ind++)
  299. {
  300. y2[ind] = (TYPE)0.0;
  301. }
  302. }
  303. STARPU_PLU(compute_lux)(size, x, y2, nblocks, rank);
  304. if (rank == 0)
  305. for (ind = 0; ind < STARPU_MIN(10, size); ind++)
  306. {
  307. fprintf(stderr, "y[%d] = %f\n", ind, (float)y2[ind]);
  308. }
  309. }
  310. /*
  311. * Termination
  312. */
  313. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  314. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  315. starpu_helper_shutdown_cublas();
  316. starpu_mpi_shutdown();
  317. starpu_shutdown();
  318. MPI_Finalize();
  319. return 0;
  320. }