cg_kernels.c 14 KB

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