heat.c 19 KB

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