plu_example.c 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428
  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. static starpu_data_handle *tmp_12_block_handles;
  42. static TYPE **tmp_12_block;
  43. static starpu_data_handle *tmp_21_block_handles;
  44. static TYPE **tmp_21_block;
  45. static void parse_args(int argc, char **argv)
  46. {
  47. int i;
  48. for (i = 1; i < argc; i++) {
  49. if (strcmp(argv[i], "-size") == 0) {
  50. char *argptr;
  51. size = strtol(argv[++i], &argptr, 10);
  52. }
  53. if (strcmp(argv[i], "-nblocks") == 0) {
  54. char *argptr;
  55. nblocks = strtol(argv[++i], &argptr, 10);
  56. }
  57. if (strcmp(argv[i], "-check") == 0) {
  58. check = 1;
  59. }
  60. if (strcmp(argv[i], "-display") == 0) {
  61. display = 1;
  62. }
  63. if (strcmp(argv[i], "-p") == 0) {
  64. char *argptr;
  65. p = strtol(argv[++i], &argptr, 10);
  66. }
  67. if (strcmp(argv[i], "-q") == 0) {
  68. char *argptr;
  69. q = strtol(argv[++i], &argptr, 10);
  70. }
  71. }
  72. }
  73. unsigned STARPU_PLU(display_flag)(void)
  74. {
  75. return display;
  76. }
  77. static void fill_block_with_random(TYPE *blockptr, unsigned size, unsigned nblocks)
  78. {
  79. const unsigned block_size = (size/nblocks);
  80. unsigned i, j;
  81. for (i = 0; i < block_size; i++)
  82. for (j = 0; j < block_size; j++)
  83. {
  84. blockptr[j+i*block_size] = (TYPE)drand48();
  85. }
  86. }
  87. #ifdef SINGLE_TMP11
  88. starpu_data_handle STARPU_PLU(get_tmp_11_block_handle)(void)
  89. {
  90. return tmp_11_block_handle;
  91. }
  92. #else
  93. starpu_data_handle STARPU_PLU(get_tmp_11_block_handle)(unsigned k)
  94. {
  95. return tmp_11_block_handles[k];
  96. }
  97. #endif
  98. starpu_data_handle STARPU_PLU(get_tmp_12_block_handle)(unsigned j)
  99. {
  100. return tmp_12_block_handles[j];
  101. }
  102. starpu_data_handle STARPU_PLU(get_tmp_21_block_handle)(unsigned i)
  103. {
  104. return tmp_21_block_handles[i];
  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. TYPE **blockptr = &dataA[j+i*nblocks];
  119. // starpu_data_handle *handleptr = &dataA_handles[j+nblocks*i];
  120. starpu_data_handle *handleptr = &dataA_handles[j+nblocks*i];
  121. if (get_block_rank(i, j) == rank)
  122. {
  123. /* This blocks should be treated by the current MPI process */
  124. /* Allocate and fill it */
  125. starpu_malloc_pinned_if_possible((void **)blockptr, blocksize);
  126. //fprintf(stderr, "Rank %d : fill block (i = %d, j = %d)\n", rank, i, j);
  127. fill_block_with_random(*blockptr, size, nblocks);
  128. //fprintf(stderr, "Rank %d : fill block (i = %d, j = %d)\n", rank, i, j);
  129. if (i == j)
  130. {
  131. unsigned tmp;
  132. for (tmp = 0; tmp < size/nblocks; tmp++)
  133. {
  134. (*blockptr)[tmp*((size/nblocks)+1)] += (TYPE)10*nblocks;
  135. }
  136. }
  137. /* Register it to StarPU */
  138. starpu_register_blas_data(handleptr, 0,
  139. (uintptr_t)*blockptr, size/nblocks,
  140. size/nblocks, size/nblocks, sizeof(TYPE));
  141. }
  142. else {
  143. *blockptr = STARPU_POISON_PTR;
  144. *handleptr = STARPU_POISON_PTR;
  145. }
  146. }
  147. }
  148. /* Allocate the temporary buffers required for the distributed algorithm */
  149. unsigned k;
  150. /* tmp buffer 11 */
  151. #ifdef SINGLE_TMP11
  152. starpu_malloc_pinned_if_possible((void **)&tmp_11_block, blocksize);
  153. starpu_register_blas_data(&tmp_11_block_handle, 0, (uintptr_t)tmp_11_block,
  154. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  155. #else
  156. tmp_11_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  157. tmp_11_block = calloc(nblocks, sizeof(TYPE *));
  158. for (k = 0; k < nblocks; k++)
  159. {
  160. starpu_malloc_pinned_if_possible((void **)&tmp_11_block[k], blocksize);
  161. STARPU_ASSERT(tmp_11_block[k]);
  162. starpu_register_blas_data(&tmp_11_block_handles[k], 0,
  163. (uintptr_t)tmp_11_block[k],
  164. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  165. }
  166. #endif
  167. /* tmp buffers 12 and 21 */
  168. tmp_12_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  169. tmp_21_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  170. tmp_12_block = calloc(nblocks, sizeof(TYPE *));
  171. tmp_21_block = calloc(nblocks, sizeof(TYPE *));
  172. for (k = 0; k < nblocks; k++)
  173. {
  174. starpu_malloc_pinned_if_possible((void **)&tmp_12_block[k], blocksize);
  175. STARPU_ASSERT(tmp_12_block[k]);
  176. starpu_register_blas_data(&tmp_12_block_handles[k], 0,
  177. (uintptr_t)tmp_12_block[k],
  178. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  179. starpu_malloc_pinned_if_possible((void **)&tmp_21_block[k], blocksize);
  180. STARPU_ASSERT(tmp_21_block[k]);
  181. starpu_register_blas_data(&tmp_21_block_handles[k], 0,
  182. (uintptr_t)tmp_21_block[k],
  183. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  184. }
  185. //display_all_blocks(nblocks, size/nblocks);
  186. }
  187. TYPE *STARPU_PLU(get_block)(unsigned i, unsigned j)
  188. {
  189. return dataA[j+i*nblocks];
  190. }
  191. int get_block_rank(unsigned i, unsigned j)
  192. {
  193. /* Take a 2D block cyclic distribution */
  194. /* NB: p (resp. q) is for "direction" i (resp. j) */
  195. return (j % q) * p + (i % p);
  196. }
  197. starpu_data_handle STARPU_PLU(get_block_handle)(unsigned i, unsigned j)
  198. {
  199. return dataA_handles[j+i*nblocks];
  200. }
  201. static void display_grid(int rank, unsigned nblocks)
  202. {
  203. if (!display)
  204. return;
  205. //if (rank == 0)
  206. {
  207. fprintf(stderr, "2D grid layout (Rank %d): \n", rank);
  208. unsigned i, j;
  209. for (j = 0; j < nblocks; j++)
  210. {
  211. for (i = 0; i < nblocks; i++)
  212. {
  213. TYPE *blockptr = STARPU_PLU(get_block)(i, j);
  214. starpu_data_handle handle = STARPU_PLU(get_block_handle)(i, j);
  215. fprintf(stderr, "%d (data %p handle %p)", get_block_rank(i, j), blockptr, handle);
  216. }
  217. fprintf(stderr, "\n");
  218. }
  219. }
  220. }
  221. int main(int argc, char **argv)
  222. {
  223. int rank;
  224. int world_size;
  225. int thread_support;
  226. /*
  227. * Initialization
  228. */
  229. if (MPI_Init_thread(&argc, &argv, MPI_THREAD_SERIALIZED, &thread_support) != MPI_SUCCESS) {
  230. fprintf(stderr,"MPI_Init_thread failed\n");
  231. exit(1);
  232. }
  233. if (thread_support == MPI_THREAD_FUNNELED)
  234. fprintf(stderr,"Warning: MPI only has funneled thread support, not serialized, hoping this will work\n");
  235. if (thread_support < MPI_THREAD_FUNNELED)
  236. fprintf(stderr,"Warning: MPI does not have thread support!\n");
  237. MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  238. MPI_Comm_size(MPI_COMM_WORLD, &world_size);
  239. srand48((long int)time(NULL));
  240. parse_args(argc, argv);
  241. STARPU_ASSERT(p*q == world_size);
  242. starpu_init(NULL);
  243. starpu_mpi_initialize();
  244. starpu_helper_init_cublas();
  245. int barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  246. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  247. /*
  248. * Problem Init
  249. */
  250. init_matrix(rank);
  251. display_grid(rank, nblocks);
  252. TYPE *a_r;
  253. // STARPU_PLU(display_data_content)(a_r, size);
  254. TYPE *x, *y;
  255. if (check)
  256. {
  257. x = calloc(size, sizeof(TYPE));
  258. STARPU_ASSERT(x);
  259. y = calloc(size, sizeof(TYPE));
  260. STARPU_ASSERT(y);
  261. if (rank == 0)
  262. {
  263. unsigned ind;
  264. for (ind = 0; ind < size; ind++)
  265. x[ind] = (TYPE)drand48();
  266. }
  267. a_r = STARPU_PLU(reconstruct_matrix)(size, nblocks);
  268. if (rank == 0)
  269. STARPU_PLU(display_data_content)(a_r, size);
  270. // STARPU_PLU(compute_ax)(size, x, y, nblocks, rank);
  271. }
  272. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  273. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  274. double timing = STARPU_PLU(plu_main)(nblocks, rank, world_size);
  275. /*
  276. * Report performance
  277. */
  278. int reduce_ret;
  279. double min_timing = timing;
  280. double max_timing = timing;
  281. double sum_timing = timing;
  282. reduce_ret = MPI_Reduce(&timing, &min_timing, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD);
  283. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  284. reduce_ret = MPI_Reduce(&timing, &max_timing, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
  285. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  286. reduce_ret = MPI_Reduce(&timing, &sum_timing, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
  287. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  288. if (rank == 0)
  289. {
  290. fprintf(stderr, "Computation took: %lf ms\n", max_timing/1000);
  291. fprintf(stderr, "\tMIN : %lf ms\n", min_timing/1000);
  292. fprintf(stderr, "\tMAX : %lf ms\n", max_timing/1000);
  293. fprintf(stderr, "\tAVG : %lf ms\n", sum_timing/(world_size*1000));
  294. unsigned n = size;
  295. double flop = (2.0f*n*n*n)/3.0f;
  296. fprintf(stderr, "Synthetic GFlops : %2.2f\n", (flop/max_timing/1000.0f));
  297. }
  298. /*
  299. * Test Result Correctness
  300. */
  301. TYPE *y2;
  302. if (check)
  303. {
  304. /*
  305. * Compute || A - LU ||
  306. */
  307. STARPU_PLU(compute_lu_matrix)(size, nblocks, a_r);
  308. #if 0
  309. /*
  310. * Compute || Ax - LUx ||
  311. */
  312. unsigned ind;
  313. y2 = calloc(size, sizeof(TYPE));
  314. STARPU_ASSERT(y);
  315. if (rank == 0)
  316. {
  317. for (ind = 0; ind < size; ind++)
  318. {
  319. y2[ind] = (TYPE)0.0;
  320. }
  321. }
  322. STARPU_PLU(compute_lux)(size, x, y2, nblocks, rank);
  323. /* Compute y2 = y2 - y */
  324. CPU_AXPY(size, -1.0, y, 1, y2, 1);
  325. TYPE err = CPU_ASUM(size, y2, 1);
  326. int max = CPU_IAMAX(size, y2, 1);
  327. fprintf(stderr, "(A - LU)X Avg error : %e\n", err/(size*size));
  328. fprintf(stderr, "(A - LU)X Max error : %e\n", y2[max]);
  329. #endif
  330. }
  331. /*
  332. * Termination
  333. */
  334. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  335. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  336. starpu_helper_shutdown_cublas();
  337. starpu_mpi_shutdown();
  338. starpu_shutdown();
  339. MPI_Finalize();
  340. return 0;
  341. }