1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
|
/*
* Copyright (c) 1997-1999, 2003 Massachusetts Institute of Technology
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 2 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
*
*/
/*
* executor.c -- execute the fft
*/
/* $Id: executor.c,v 1.1 2008/10/17 06:13:18 scuri Exp $ */
#include "fftw-int.h"
#include <stdio.h>
#include <stdlib.h>
const char *fftw_version = "FFTW V" FFTW_VERSION " ($Id: executor.c,v 1.1 2008/10/17 06:13:18 scuri Exp $)";
/*
* This function is called in other files, so we cannot declare
* it static.
*/
void fftw_strided_copy(int n, fftw_complex *in, int ostride,
fftw_complex *out)
{
int i;
fftw_real r0, r1, i0, i1;
fftw_real r2, r3, i2, i3;
i = 0;
for (; i < (n & 3); ++i) {
out[i * ostride] = in[i];
}
for (; i < n; i += 4) {
r0 = c_re(in[i]);
i0 = c_im(in[i]);
r1 = c_re(in[i + 1]);
i1 = c_im(in[i + 1]);
r2 = c_re(in[i + 2]);
i2 = c_im(in[i + 2]);
r3 = c_re(in[i + 3]);
i3 = c_im(in[i + 3]);
c_re(out[i * ostride]) = r0;
c_im(out[i * ostride]) = i0;
c_re(out[(i + 1) * ostride]) = r1;
c_im(out[(i + 1) * ostride]) = i1;
c_re(out[(i + 2) * ostride]) = r2;
c_im(out[(i + 2) * ostride]) = i2;
c_re(out[(i + 3) * ostride]) = r3;
c_im(out[(i + 3) * ostride]) = i3;
}
}
static void executor_many(int n, const fftw_complex *in,
fftw_complex *out,
fftw_plan_node *p,
int istride,
int ostride,
int howmany, int idist, int odist,
fftw_recurse_kind recurse_kind)
{
int s;
switch (p->type) {
case FFTW_NOTW:
{
fftw_notw_codelet *codelet = p->nodeu.notw.codelet;
HACK_ALIGN_STACK_ODD;
for (s = 0; s < howmany; ++s)
codelet(in + s * idist,
out + s * odist,
istride, ostride);
break;
}
default:
for (s = 0; s < howmany; ++s)
fftw_executor_simple(n, in + s * idist,
out + s * odist,
p, istride, ostride,
recurse_kind);
}
}
#ifdef FFTW_ENABLE_VECTOR_RECURSE
/* executor_many_vector is like executor_many, but it pushes the
howmany loop down to the leaves of the transform: */
static void executor_many_vector(int n, const fftw_complex *in,
fftw_complex *out,
fftw_plan_node *p,
int istride,
int ostride,
int howmany, int idist, int odist)
{
int s;
switch (p->type) {
case FFTW_NOTW:
{
fftw_notw_codelet *codelet = p->nodeu.notw.codelet;
HACK_ALIGN_STACK_ODD;
for (s = 0; s < howmany; ++s)
codelet(in + s * idist,
out + s * odist,
istride, ostride);
break;
}
case FFTW_TWIDDLE:
{
int r = p->nodeu.twiddle.size;
int m = n / r;
fftw_twiddle_codelet *codelet;
fftw_complex *W;
for (s = 0; s < r; ++s)
executor_many_vector(m, in + s * istride,
out + s * (m * ostride),
p->nodeu.twiddle.recurse,
istride * r, ostride,
howmany, idist, odist);
codelet = p->nodeu.twiddle.codelet;
W = p->nodeu.twiddle.tw->twarray;
/* This may not be the right thing. We maybe should have
the howmany loop for the twiddle codelets at the
topmost level of the recursion, since odist is big;
i.e. separate recursions for twiddle and notwiddle. */
HACK_ALIGN_STACK_EVEN;
for (s = 0; s < howmany; ++s)
codelet(out + s * odist, W, m * ostride, m, ostride);
break;
}
case FFTW_GENERIC:
{
int r = p->nodeu.generic.size;
int m = n / r;
fftw_generic_codelet *codelet;
fftw_complex *W;
for (s = 0; s < r; ++s)
executor_many_vector(m, in + s * istride,
out + s * (m * ostride),
p->nodeu.generic.recurse,
istride * r, ostride,
howmany, idist, odist);
codelet = p->nodeu.generic.codelet;
W = p->nodeu.generic.tw->twarray;
for (s = 0; s < howmany; ++s)
codelet(out + s * odist, W, m, r, n, ostride);
break;
}
case FFTW_RADER:
{
int r = p->nodeu.rader.size;
int m = n / r;
fftw_rader_codelet *codelet;
fftw_complex *W;
for (s = 0; s < r; ++s)
executor_many_vector(m, in + s * istride,
out + s * (m * ostride),
p->nodeu.rader.recurse,
istride * r, ostride,
howmany, idist, odist);
codelet = p->nodeu.rader.codelet;
W = p->nodeu.rader.tw->twarray;
for (s = 0; s < howmany; ++s)
codelet(out + s * odist, W, m, r, ostride,
p->nodeu.rader.rader_data);
break;
}
default:
fftw_die("BUG in executor: invalid plan\n");
break;
}
}
#endif /* FFTW_ENABLE_VECTOR_RECURSE */
/*
* Do *not* declare simple executor static--we need to call it
* from other files...also, preface its name with "fftw_"
* to avoid any possible name collisions.
*/
void fftw_executor_simple(int n, const fftw_complex *in,
fftw_complex *out,
fftw_plan_node *p,
int istride,
int ostride,
fftw_recurse_kind recurse_kind)
{
switch (p->type) {
case FFTW_NOTW:
HACK_ALIGN_STACK_ODD;
(p->nodeu.notw.codelet)(in, out, istride, ostride);
break;
case FFTW_TWIDDLE:
{
int r = p->nodeu.twiddle.size;
int m = n / r;
fftw_twiddle_codelet *codelet;
fftw_complex *W;
#ifdef FFTW_ENABLE_VECTOR_RECURSE
if (recurse_kind == FFTW_NORMAL_RECURSE)
#endif
executor_many(m, in, out,
p->nodeu.twiddle.recurse,
istride * r, ostride,
r, istride, m * ostride,
FFTW_NORMAL_RECURSE);
#ifdef FFTW_ENABLE_VECTOR_RECURSE
else
executor_many_vector(m, in, out,
p->nodeu.twiddle.recurse,
istride * r, ostride,
r, istride, m * ostride);
#endif
codelet = p->nodeu.twiddle.codelet;
W = p->nodeu.twiddle.tw->twarray;
HACK_ALIGN_STACK_EVEN;
codelet(out, W, m * ostride, m, ostride);
break;
}
case FFTW_GENERIC:
{
int r = p->nodeu.generic.size;
int m = n / r;
fftw_generic_codelet *codelet;
fftw_complex *W;
#ifdef FFTW_ENABLE_VECTOR_RECURSE
if (recurse_kind == FFTW_NORMAL_RECURSE)
#endif
executor_many(m, in, out,
p->nodeu.generic.recurse,
istride * r, ostride,
r, istride, m * ostride,
FFTW_NORMAL_RECURSE);
#ifdef FFTW_ENABLE_VECTOR_RECURSE
else
executor_many_vector(m, in, out,
p->nodeu.generic.recurse,
istride * r, ostride,
r, istride, m * ostride);
#endif
codelet = p->nodeu.generic.codelet;
W = p->nodeu.generic.tw->twarray;
codelet(out, W, m, r, n, ostride);
break;
}
case FFTW_RADER:
{
int r = p->nodeu.rader.size;
int m = n / r;
fftw_rader_codelet *codelet;
fftw_complex *W;
#ifdef FFTW_ENABLE_VECTOR_RECURSE
if (recurse_kind == FFTW_NORMAL_RECURSE)
#endif
executor_many(m, in, out,
p->nodeu.rader.recurse,
istride * r, ostride,
r, istride, m * ostride,
FFTW_NORMAL_RECURSE);
#ifdef FFTW_ENABLE_VECTOR_RECURSE
else
executor_many_vector(m, in, out,
p->nodeu.rader.recurse,
istride * r, ostride,
r, istride, m * ostride);
#endif
codelet = p->nodeu.rader.codelet;
W = p->nodeu.rader.tw->twarray;
codelet(out, W, m, r, ostride,
p->nodeu.rader.rader_data);
break;
}
default:
fftw_die("BUG in executor: invalid plan\n");
break;
}
}
static void executor_simple_inplace(int n, fftw_complex *in,
fftw_complex *out,
fftw_plan_node *p,
int istride,
fftw_recurse_kind recurse_kind)
{
switch (p->type) {
case FFTW_NOTW:
HACK_ALIGN_STACK_ODD;
(p->nodeu.notw.codelet)(in, in, istride, istride);
break;
default:
{
fftw_complex *tmp;
if (out)
tmp = out;
else
tmp = (fftw_complex *)
fftw_malloc(n * sizeof(fftw_complex));
fftw_executor_simple(n, in, tmp, p, istride, 1,
recurse_kind);
fftw_strided_copy(n, tmp, istride, in);
if (!out)
fftw_free(tmp);
}
}
}
static void executor_many_inplace(int n, fftw_complex *in,
fftw_complex *out,
fftw_plan_node *p,
int istride,
int howmany, int idist,
fftw_recurse_kind recurse_kind)
{
switch (p->type) {
case FFTW_NOTW:
{
fftw_notw_codelet *codelet = p->nodeu.notw.codelet;
int s;
HACK_ALIGN_STACK_ODD;
for (s = 0; s < howmany; ++s)
codelet(in + s * idist,
in + s * idist,
istride, istride);
break;
}
default:
{
int s;
fftw_complex *tmp;
if (out)
tmp = out;
else
tmp = (fftw_complex *)
fftw_malloc(n * sizeof(fftw_complex));
for (s = 0; s < howmany; ++s) {
fftw_executor_simple(n,
in + s * idist,
tmp,
p, istride, 1, recurse_kind);
fftw_strided_copy(n, tmp, istride, in + s * idist);
}
if (!out)
fftw_free(tmp);
}
}
}
/* user interface */
void fftw(fftw_plan plan, int howmany, fftw_complex *in, int istride,
int idist, fftw_complex *out, int ostride, int odist)
{
int n = plan->n;
if (plan->flags & FFTW_IN_PLACE) {
if (howmany == 1) {
executor_simple_inplace(n, in, out, plan->root, istride,
plan->recurse_kind);
} else {
executor_many_inplace(n, in, out, plan->root, istride, howmany,
idist, plan->recurse_kind);
}
} else {
if (howmany == 1) {
fftw_executor_simple(n, in, out, plan->root, istride, ostride,
plan->recurse_kind);
} else {
#ifdef FFTW_ENABLE_VECTOR_RECURSE
int vector_size = plan->vector_size;
if (vector_size <= 1)
#endif
executor_many(n, in, out, plan->root, istride, ostride,
howmany, idist, odist, plan->recurse_kind);
#ifdef FFTW_ENABLE_VECTOR_RECURSE
else {
int s;
int num_vects = howmany / vector_size;
fftw_plan_node *root = plan->root;
for (s = 0; s < num_vects; ++s)
executor_many_vector(n,
in + s * (vector_size * idist),
out + s * (vector_size * odist),
root,
istride, ostride,
vector_size, idist, odist);
s = howmany % vector_size;
if (s > 0)
executor_many(n,
in + num_vects * (vector_size * idist),
out + num_vects * (vector_size * odist),
root,
istride, ostride,
s, idist, odist,
FFTW_NORMAL_RECURSE);
}
#endif
}
}
}
void fftw_one(fftw_plan plan, fftw_complex *in, fftw_complex *out)
{
int n = plan->n;
if (plan->flags & FFTW_IN_PLACE)
executor_simple_inplace(n, in, out, plan->root, 1,
plan->recurse_kind);
else
fftw_executor_simple(n, in, out, plan->root, 1, 1,
plan->recurse_kind);
}
|