forked from cyclops-community/ctf
-
Notifications
You must be signed in to change notification settings - Fork 0
/
fft.cxx
217 lines (181 loc) · 5.61 KB
/
fft.cxx
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
/** \addtogroup examples
* @{
* \defgroup fft fft
* @{
* \brief FFT iterative method using gemv and spmv
*/
#include <ctf.hpp>
using namespace CTF;
namespace CTF_int{
// factorizes n into nfactor factors
void factorize(int n, int *nfactor, int **factor);
}
// gets ith power of nth imaginary root of the identity, w^i_n
std::complex<double> omega(int i, int n){
std::complex<double> imag(0,1);
return exp(-2.*i*(M_PI/n)*imag);
}
// gets DFT matrix D(i,j) = w_n^(ij)
Matrix< std::complex<double> > DFT_matrix(int n, World & wrld){
int64_t np;
int64_t * idx;
std::complex<double> * data;
Matrix < std::complex<double> >DFT(n, n, NS, wrld, "DFT");
DFT.get_local_data(&np, &idx, &data);
for (int64_t i=0; i<np; i++){
data[i] = omega((idx[i]/n)*(idx[i]%n), n);
}
DFT.write(np, idx, data);
free(idx);
delete [] data;
return DFT;
}
// gets twiddle factor matirx T = [1, 1; 1, w^1_n]
Matrix< std::complex<double> > twiddle_matrix(int n, World & wrld){
int64_t np;
int64_t * idx;
std::complex<double> * data;
Matrix < std::complex<double> >T(2, 2, NS, wrld, "T");
T.get_local_data(&np, &idx, &data);
for (int64_t i=0; i<np; i++){
if (idx[i]<3) data[i] = std::complex<double>(1,0);
else data[i] = omega(1,n);
}
T.write(np, idx, data);
free(idx);
delete [] data;
return T;
}
void fft(Vector< std::complex<double> > & v, int n){
// assert that n is a power of two
int logn=0;
while (1<<logn < n){
logn++;
}
assert(1<<logn == n);
// factorize n, factors will all just be 2 for now
int nfact;
int * factors;
CTF_int::factorize(n, &nfact, &factors);
assert(nfact == logn);
// Fold v into log_2(n) tensor V
Tensor< std::complex<double> > V = v.reshape(nfact, factors);
// define range of indices [a, b, c, ...]
char inds[nfact+1];
char rev_inds[nfact];
for (int i=0; i<nfact; i++){
inds[i]='a'+i;
rev_inds[i]='a'+nfact-i-1;
}
inds[nfact] = 'a';
// reverse tensor index order
// now that each FFT step transforms the lowest level modes, which should be faster
V[rev_inds] = V[inds];
for (int i=0; i<nfact; i++){
inds[nfact+i]='a'+i;
}
Matrix< std::complex<double> > DFT = DFT_matrix(2,*v.wrld);
Tensor< std::complex<double> > S_old(1, factors, *v.wrld, *v.sr, "S");
S_old["i"] = 1;
char inds_in[nfact];
char inds_DFT[2];
// m = 2^{i+1}
int m = 1;
// execute FFT from the bottom level of recursion upwards
for (int i=0; i<nfact; i++){
m*=2;
// S(m) when unfolded looks like [ 1 1 1 ... 1 ]
// [ w_m w_m^2 w_m^3 ... w_m^{m/2} ]
// it is used to scale the odd elements of the FFT, i.e., we compute recursively FFT(v) via
// for i=1 to m/2, A_i = sum_{j=0}^{m/2-1} v_{2j} *w_{m/2}^{ij}
// B_i = w_m^i sum_{j=0}^{m/2-1} v_{2j+1}*w_{m/2}^{ij}
// multiplying V by S will do the scaling of the elements of B_i by w_m^i
// we do it first, since we want to apply it to the sequences obtain recursively
// note that for i=0, m=2 the application of S(2) = [1; 1] is trivial
Tensor< std::complex<double> > S;
if (i==0){
S = S_old;
} else {
// we construct S(m) by the Kharti-Rao product of S(m/2) and the twiddle factor [1, 1; 1, w_m]
S = Tensor< std::complex<double> >(i+1, factors, *v.wrld, *v.sr, "S");
Matrix< std::complex<double> > T = twiddle_matrix(m, *v.wrld);
char inds_T[2] = {'a', inds[i]};
S[inds] = T[inds_T]*S_old[inds+1];
S_old = S;
}
V[inds] *= S[inds];
// then we multiply by the 2-by-2 DFT matrix, which is just [1, 1; 1, -1]
// this gives FFT(v) = [A+B; A-B]
// this combines the recursive sequences or does the DFT at the base case when i=0
memcpy(inds_in, inds, nfact*sizeof(char));
inds_in[i] = 'a'-1;
inds_DFT[0] = inds[i];
inds_DFT[1] = 'a'-1;
V[inds] = DFT[inds_DFT]*V[inds_in];
}
// we now unfold the tensor back into vector form
v.reshape(V);
}
int fft(int n,
World & dw){
Vector< std::complex<double> > a(n, dw, "a");
Vector< std::complex<double> > da(n, dw, "da");
srand48(2*dw.rank);
Transform< std::complex<double> > set_rand([](std::complex<double> & d){ d=std::complex<double>(drand48(),drand48()); } );
set_rand(a["i"]);
Matrix< std::complex<double> > DFT = DFT_matrix(n,dw);
da["i"] = DFT["ij"]*a["j"];
Vector< std::complex<double> > fa(n, dw, "fa");
fa["i"] = a["i"];
fft(fa, n);
da["i"] -= fa["i"];
std::complex<double> dnrm;
Scalar< std::complex<double> > s(dnrm, dw);
s[""] += da["i"]*da["i"];
dnrm = s;
bool pass = fabs(dnrm.real()) <= 1.E-6 && fabs(dnrm.imag()) <= 1.E-6;
if (dw.rank == 0){
if (pass)
printf("{ FFT = DFT } passed \n");
else
printf("{ FFT = DFT } failed \n");
}
return pass;
}
#ifndef TEST_SUITE
char* getCmdOption(char ** begin,
char ** end,
const std::string & option){
char ** itr = std::find(begin, end, option);
if (itr != end && ++itr != end){
return *itr;
}
return 0;
}
int main(int argc, char ** argv){
int rank, np, n, pass;
int const in_num = argc;
char ** input_str = argv;
MPI_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &np);
if (getCmdOption(input_str, input_str+in_num, "-n")){
n = atoi(getCmdOption(input_str, input_str+in_num, "-n"));
if (n < 0) n = 16;
} else n = 16;
{
World dw(argc, argv);
if (rank == 0){
printf("Running FFT on random dimension %d vector\n",n);
}
pass = fft(n, dw);
assert(pass);
}
MPI_Finalize();
return 0;
}
/**
* @}
* @}
*/
#endif