You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

265 lines
8.9 KiB

  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. #include "main.h"
  10. #include <unsupported/Eigen/FFT>
  11. template <typename T>
  12. std::complex<T> RandomCpx() { return std::complex<T>( (T)(rand()/(T)RAND_MAX - .5), (T)(rand()/(T)RAND_MAX - .5) ); }
  13. using namespace std;
  14. using namespace Eigen;
  15. float norm(float x) {return x*x;}
  16. double norm(double x) {return x*x;}
  17. long double norm(long double x) {return x*x;}
  18. template < typename T>
  19. complex<long double> promote(complex<T> x) { return complex<long double>(x.real(),x.imag()); }
  20. complex<long double> promote(float x) { return complex<long double>( x); }
  21. complex<long double> promote(double x) { return complex<long double>( x); }
  22. complex<long double> promote(long double x) { return complex<long double>( x); }
  23. template <typename VT1,typename VT2>
  24. long double fft_rmse( const VT1 & fftbuf,const VT2 & timebuf)
  25. {
  26. long double totalpower=0;
  27. long double difpower=0;
  28. long double pi = acos((long double)-1 );
  29. for (size_t k0=0;k0<(size_t)fftbuf.size();++k0) {
  30. complex<long double> acc = 0;
  31. long double phinc = -2.*k0* pi / timebuf.size();
  32. for (size_t k1=0;k1<(size_t)timebuf.size();++k1) {
  33. acc += promote( timebuf[k1] ) * exp( complex<long double>(0,k1*phinc) );
  34. }
  35. totalpower += norm(acc);
  36. complex<long double> x = promote(fftbuf[k0]);
  37. complex<long double> dif = acc - x;
  38. difpower += norm(dif);
  39. //cerr << k0 << "\t" << acc << "\t" << x << "\t" << sqrt(norm(dif)) << endl;
  40. }
  41. cerr << "rmse:" << sqrt(difpower/totalpower) << endl;
  42. return sqrt(difpower/totalpower);
  43. }
  44. template <typename VT1,typename VT2>
  45. long double dif_rmse( const VT1 buf1,const VT2 buf2)
  46. {
  47. long double totalpower=0;
  48. long double difpower=0;
  49. size_t n = (min)( buf1.size(),buf2.size() );
  50. for (size_t k=0;k<n;++k) {
  51. totalpower += (norm( buf1[k] ) + norm(buf2[k]) )/2.;
  52. difpower += norm(buf1[k] - buf2[k]);
  53. }
  54. return sqrt(difpower/totalpower);
  55. }
  56. enum { StdVectorContainer, EigenVectorContainer };
  57. template<int Container, typename Scalar> struct VectorType;
  58. template<typename Scalar> struct VectorType<StdVectorContainer,Scalar>
  59. {
  60. typedef vector<Scalar> type;
  61. };
  62. template<typename Scalar> struct VectorType<EigenVectorContainer,Scalar>
  63. {
  64. typedef Matrix<Scalar,Dynamic,1> type;
  65. };
  66. template <int Container, typename T>
  67. void test_scalar_generic(int nfft)
  68. {
  69. typedef typename FFT<T>::Complex Complex;
  70. typedef typename FFT<T>::Scalar Scalar;
  71. typedef typename VectorType<Container,Scalar>::type ScalarVector;
  72. typedef typename VectorType<Container,Complex>::type ComplexVector;
  73. FFT<T> fft;
  74. ScalarVector tbuf(nfft);
  75. ComplexVector freqBuf;
  76. for (int k=0;k<nfft;++k)
  77. tbuf[k]= (T)( rand()/(double)RAND_MAX - .5);
  78. // make sure it DOESN'T give the right full spectrum answer
  79. // if we've asked for half-spectrum
  80. fft.SetFlag(fft.HalfSpectrum );
  81. fft.fwd( freqBuf,tbuf);
  82. VERIFY((size_t)freqBuf.size() == (size_t)( (nfft>>1)+1) );
  83. VERIFY( fft_rmse(freqBuf,tbuf) < test_precision<T>() );// gross check
  84. fft.ClearFlag(fft.HalfSpectrum );
  85. fft.fwd( freqBuf,tbuf);
  86. VERIFY( (size_t)freqBuf.size() == (size_t)nfft);
  87. VERIFY( fft_rmse(freqBuf,tbuf) < test_precision<T>() );// gross check
  88. if (nfft&1)
  89. return; // odd FFTs get the wrong size inverse FFT
  90. ScalarVector tbuf2;
  91. fft.inv( tbuf2 , freqBuf);
  92. VERIFY( dif_rmse(tbuf,tbuf2) < test_precision<T>() );// gross check
  93. // verify that the Unscaled flag takes effect
  94. ScalarVector tbuf3;
  95. fft.SetFlag(fft.Unscaled);
  96. fft.inv( tbuf3 , freqBuf);
  97. for (int k=0;k<nfft;++k)
  98. tbuf3[k] *= T(1./nfft);
  99. //for (size_t i=0;i<(size_t) tbuf.size();++i)
  100. // cout << "freqBuf=" << freqBuf[i] << " in2=" << tbuf3[i] << " - in=" << tbuf[i] << " => " << (tbuf3[i] - tbuf[i] ) << endl;
  101. VERIFY( dif_rmse(tbuf,tbuf3) < test_precision<T>() );// gross check
  102. // verify that ClearFlag works
  103. fft.ClearFlag(fft.Unscaled);
  104. fft.inv( tbuf2 , freqBuf);
  105. VERIFY( dif_rmse(tbuf,tbuf2) < test_precision<T>() );// gross check
  106. }
  107. template <typename T>
  108. void test_scalar(int nfft)
  109. {
  110. test_scalar_generic<StdVectorContainer,T>(nfft);
  111. //test_scalar_generic<EigenVectorContainer,T>(nfft);
  112. }
  113. template <int Container, typename T>
  114. void test_complex_generic(int nfft)
  115. {
  116. typedef typename FFT<T>::Complex Complex;
  117. typedef typename VectorType<Container,Complex>::type ComplexVector;
  118. FFT<T> fft;
  119. ComplexVector inbuf(nfft);
  120. ComplexVector outbuf;
  121. ComplexVector buf3;
  122. for (int k=0;k<nfft;++k)
  123. inbuf[k]= Complex( (T)(rand()/(double)RAND_MAX - .5), (T)(rand()/(double)RAND_MAX - .5) );
  124. fft.fwd( outbuf , inbuf);
  125. VERIFY( fft_rmse(outbuf,inbuf) < test_precision<T>() );// gross check
  126. fft.inv( buf3 , outbuf);
  127. VERIFY( dif_rmse(inbuf,buf3) < test_precision<T>() );// gross check
  128. // verify that the Unscaled flag takes effect
  129. ComplexVector buf4;
  130. fft.SetFlag(fft.Unscaled);
  131. fft.inv( buf4 , outbuf);
  132. for (int k=0;k<nfft;++k)
  133. buf4[k] *= T(1./nfft);
  134. VERIFY( dif_rmse(inbuf,buf4) < test_precision<T>() );// gross check
  135. // verify that ClearFlag works
  136. fft.ClearFlag(fft.Unscaled);
  137. fft.inv( buf3 , outbuf);
  138. VERIFY( dif_rmse(inbuf,buf3) < test_precision<T>() );// gross check
  139. }
  140. template <typename T>
  141. void test_complex(int nfft)
  142. {
  143. test_complex_generic<StdVectorContainer,T>(nfft);
  144. test_complex_generic<EigenVectorContainer,T>(nfft);
  145. }
  146. /*
  147. template <typename T,int nrows,int ncols>
  148. void test_complex2d()
  149. {
  150. typedef typename Eigen::FFT<T>::Complex Complex;
  151. FFT<T> fft;
  152. Eigen::Matrix<Complex,nrows,ncols> src,src2,dst,dst2;
  153. src = Eigen::Matrix<Complex,nrows,ncols>::Random();
  154. //src = Eigen::Matrix<Complex,nrows,ncols>::Identity();
  155. for (int k=0;k<ncols;k++) {
  156. Eigen::Matrix<Complex,nrows,1> tmpOut;
  157. fft.fwd( tmpOut,src.col(k) );
  158. dst2.col(k) = tmpOut;
  159. }
  160. for (int k=0;k<nrows;k++) {
  161. Eigen::Matrix<Complex,1,ncols> tmpOut;
  162. fft.fwd( tmpOut, dst2.row(k) );
  163. dst2.row(k) = tmpOut;
  164. }
  165. fft.fwd2(dst.data(),src.data(),ncols,nrows);
  166. fft.inv2(src2.data(),dst.data(),ncols,nrows);
  167. VERIFY( (src-src2).norm() < test_precision<T>() );
  168. VERIFY( (dst-dst2).norm() < test_precision<T>() );
  169. }
  170. */
  171. void test_return_by_value(int len)
  172. {
  173. VectorXf in;
  174. VectorXf in1;
  175. in.setRandom( len );
  176. VectorXcf out1,out2;
  177. FFT<float> fft;
  178. fft.SetFlag(fft.HalfSpectrum );
  179. fft.fwd(out1,in);
  180. out2 = fft.fwd(in);
  181. VERIFY( (out1-out2).norm() < test_precision<float>() );
  182. in1 = fft.inv(out1);
  183. VERIFY( (in1-in).norm() < test_precision<float>() );
  184. }
  185. void test_FFTW()
  186. {
  187. CALL_SUBTEST( test_return_by_value(32) );
  188. //CALL_SUBTEST( ( test_complex2d<float,4,8> () ) ); CALL_SUBTEST( ( test_complex2d<double,4,8> () ) );
  189. //CALL_SUBTEST( ( test_complex2d<long double,4,8> () ) );
  190. CALL_SUBTEST( test_complex<float>(32) ); CALL_SUBTEST( test_complex<double>(32) );
  191. CALL_SUBTEST( test_complex<float>(256) ); CALL_SUBTEST( test_complex<double>(256) );
  192. CALL_SUBTEST( test_complex<float>(3*8) ); CALL_SUBTEST( test_complex<double>(3*8) );
  193. CALL_SUBTEST( test_complex<float>(5*32) ); CALL_SUBTEST( test_complex<double>(5*32) );
  194. CALL_SUBTEST( test_complex<float>(2*3*4) ); CALL_SUBTEST( test_complex<double>(2*3*4) );
  195. CALL_SUBTEST( test_complex<float>(2*3*4*5) ); CALL_SUBTEST( test_complex<double>(2*3*4*5) );
  196. CALL_SUBTEST( test_complex<float>(2*3*4*5*7) ); CALL_SUBTEST( test_complex<double>(2*3*4*5*7) );
  197. CALL_SUBTEST( test_scalar<float>(32) ); CALL_SUBTEST( test_scalar<double>(32) );
  198. CALL_SUBTEST( test_scalar<float>(45) ); CALL_SUBTEST( test_scalar<double>(45) );
  199. CALL_SUBTEST( test_scalar<float>(50) ); CALL_SUBTEST( test_scalar<double>(50) );
  200. CALL_SUBTEST( test_scalar<float>(256) ); CALL_SUBTEST( test_scalar<double>(256) );
  201. CALL_SUBTEST( test_scalar<float>(2*3*4*5*7) ); CALL_SUBTEST( test_scalar<double>(2*3*4*5*7) );
  202. #ifdef EIGEN_HAS_FFTWL
  203. CALL_SUBTEST( test_complex<long double>(32) );
  204. CALL_SUBTEST( test_complex<long double>(256) );
  205. CALL_SUBTEST( test_complex<long double>(3*8) );
  206. CALL_SUBTEST( test_complex<long double>(5*32) );
  207. CALL_SUBTEST( test_complex<long double>(2*3*4) );
  208. CALL_SUBTEST( test_complex<long double>(2*3*4*5) );
  209. CALL_SUBTEST( test_complex<long double>(2*3*4*5*7) );
  210. CALL_SUBTEST( test_scalar<long double>(32) );
  211. CALL_SUBTEST( test_scalar<long double>(45) );
  212. CALL_SUBTEST( test_scalar<long double>(50) );
  213. CALL_SUBTEST( test_scalar<long double>(256) );
  214. CALL_SUBTEST( test_scalar<long double>(2*3*4*5*7) );
  215. #endif
  216. }