starpufftx1d.c 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009-2011 Université de Bordeaux 1
  4. * Copyright (C) 2010 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. #define PARALLEL
  18. #ifdef PARALLEL
  19. /* Dumb parallel version, disabled */
  20. #define DIV_1D 64
  21. /*
  22. * Overall strategy for an fft of size n:
  23. * - perform n1 ffts of size n2
  24. * - twiddle
  25. * - perform n2 ffts of size n1
  26. *
  27. * - n1 defaults to DIV_1D, thus n2 defaults to n / DIV_1D.
  28. *
  29. * Precise tasks:
  30. *
  31. * - twist1: twist the whole n-element input (called "in") into n1 chunks of
  32. * size n2, by using n1 tasks taking the whole n-element input as a
  33. * R parameter and one n2 output as a W parameter. The result is
  34. * called twisted1.
  35. * - fft1: perform n1 (n2) ffts, by using n1 tasks doing one fft each. Also
  36. * twiddle the result to prepare for the fft2. The result is called
  37. * fft1.
  38. * - join: depends on all the fft1s, to gather the n1 results of size n2 in
  39. * the fft1 vector.
  40. * - twist2: twist the fft1 vector into n2 chunks of size n1, called twisted2.
  41. * since n2 is typically very large, this step is divided in DIV_1D
  42. * tasks, each of them performing n2/DIV_1D of them
  43. * - fft2: perform n2 ffts of size n1. This is divided in DIV_1D tasks of
  44. * n2/DIV_1D ffts, to be performed in batches. The result is called
  45. * fft2.
  46. * - twist3: twist back the result of the fft2s above into the output buffer.
  47. * Only implemented on CPUs for simplicity of the gathering.
  48. *
  49. * The tag space thus uses 3 dimensions:
  50. * - the number of the plan.
  51. * - the step (TWIST1, FFT1, JOIN, TWIST2, FFT2, TWIST3, END)
  52. * - an index i between 0 and DIV_1D-1.
  53. */
  54. #define STEP_TAG_1D(plan, step, i) _STEP_TAG(plan, step, i)
  55. #ifdef STARPU_USE_CUDA
  56. /* twist1:
  57. *
  58. * Twist the full input vector (first parameter) into one chunk of size n2
  59. * (second parameter) */
  60. static void
  61. STARPUFFT(twist1_1d_kernel_gpu)(void *descr[], void *_args)
  62. {
  63. struct STARPUFFT(args) *args = _args;
  64. STARPUFFT(plan) plan = args->plan;
  65. int i = args->i;
  66. int n1 = plan->n1[0];
  67. int n2 = plan->n2[0];
  68. _cufftComplex * restrict in = (_cufftComplex *)STARPU_VECTOR_GET_PTR(descr[0]);
  69. _cufftComplex * restrict twisted1 = (_cufftComplex *)STARPU_VECTOR_GET_PTR(descr[1]);
  70. STARPUFFT(cuda_twist1_1d_host)(in, twisted1, i, n1, n2);
  71. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  72. }
  73. /* fft1:
  74. *
  75. * Perform one fft of size n2 */
  76. static void
  77. STARPUFFT(fft1_1d_plan_gpu)(void *args)
  78. {
  79. STARPUFFT(plan) plan = args;
  80. cufftResult cures;
  81. int n2 = plan->n2[0];
  82. int workerid = starpu_worker_get_id();
  83. cures = cufftPlan1d(&plan->plans[workerid].plan1_cuda, n2, _CUFFT_C2C, 1);
  84. STARPU_ASSERT(cures == CUFFT_SUCCESS);
  85. cufftSetStream(plan->plans[workerid].plan1_cuda, starpu_cuda_get_local_stream());
  86. STARPU_ASSERT(cures == CUFFT_SUCCESS);
  87. }
  88. static void
  89. STARPUFFT(fft1_1d_kernel_gpu)(void *descr[], void *_args)
  90. {
  91. struct STARPUFFT(args) *args = _args;
  92. STARPUFFT(plan) plan = args->plan;
  93. int i = args->i;
  94. int n2 = plan->n2[0];
  95. cufftResult cures;
  96. _cufftComplex * restrict in = (_cufftComplex *)STARPU_VECTOR_GET_PTR(descr[0]);
  97. _cufftComplex * restrict out = (_cufftComplex *)STARPU_VECTOR_GET_PTR(descr[1]);
  98. const _cufftComplex * restrict roots = (_cufftComplex *)STARPU_VECTOR_GET_PTR(descr[2]);
  99. int workerid = starpu_worker_get_id();
  100. task_per_worker[workerid]++;
  101. cures = _cufftExecC2C(plan->plans[workerid].plan1_cuda, in, out, plan->sign == -1 ? CUFFT_FORWARD : CUFFT_INVERSE);
  102. STARPU_ASSERT(cures == CUFFT_SUCCESS);
  103. STARPUFFT(cuda_twiddle_1d_host)(out, roots, n2, i);
  104. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  105. }
  106. /* fft2:
  107. *
  108. * Perform n3 = n2/DIV_1D ffts of size n1 */
  109. static void
  110. STARPUFFT(fft2_1d_plan_gpu)(void *args)
  111. {
  112. STARPUFFT(plan) plan = args;
  113. cufftResult cures;
  114. int n1 = plan->n1[0];
  115. int n2 = plan->n2[0];
  116. int n3 = n2/DIV_1D;
  117. int workerid = starpu_worker_get_id();
  118. cures = cufftPlan1d(&plan->plans[workerid].plan2_cuda, n1, _CUFFT_C2C, n3);
  119. STARPU_ASSERT(cures == CUFFT_SUCCESS);
  120. cufftSetStream(plan->plans[workerid].plan2_cuda, starpu_cuda_get_local_stream());
  121. STARPU_ASSERT(cures == CUFFT_SUCCESS);
  122. }
  123. static void
  124. STARPUFFT(fft2_1d_kernel_gpu)(void *descr[], void *_args)
  125. {
  126. struct STARPUFFT(args) *args = _args;
  127. STARPUFFT(plan) plan = args->plan;
  128. cufftResult cures;
  129. _cufftComplex * restrict in = (_cufftComplex *)STARPU_VECTOR_GET_PTR(descr[0]);
  130. _cufftComplex * restrict out = (_cufftComplex *)STARPU_VECTOR_GET_PTR(descr[1]);
  131. int workerid = starpu_worker_get_id();
  132. task_per_worker[workerid]++;
  133. /* NOTE using batch support */
  134. cures = _cufftExecC2C(plan->plans[workerid].plan2_cuda, in, out, plan->sign == -1 ? CUFFT_FORWARD : CUFFT_INVERSE);
  135. STARPU_ASSERT(cures == CUFFT_SUCCESS);
  136. cudaStreamSynchronize(starpu_cuda_get_local_stream());
  137. }
  138. #endif
  139. /* twist1:
  140. *
  141. * Twist the full input vector (first parameter) into one chunk of size n2
  142. * (second parameter) */
  143. static void
  144. STARPUFFT(twist1_1d_kernel_cpu)(void *descr[], void *_args)
  145. {
  146. struct STARPUFFT(args) *args = _args;
  147. STARPUFFT(plan) plan = args->plan;
  148. int i = args->i;
  149. int j;
  150. int n1 = plan->n1[0];
  151. int n2 = plan->n2[0];
  152. STARPUFFT(complex) * restrict in = (STARPUFFT(complex) *)STARPU_VECTOR_GET_PTR(descr[0]);
  153. STARPUFFT(complex) * restrict twisted1 = (STARPUFFT(complex) *)STARPU_VECTOR_GET_PTR(descr[1]);
  154. /* printf("twist1 %d %g\n", i, (double) cabs(plan->in[i])); */
  155. for (j = 0; j < n2; j++)
  156. twisted1[j] = in[i+j*n1];
  157. }
  158. #ifdef STARPU_HAVE_FFTW
  159. /* fft1:
  160. *
  161. * Perform one fft of size n2 */
  162. static void
  163. STARPUFFT(fft1_1d_kernel_cpu)(void *descr[], void *_args)
  164. {
  165. struct STARPUFFT(args) *args = _args;
  166. STARPUFFT(plan) plan = args->plan;
  167. int i = args->i;
  168. int j;
  169. int n2 = plan->n2[0];
  170. int workerid = starpu_worker_get_id();
  171. task_per_worker[workerid]++;
  172. STARPUFFT(complex) * restrict twisted1 = (STARPUFFT(complex) *)STARPU_VECTOR_GET_PTR(descr[0]);
  173. STARPUFFT(complex) * restrict fft1 = (STARPUFFT(complex) *)STARPU_VECTOR_GET_PTR(descr[1]);
  174. /* printf("fft1 %d %g\n", i, (double) cabs(twisted1[0])); */
  175. _FFTW(execute_dft)(plan->plans[workerid].plan1_cpu, twisted1, fft1);
  176. /* twiddle fft1 buffer */
  177. for (j = 0; j < n2; j++)
  178. fft1[j] = fft1[j] * plan->roots[0][i*j];
  179. }
  180. #endif
  181. /* twist2:
  182. *
  183. * Twist the full vector (results of the fft1s) into one package of n2/DIV_1D
  184. * chunks of size n1 */
  185. static void
  186. STARPUFFT(twist2_1d_kernel_cpu)(void *descr[], void *_args)
  187. {
  188. struct STARPUFFT(args) *args = _args;
  189. STARPUFFT(plan) plan = args->plan;
  190. int jj = args->jj; /* between 0 and DIV_1D */
  191. int jjj; /* beetween 0 and n3 */
  192. int i;
  193. int n1 = plan->n1[0];
  194. int n2 = plan->n2[0];
  195. int n3 = n2/DIV_1D;
  196. STARPUFFT(complex) * restrict twisted2 = (STARPUFFT(complex) *)STARPU_VECTOR_GET_PTR(descr[0]);
  197. /* printf("twist2 %d %g\n", jj, (double) cabs(plan->fft1[jj])); */
  198. for (jjj = 0; jjj < n3; jjj++) {
  199. int j = jj * n3 + jjj;
  200. for (i = 0; i < n1; i++)
  201. twisted2[jjj*n1+i] = plan->fft1[i*n2+j];
  202. }
  203. }
  204. #ifdef STARPU_HAVE_FFTW
  205. /* fft2:
  206. *
  207. * Perform n3 = n2/DIV_1D ffts of size n1 */
  208. static void
  209. STARPUFFT(fft2_1d_kernel_cpu)(void *descr[], void *_args)
  210. {
  211. struct STARPUFFT(args) *args = _args;
  212. STARPUFFT(plan) plan = args->plan;
  213. /* int jj = args->jj; */
  214. int workerid = starpu_worker_get_id();
  215. task_per_worker[workerid]++;
  216. STARPUFFT(complex) * restrict twisted2 = (STARPUFFT(complex) *)STARPU_VECTOR_GET_PTR(descr[0]);
  217. STARPUFFT(complex) * restrict fft2 = (STARPUFFT(complex) *)STARPU_VECTOR_GET_PTR(descr[1]);
  218. /* printf("fft2 %d %g\n", jj, (double) cabs(twisted2[plan->totsize4-1])); */
  219. _FFTW(execute_dft)(plan->plans[workerid].plan2_cpu, twisted2, fft2);
  220. }
  221. #endif
  222. /* twist3:
  223. *
  224. * Spread the package of n2/DIV_1D chunks of size n1 into the output vector */
  225. static void
  226. STARPUFFT(twist3_1d_kernel_cpu)(void *descr[], void *_args)
  227. {
  228. struct STARPUFFT(args) *args = _args;
  229. STARPUFFT(plan) plan = args->plan;
  230. int jj = args->jj; /* between 0 and DIV_1D */
  231. int jjj; /* beetween 0 and n3 */
  232. int i;
  233. int n1 = plan->n1[0];
  234. int n2 = plan->n2[0];
  235. int n3 = n2/DIV_1D;
  236. const STARPUFFT(complex) * restrict fft2 = (STARPUFFT(complex) *)STARPU_VECTOR_GET_PTR(descr[0]);
  237. /* printf("twist3 %d %g\n", jj, (double) cabs(fft2[0])); */
  238. for (jjj = 0; jjj < n3; jjj++) {
  239. int j = jj * n3 + jjj;
  240. for (i = 0; i < n1; i++)
  241. plan->out[i*n2+j] = fft2[jjj*n1+i];
  242. }
  243. }
  244. /* Performance models for the 5 kinds of tasks */
  245. static struct starpu_perfmodel STARPUFFT(twist1_1d_model) = {
  246. .type = STARPU_HISTORY_BASED,
  247. .symbol = TYPE"twist1_1d"
  248. };
  249. static struct starpu_perfmodel STARPUFFT(fft1_1d_model) = {
  250. .type = STARPU_HISTORY_BASED,
  251. .symbol = TYPE"fft1_1d"
  252. };
  253. static struct starpu_perfmodel STARPUFFT(twist2_1d_model) = {
  254. .type = STARPU_HISTORY_BASED,
  255. .symbol = TYPE"twist2_1d"
  256. };
  257. static struct starpu_perfmodel STARPUFFT(fft2_1d_model) = {
  258. .type = STARPU_HISTORY_BASED,
  259. .symbol = TYPE"fft2_1d"
  260. };
  261. static struct starpu_perfmodel STARPUFFT(twist3_1d_model) = {
  262. .type = STARPU_HISTORY_BASED,
  263. .symbol = TYPE"twist3_1d"
  264. };
  265. /* codelet pointers for the 5 kinds of tasks */
  266. static struct starpu_codelet STARPUFFT(twist1_1d_codelet) = {
  267. .where =
  268. #ifdef STARPU_USE_CUDA
  269. STARPU_CUDA|
  270. #endif
  271. STARPU_CPU,
  272. #ifdef STARPU_USE_CUDA
  273. .cuda_func = STARPUFFT(twist1_1d_kernel_gpu),
  274. #endif
  275. .cpu_func = STARPUFFT(twist1_1d_kernel_cpu),
  276. .model = &STARPUFFT(twist1_1d_model),
  277. .nbuffers = 2
  278. };
  279. static struct starpu_codelet STARPUFFT(fft1_1d_codelet) = {
  280. .where =
  281. #ifdef STARPU_USE_CUDA
  282. STARPU_CUDA|
  283. #endif
  284. #ifdef STARPU_HAVE_FFTW
  285. STARPU_CPU|
  286. #endif
  287. 0,
  288. #ifdef STARPU_USE_CUDA
  289. .cuda_func = STARPUFFT(fft1_1d_kernel_gpu),
  290. #endif
  291. #ifdef STARPU_HAVE_FFTW
  292. .cpu_func = STARPUFFT(fft1_1d_kernel_cpu),
  293. #endif
  294. .model = &STARPUFFT(fft1_1d_model),
  295. .nbuffers = 3
  296. };
  297. static struct starpu_codelet STARPUFFT(twist2_1d_codelet) = {
  298. .where = STARPU_CPU,
  299. .cpu_func = STARPUFFT(twist2_1d_kernel_cpu),
  300. .model = &STARPUFFT(twist2_1d_model),
  301. .nbuffers = 1
  302. };
  303. static struct starpu_codelet STARPUFFT(fft2_1d_codelet) = {
  304. .where =
  305. #ifdef STARPU_USE_CUDA
  306. STARPU_CUDA|
  307. #endif
  308. #ifdef STARPU_HAVE_FFTW
  309. STARPU_CPU|
  310. #endif
  311. 0,
  312. #ifdef STARPU_USE_CUDA
  313. .cuda_func = STARPUFFT(fft2_1d_kernel_gpu),
  314. #endif
  315. #ifdef STARPU_HAVE_FFTW
  316. .cpu_func = STARPUFFT(fft2_1d_kernel_cpu),
  317. #endif
  318. .model = &STARPUFFT(fft2_1d_model),
  319. .nbuffers = 2
  320. };
  321. static struct starpu_codelet STARPUFFT(twist3_1d_codelet) = {
  322. .where = STARPU_CPU,
  323. .cpu_func = STARPUFFT(twist3_1d_kernel_cpu),
  324. .model = &STARPUFFT(twist3_1d_model),
  325. .nbuffers = 1
  326. };
  327. #endif /* PARALLEL */
  328. /* Planning:
  329. *
  330. * - For each CPU worker, we need to plan the two fftw stages.
  331. * - For GPU workers, we need to do the planning in the CUDA context, so we do
  332. * this lazily through the initialised1 and initialised2 flags ; TODO: use
  333. * starpu_execute_on_each_worker instead (done in the omp branch).
  334. * - We allocate all the temporary buffers and register them to starpu.
  335. * - We create all the tasks, but do not submit them yet. It will be possible
  336. * to reuse them at will to perform several ffts with the same planning.
  337. */
  338. STARPUFFT(plan)
  339. STARPUFFT(plan_dft_1d)(int n, int sign, unsigned flags)
  340. {
  341. int workerid;
  342. int n1 = DIV_1D;
  343. int n2 = n / n1;
  344. int n3;
  345. int z;
  346. struct starpu_task *task;
  347. #ifdef STARPU_USE_CUDA
  348. /* cufft 1D limited to 8M elements */
  349. while (n2 > 8 << 20) {
  350. n1 *= 2;
  351. n2 /= 2;
  352. }
  353. #endif
  354. STARPU_ASSERT(n == n1*n2);
  355. STARPU_ASSERT(n1 < (1ULL << I_BITS));
  356. /* distribute the n2 second ffts into DIV_1D packages */
  357. n3 = n2 / DIV_1D;
  358. STARPU_ASSERT(n2 == n3*DIV_1D);
  359. /* TODO: flags? Automatically set FFTW_MEASURE on calibration? */
  360. STARPU_ASSERT(flags == 0);
  361. STARPUFFT(plan) plan = malloc(sizeof(*plan));
  362. memset(plan, 0, sizeof(*plan));
  363. plan->number = STARPU_ATOMIC_ADD(&starpufft_last_plan_number, 1) - 1;
  364. /* The plan number has a limited size */
  365. STARPU_ASSERT(plan->number < (1ULL << NUMBER_BITS));
  366. /* Just one dimension */
  367. plan->dim = 1;
  368. plan->n = malloc(plan->dim * sizeof(*plan->n));
  369. plan->n[0] = n;
  370. check_dims(plan);
  371. plan->n1 = malloc(plan->dim * sizeof(*plan->n1));
  372. plan->n1[0] = n1;
  373. plan->n2 = malloc(plan->dim * sizeof(*plan->n2));
  374. plan->n2[0] = n2;
  375. /* Note: this is for coherency with the 2D case */
  376. plan->totsize = n;
  377. plan->totsize1 = n1;
  378. plan->totsize2 = n2;
  379. plan->totsize3 = DIV_1D;
  380. plan->totsize4 = plan->totsize / plan->totsize3;
  381. plan->type = C2C;
  382. plan->sign = sign;
  383. /* Compute the w^k just once. */
  384. compute_roots(plan);
  385. /* Initialize per-worker working set */
  386. for (workerid = 0; workerid < starpu_worker_get_count(); workerid++) {
  387. switch (starpu_worker_get_type(workerid)) {
  388. case STARPU_CPU_WORKER:
  389. #ifdef STARPU_HAVE_FFTW
  390. /* first fft plan: one fft of size n2.
  391. * FFTW imposes that buffer pointers are known at
  392. * planning time. */
  393. plan->plans[workerid].plan1_cpu = _FFTW(plan_dft_1d)(n2, NULL, NULL, sign, _FFTW_FLAGS);
  394. STARPU_ASSERT(plan->plans[workerid].plan1_cpu);
  395. /* second fft plan: n3 ffts of size n1 */
  396. plan->plans[workerid].plan2_cpu = _FFTW(plan_many_dft)(plan->dim,
  397. plan->n1, n3,
  398. NULL, NULL, 1, plan->totsize1,
  399. NULL, NULL, 1, plan->totsize1,
  400. sign, _FFTW_FLAGS);
  401. STARPU_ASSERT(plan->plans[workerid].plan2_cpu);
  402. #else
  403. #warning libstarpufft can not work correctly if libfftw3 is not installed
  404. #endif
  405. break;
  406. case STARPU_CUDA_WORKER:
  407. break;
  408. default:
  409. STARPU_ABORT();
  410. break;
  411. }
  412. }
  413. #ifdef STARPU_USE_CUDA
  414. starpu_execute_on_each_worker(STARPUFFT(fft1_1d_plan_gpu), plan, STARPU_CUDA);
  415. starpu_execute_on_each_worker(STARPUFFT(fft2_1d_plan_gpu), plan, STARPU_CUDA);
  416. #endif
  417. /* Allocate buffers. */
  418. plan->twisted1 = STARPUFFT(malloc)(plan->totsize * sizeof(*plan->twisted1));
  419. memset(plan->twisted1, 0, plan->totsize * sizeof(*plan->twisted1));
  420. plan->fft1 = STARPUFFT(malloc)(plan->totsize * sizeof(*plan->fft1));
  421. memset(plan->fft1, 0, plan->totsize * sizeof(*plan->fft1));
  422. plan->twisted2 = STARPUFFT(malloc)(plan->totsize * sizeof(*plan->twisted2));
  423. memset(plan->twisted2, 0, plan->totsize * sizeof(*plan->twisted2));
  424. plan->fft2 = STARPUFFT(malloc)(plan->totsize * sizeof(*plan->fft2));
  425. memset(plan->fft2, 0, plan->totsize * sizeof(*plan->fft2));
  426. /* Allocate handle arrays */
  427. plan->twisted1_handle = malloc(plan->totsize1 * sizeof(*plan->twisted1_handle));
  428. plan->fft1_handle = malloc(plan->totsize1 * sizeof(*plan->fft1_handle));
  429. plan->twisted2_handle = malloc(plan->totsize3 * sizeof(*plan->twisted2_handle));
  430. plan->fft2_handle = malloc(plan->totsize3 * sizeof(*plan->fft2_handle));
  431. /* Allocate task arrays */
  432. plan->twist1_tasks = malloc(plan->totsize1 * sizeof(*plan->twist1_tasks));
  433. plan->fft1_tasks = malloc(plan->totsize1 * sizeof(*plan->fft1_tasks));
  434. plan->twist2_tasks = malloc(plan->totsize3 * sizeof(*plan->twist2_tasks));
  435. plan->fft2_tasks = malloc(plan->totsize3 * sizeof(*plan->fft2_tasks));
  436. plan->twist3_tasks = malloc(plan->totsize3 * sizeof(*plan->twist3_tasks));
  437. /* Allocate codelet argument arrays */
  438. plan->fft1_args = malloc(plan->totsize1 * sizeof(*plan->fft1_args));
  439. plan->fft2_args = malloc(plan->totsize3 * sizeof(*plan->fft2_args));
  440. /* Create first-round tasks: DIV_1D tasks of type twist1 and fft1 */
  441. for (z = 0; z < plan->totsize1; z++) {
  442. int i = z;
  443. #define STEP_TAG(step) STEP_TAG_1D(plan, step, i)
  444. plan->fft1_args[z].plan = plan;
  445. plan->fft1_args[z].i = i;
  446. /* Register the twisted1 buffer of size n2. */
  447. starpu_vector_data_register(&plan->twisted1_handle[z], 0, (uintptr_t) &plan->twisted1[z*plan->totsize2], plan->totsize2, sizeof(*plan->twisted1));
  448. /* Register the fft1 buffer of size n2. */
  449. starpu_vector_data_register(&plan->fft1_handle[z], 0, (uintptr_t) &plan->fft1[z*plan->totsize2], plan->totsize2, sizeof(*plan->fft1));
  450. /* We'll need the result of fft1 on the CPU for the second
  451. * twist anyway, so tell starpu to not keep the fft1 buffer in
  452. * the GPU. */
  453. starpu_data_set_wt_mask(plan->fft1_handle[z], 1<<0);
  454. /* Create twist1 task */
  455. plan->twist1_tasks[z] = task = starpu_task_create();
  456. task->cl = &STARPUFFT(twist1_1d_codelet);
  457. /* task->buffers[0].handle = to be filled at execution to point
  458. to the application input. */
  459. task->buffers[0].mode = STARPU_R;
  460. task->buffers[1].handle = plan->twisted1_handle[z];
  461. task->buffers[1].mode = STARPU_W;
  462. task->cl_arg = &plan->fft1_args[z];
  463. task->tag_id = STEP_TAG(TWIST1);
  464. task->use_tag = 1;
  465. task->detach = 1;
  466. task->destroy = 0;
  467. /* Tell that fft1 depends on twisted1 */
  468. starpu_tag_declare_deps(STEP_TAG(FFT1),
  469. 1, STEP_TAG(TWIST1));
  470. /* Create FFT1 task */
  471. plan->fft1_tasks[z] = task = starpu_task_create();
  472. task->cl = &STARPUFFT(fft1_1d_codelet);
  473. task->buffers[0].handle = plan->twisted1_handle[z];
  474. task->buffers[0].mode = STARPU_R;
  475. task->buffers[1].handle = plan->fft1_handle[z];
  476. task->buffers[1].mode = STARPU_W;
  477. task->buffers[2].handle = plan->roots_handle[0];
  478. task->buffers[2].mode = STARPU_R;
  479. task->cl_arg = &plan->fft1_args[z];
  480. task->tag_id = STEP_TAG(FFT1);
  481. task->use_tag = 1;
  482. task->detach = 1;
  483. task->destroy = 0;
  484. /* Tell that the join task will depend on the fft1 task. */
  485. starpu_tag_declare_deps(STEP_TAG_1D(plan, JOIN, 0),
  486. 1, STEP_TAG(FFT1));
  487. #undef STEP_TAG
  488. }
  489. /* Create the join task, only serving as a dependency point between
  490. * fft1 and twist2 tasks */
  491. plan->join_task = task = starpu_task_create();
  492. task->cl = NULL;
  493. task->tag_id = STEP_TAG_1D(plan, JOIN, 0);
  494. task->use_tag = 1;
  495. task->detach = 1;
  496. task->destroy = 0;
  497. /* Create second-round tasks: DIV_1D batches of n2/DIV_1D twist2, fft2,
  498. * and twist3 */
  499. for (z = 0; z < plan->totsize3; z++) {
  500. int jj = z;
  501. #define STEP_TAG(step) STEP_TAG_1D(plan, step, jj)
  502. plan->fft2_args[z].plan = plan;
  503. plan->fft2_args[z].jj = jj;
  504. /* Register n3 twisted2 buffers of size n1 */
  505. starpu_vector_data_register(&plan->twisted2_handle[z], 0, (uintptr_t) &plan->twisted2[z*plan->totsize4], plan->totsize4, sizeof(*plan->twisted2));
  506. starpu_vector_data_register(&plan->fft2_handle[z], 0, (uintptr_t) &plan->fft2[z*plan->totsize4], plan->totsize4, sizeof(*plan->fft2));
  507. /* We'll need the result of fft2 on the CPU for the third
  508. * twist anyway, so tell starpu to not keep the fft2 buffer in
  509. * the GPU. */
  510. starpu_data_set_wt_mask(plan->fft2_handle[z], 1<<0);
  511. /* Tell that twisted2 depends on the join task */
  512. starpu_tag_declare_deps(STEP_TAG(TWIST2),
  513. 1, STEP_TAG_1D(plan, JOIN, 0));
  514. /* Create twist2 task */
  515. plan->twist2_tasks[z] = task = starpu_task_create();
  516. task->cl = &STARPUFFT(twist2_1d_codelet);
  517. task->buffers[0].handle = plan->twisted2_handle[z];
  518. task->buffers[0].mode = STARPU_W;
  519. task->cl_arg = &plan->fft2_args[z];
  520. task->tag_id = STEP_TAG(TWIST2);
  521. task->use_tag = 1;
  522. task->detach = 1;
  523. task->destroy = 0;
  524. /* Tell that fft2 depends on twisted2 */
  525. starpu_tag_declare_deps(STEP_TAG(FFT2),
  526. 1, STEP_TAG(TWIST2));
  527. /* Create FFT2 task */
  528. plan->fft2_tasks[z] = task = starpu_task_create();
  529. task->cl = &STARPUFFT(fft2_1d_codelet);
  530. task->buffers[0].handle = plan->twisted2_handle[z];
  531. task->buffers[0].mode = STARPU_R;
  532. task->buffers[1].handle = plan->fft2_handle[z];
  533. task->buffers[1].mode = STARPU_W;
  534. task->cl_arg = &plan->fft2_args[z];
  535. task->tag_id = STEP_TAG(FFT2);
  536. task->use_tag = 1;
  537. task->detach = 1;
  538. task->destroy = 0;
  539. /* Tell that twist3 depends on fft2 */
  540. starpu_tag_declare_deps(STEP_TAG(TWIST3),
  541. 1, STEP_TAG(FFT2));
  542. /* Create twist3 tasks */
  543. /* These run only on CPUs and thus write directly into the
  544. * application output buffer. */
  545. plan->twist3_tasks[z] = task = starpu_task_create();
  546. task->cl = &STARPUFFT(twist3_1d_codelet);
  547. task->buffers[0].handle = plan->fft2_handle[z];
  548. task->buffers[0].mode = STARPU_R;
  549. task->cl_arg = &plan->fft2_args[z];
  550. task->tag_id = STEP_TAG(TWIST3);
  551. task->use_tag = 1;
  552. task->detach = 1;
  553. task->destroy = 0;
  554. /* Tell that to be completely finished we need to have finished
  555. * this twisted3 */
  556. starpu_tag_declare_deps(STEP_TAG_1D(plan, END, 0),
  557. 1, STEP_TAG(TWIST3));
  558. #undef STEP_TAG
  559. }
  560. /* Create end task, only serving as a join point. */
  561. plan->end_task = task = starpu_task_create();
  562. task->cl = NULL;
  563. task->tag_id = STEP_TAG_1D(plan, END, 0);
  564. task->use_tag = 1;
  565. task->detach = 1;
  566. task->destroy = 0;
  567. return plan;
  568. }
  569. /* Actually submit all the tasks. */
  570. static starpu_tag_t
  571. STARPUFFT(start1dC2C)(STARPUFFT(plan) plan)
  572. {
  573. STARPU_ASSERT(plan->type == C2C);
  574. int z;
  575. for (z=0; z < plan->totsize1; z++) {
  576. starpu_task_submit(plan->twist1_tasks[z]);
  577. starpu_task_submit(plan->fft1_tasks[z]);
  578. }
  579. starpu_task_submit(plan->join_task);
  580. for (z=0; z < plan->totsize3; z++) {
  581. starpu_task_submit(plan->twist2_tasks[z]);
  582. starpu_task_submit(plan->fft2_tasks[z]);
  583. starpu_task_submit(plan->twist3_tasks[z]);
  584. }
  585. starpu_task_submit(plan->end_task);
  586. return STEP_TAG_1D(plan, END, 0);
  587. }
  588. /* Free all the tags. The generic code handles freeing the buffers. */
  589. static void
  590. STARPUFFT(free_1d_tags)(STARPUFFT(plan) plan)
  591. {
  592. unsigned i;
  593. int n1 = plan->n1[0];
  594. for (i = 0; i < n1; i++) {
  595. starpu_tag_remove(STEP_TAG_1D(plan, TWIST1, i));
  596. starpu_tag_remove(STEP_TAG_1D(plan, FFT1, i));
  597. }
  598. starpu_tag_remove(STEP_TAG_1D(plan, JOIN, 0));
  599. for (i = 0; i < DIV_1D; i++) {
  600. starpu_tag_remove(STEP_TAG_1D(plan, TWIST2, i));
  601. starpu_tag_remove(STEP_TAG_1D(plan, FFT2, i));
  602. starpu_tag_remove(STEP_TAG_1D(plan, TWIST3, i));
  603. }
  604. starpu_tag_remove(STEP_TAG_1D(plan, END, 0));
  605. }