plu_solve.c 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407
  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 <starpu.h>
  17. #include <math.h>
  18. #include "pxlu.h"
  19. static double frobenius_norm(TYPE *v, unsigned n)
  20. {
  21. double sum2 = 0.0;
  22. /* compute sqrt(Sum(|x|^2)) */
  23. unsigned i,j;
  24. for (j = 0; j < n; j++)
  25. for (i = 0; i < n; i++)
  26. {
  27. double a = fabsl((double)v[i+n*j]);
  28. sum2 += a*a;
  29. }
  30. return sqrt(sum2);
  31. }
  32. void STARPU_PLU(display_data_content)(TYPE *data, unsigned blocksize)
  33. {
  34. fprintf(stderr, "DISPLAY BLOCK\n");
  35. unsigned i, j;
  36. for (j = 0; j < blocksize; j++)
  37. {
  38. for (i = 0; i < blocksize; i++)
  39. {
  40. fprintf(stderr, "%f ", data[j+i*blocksize]);
  41. }
  42. fprintf(stderr, "\n");
  43. }
  44. fprintf(stderr, "****\n");
  45. }
  46. static STARPU_PLU(compute_ax_block)(unsigned block_size, TYPE *block_data, TYPE *sub_x, TYPE *sub_y)
  47. {
  48. fprintf(stderr, "block data %p sub x %p sub y %p\n", block_data, sub_x, sub_y);
  49. CPU_GEMV("N", block_size, block_size, 1.0, block_data, block_size, sub_x, 1, 1.0, sub_y, 1);
  50. }
  51. void STARPU_PLU(extract_upper)(unsigned block_size, TYPE *inblock, TYPE *outblock)
  52. {
  53. unsigned li, lj;
  54. for (lj = 0; lj < block_size; lj++)
  55. {
  56. /* Upper block diag is 1 */
  57. outblock[lj*(block_size + 1)] = (TYPE)1.0;
  58. for (li = lj + 1; li < block_size; li++)
  59. {
  60. outblock[lj + li*block_size] = inblock[lj + li*block_size];
  61. }
  62. }
  63. }
  64. void STARPU_PLU(extract_lower)(unsigned block_size, TYPE *inblock, TYPE *outblock)
  65. {
  66. unsigned li, lj;
  67. for (lj = 0; lj < block_size; lj++)
  68. {
  69. for (li = 0; li <= lj; li++)
  70. {
  71. outblock[lj + li*block_size] = inblock[lj + li*block_size];
  72. }
  73. }
  74. }
  75. static STARPU_PLU(compute_ax_block_upper)(unsigned size, unsigned nblocks,
  76. TYPE *block_data, TYPE *sub_x, TYPE *sub_y)
  77. {
  78. unsigned block_size = size/nblocks;
  79. // fprintf(stderr, "KEEP UPPER\n");
  80. // STARPU_PLU(display_data_content)(block_data, block_size);
  81. /* Take a copy of the upper part of the diagonal block */
  82. TYPE *upper_block_copy = calloc((block_size)*(block_size), sizeof(TYPE));
  83. STARPU_PLU(extract_upper)(block_size, block_data, upper_block_copy);
  84. // STARPU_PLU(display_data_content)(upper_block_copy, block_size);
  85. STARPU_PLU(compute_ax_block)(block_size, upper_block_copy, sub_x, sub_y);
  86. free(upper_block_copy);
  87. }
  88. TYPE *STARPU_PLU(reconstruct_matrix)(unsigned size, unsigned nblocks)
  89. {
  90. // fprintf(stderr, "RECONSTRUCT MATRIX size %d nblocks %d\n", size, nblocks);
  91. TYPE *bigmatrix = calloc(size*size, sizeof(TYPE));
  92. unsigned block_size = size/nblocks;
  93. int rank;
  94. MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  95. unsigned bi, bj;
  96. for (bj = 0; bj < nblocks; bj++)
  97. for (bi = 0; bi < nblocks; bi++)
  98. {
  99. TYPE *block;
  100. if (get_block_rank(bi, bj) == 0)
  101. {
  102. // if (rank == 0)
  103. // fprintf(stderr, "RECONSTRUCT bi %d bj %d\n", bi, bj);
  104. block = STARPU_PLU(get_block)(bi, bj);
  105. }
  106. else {
  107. MPI_Status status;
  108. if (rank == 0)
  109. {
  110. block = calloc(block_size*block_size, sizeof(TYPE));
  111. // fprintf(stderr, "RECONSTRUCT bi %d bj %d - receive in %p from node %d\n", bi, bj, block, get_block_rank(bi, bj));
  112. int ret = MPI_Recv(block, block_size*block_size, MPI_TYPE, get_block_rank(bi, bj), 0, MPI_COMM_WORLD, &status);
  113. STARPU_ASSERT(ret == MPI_SUCCESS);
  114. }
  115. else {
  116. block = STARPU_PLU(get_block)(bi, bj);
  117. // fprintf(stderr, "RECONSTRUCT bi %d bj %d - send %p from node %d\n", bi, bj, block, get_block_rank(bi, bj));
  118. int ret = MPI_Send(block, block_size*block_size, MPI_TYPE, 0, 0, MPI_COMM_WORLD);
  119. STARPU_ASSERT(ret == MPI_SUCCESS);
  120. }
  121. }
  122. if (rank == 0)
  123. {
  124. unsigned j, i;
  125. for (j = 0; j < block_size; j++)
  126. for (i = 0; i < block_size; i++)
  127. {
  128. bigmatrix[(j + bj*block_size)+(i+bi*block_size)*size] =
  129. block[j+i*block_size];
  130. }
  131. if (get_block_rank(bi, bj) != 0)
  132. free(block);
  133. }
  134. }
  135. return bigmatrix;
  136. }
  137. static TYPE *reconstruct_lower(unsigned size, unsigned nblocks)
  138. {
  139. TYPE *lower = calloc(size*size, sizeof(TYPE));
  140. TYPE *bigmatrix = STARPU_PLU(reconstruct_matrix)(size, nblocks);
  141. STARPU_PLU(extract_lower)(size, bigmatrix, lower);
  142. return lower;
  143. }
  144. static TYPE *reconstruct_upper(unsigned size, unsigned nblocks)
  145. {
  146. TYPE *upper = calloc(size*size, sizeof(TYPE));
  147. TYPE *bigmatrix = STARPU_PLU(reconstruct_matrix)(size, nblocks);
  148. STARPU_PLU(extract_upper)(size, bigmatrix, upper);
  149. return upper;
  150. }
  151. void STARPU_PLU(compute_lu_matrix)(unsigned size, unsigned nblocks, TYPE *Asaved)
  152. {
  153. TYPE *all_r = STARPU_PLU(reconstruct_matrix)(size, nblocks);
  154. int rank;
  155. MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  156. if (rank == 0)
  157. {
  158. TYPE *L = malloc((size_t)size*size*sizeof(TYPE));
  159. TYPE *U = malloc((size_t)size*size*sizeof(TYPE));
  160. memset(L, 0, size*size*sizeof(TYPE));
  161. memset(U, 0, size*size*sizeof(TYPE));
  162. /* only keep the lower part */
  163. unsigned i, j;
  164. for (j = 0; j < size; j++)
  165. {
  166. for (i = 0; i < j; i++)
  167. {
  168. L[j+i*size] = all_r[j+i*size];
  169. }
  170. /* diag i = j */
  171. L[j+j*size] = all_r[j+j*size];
  172. U[j+j*size] = 1.0;
  173. for (i = j+1; i < size; i++)
  174. {
  175. U[j+i*size] = all_r[j+i*size];
  176. }
  177. }
  178. STARPU_PLU(display_data_content)(L, size);
  179. STARPU_PLU(display_data_content)(U, size);
  180. /* now A_err = L, compute L*U */
  181. CPU_TRMM("R", "U", "N", "U", size, size, 1.0f, U, size, L, size);
  182. fprintf(stderr, "\nLU\n");
  183. STARPU_PLU(display_data_content)(L, size);
  184. /* compute "LU - A" in L*/
  185. CPU_AXPY(size*size, -1.0, Asaved, 1, L, 1);
  186. TYPE err = CPU_ASUM(size*size, L, 1);
  187. int max = CPU_IAMAX(size*size, L, 1);
  188. fprintf(stderr, "DISPLAY ERROR\n");
  189. STARPU_PLU(display_data_content)(L, size);
  190. fprintf(stderr, "(A - LU) Avg error : %e\n", err/(size*size));
  191. fprintf(stderr, "(A - LU) Max error : %e\n", L[max]);
  192. double residual = frobenius_norm(L, size);
  193. double matnorm = frobenius_norm(Asaved, size);
  194. fprintf(stderr, "||A-LU|| / (||A||*N) : %e\n", residual/(matnorm*size));
  195. }
  196. }
  197. static STARPU_PLU(compute_ax_block_lower)(unsigned size, unsigned nblocks,
  198. TYPE *block_data, TYPE *sub_x, TYPE *sub_y)
  199. {
  200. unsigned block_size = size/nblocks;
  201. // fprintf(stderr, "KEEP LOWER\n");
  202. // STARPU_PLU(display_data_content)(block_data, block_size);
  203. /* Take a copy of the upper part of the diagonal block */
  204. TYPE *lower_block_copy = calloc((block_size)*(block_size), sizeof(TYPE));
  205. STARPU_PLU(extract_lower)(block_size, block_data, lower_block_copy);
  206. // STARPU_PLU(display_data_content)(lower_block_copy, block_size);
  207. STARPU_PLU(compute_ax_block)(size/nblocks, lower_block_copy, sub_x, sub_y);
  208. free(lower_block_copy);
  209. }
  210. void STARPU_PLU(compute_lux)(unsigned size, TYPE *x, TYPE *y, unsigned nblocks, int rank)
  211. {
  212. /* Create temporary buffers where all MPI processes are going to
  213. * compute Ui x = yi where Ai is the matrix containing the blocks of U
  214. * affected to process i, and 0 everywhere else. We then have y as the
  215. * sum of all yi. */
  216. TYPE *yi = calloc(size, sizeof(TYPE));
  217. fprintf(stderr, "Compute LU\n");
  218. unsigned block_size = size/nblocks;
  219. /* Compute UiX = Yi */
  220. unsigned long i,j;
  221. for (j = 0; j < nblocks; j++)
  222. {
  223. if (get_block_rank(j, j) == rank)
  224. {
  225. TYPE *block_data = STARPU_PLU(get_block)(j, j);
  226. TYPE *sub_x = &x[j*(block_size)];
  227. TYPE *sub_yi = &yi[j*(block_size)];
  228. STARPU_PLU(compute_ax_block_upper)(size, nblocks, block_data, sub_x, sub_yi);
  229. }
  230. for (i = j + 1; i < nblocks; i++)
  231. {
  232. if (get_block_rank(i, j) == rank)
  233. {
  234. /* That block belongs to the current MPI process */
  235. TYPE *block_data = STARPU_PLU(get_block)(i, j);
  236. TYPE *sub_x = &x[i*(block_size)];
  237. TYPE *sub_yi = &yi[j*(block_size)];
  238. STARPU_PLU(compute_ax_block)(size/nblocks, block_data, sub_x, sub_yi);
  239. }
  240. }
  241. }
  242. /* Grab Sum Yi in X */
  243. MPI_Reduce(yi, x, size, MPI_TYPE, MPI_SUM, 0, MPI_COMM_WORLD);
  244. memset(yi, 0, size*sizeof(TYPE));
  245. unsigned ind;
  246. // if (rank == 0)
  247. // {
  248. // fprintf(stderr, "INTERMEDIATE\n");
  249. // for (ind = 0; ind < STARPU_MIN(10, size); ind++)
  250. // {
  251. // fprintf(stderr, "x[%d] = %f\n", ind, (float)x[ind]);
  252. // }
  253. // fprintf(stderr, "****\n");
  254. // }
  255. /* Everyone needs x */
  256. int bcst_ret;
  257. bcst_ret = MPI_Bcast(&x, size, MPI_TYPE, 0, MPI_COMM_WORLD);
  258. STARPU_ASSERT(bcst_ret == MPI_SUCCESS);
  259. /* Compute LiX = Yi (with X = UX) */
  260. for (j = 0; j < nblocks; j++)
  261. {
  262. if (j > 0)
  263. for (i = 0; i < j; i++)
  264. {
  265. if (get_block_rank(i, j) == rank)
  266. {
  267. /* That block belongs to the current MPI process */
  268. TYPE *block_data = STARPU_PLU(get_block)(i, j);
  269. TYPE *sub_x = &x[i*(block_size)];
  270. TYPE *sub_yi = &yi[j*(block_size)];
  271. STARPU_PLU(compute_ax_block)(size/nblocks, block_data, sub_x, sub_yi);
  272. }
  273. }
  274. if (get_block_rank(j, j) == rank)
  275. {
  276. TYPE *block_data = STARPU_PLU(get_block)(j, j);
  277. TYPE *sub_x = &x[j*(block_size)];
  278. TYPE *sub_yi = &yi[j*(block_size)];
  279. STARPU_PLU(compute_ax_block_lower)(size, nblocks, block_data, sub_x, sub_yi);
  280. }
  281. }
  282. /* Grab Sum Yi in Y */
  283. MPI_Reduce(yi, y, size, MPI_TYPE, MPI_SUM, 0, MPI_COMM_WORLD);
  284. free(yi);
  285. }
  286. /* x and y must be valid (at least) on 0 */
  287. void STARPU_PLU(compute_ax)(unsigned size, TYPE *x, TYPE *y, unsigned nblocks, int rank)
  288. {
  289. unsigned block_size = size/nblocks;
  290. /* Send x to everyone */
  291. int bcst_ret;
  292. bcst_ret = MPI_Bcast(&x, size, MPI_TYPE, 0, MPI_COMM_WORLD);
  293. STARPU_ASSERT(bcst_ret == MPI_SUCCESS);
  294. /* Create temporary buffers where all MPI processes are going to
  295. * compute Ai x = yi where Ai is the matrix containing the blocks of A
  296. * affected to process i, and 0 everywhere else. We then have y as the
  297. * sum of all yi. */
  298. TYPE *yi = calloc(size, sizeof(TYPE));
  299. /* Compute Aix = yi */
  300. unsigned long i,j;
  301. for (j = 0; j < nblocks; j++)
  302. {
  303. for (i = 0; i < nblocks; i++)
  304. {
  305. if (get_block_rank(i, j) == rank)
  306. {
  307. /* That block belongs to the current MPI process */
  308. TYPE *block_data = STARPU_PLU(get_block)(i, j);
  309. TYPE *sub_x = &x[i*block_size];
  310. TYPE *sub_yi = &yi[j*block_size];
  311. STARPU_PLU(compute_ax_block)(block_size, block_data, sub_x, sub_yi);
  312. }
  313. }
  314. }
  315. /* Compute the Sum of all yi = y */
  316. MPI_Reduce(yi, y, size, MPI_TYPE, MPI_SUM, 0, MPI_COMM_WORLD);
  317. fprintf(stderr, "RANK %d - FOO 1 y[0] %f\n", rank, y[0]);
  318. free(yi);
  319. }