plu_example.c 7.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330
  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 (j = 0; j < block_size; j++)
  69. for (i = 0; i < block_size; i++)
  70. {
  71. blockptr[i+j*block_size] = (TYPE)drand48();
  72. // blockptr[i+j*block_size] = (i == j)?2.0:0.0;
  73. }
  74. }
  75. starpu_data_handle STARPU_PLU(get_block_handle)(unsigned j, unsigned i)
  76. {
  77. return dataA_handles[i+j*nblocks];
  78. }
  79. starpu_data_handle STARPU_PLU(get_tmp_11_block_handle)(void)
  80. {
  81. return tmp_11_block_handle;
  82. }
  83. starpu_data_handle STARPU_PLU(get_tmp_12_block_handle)(unsigned j)
  84. {
  85. return tmp_12_block_handles[j];
  86. }
  87. starpu_data_handle STARPU_PLU(get_tmp_21_block_handle)(unsigned i)
  88. {
  89. return tmp_21_block_handles[i];
  90. }
  91. TYPE *STARPU_PLU(get_block)(unsigned j, unsigned i)
  92. {
  93. return dataA[i+j*nblocks];
  94. }
  95. static void init_matrix(int rank)
  96. {
  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. if (get_block_rank(i, j) == rank)
  108. {
  109. /* This blocks should be treated by the current MPI process */
  110. /* Allocate and fill it */
  111. starpu_malloc_pinned_if_possible((void **)&dataA[i+j*nblocks], blocksize);
  112. fill_block_with_random(STARPU_PLU(get_block)(j, i), size, nblocks);
  113. if (i == j)
  114. {
  115. TYPE *b = STARPU_PLU(get_block)(j, i);
  116. unsigned tmp;
  117. for (tmp = 0; tmp < size/nblocks; tmp++)
  118. {
  119. b[tmp*((size/nblocks)+1)] += (TYPE)10*nblocks;
  120. }
  121. }
  122. /* Register it to StarPU */
  123. starpu_register_blas_data(&dataA_handles[i+nblocks*j], 0,
  124. (uintptr_t)dataA[i+nblocks*j], size/nblocks,
  125. size/nblocks, size/nblocks, sizeof(TYPE));
  126. }
  127. else {
  128. dataA[i+j*nblocks] = STARPU_POISON_PTR;
  129. dataA_handles[i+j*nblocks] = STARPU_POISON_PTR;
  130. }
  131. }
  132. }
  133. /* Allocate the temporary buffers required for the distributed algorithm */
  134. /* tmp buffer 11 */
  135. starpu_malloc_pinned_if_possible((void **)&tmp_11_block, blocksize);
  136. starpu_register_blas_data(&tmp_11_block_handle, 0, (uintptr_t)tmp_11_block,
  137. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  138. /* tmp buffers 12 and 21 */
  139. tmp_12_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  140. tmp_21_block_handles = calloc(nblocks, sizeof(starpu_data_handle));
  141. tmp_12_block = calloc(nblocks, sizeof(TYPE *));
  142. tmp_21_block = calloc(nblocks, sizeof(TYPE *));
  143. unsigned k;
  144. for (k = 0; k < nblocks; k++)
  145. {
  146. starpu_malloc_pinned_if_possible((void **)&tmp_12_block[k], blocksize);
  147. STARPU_ASSERT(tmp_12_block[k]);
  148. starpu_register_blas_data(&tmp_12_block_handles[k], 0,
  149. (uintptr_t)tmp_12_block[k],
  150. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  151. starpu_malloc_pinned_if_possible((void **)&tmp_21_block[k], blocksize);
  152. STARPU_ASSERT(tmp_21_block[k]);
  153. starpu_register_blas_data(&tmp_21_block_handles[k], 0,
  154. (uintptr_t)tmp_21_block[k],
  155. size/nblocks, size/nblocks, size/nblocks, sizeof(TYPE));
  156. }
  157. }
  158. int get_block_rank(unsigned i, unsigned j)
  159. {
  160. /* Take a 2D block cyclic distribution */
  161. /* NB: p (resp. q) is for "direction" i (resp. j) */
  162. return (j % q) * p + (i % p);
  163. }
  164. static void display_grid(int rank, unsigned nblocks)
  165. {
  166. if (rank == 0)
  167. {
  168. fprintf(stderr, "2D grid layout: \n");
  169. unsigned i, j;
  170. for (j = 0; j < nblocks; j++)
  171. {
  172. for (i = 0; i < nblocks; i++)
  173. {
  174. fprintf(stderr, "%d ", get_block_rank(i, j));
  175. }
  176. fprintf(stderr, "\n");
  177. }
  178. }
  179. }
  180. int main(int argc, char **argv)
  181. {
  182. int rank;
  183. int world_size;
  184. /*
  185. * Initialization
  186. */
  187. MPI_Init(&argc, &argv);
  188. MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  189. MPI_Comm_size(MPI_COMM_WORLD, &world_size);
  190. srand48((long int)time(NULL));
  191. parse_args(argc, argv);
  192. STARPU_ASSERT(p*q == world_size);
  193. //display_grid(rank, nblocks);
  194. starpu_init(NULL);
  195. starpu_mpi_initialize();
  196. starpu_helper_init_cublas();
  197. int barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  198. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  199. /*
  200. * Problem Init
  201. */
  202. init_matrix(rank);
  203. TYPE *x, *y;
  204. if (check)
  205. {
  206. if (rank == 0)
  207. fprintf(stderr, "Compute AX = B\n");
  208. x = calloc(size, sizeof(TYPE));
  209. STARPU_ASSERT(x);
  210. y = calloc(size, sizeof(TYPE));
  211. STARPU_ASSERT(y);
  212. unsigned ind;
  213. for (ind = 0; ind < size; ind++)
  214. {
  215. //x[ind] = (TYPE)1.0;
  216. x[ind] = (TYPE)drand48();
  217. y[ind] = (TYPE)0.0;
  218. }
  219. STARPU_PLU(compute_ax)(size, x, y, nblocks, rank);
  220. if (rank == 0)
  221. for (ind = 0; ind < 10; ind++)
  222. {
  223. fprintf(stderr, "y[%d] = %f\n", ind, (float)y[ind]);
  224. }
  225. }
  226. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  227. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  228. double timing = STARPU_PLU(plu_main)(nblocks, rank, world_size);
  229. /*
  230. * Report performance
  231. */
  232. int reduce_ret;
  233. double min_timing = timing;
  234. double max_timing = timing;
  235. double sum_timing = timing;
  236. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  237. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  238. reduce_ret = MPI_Reduce(&timing, &min_timing, 1, MPI_DOUBLE, MPI_MIN, 0, MPI_COMM_WORLD);
  239. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  240. reduce_ret = MPI_Reduce(&timing, &max_timing, 1, MPI_DOUBLE, MPI_MAX, 0, MPI_COMM_WORLD);
  241. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  242. reduce_ret = MPI_Reduce(&timing, &sum_timing, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD);
  243. STARPU_ASSERT(reduce_ret == MPI_SUCCESS);
  244. if (rank == 0)
  245. {
  246. fprintf(stderr, "Computation took: %lf ms\n", max_timing/1000);
  247. fprintf(stderr, "\tMIN : %lf ms\n", min_timing/1000);
  248. fprintf(stderr, "\tMAX : %lf ms\n", max_timing/1000);
  249. fprintf(stderr, "\tAVG : %lf ms\n", sum_timing/(world_size*1000));
  250. unsigned n = size;
  251. double flop = (2.0f*n*n*n)/3.0f;
  252. fprintf(stderr, "Synthetic GFlops : %2.2f\n", (flop/max_timing/1000.0f));
  253. }
  254. /*
  255. * Termination
  256. */
  257. barrier_ret = MPI_Barrier(MPI_COMM_WORLD);
  258. STARPU_ASSERT(barrier_ret == MPI_SUCCESS);
  259. starpu_helper_shutdown_cublas();
  260. starpu_mpi_shutdown();
  261. starpu_shutdown();
  262. MPI_Finalize();
  263. return 0;
  264. }