starpufftx1d.c 23 KB

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