heat.c 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009, 2010, 2012, 2015 Université de Bordeaux
  4. * Copyright (C) 2010, 2011, 2012, 2016 CNRS
  5. *
  6. * StarPU is free software; you can redistribute it and/or modify
  7. * it under the terms of the GNU Lesser General Public License as published by
  8. * the Free Software Foundation; either version 2.1 of the License, or (at
  9. * your option) any later version.
  10. *
  11. * StarPU is distributed in the hope that it will be useful, but
  12. * WITHOUT ANY WARRANTY; without even the implied warranty of
  13. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  14. *
  15. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  16. */
  17. /*
  18. * heat propagation simulation through either direct LU factorization or
  19. * iterative conjugate gradient.
  20. */
  21. #include "heat.h"
  22. /* default values */
  23. static unsigned ntheta = 32+2;
  24. static unsigned nthick = 32+2;
  25. static unsigned nblocks = 16;
  26. static unsigned nbigblocks = 8;
  27. static unsigned shape = 0;
  28. static unsigned pinned = 0;
  29. static unsigned check = 0;
  30. static unsigned version = 2;
  31. static unsigned use_cg = 0; /* use a LU decomposition of CG ? */
  32. static unsigned no_prio = 0;
  33. extern void do_conjugate_gradient(float *nzvalA, float *vecb, float *vecx, uint32_t nnz,
  34. unsigned nrow, uint32_t *colind, uint32_t *rowptr);
  35. static void parse_args(int argc, char **argv)
  36. {
  37. int i;
  38. for (i = 1; i < argc; i++)
  39. {
  40. if (strcmp(argv[i], "-cg") == 0)
  41. {
  42. use_cg = 1;
  43. }
  44. if (strcmp(argv[i], "-shape") == 0)
  45. {
  46. char *argptr;
  47. shape = strtol(argv[++i], &argptr, 10);
  48. }
  49. if (strcmp(argv[i], "-nthick") == 0)
  50. {
  51. char *argptr;
  52. nthick = strtol(argv[++i], &argptr, 10);
  53. }
  54. if (strcmp(argv[i], "-ntheta") == 0)
  55. {
  56. char *argptr;
  57. ntheta = strtol(argv[++i], &argptr, 10);
  58. }
  59. if (strcmp(argv[i], "-nblocks") == 0)
  60. {
  61. char *argptr;
  62. nblocks = strtol(argv[++i], &argptr, 10);
  63. }
  64. if (strcmp(argv[i], "-nbigblocks") == 0)
  65. {
  66. char *argptr;
  67. nbigblocks = strtol(argv[++i], &argptr, 10);
  68. }
  69. if (strcmp(argv[i], "-v1") == 0)
  70. {
  71. version = 1;
  72. }
  73. if (strcmp(argv[i], "-v2") == 0)
  74. {
  75. version = 2;
  76. }
  77. if (strcmp(argv[i], "-v3") == 0)
  78. {
  79. version = 3;
  80. }
  81. if (strcmp(argv[i], "-v4") == 0)
  82. {
  83. version = 4;
  84. }
  85. if (strcmp(argv[i], "-pin") == 0)
  86. {
  87. pinned = 1;
  88. }
  89. if (strcmp(argv[i], "-check") == 0)
  90. {
  91. check = 1;
  92. }
  93. if (strcmp(argv[i], "-no-prio") == 0)
  94. {
  95. no_prio = 1;
  96. }
  97. if (strcmp(argv[i], "-size") == 0)
  98. {
  99. char *argptr;
  100. unsigned size = strtol(argv[++i], &argptr, 10);
  101. nthick = 130;
  102. ntheta = (size/128) + 2;
  103. STARPU_ASSERT((nthick - 2)*(ntheta - 2) == size);
  104. }
  105. if (strcmp(argv[i], "-h") == 0)
  106. {
  107. printf("usage : %s [-v1|-v2|-v3] [-pin] [-nthick number] [-ntheta number] [-shape [0|1|2]] [-cg] [-size number] [-no-prio]\n", argv[0]);
  108. }
  109. }
  110. }
  111. /*
  112. * The Finite element method code
  113. *
  114. * B C
  115. * **********
  116. * * 0 * *
  117. * * * *
  118. * * * 1 *
  119. * **********
  120. * A D
  121. */
  122. static inline float diff_psi(unsigned theta_tr, unsigned thick_tr, unsigned side_tr,
  123. unsigned theta_psi, unsigned thick_psi, unsigned xy, point *pmesh)
  124. {
  125. float xa,ya,xb,yb,xc,yc;
  126. float tmp;
  127. assert(theta_tr + 2 <= ntheta);
  128. assert(thick_tr + 2 <= nthick);
  129. /* A */
  130. xa = pmesh[NODE_NUMBER(theta_tr, thick_tr)].x;
  131. ya = pmesh[NODE_NUMBER(theta_tr, thick_tr)].y;
  132. /* B */
  133. if (side_tr)
  134. {
  135. /* lower D is actually B here */
  136. xb = pmesh[NODE_NUMBER(theta_tr+1, thick_tr)].x;
  137. yb = pmesh[NODE_NUMBER(theta_tr+1, thick_tr)].y;
  138. }
  139. else
  140. {
  141. /* upper */
  142. xb = pmesh[NODE_NUMBER(theta_tr, thick_tr+1)].x;
  143. yb = pmesh[NODE_NUMBER(theta_tr, thick_tr+1)].y;
  144. }
  145. xc = pmesh[NODE_NUMBER(theta_tr+1, thick_tr+1)].x;
  146. yc = pmesh[NODE_NUMBER(theta_tr+1, thick_tr+1)].y;
  147. /* now look for the actual psi node */
  148. if (NODE_NUMBER(theta_tr, thick_tr) == NODE_NUMBER(theta_psi, thick_psi))
  149. {
  150. /* A nothing to do */
  151. }
  152. else if (NODE_NUMBER(theta_tr+1, thick_tr+1) == NODE_NUMBER(theta_psi, thick_psi))
  153. {
  154. /* psi matches C */
  155. /* swap A and C coordinates */
  156. tmp = xa; xa = xc; xc = tmp;
  157. tmp = ya; ya = yc; yc = tmp;
  158. }
  159. else if (side_tr && (NODE_NUMBER(theta_tr+1, thick_tr) == NODE_NUMBER(theta_psi, thick_psi)))
  160. {
  161. /* psi is D (that was stored in C) XXX */
  162. tmp = xa; xa = xb; xb = tmp;
  163. tmp = ya; ya = yb; yb = tmp;
  164. }
  165. else if (!side_tr && (NODE_NUMBER(theta_tr, thick_tr+1) == NODE_NUMBER(theta_psi, thick_psi)))
  166. {
  167. /* psi is C */
  168. tmp = xa; xa = xb; xb = tmp;
  169. tmp = ya; ya = yb; yb = tmp;
  170. }
  171. else
  172. {
  173. /* the psi node is not a node of the current triangle */
  174. return 0.0f;
  175. }
  176. /* now the triangle should have A as the psi node */
  177. float denom;
  178. float value;
  179. denom = (xa - xb)*(yc - ya) - (xc - xb)*(ya - yb);
  180. switch (xy)
  181. {
  182. case X:
  183. value = (yc - yb)/denom;
  184. break;
  185. case Y:
  186. value = -(xc - xb)/denom;
  187. break;
  188. default:
  189. assert(0);
  190. }
  191. return value;
  192. }
  193. static inline float diff_y_psi(unsigned theta_tr, unsigned thick_tr, unsigned side_tr,
  194. unsigned theta_psi, unsigned thick_psi, point *pmesh)
  195. {
  196. return diff_psi(theta_tr, thick_tr, side_tr, theta_psi, thick_psi, Y, pmesh);
  197. }
  198. static inline float diff_x_psi(unsigned theta_tr, unsigned thick_tr, unsigned side_tr,
  199. unsigned theta_psi, unsigned thick_psi, point *pmesh)
  200. {
  201. return diff_psi(theta_tr, thick_tr, side_tr, theta_psi, thick_psi, X, pmesh);
  202. }
  203. static inline float surface_triangle(unsigned theta_tr, unsigned thick_tr, unsigned side_tr, point *pmesh)
  204. {
  205. float surface;
  206. float tmp;
  207. float xi, xj, xk, yi, yj, yk;
  208. STARPU_ASSERT(theta_tr + 2 <= ntheta);
  209. STARPU_ASSERT(thick_tr + 2 <= nthick);
  210. xi = pmesh[NODE_NUMBER(theta_tr, thick_tr)].x;
  211. yi = pmesh[NODE_NUMBER(theta_tr, thick_tr)].y;
  212. xj = pmesh[NODE_NUMBER(theta_tr+1, thick_tr+1)].x;
  213. yj = pmesh[NODE_NUMBER(theta_tr+1, thick_tr+1)].y;
  214. if (side_tr)
  215. {
  216. /* lower */
  217. xk = pmesh[NODE_NUMBER(theta_tr+1, thick_tr)].x;
  218. yk = pmesh[NODE_NUMBER(theta_tr+1, thick_tr)].y;
  219. }
  220. else
  221. {
  222. xk = pmesh[NODE_NUMBER(theta_tr, thick_tr+1)].x;
  223. yk = pmesh[NODE_NUMBER(theta_tr, thick_tr+1)].y;
  224. }
  225. tmp = (xi - xj)*(yk -yj) - (xk - xj)*(yi -yj);
  226. surface = 0.5*fabs(tmp);
  227. return surface;
  228. }
  229. static inline float integral_triangle(int theta_tr, int thick_tr, unsigned side_tr,
  230. unsigned theta_i, unsigned thick_i, unsigned theta_j, unsigned thick_j, point *pmesh)
  231. {
  232. float surface;
  233. float value;
  234. float dxi, dxj, dyi, dyj;
  235. if (theta_tr < 0) return 0.0f;
  236. if (theta_tr + 2 > (int)ntheta) return 0.0f;
  237. if (thick_tr < 0) return 0.0f;
  238. if (thick_tr + 2 > (int)nthick) return 0.0f;
  239. dxi = diff_x_psi(theta_tr, thick_tr, side_tr, theta_i, thick_i, pmesh);
  240. dyi = diff_y_psi(theta_tr, thick_tr, side_tr, theta_i, thick_i, pmesh);
  241. dxj = diff_x_psi(theta_tr, thick_tr, side_tr, theta_j, thick_j, pmesh);
  242. dyj = diff_y_psi(theta_tr, thick_tr, side_tr, theta_j, thick_j, pmesh);
  243. surface = surface_triangle(theta_tr, thick_tr, side_tr, pmesh);
  244. value = (dxi*dxj + dyi*dyj)*surface;
  245. return value;
  246. }
  247. static inline float integrale_sum(unsigned theta_i, unsigned thick_i, unsigned theta_j, unsigned thick_j, point *pmesh)
  248. {
  249. float integral = 0.0f;
  250. integral += integral_triangle(theta_i - 1, thick_i - 1, 1, theta_i, thick_i, theta_j, thick_j, pmesh);
  251. integral += integral_triangle(theta_i - 1, thick_i - 1, 0, theta_i, thick_i, theta_j, thick_j, pmesh);
  252. integral += integral_triangle(theta_i - 1, thick_i, 1, theta_i, thick_i, theta_j, thick_j, pmesh);
  253. integral += integral_triangle(theta_i, thick_i, 0, theta_i, thick_i, theta_j, thick_j, pmesh);
  254. integral += integral_triangle(theta_i, thick_i, 1, theta_i, thick_i, theta_j, thick_j, pmesh);
  255. integral += integral_triangle(theta_i, thick_i - 1, 0, theta_i, thick_i, theta_j, thick_j, pmesh);
  256. return integral;
  257. }
  258. static float compute_A_value(unsigned i, unsigned j, point *pmesh)
  259. {
  260. float value = 0.0f;
  261. unsigned thick_i, thick_j;
  262. unsigned theta_i, theta_j;
  263. /* add all contributions from all connex triangles */
  264. thick_i = NODE_TO_THICK(i);
  265. thick_j = NODE_TO_THICK(j);
  266. theta_i = NODE_TO_THETA(i);
  267. theta_j = NODE_TO_THETA(j);
  268. /* Compute the Sum of all the integral over all triangles */
  269. if ( (abs(thick_i - thick_j) <= 1) && (abs(theta_i - theta_j) <= 1) )
  270. {
  271. if ( (theta_j == theta_i -1) && (thick_j == thick_i +1))
  272. goto done;
  273. if ( (theta_j == theta_i + 1) && (thick_j == thick_i - 1))
  274. goto done;
  275. /* this may not be a null entry */
  276. value += integrale_sum(theta_i, thick_i, theta_j, thick_j, pmesh);
  277. }
  278. done:
  279. return value;
  280. }
  281. #define TRANSLATE(k) (RefArray[(k)])
  282. #define TRANSLATEBACK(k) (RefArrayBack[(k)])
  283. static void solve_system(unsigned size, unsigned subsize, float *result, int *RefArray, float *Bformer, float *A, float *B)
  284. {
  285. unsigned i;
  286. /* solve the actual problem LU X = B */
  287. /* solve LX' = Y with X' = UX */
  288. /* solve UX = X' */
  289. FPRINTF(stderr, "Solving the problem ...\n");
  290. float *savedB = NULL;
  291. float *LUB = NULL;
  292. if (check)
  293. {
  294. savedB = malloc(subsize*sizeof(float));
  295. memcpy(savedB, B, subsize*sizeof(float));
  296. LUB = malloc(subsize*sizeof(float));
  297. }
  298. /* L */
  299. STARPU_STRSV("L", "N", "N", subsize, A, subsize, B, 1);
  300. /* U */
  301. STARPU_STRSV("U", "N", "U", subsize, A, subsize, B, 1);
  302. STARPU_ASSERT(DIM == size);
  303. if (check)
  304. {
  305. /* compute the error on (LUB - savedB) which should be 0 */
  306. /* LUB = B */
  307. memcpy(LUB, B, subsize*sizeof(float));
  308. /* LUB = U * LUB */
  309. STARPU_STRMV("U", "N", "U", subsize, A, subsize, LUB, 1);
  310. /* LUB = L * LUB */
  311. STARPU_STRMV("L", "N", "N", subsize, A, subsize, LUB, 1);
  312. /* LUB -= B */
  313. STARPU_SAXPY(subsize, -1.0f, savedB, 1, LUB, 1);
  314. /* check if LUB is close to the 0 vector */
  315. int maxind = STARPU_ISAMAX(subsize, LUB, 1);
  316. FPRINTF(stderr, "max error (LUX - B) = %e\n",LUB[maxind - 1]);
  317. float sum = STARPU_SASUM(subsize, LUB, 1);
  318. FPRINTF(stderr,"avg. error %e\n", sum/subsize);
  319. free(LUB);
  320. free(savedB);
  321. }
  322. /* now display back the ACTUAL result */
  323. for (i = 0; i < subsize; i++)
  324. {
  325. result[TRANSLATE(i)] = B[i];
  326. }
  327. for (i = subsize ; i < size; i++)
  328. {
  329. result[TRANSLATE(i)] = Bformer[TRANSLATE(i)];
  330. }
  331. }
  332. unsigned compute_pivot_array(int *RefArray, int *RefArrayBack, unsigned size)
  333. {
  334. unsigned k;
  335. unsigned index = 0;
  336. unsigned theta, thick;
  337. unsigned newsize;
  338. for (k = 0; k < size; k++)
  339. {
  340. RefArray[k] = k;
  341. RefArrayBack[k] = k;
  342. }
  343. /* first inner nodes */
  344. for (theta = 1; theta < ntheta - 1 ; theta++)
  345. {
  346. for (thick = 1; thick < nthick - 1; thick++)
  347. {
  348. /* inner nodes are unknown */
  349. RefArrayBack[NODE_NUMBER(theta, thick)] = index;
  350. RefArray[index] = NODE_NUMBER(theta, thick);
  351. index++;
  352. }
  353. }
  354. newsize = index;
  355. for (theta=0; theta < ntheta; theta++)
  356. {
  357. /* Lower boundary "South" */
  358. RefArrayBack[NODE_NUMBER(theta, 0)] = index;
  359. RefArray[index++] = NODE_NUMBER(theta, 0);
  360. /* Upper boundary "North" */
  361. RefArrayBack[NODE_NUMBER(theta, nthick-1)] = index;
  362. RefArray[index++] = NODE_NUMBER(theta, nthick-1);
  363. }
  364. for (thick = 1; thick < nthick -1; thick++)
  365. {
  366. /* "West "*/
  367. RefArrayBack[NODE_NUMBER(0, thick)] = index;
  368. RefArray[index++] = NODE_NUMBER(0, thick);
  369. /* "East" */
  370. RefArrayBack[NODE_NUMBER(ntheta-1, thick)] = index;
  371. RefArray[index++] = NODE_NUMBER(ntheta-1, thick);
  372. }
  373. assert(index == size);
  374. return newsize;
  375. }
  376. void build_mesh(point *mesh)
  377. {
  378. unsigned theta, thick;
  379. /* first build the mesh by determining all points positions */
  380. for (theta = 0; theta < ntheta; theta++)
  381. {
  382. float angle;
  383. angle = (ntheta - 1 - theta) * Pi/(ntheta-1);
  384. for (thick = 0; thick < nthick; thick++)
  385. {
  386. float r;
  387. r = thick * (RMAX - RMIN)/(nthick - 1) + RMIN;
  388. switch (shape)
  389. {
  390. default:
  391. case 0:
  392. mesh[NODE_NUMBER(theta,thick)].x = r*cosf(angle);
  393. mesh[NODE_NUMBER(theta,thick)].y = r*sinf(angle);
  394. break;
  395. case 1:
  396. mesh[NODE_NUMBER(theta,thick)].x =
  397. -100 + RMIN+((RMAX-RMIN)*theta)/(ntheta - 1);
  398. mesh[NODE_NUMBER(theta,thick)].y =
  399. RMIN+((RMAX-RMIN)*thick)/(nthick - 1);
  400. break;
  401. case 2:
  402. mesh[NODE_NUMBER(theta,thick)].x = r*(2.0f*theta/(ntheta - 1)- 1.0f);
  403. mesh[NODE_NUMBER(theta,thick)].y = r*(2.0f*thick/(nthick - 1)- 1.0f);
  404. break;
  405. }
  406. }
  407. }
  408. }
  409. static unsigned long build_neighbour_vector(unsigned long*neighbours, unsigned node, int *RefArray, int *RefArrayBack)
  410. {
  411. /* where is that point in the former space ? */
  412. int former = TRANSLATE(node);
  413. int former_thick, former_theta;
  414. former_thick= (int)NODE_TO_THICK(former);
  415. former_theta = (int)NODE_TO_THETA(former);
  416. /* do a list of all the possible neighbours */
  417. unsigned nneighbours = 0;
  418. int dtheta, dthick;
  419. for (dthick = -1; dthick <= 1; dthick++)
  420. {
  421. if ((former_thick + dthick) >= 0 && (former_thick + dthick) <= (int)nthick )
  422. {
  423. for (dtheta = -1; dtheta <= 1; dtheta++)
  424. {
  425. if ((former_theta + dtheta) >= 0 && (former_theta + dtheta) <= (int)ntheta )
  426. {
  427. /* we got a possible neighbour */
  428. unsigned pnode =
  429. NODE_NUMBER((former_theta + dtheta), (former_thick + dthick));
  430. neighbours[nneighbours++] = TRANSLATEBACK(pnode);
  431. }
  432. }
  433. }
  434. }
  435. unsigned i;
  436. /* order that list */
  437. for (i = 0; i < nneighbours; i++)
  438. {
  439. /* find the i^th smallest entry for position i */
  440. unsigned index;
  441. unsigned min , min_index;
  442. min = neighbours[i];
  443. min_index = i;
  444. for (index = i+1; index < nneighbours; index++)
  445. {
  446. STARPU_ASSERT(neighbours[i] != neighbours[index]);
  447. if (neighbours[index] < min)
  448. {
  449. min = neighbours[index];
  450. min_index = index;
  451. }
  452. }
  453. /* swap values */
  454. neighbours[min_index] = neighbours[i];
  455. neighbours[i] = min;
  456. }
  457. return nneighbours;
  458. }
  459. static void build_sparse_stiffness_matrix_B(point *pmesh, float *B, float *Bformer, unsigned size, unsigned newsize, int *RefArray, int *RefArrayBack)
  460. {
  461. unsigned i,j;
  462. /* first give the value of known nodes (at boundaries) */
  463. for (i = 0; i < size; i++)
  464. {
  465. Bformer[i] = 0.0f;
  466. }
  467. for (i = 0; i < nthick; i++)
  468. {
  469. Bformer[i] = 200.0f;
  470. Bformer[size-1-i] = 200.0f;
  471. }
  472. for (i = 1; i < ntheta-1; i++)
  473. {
  474. Bformer[i*nthick] = 200.0f;
  475. Bformer[(i+1)*nthick-1] = 100.0f;
  476. }
  477. /* now the actual stiffness (reordered) matrix*/
  478. for (j = 0 ; j < newsize ; j++)
  479. {
  480. unsigned long neighbour;
  481. unsigned long nneighbours;
  482. unsigned long neighbours[9];
  483. nneighbours = build_neighbour_vector(&neighbours[0], j, RefArray, RefArrayBack);
  484. B[j] = Bformer[TRANSLATE(j)];
  485. for (neighbour = 0; neighbour < nneighbours; neighbour++)
  486. {
  487. unsigned n = neighbours[neighbour];
  488. if (n >= newsize)
  489. {
  490. B[j] -= compute_A_value(TRANSLATE(n), TRANSLATE(j), pmesh)*Bformer[TRANSLATE(n)];
  491. }
  492. }
  493. }
  494. }
  495. static unsigned build_sparse_stiffness_matrix_A(point *pmesh, float **nzval, uint32_t **colind,
  496. uint32_t *rowptr, unsigned newsize, int *RefArray, int *RefArrayBack)
  497. {
  498. unsigned j;
  499. unsigned pos = 0;
  500. *nzval = NULL;
  501. *colind = NULL;
  502. /* now the actual stiffness (reordered) matrix*/
  503. for (j = 0 ; j < newsize ; j++)
  504. {
  505. rowptr[j] = pos;
  506. unsigned long neighbour;
  507. unsigned long nneighbours;
  508. unsigned long neighbours[9];
  509. nneighbours = build_neighbour_vector(&neighbours[0], j, RefArray, RefArrayBack);
  510. for (neighbour = 0; neighbour < nneighbours; neighbour++)
  511. {
  512. unsigned nodeneighbour = neighbours[neighbour];
  513. if (nodeneighbour < newsize)
  514. {
  515. float val;
  516. val = compute_A_value(TRANSLATE(j), TRANSLATE(nodeneighbour), pmesh);
  517. if (val != 0.0f)
  518. {
  519. *nzval = realloc(*nzval, (pos+1)*sizeof(float));
  520. *colind = realloc(*colind, (pos+1)*sizeof(uint32_t));
  521. (*nzval)[pos] = val;
  522. (*colind)[pos] = nodeneighbour;
  523. pos++;
  524. }
  525. }
  526. }
  527. }
  528. rowptr[newsize] = pos;
  529. return pos;
  530. }
  531. static void build_dense_stiffness_matrix_A(point *pmesh, float *A, unsigned newsize, int *RefArray, int *RefArrayBack)
  532. {
  533. unsigned long j;
  534. /* touch all the memory */
  535. memset(A, 0, newsize*newsize*sizeof(float));
  536. /* now the actual stiffness (reordered) matrix*/
  537. for (j = 0 ; j < newsize ; j++)
  538. {
  539. unsigned long neighbour;
  540. unsigned long nneighbours;
  541. unsigned long neighbours[9];
  542. nneighbours = build_neighbour_vector(&neighbours[0], j, RefArray, RefArrayBack);
  543. for (neighbour = 0; neighbour < nneighbours; neighbour++)
  544. {
  545. unsigned long nodeneighbour = neighbours[neighbour];
  546. if (nodeneighbour < newsize)
  547. {
  548. float val;
  549. val = compute_A_value(TRANSLATE(j), TRANSLATE(nodeneighbour), pmesh);
  550. A[j+ (unsigned long)newsize*nodeneighbour] = val;
  551. }
  552. }
  553. }
  554. }
  555. int main(int argc, char **argv)
  556. {
  557. float *A;
  558. float *B;
  559. unsigned newsize;
  560. float *result;
  561. int *RefArray, *RefArrayBack;
  562. point *pmesh;
  563. float *Bformer;
  564. parse_args(argc, argv);
  565. pmesh = malloc(DIM*sizeof(point));
  566. RefArray = malloc(DIM*sizeof(int));
  567. RefArrayBack = malloc(DIM*sizeof(int));
  568. Bformer = malloc(DIM*sizeof(float));
  569. result = malloc(DIM*sizeof(float));
  570. build_mesh(pmesh);
  571. /* now simplify that problem given the boundary conditions
  572. * to do so, we remove the already known variables from the system
  573. * by pivoting the various know variable, RefArray keep track of that
  574. * pivoting */
  575. newsize = compute_pivot_array(RefArray, RefArrayBack, DIM);
  576. /* we can either use a direct method (LU decomposition here) or an
  577. * iterative method (conjugate gradient here) */
  578. if (use_cg)
  579. {
  580. unsigned nnz;
  581. float *nzval;
  582. uint32_t *colind;
  583. uint32_t *rowptr;
  584. rowptr = malloc((newsize+1)*sizeof(uint32_t));
  585. B = malloc(newsize*sizeof(float));
  586. build_sparse_stiffness_matrix_B(pmesh, B, Bformer, DIM, newsize, RefArray, RefArrayBack);
  587. nnz = build_sparse_stiffness_matrix_A(pmesh, &nzval, &colind, rowptr, newsize, RefArray, RefArrayBack);
  588. do_conjugate_gradient(nzval, B, result, nnz, newsize, colind, rowptr);
  589. /* XXX */
  590. memcpy(B, result, newsize*sizeof(float));
  591. /* now display back the ACTUAL result */
  592. unsigned i;
  593. for (i = 0; i < newsize; i++)
  594. {
  595. result[TRANSLATE(i)] = B[i];
  596. }
  597. for (i = newsize ; i < DIM; i++)
  598. {
  599. result[TRANSLATE(i)] = Bformer[TRANSLATE(i)];
  600. }
  601. }
  602. else
  603. {
  604. /* unfortunately CUDA does not allow late memory registration,
  605. * we need to do the malloc using CUDA itself ... */
  606. initialize_system(&A, &B, newsize, pinned);
  607. /* then build the stiffness matrix A */
  608. build_sparse_stiffness_matrix_B(pmesh, B, Bformer, DIM, newsize, RefArray, RefArrayBack);
  609. build_dense_stiffness_matrix_A(pmesh, A, newsize, RefArray, RefArrayBack);
  610. FPRINTF(stderr, "Problem size : %ux%u (%ux%u) (%lu MB)\n", newsize, newsize, DIM, DIM, ((unsigned long)newsize*newsize*4UL)/(1024*1024));
  611. STARPU_ASSERT(newsize % nblocks == 0);
  612. switch (version)
  613. {
  614. case 1:
  615. case 2:
  616. dw_factoLU(A, newsize, newsize, nblocks, version, no_prio);
  617. break;
  618. case 3:
  619. dw_factoLU_tag(A, newsize, newsize, nblocks, no_prio);
  620. break;
  621. case 4:
  622. dw_factoLU_grain(A, newsize, newsize, nblocks, nbigblocks);
  623. break;
  624. default:
  625. STARPU_ABORT();
  626. }
  627. display_stat_heat();
  628. if (check)
  629. solve_system(DIM, newsize, result, RefArray, Bformer, A, B);
  630. starpu_cublas_shutdown();
  631. starpu_shutdown();
  632. free_system(A, B, newsize, pinned);
  633. }
  634. #ifdef STARPU_OPENGL_RENDER
  635. opengl_render(ntheta, nthick, result, pmesh, argc, argv);
  636. #endif
  637. free(pmesh);
  638. free(RefArray);
  639. free(RefArrayBack);
  640. free(Bformer);
  641. free(result);
  642. return 0;
  643. }