1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Mark Borgerding mark a borgerding net
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10 #include "main.h"
11 #include <unsupported/Eigen/FFT>
12
13 template <typename T>
RandomCpx()14 std::complex<T> RandomCpx() { return std::complex<T>( (T)(rand()/(T)RAND_MAX - .5), (T)(rand()/(T)RAND_MAX - .5) ); }
15
16 using namespace std;
17 using namespace Eigen;
18
norm(float x)19 float norm(float x) {return x*x;}
norm(double x)20 double norm(double x) {return x*x;}
norm(long double x)21 long double norm(long double x) {return x*x;}
22
23 template < typename T>
promote(complex<T> x)24 complex<long double> promote(complex<T> x) { return complex<long double>(x.real(),x.imag()); }
25
promote(float x)26 complex<long double> promote(float x) { return complex<long double>( x); }
promote(double x)27 complex<long double> promote(double x) { return complex<long double>( x); }
promote(long double x)28 complex<long double> promote(long double x) { return complex<long double>( x); }
29
30
31 template <typename VT1,typename VT2>
fft_rmse(const VT1 & fftbuf,const VT2 & timebuf)32 long double fft_rmse( const VT1 & fftbuf,const VT2 & timebuf)
33 {
34 long double totalpower=0;
35 long double difpower=0;
36 long double pi = acos((long double)-1 );
37 for (size_t k0=0;k0<(size_t)fftbuf.size();++k0) {
38 complex<long double> acc = 0;
39 long double phinc = -2.*k0* pi / timebuf.size();
40 for (size_t k1=0;k1<(size_t)timebuf.size();++k1) {
41 acc += promote( timebuf[k1] ) * exp( complex<long double>(0,k1*phinc) );
42 }
43 totalpower += norm(acc);
44 complex<long double> x = promote(fftbuf[k0]);
45 complex<long double> dif = acc - x;
46 difpower += norm(dif);
47 //cerr << k0 << "\t" << acc << "\t" << x << "\t" << sqrt(norm(dif)) << endl;
48 }
49 cerr << "rmse:" << sqrt(difpower/totalpower) << endl;
50 return sqrt(difpower/totalpower);
51 }
52
53 template <typename VT1,typename VT2>
dif_rmse(const VT1 buf1,const VT2 buf2)54 long double dif_rmse( const VT1 buf1,const VT2 buf2)
55 {
56 long double totalpower=0;
57 long double difpower=0;
58 size_t n = (min)( buf1.size(),buf2.size() );
59 for (size_t k=0;k<n;++k) {
60 totalpower += (norm( buf1[k] ) + norm(buf2[k]) )/2.;
61 difpower += norm(buf1[k] - buf2[k]);
62 }
63 return sqrt(difpower/totalpower);
64 }
65
66 enum { StdVectorContainer, EigenVectorContainer };
67
68 template<int Container, typename Scalar> struct VectorType;
69
70 template<typename Scalar> struct VectorType<StdVectorContainer,Scalar>
71 {
72 typedef vector<Scalar> type;
73 };
74
75 template<typename Scalar> struct VectorType<EigenVectorContainer,Scalar>
76 {
77 typedef Matrix<Scalar,Dynamic,1> type;
78 };
79
80 template <int Container, typename T>
test_scalar_generic(int nfft)81 void test_scalar_generic(int nfft)
82 {
83 typedef typename FFT<T>::Complex Complex;
84 typedef typename FFT<T>::Scalar Scalar;
85 typedef typename VectorType<Container,Scalar>::type ScalarVector;
86 typedef typename VectorType<Container,Complex>::type ComplexVector;
87
88 FFT<T> fft;
89 ScalarVector tbuf(nfft);
90 ComplexVector freqBuf;
91 for (int k=0;k<nfft;++k)
92 tbuf[k]= (T)( rand()/(double)RAND_MAX - .5);
93
94 // make sure it DOESN'T give the right full spectrum answer
95 // if we've asked for half-spectrum
96 fft.SetFlag(fft.HalfSpectrum );
97 fft.fwd( freqBuf,tbuf);
98 VERIFY((size_t)freqBuf.size() == (size_t)( (nfft>>1)+1) );
99 VERIFY( fft_rmse(freqBuf,tbuf) < test_precision<T>() );// gross check
100
101 fft.ClearFlag(fft.HalfSpectrum );
102 fft.fwd( freqBuf,tbuf);
103 VERIFY( (size_t)freqBuf.size() == (size_t)nfft);
104 VERIFY( fft_rmse(freqBuf,tbuf) < test_precision<T>() );// gross check
105
106 if (nfft&1)
107 return; // odd FFTs get the wrong size inverse FFT
108
109 ScalarVector tbuf2;
110 fft.inv( tbuf2 , freqBuf);
111 VERIFY( dif_rmse(tbuf,tbuf2) < test_precision<T>() );// gross check
112
113
114 // verify that the Unscaled flag takes effect
115 ScalarVector tbuf3;
116 fft.SetFlag(fft.Unscaled);
117
118 fft.inv( tbuf3 , freqBuf);
119
120 for (int k=0;k<nfft;++k)
121 tbuf3[k] *= T(1./nfft);
122
123
124 //for (size_t i=0;i<(size_t) tbuf.size();++i)
125 // cout << "freqBuf=" << freqBuf[i] << " in2=" << tbuf3[i] << " - in=" << tbuf[i] << " => " << (tbuf3[i] - tbuf[i] ) << endl;
126
127 VERIFY( dif_rmse(tbuf,tbuf3) < test_precision<T>() );// gross check
128
129 // verify that ClearFlag works
130 fft.ClearFlag(fft.Unscaled);
131 fft.inv( tbuf2 , freqBuf);
132 VERIFY( dif_rmse(tbuf,tbuf2) < test_precision<T>() );// gross check
133 }
134
135 template <typename T>
test_scalar(int nfft)136 void test_scalar(int nfft)
137 {
138 test_scalar_generic<StdVectorContainer,T>(nfft);
139 //test_scalar_generic<EigenVectorContainer,T>(nfft);
140 }
141
142
143 template <int Container, typename T>
test_complex_generic(int nfft)144 void test_complex_generic(int nfft)
145 {
146 typedef typename FFT<T>::Complex Complex;
147 typedef typename VectorType<Container,Complex>::type ComplexVector;
148
149 FFT<T> fft;
150
151 ComplexVector inbuf(nfft);
152 ComplexVector outbuf;
153 ComplexVector buf3;
154 for (int k=0;k<nfft;++k)
155 inbuf[k]= Complex( (T)(rand()/(double)RAND_MAX - .5), (T)(rand()/(double)RAND_MAX - .5) );
156 fft.fwd( outbuf , inbuf);
157
158 VERIFY( fft_rmse(outbuf,inbuf) < test_precision<T>() );// gross check
159 fft.inv( buf3 , outbuf);
160
161 VERIFY( dif_rmse(inbuf,buf3) < test_precision<T>() );// gross check
162
163 // verify that the Unscaled flag takes effect
164 ComplexVector buf4;
165 fft.SetFlag(fft.Unscaled);
166 fft.inv( buf4 , outbuf);
167 for (int k=0;k<nfft;++k)
168 buf4[k] *= T(1./nfft);
169 VERIFY( dif_rmse(inbuf,buf4) < test_precision<T>() );// gross check
170
171 // verify that ClearFlag works
172 fft.ClearFlag(fft.Unscaled);
173 fft.inv( buf3 , outbuf);
174 VERIFY( dif_rmse(inbuf,buf3) < test_precision<T>() );// gross check
175 }
176
177 template <typename T>
test_complex(int nfft)178 void test_complex(int nfft)
179 {
180 test_complex_generic<StdVectorContainer,T>(nfft);
181 test_complex_generic<EigenVectorContainer,T>(nfft);
182 }
183 /*
184 template <typename T,int nrows,int ncols>
185 void test_complex2d()
186 {
187 typedef typename Eigen::FFT<T>::Complex Complex;
188 FFT<T> fft;
189 Eigen::Matrix<Complex,nrows,ncols> src,src2,dst,dst2;
190
191 src = Eigen::Matrix<Complex,nrows,ncols>::Random();
192 //src = Eigen::Matrix<Complex,nrows,ncols>::Identity();
193
194 for (int k=0;k<ncols;k++) {
195 Eigen::Matrix<Complex,nrows,1> tmpOut;
196 fft.fwd( tmpOut,src.col(k) );
197 dst2.col(k) = tmpOut;
198 }
199
200 for (int k=0;k<nrows;k++) {
201 Eigen::Matrix<Complex,1,ncols> tmpOut;
202 fft.fwd( tmpOut, dst2.row(k) );
203 dst2.row(k) = tmpOut;
204 }
205
206 fft.fwd2(dst.data(),src.data(),ncols,nrows);
207 fft.inv2(src2.data(),dst.data(),ncols,nrows);
208 VERIFY( (src-src2).norm() < test_precision<T>() );
209 VERIFY( (dst-dst2).norm() < test_precision<T>() );
210 }
211 */
212
213
test_return_by_value(int len)214 void test_return_by_value(int len)
215 {
216 VectorXf in;
217 VectorXf in1;
218 in.setRandom( len );
219 VectorXcf out1,out2;
220 FFT<float> fft;
221
222 fft.SetFlag(fft.HalfSpectrum );
223
224 fft.fwd(out1,in);
225 out2 = fft.fwd(in);
226 VERIFY( (out1-out2).norm() < test_precision<float>() );
227 in1 = fft.inv(out1);
228 VERIFY( (in1-in).norm() < test_precision<float>() );
229 }
230
test_FFTW()231 void test_FFTW()
232 {
233 CALL_SUBTEST( test_return_by_value(32) );
234 //CALL_SUBTEST( ( test_complex2d<float,4,8> () ) ); CALL_SUBTEST( ( test_complex2d<double,4,8> () ) );
235 //CALL_SUBTEST( ( test_complex2d<long double,4,8> () ) );
236 CALL_SUBTEST( test_complex<float>(32) ); CALL_SUBTEST( test_complex<double>(32) );
237 CALL_SUBTEST( test_complex<float>(256) ); CALL_SUBTEST( test_complex<double>(256) );
238 CALL_SUBTEST( test_complex<float>(3*8) ); CALL_SUBTEST( test_complex<double>(3*8) );
239 CALL_SUBTEST( test_complex<float>(5*32) ); CALL_SUBTEST( test_complex<double>(5*32) );
240 CALL_SUBTEST( test_complex<float>(2*3*4) ); CALL_SUBTEST( test_complex<double>(2*3*4) );
241 CALL_SUBTEST( test_complex<float>(2*3*4*5) ); CALL_SUBTEST( test_complex<double>(2*3*4*5) );
242 CALL_SUBTEST( test_complex<float>(2*3*4*5*7) ); CALL_SUBTEST( test_complex<double>(2*3*4*5*7) );
243
244 CALL_SUBTEST( test_scalar<float>(32) ); CALL_SUBTEST( test_scalar<double>(32) );
245 CALL_SUBTEST( test_scalar<float>(45) ); CALL_SUBTEST( test_scalar<double>(45) );
246 CALL_SUBTEST( test_scalar<float>(50) ); CALL_SUBTEST( test_scalar<double>(50) );
247 CALL_SUBTEST( test_scalar<float>(256) ); CALL_SUBTEST( test_scalar<double>(256) );
248 CALL_SUBTEST( test_scalar<float>(2*3*4*5*7) ); CALL_SUBTEST( test_scalar<double>(2*3*4*5*7) );
249
250 #ifdef EIGEN_HAS_FFTWL
251 CALL_SUBTEST( test_complex<long double>(32) );
252 CALL_SUBTEST( test_complex<long double>(256) );
253 CALL_SUBTEST( test_complex<long double>(3*8) );
254 CALL_SUBTEST( test_complex<long double>(5*32) );
255 CALL_SUBTEST( test_complex<long double>(2*3*4) );
256 CALL_SUBTEST( test_complex<long double>(2*3*4*5) );
257 CALL_SUBTEST( test_complex<long double>(2*3*4*5*7) );
258
259 CALL_SUBTEST( test_scalar<long double>(32) );
260 CALL_SUBTEST( test_scalar<long double>(45) );
261 CALL_SUBTEST( test_scalar<long double>(50) );
262 CALL_SUBTEST( test_scalar<long double>(256) );
263 CALL_SUBTEST( test_scalar<long double>(2*3*4*5*7) );
264 #endif
265 }
266