plu_example.c 9.3 KB

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