heat.c 19 KB

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