starpufftx1d.c 22 KB

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