cg_kernels.c 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010 Université de Bordeaux 1
  4. *
  5. * StarPU 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. * StarPU 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 "cg.h"
  17. #include <math.h>
  18. #include <limits.h>
  19. #if 0
  20. static void print_vector_from_descr(unsigned nx, TYPE *v)
  21. {
  22. unsigned i;
  23. for (i = 0; i < nx; i++)
  24. {
  25. fprintf(stderr, "%2.2e ", v[i]);
  26. }
  27. fprintf(stderr, "\n");
  28. }
  29. static void print_matrix_from_descr(unsigned nx, unsigned ny, unsigned ld, TYPE *mat)
  30. {
  31. unsigned i, j;
  32. for (j = 0; j < nx; j++)
  33. {
  34. for (i = 0; i < ny; i++)
  35. {
  36. fprintf(stderr, "%2.2e ", mat[j+i*ld]);
  37. }
  38. fprintf(stderr, "\n");
  39. }
  40. }
  41. #endif
  42. /*
  43. * Reduction accumulation methods
  44. */
  45. #ifdef STARPU_USE_CUDA
  46. static void accumulate_variable_cuda(void *descr[], void *cl_arg)
  47. {
  48. TYPE *v_dst = (TYPE *)STARPU_VARIABLE_GET_PTR(descr[0]);
  49. TYPE *v_src = (TYPE *)STARPU_VARIABLE_GET_PTR(descr[1]);
  50. cublasaxpy(1, (TYPE)1.0, v_src, 1, v_dst, 1);
  51. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  52. }
  53. #endif
  54. static void accumulate_variable_cpu(void *descr[], void *cl_arg)
  55. {
  56. TYPE *v_dst = (TYPE *)STARPU_VARIABLE_GET_PTR(descr[0]);
  57. TYPE *v_src = (TYPE *)STARPU_VARIABLE_GET_PTR(descr[1]);
  58. *v_dst = *v_dst + *v_src;
  59. }
  60. static struct starpu_perfmodel accumulate_variable_model =
  61. {
  62. .type = STARPU_HISTORY_BASED,
  63. .symbol = "accumulate_variable"
  64. };
  65. struct starpu_codelet accumulate_variable_cl =
  66. {
  67. .where = STARPU_CPU|STARPU_CUDA,
  68. .cpu_funcs = {accumulate_variable_cpu, NULL},
  69. #ifdef STARPU_USE_CUDA
  70. .cuda_funcs = {accumulate_variable_cuda, NULL},
  71. #endif
  72. .nbuffers = 2,
  73. .model = &accumulate_variable_model
  74. };
  75. #ifdef STARPU_USE_CUDA
  76. static void accumulate_vector_cuda(void *descr[], void *cl_arg)
  77. {
  78. TYPE *v_dst = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  79. TYPE *v_src = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  80. unsigned n = STARPU_VECTOR_GET_NX(descr[0]);
  81. cublasaxpy(n, (TYPE)1.0, v_src, 1, v_dst, 1);
  82. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  83. }
  84. #endif
  85. static void accumulate_vector_cpu(void *descr[], void *cl_arg)
  86. {
  87. TYPE *v_dst = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  88. TYPE *v_src = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  89. unsigned n = STARPU_VECTOR_GET_NX(descr[0]);
  90. AXPY(n, (TYPE)1.0, v_src, 1, v_dst, 1);
  91. }
  92. static struct starpu_perfmodel accumulate_vector_model =
  93. {
  94. .type = STARPU_HISTORY_BASED,
  95. .symbol = "accumulate_vector"
  96. };
  97. struct starpu_codelet accumulate_vector_cl =
  98. {
  99. .where = STARPU_CPU|STARPU_CUDA,
  100. .cpu_funcs = {accumulate_vector_cpu, NULL},
  101. #ifdef STARPU_USE_CUDA
  102. .cuda_funcs = {accumulate_vector_cuda, NULL},
  103. #endif
  104. .nbuffers = 2,
  105. .model = &accumulate_vector_model
  106. };
  107. /*
  108. * Reduction initialization methods
  109. */
  110. #ifdef STARPU_USE_CUDA
  111. extern void zero_vector(TYPE *x, unsigned nelems);
  112. static void bzero_variable_cuda(void *descr[], void *cl_arg)
  113. {
  114. TYPE *v = (TYPE *)STARPU_VARIABLE_GET_PTR(descr[0]);
  115. zero_vector(v, 1);
  116. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  117. }
  118. #endif
  119. static void bzero_variable_cpu(void *descr[], void *cl_arg)
  120. {
  121. TYPE *v = (TYPE *)STARPU_VARIABLE_GET_PTR(descr[0]);
  122. *v = (TYPE)0.0;
  123. }
  124. static struct starpu_perfmodel bzero_variable_model =
  125. {
  126. .type = STARPU_HISTORY_BASED,
  127. .symbol = "bzero_variable"
  128. };
  129. struct starpu_codelet bzero_variable_cl =
  130. {
  131. .where = STARPU_CPU|STARPU_CUDA,
  132. .cpu_funcs = {bzero_variable_cpu, NULL},
  133. #ifdef STARPU_USE_CUDA
  134. .cuda_funcs = {bzero_variable_cuda, NULL},
  135. #endif
  136. .nbuffers = 1,
  137. .model = &bzero_variable_model
  138. };
  139. #ifdef STARPU_USE_CUDA
  140. static void bzero_vector_cuda(void *descr[], void *cl_arg)
  141. {
  142. TYPE *v = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  143. unsigned n = STARPU_VECTOR_GET_NX(descr[0]);
  144. zero_vector(v, n);
  145. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  146. }
  147. #endif
  148. static void bzero_vector_cpu(void *descr[], void *cl_arg)
  149. {
  150. TYPE *v = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  151. unsigned n = STARPU_VECTOR_GET_NX(descr[0]);
  152. memset(v, 0, n*sizeof(TYPE));
  153. }
  154. static struct starpu_perfmodel bzero_vector_model =
  155. {
  156. .type = STARPU_HISTORY_BASED,
  157. .symbol = "bzero_vector"
  158. };
  159. struct starpu_codelet bzero_vector_cl =
  160. {
  161. .where = STARPU_CPU|STARPU_CUDA,
  162. .cpu_funcs = {bzero_vector_cpu, NULL},
  163. #ifdef STARPU_USE_CUDA
  164. .cuda_funcs = {bzero_vector_cuda, NULL},
  165. #endif
  166. .nbuffers = 1,
  167. .model = &bzero_vector_model
  168. };
  169. /*
  170. * DOT kernel : s = dot(v1, v2)
  171. */
  172. #ifdef STARPU_USE_CUDA
  173. extern void dot_host(TYPE *x, TYPE *y, unsigned nelems, TYPE *dot);
  174. static void dot_kernel_cuda(void *descr[], void *cl_arg)
  175. {
  176. TYPE *dot = (TYPE *)STARPU_VARIABLE_GET_PTR(descr[0]);
  177. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  178. TYPE *v2 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[2]);
  179. unsigned n = STARPU_VECTOR_GET_NX(descr[1]);
  180. /* Contrary to cublasSdot, this function puts its result directly in
  181. * device memory, so that we don't have to transfer that value back and
  182. * forth. */
  183. dot_host(v1, v2, n, dot);
  184. }
  185. #endif
  186. static void dot_kernel_cpu(void *descr[], void *cl_arg)
  187. {
  188. TYPE *dot = (TYPE *)STARPU_VARIABLE_GET_PTR(descr[0]);
  189. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  190. TYPE *v2 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[2]);
  191. unsigned n = STARPU_VECTOR_GET_NX(descr[1]);
  192. TYPE local_dot = 0.0;
  193. /* Note that we explicitely cast the result of the DOT kernel because
  194. * some BLAS library will return a double for sdot for instance. */
  195. local_dot = (TYPE)DOT(n, v1, 1, v2, 1);
  196. *dot = *dot + local_dot;
  197. }
  198. static struct starpu_perfmodel dot_kernel_model =
  199. {
  200. .type = STARPU_HISTORY_BASED,
  201. .symbol = "dot_kernel"
  202. };
  203. static struct starpu_codelet dot_kernel_cl =
  204. {
  205. .where = STARPU_CPU|STARPU_CUDA,
  206. .cpu_funcs = {dot_kernel_cpu, NULL},
  207. #ifdef STARPU_USE_CUDA
  208. .cuda_funcs = {dot_kernel_cuda, NULL},
  209. #endif
  210. .nbuffers = 3,
  211. .model = &dot_kernel_model
  212. };
  213. void dot_kernel(starpu_data_handle_t v1,
  214. starpu_data_handle_t v2,
  215. starpu_data_handle_t s,
  216. unsigned nblocks,
  217. int use_reduction)
  218. {
  219. /* Blank the accumulation variable */
  220. starpu_insert_task(&bzero_variable_cl, STARPU_W, s, 0);
  221. unsigned b;
  222. for (b = 0; b < nblocks; b++)
  223. {
  224. starpu_insert_task(&dot_kernel_cl,
  225. use_reduction?STARPU_REDUX:STARPU_RW, s,
  226. STARPU_R, starpu_data_get_sub_data(v1, 1, b),
  227. STARPU_R, starpu_data_get_sub_data(v2, 1, b),
  228. 0);
  229. }
  230. }
  231. /*
  232. * SCAL kernel : v1 = p1 v1
  233. */
  234. #ifdef STARPU_USE_CUDA
  235. static void scal_kernel_cuda(void *descr[], void *cl_arg)
  236. {
  237. TYPE p1;
  238. starpu_unpack_cl_args(cl_arg, &p1);
  239. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  240. unsigned n = STARPU_VECTOR_GET_NX(descr[0]);
  241. /* v1 = p1 v1 */
  242. TYPE alpha = p1;
  243. cublasscal(n, alpha, v1, 1);
  244. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  245. }
  246. #endif
  247. static void scal_kernel_cpu(void *descr[], void *cl_arg)
  248. {
  249. TYPE alpha;
  250. starpu_unpack_cl_args(cl_arg, &alpha);
  251. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  252. unsigned n = STARPU_VECTOR_GET_NX(descr[0]);
  253. /* v1 = alpha v1 */
  254. SCAL(n, alpha, v1, 1);
  255. }
  256. static struct starpu_perfmodel scal_kernel_model =
  257. {
  258. .type = STARPU_HISTORY_BASED,
  259. .symbol = "scal_kernel"
  260. };
  261. static struct starpu_codelet scal_kernel_cl =
  262. {
  263. .where = STARPU_CPU|STARPU_CUDA,
  264. .cpu_funcs = {scal_kernel_cpu, NULL},
  265. #ifdef STARPU_USE_CUDA
  266. .cuda_funcs = {scal_kernel_cuda, NULL},
  267. #endif
  268. .nbuffers = 1,
  269. .model = &scal_kernel_model
  270. };
  271. /*
  272. * GEMV kernel : v1 = p1 * v1 + p2 * M v2
  273. */
  274. #ifdef STARPU_USE_CUDA
  275. static void gemv_kernel_cuda(void *descr[], void *cl_arg)
  276. {
  277. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  278. TYPE *v2 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[2]);
  279. TYPE *M = (TYPE *)STARPU_MATRIX_GET_PTR(descr[1]);
  280. unsigned ld = STARPU_MATRIX_GET_LD(descr[1]);
  281. unsigned nx = STARPU_MATRIX_GET_NX(descr[1]);
  282. unsigned ny = STARPU_MATRIX_GET_NY(descr[1]);
  283. TYPE alpha, beta;
  284. starpu_unpack_cl_args(cl_arg, &beta, &alpha);
  285. /* Compute v1 = alpha M v2 + beta v1 */
  286. cublasgemv('N', nx, ny, alpha, M, ld, v2, 1, beta, v1, 1);
  287. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  288. }
  289. #endif
  290. static void gemv_kernel_cpu(void *descr[], void *cl_arg)
  291. {
  292. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  293. TYPE *v2 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[2]);
  294. TYPE *M = (TYPE *)STARPU_MATRIX_GET_PTR(descr[1]);
  295. unsigned ld = STARPU_MATRIX_GET_LD(descr[1]);
  296. unsigned nx = STARPU_MATRIX_GET_NX(descr[1]);
  297. unsigned ny = STARPU_MATRIX_GET_NY(descr[1]);
  298. TYPE alpha, beta;
  299. starpu_unpack_cl_args(cl_arg, &beta, &alpha);
  300. int worker_size = starpu_combined_worker_get_size();
  301. if (worker_size > 1)
  302. {
  303. /* Parallel CPU task */
  304. int rank = starpu_combined_worker_get_rank();
  305. int block_size = (ny + worker_size - 1)/worker_size;
  306. int new_nx = STARPU_MIN(nx, block_size*(rank+1)) - block_size*rank;
  307. nx = new_nx;
  308. v1 = &v1[block_size*rank];
  309. M = &M[block_size*rank];
  310. }
  311. /* Compute v1 = alpha M v2 + beta v1 */
  312. GEMV("N", nx, ny, alpha, M, ld, v2, 1, beta, v1, 1);
  313. }
  314. static struct starpu_perfmodel gemv_kernel_model =
  315. {
  316. .type = STARPU_HISTORY_BASED,
  317. .symbol = "gemv_kernel"
  318. };
  319. static struct starpu_codelet gemv_kernel_cl =
  320. {
  321. .where = STARPU_CPU|STARPU_CUDA,
  322. .type = STARPU_SPMD,
  323. .max_parallelism = INT_MAX,
  324. .cpu_funcs = {gemv_kernel_cpu, NULL},
  325. #ifdef STARPU_USE_CUDA
  326. .cuda_funcs = {gemv_kernel_cuda, NULL},
  327. #endif
  328. .nbuffers = 3,
  329. .model = &gemv_kernel_model
  330. };
  331. void gemv_kernel(starpu_data_handle_t v1,
  332. starpu_data_handle_t matrix,
  333. starpu_data_handle_t v2,
  334. TYPE p1, TYPE p2,
  335. unsigned nblocks,
  336. int use_reduction)
  337. {
  338. unsigned b1, b2;
  339. for (b2 = 0; b2 < nblocks; b2++)
  340. {
  341. starpu_insert_task(&scal_kernel_cl,
  342. STARPU_RW, starpu_data_get_sub_data(v1, 1, b2),
  343. STARPU_VALUE, &p1, sizeof(p1),
  344. 0);
  345. }
  346. for (b2 = 0; b2 < nblocks; b2++)
  347. {
  348. for (b1 = 0; b1 < nblocks; b1++)
  349. {
  350. TYPE one = 1.0;
  351. starpu_insert_task(&gemv_kernel_cl,
  352. use_reduction?STARPU_REDUX:STARPU_RW, starpu_data_get_sub_data(v1, 1, b2),
  353. STARPU_R, starpu_data_get_sub_data(matrix, 2, b2, b1),
  354. STARPU_R, starpu_data_get_sub_data(v2, 1, b1),
  355. STARPU_VALUE, &one, sizeof(one),
  356. STARPU_VALUE, &p2, sizeof(p2),
  357. 0);
  358. }
  359. }
  360. }
  361. /*
  362. * AXPY + SCAL kernel : v1 = p1 * v1 + p2 * v2
  363. */
  364. #ifdef STARPU_USE_CUDA
  365. static void scal_axpy_kernel_cuda(void *descr[], void *cl_arg)
  366. {
  367. TYPE p1, p2;
  368. starpu_unpack_cl_args(cl_arg, &p1, &p2);
  369. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  370. TYPE *v2 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  371. unsigned n = STARPU_VECTOR_GET_NX(descr[0]);
  372. /* Compute v1 = p1 * v1 + p2 * v2.
  373. * v1 = p1 v1
  374. * v1 = v1 + p2 v2
  375. */
  376. cublasscal(n, p1, v1, 1);
  377. cublasaxpy(n, p2, v2, 1, v1, 1);
  378. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  379. }
  380. #endif
  381. static void scal_axpy_kernel_cpu(void *descr[], void *cl_arg)
  382. {
  383. TYPE p1, p2;
  384. starpu_unpack_cl_args(cl_arg, &p1, &p2);
  385. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  386. TYPE *v2 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  387. unsigned nx = STARPU_VECTOR_GET_NX(descr[0]);
  388. /* Compute v1 = p1 * v1 + p2 * v2.
  389. * v1 = p1 v1
  390. * v1 = v1 + p2 v2
  391. */
  392. SCAL(nx, p1, v1, 1);
  393. AXPY(nx, p2, v2, 1, v1, 1);
  394. }
  395. static struct starpu_perfmodel scal_axpy_kernel_model =
  396. {
  397. .type = STARPU_HISTORY_BASED,
  398. .symbol = "scal_axpy_kernel"
  399. };
  400. static struct starpu_codelet scal_axpy_kernel_cl =
  401. {
  402. .where = STARPU_CPU|STARPU_CUDA,
  403. .cpu_funcs = {scal_axpy_kernel_cpu, NULL},
  404. #ifdef STARPU_USE_CUDA
  405. .cuda_funcs = {scal_axpy_kernel_cuda, NULL},
  406. #endif
  407. .nbuffers = 2,
  408. .model = &scal_axpy_kernel_model
  409. };
  410. void scal_axpy_kernel(starpu_data_handle_t v1, TYPE p1,
  411. starpu_data_handle_t v2, TYPE p2,
  412. unsigned nblocks)
  413. {
  414. unsigned b;
  415. for (b = 0; b < nblocks; b++)
  416. {
  417. starpu_insert_task(&scal_axpy_kernel_cl,
  418. STARPU_RW, starpu_data_get_sub_data(v1, 1, b),
  419. STARPU_R, starpu_data_get_sub_data(v2, 1, b),
  420. STARPU_VALUE, &p1, sizeof(p1),
  421. STARPU_VALUE, &p2, sizeof(p2),
  422. 0);
  423. }
  424. }
  425. /*
  426. * AXPY kernel : v1 = v1 + p1 * v2
  427. */
  428. #ifdef STARPU_USE_CUDA
  429. static void axpy_kernel_cuda(void *descr[], void *cl_arg)
  430. {
  431. TYPE p1;
  432. starpu_unpack_cl_args(cl_arg, &p1);
  433. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  434. TYPE *v2 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  435. unsigned n = STARPU_VECTOR_GET_NX(descr[0]);
  436. /* Compute v1 = v1 + p1 * v2.
  437. */
  438. cublasaxpy(n, p1, v2, 1, v1, 1);
  439. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  440. }
  441. #endif
  442. static void axpy_kernel_cpu(void *descr[], void *cl_arg)
  443. {
  444. TYPE p1;
  445. starpu_unpack_cl_args(cl_arg, &p1);
  446. TYPE *v1 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  447. TYPE *v2 = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  448. unsigned nx = STARPU_VECTOR_GET_NX(descr[0]);
  449. /* Compute v1 = p1 * v1 + p2 * v2.
  450. */
  451. AXPY(nx, p1, v2, 1, v1, 1);
  452. }
  453. static struct starpu_perfmodel axpy_kernel_model =
  454. {
  455. .type = STARPU_HISTORY_BASED,
  456. .symbol = "axpy_kernel"
  457. };
  458. static struct starpu_codelet axpy_kernel_cl =
  459. {
  460. .where = STARPU_CPU|STARPU_CUDA,
  461. .cpu_funcs = {axpy_kernel_cpu, NULL},
  462. #ifdef STARPU_USE_CUDA
  463. .cuda_funcs = {axpy_kernel_cuda, NULL},
  464. #endif
  465. .nbuffers = 2,
  466. .model = &axpy_kernel_model
  467. };
  468. void axpy_kernel(starpu_data_handle_t v1,
  469. starpu_data_handle_t v2, TYPE p1,
  470. unsigned nblocks)
  471. {
  472. unsigned b;
  473. for (b = 0; b < nblocks; b++)
  474. {
  475. starpu_insert_task(&axpy_kernel_cl,
  476. STARPU_RW, starpu_data_get_sub_data(v1, 1, b),
  477. STARPU_R, starpu_data_get_sub_data(v2, 1, b),
  478. STARPU_VALUE, &p1, sizeof(p1),
  479. 0);
  480. }
  481. }
  482. /*
  483. * COPY kernel : vector_dst <- vector_src
  484. */
  485. static void copy_handle_cpu(void *descr[], void *cl_arg)
  486. {
  487. TYPE *dst = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  488. TYPE *src = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  489. unsigned nx = STARPU_VECTOR_GET_NX(descr[0]);
  490. size_t elemsize = STARPU_VECTOR_GET_ELEMSIZE(descr[0]);
  491. memcpy(dst, src, nx*elemsize);
  492. }
  493. #ifdef STARPU_USE_CUDA
  494. static void copy_handle_cuda(void *descr[], void *cl_arg)
  495. {
  496. TYPE *dst = (TYPE *)STARPU_VECTOR_GET_PTR(descr[0]);
  497. TYPE *src = (TYPE *)STARPU_VECTOR_GET_PTR(descr[1]);
  498. unsigned nx = STARPU_VECTOR_GET_NX(descr[0]);
  499. size_t elemsize = STARPU_VECTOR_GET_ELEMSIZE(descr[0]);
  500. cudaMemcpyAsync(dst, src, nx*elemsize, cudaMemcpyDeviceToDevice, starpu_cuda_get_local_stream());
  501. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  502. }
  503. #endif
  504. static struct starpu_perfmodel copy_handle_model =
  505. {
  506. .type = STARPU_HISTORY_BASED,
  507. .symbol = "copy_handle"
  508. };
  509. static struct starpu_codelet copy_handle_cl =
  510. {
  511. .where = STARPU_CPU|STARPU_CUDA,
  512. .cpu_funcs = {copy_handle_cpu, NULL},
  513. #ifdef STARPU_USE_CUDA
  514. .cuda_funcs = {copy_handle_cuda, NULL},
  515. #endif
  516. .nbuffers = 2,
  517. .model = &copy_handle_model
  518. };
  519. void copy_handle(starpu_data_handle_t dst, starpu_data_handle_t src, unsigned nblocks)
  520. {
  521. unsigned b;
  522. for (b = 0; b < nblocks; b++)
  523. {
  524. starpu_insert_task(&copy_handle_cl,
  525. STARPU_W, starpu_data_get_sub_data(dst, 1, b),
  526. STARPU_R, starpu_data_get_sub_data(src, 1, b),
  527. 0);
  528. }
  529. }