plu_example.c 9.0 KB

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