plu_example.c 9.1 KB

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