test_pseudo_diffs.py 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380
  1. # Created by Pearu Peterson, September 2002
  2. __usage__ = """
  3. Build fftpack:
  4. python setup_fftpack.py build
  5. Run tests if scipy is installed:
  6. python -c 'import scipy;scipy.fftpack.test(<level>)'
  7. Run tests if fftpack is not installed:
  8. python tests/test_pseudo_diffs.py [<level>]
  9. """
  10. from numpy.testing import (assert_equal, assert_almost_equal,
  11. assert_array_almost_equal)
  12. from scipy.fftpack import (diff, fft, ifft, tilbert, itilbert, hilbert,
  13. ihilbert, shift, fftfreq, cs_diff, sc_diff,
  14. ss_diff, cc_diff)
  15. import numpy as np
  16. from numpy import arange, sin, cos, pi, exp, tanh, sum, sign
  17. from numpy.random import random
  18. def direct_diff(x,k=1,period=None):
  19. fx = fft(x)
  20. n = len(fx)
  21. if period is None:
  22. period = 2*pi
  23. w = fftfreq(n)*2j*pi/period*n
  24. if k < 0:
  25. w = 1 / w**k
  26. w[0] = 0.0
  27. else:
  28. w = w**k
  29. if n > 2000:
  30. w[250:n-250] = 0.0
  31. return ifft(w*fx).real
  32. def direct_tilbert(x,h=1,period=None):
  33. fx = fft(x)
  34. n = len(fx)
  35. if period is None:
  36. period = 2*pi
  37. w = fftfreq(n)*h*2*pi/period*n
  38. w[0] = 1
  39. w = 1j/tanh(w)
  40. w[0] = 0j
  41. return ifft(w*fx)
  42. def direct_itilbert(x,h=1,period=None):
  43. fx = fft(x)
  44. n = len(fx)
  45. if period is None:
  46. period = 2*pi
  47. w = fftfreq(n)*h*2*pi/period*n
  48. w = -1j*tanh(w)
  49. return ifft(w*fx)
  50. def direct_hilbert(x):
  51. fx = fft(x)
  52. n = len(fx)
  53. w = fftfreq(n)*n
  54. w = 1j*sign(w)
  55. return ifft(w*fx)
  56. def direct_ihilbert(x):
  57. return -direct_hilbert(x)
  58. def direct_shift(x,a,period=None):
  59. n = len(x)
  60. if period is None:
  61. k = fftfreq(n)*1j*n
  62. else:
  63. k = fftfreq(n)*2j*pi/period*n
  64. return ifft(fft(x)*exp(k*a)).real
  65. class TestDiff:
  66. def test_definition(self):
  67. for n in [16,17,64,127,32]:
  68. x = arange(n)*2*pi/n
  69. assert_array_almost_equal(diff(sin(x)),direct_diff(sin(x)))
  70. assert_array_almost_equal(diff(sin(x),2),direct_diff(sin(x),2))
  71. assert_array_almost_equal(diff(sin(x),3),direct_diff(sin(x),3))
  72. assert_array_almost_equal(diff(sin(x),4),direct_diff(sin(x),4))
  73. assert_array_almost_equal(diff(sin(x),5),direct_diff(sin(x),5))
  74. assert_array_almost_equal(diff(sin(2*x),3),direct_diff(sin(2*x),3))
  75. assert_array_almost_equal(diff(sin(2*x),4),direct_diff(sin(2*x),4))
  76. assert_array_almost_equal(diff(cos(x)),direct_diff(cos(x)))
  77. assert_array_almost_equal(diff(cos(x),2),direct_diff(cos(x),2))
  78. assert_array_almost_equal(diff(cos(x),3),direct_diff(cos(x),3))
  79. assert_array_almost_equal(diff(cos(x),4),direct_diff(cos(x),4))
  80. assert_array_almost_equal(diff(cos(2*x)),direct_diff(cos(2*x)))
  81. assert_array_almost_equal(diff(sin(x*n/8)),direct_diff(sin(x*n/8)))
  82. assert_array_almost_equal(diff(cos(x*n/8)),direct_diff(cos(x*n/8)))
  83. for k in range(5):
  84. assert_array_almost_equal(diff(sin(4*x),k),direct_diff(sin(4*x),k))
  85. assert_array_almost_equal(diff(cos(4*x),k),direct_diff(cos(4*x),k))
  86. def test_period(self):
  87. for n in [17,64]:
  88. x = arange(n)/float(n)
  89. assert_array_almost_equal(diff(sin(2*pi*x),period=1),
  90. 2*pi*cos(2*pi*x))
  91. assert_array_almost_equal(diff(sin(2*pi*x),3,period=1),
  92. -(2*pi)**3*cos(2*pi*x))
  93. def test_sin(self):
  94. for n in [32,64,77]:
  95. x = arange(n)*2*pi/n
  96. assert_array_almost_equal(diff(sin(x)),cos(x))
  97. assert_array_almost_equal(diff(cos(x)),-sin(x))
  98. assert_array_almost_equal(diff(sin(x),2),-sin(x))
  99. assert_array_almost_equal(diff(sin(x),4),sin(x))
  100. assert_array_almost_equal(diff(sin(4*x)),4*cos(4*x))
  101. assert_array_almost_equal(diff(sin(sin(x))),cos(x)*cos(sin(x)))
  102. def test_expr(self):
  103. for n in [64,77,100,128,256,512,1024,2048,4096,8192][:5]:
  104. x = arange(n)*2*pi/n
  105. f = sin(x)*cos(4*x)+exp(sin(3*x))
  106. df = cos(x)*cos(4*x)-4*sin(x)*sin(4*x)+3*cos(3*x)*exp(sin(3*x))
  107. ddf = -17*sin(x)*cos(4*x)-8*cos(x)*sin(4*x)\
  108. - 9*sin(3*x)*exp(sin(3*x))+9*cos(3*x)**2*exp(sin(3*x))
  109. d1 = diff(f)
  110. assert_array_almost_equal(d1,df)
  111. assert_array_almost_equal(diff(df),ddf)
  112. assert_array_almost_equal(diff(f,2),ddf)
  113. assert_array_almost_equal(diff(ddf,-1),df)
  114. def test_expr_large(self):
  115. for n in [2048,4096]:
  116. x = arange(n)*2*pi/n
  117. f = sin(x)*cos(4*x)+exp(sin(3*x))
  118. df = cos(x)*cos(4*x)-4*sin(x)*sin(4*x)+3*cos(3*x)*exp(sin(3*x))
  119. ddf = -17*sin(x)*cos(4*x)-8*cos(x)*sin(4*x)\
  120. - 9*sin(3*x)*exp(sin(3*x))+9*cos(3*x)**2*exp(sin(3*x))
  121. assert_array_almost_equal(diff(f),df)
  122. assert_array_almost_equal(diff(df),ddf)
  123. assert_array_almost_equal(diff(ddf,-1),df)
  124. assert_array_almost_equal(diff(f,2),ddf)
  125. def test_int(self):
  126. n = 64
  127. x = arange(n)*2*pi/n
  128. assert_array_almost_equal(diff(sin(x),-1),-cos(x))
  129. assert_array_almost_equal(diff(sin(x),-2),-sin(x))
  130. assert_array_almost_equal(diff(sin(x),-4),sin(x))
  131. assert_array_almost_equal(diff(2*cos(2*x),-1),sin(2*x))
  132. def test_random_even(self):
  133. for k in [0,2,4,6]:
  134. for n in [60,32,64,56,55]:
  135. f = random((n,))
  136. af = sum(f,axis=0)/n
  137. f = f-af
  138. # zeroing Nyquist mode:
  139. f = diff(diff(f,1),-1)
  140. assert_almost_equal(sum(f,axis=0),0.0)
  141. assert_array_almost_equal(diff(diff(f,k),-k),f)
  142. assert_array_almost_equal(diff(diff(f,-k),k),f)
  143. def test_random_odd(self):
  144. for k in [0,1,2,3,4,5,6]:
  145. for n in [33,65,55]:
  146. f = random((n,))
  147. af = sum(f,axis=0)/n
  148. f = f-af
  149. assert_almost_equal(sum(f,axis=0),0.0)
  150. assert_array_almost_equal(diff(diff(f,k),-k),f)
  151. assert_array_almost_equal(diff(diff(f,-k),k),f)
  152. def test_zero_nyquist(self):
  153. for k in [0,1,2,3,4,5,6]:
  154. for n in [32,33,64,56,55]:
  155. f = random((n,))
  156. af = sum(f,axis=0)/n
  157. f = f-af
  158. # zeroing Nyquist mode:
  159. f = diff(diff(f,1),-1)
  160. assert_almost_equal(sum(f,axis=0),0.0)
  161. assert_array_almost_equal(diff(diff(f,k),-k),f)
  162. assert_array_almost_equal(diff(diff(f,-k),k),f)
  163. class TestTilbert:
  164. def test_definition(self):
  165. for h in [0.1,0.5,1,5.5,10]:
  166. for n in [16,17,64,127]:
  167. x = arange(n)*2*pi/n
  168. y = tilbert(sin(x),h)
  169. y1 = direct_tilbert(sin(x),h)
  170. assert_array_almost_equal(y,y1)
  171. assert_array_almost_equal(tilbert(sin(x),h),
  172. direct_tilbert(sin(x),h))
  173. assert_array_almost_equal(tilbert(sin(2*x),h),
  174. direct_tilbert(sin(2*x),h))
  175. def test_random_even(self):
  176. for h in [0.1,0.5,1,5.5,10]:
  177. for n in [32,64,56]:
  178. f = random((n,))
  179. af = sum(f,axis=0)/n
  180. f = f-af
  181. assert_almost_equal(sum(f,axis=0),0.0)
  182. assert_array_almost_equal(direct_tilbert(direct_itilbert(f,h),h),f)
  183. def test_random_odd(self):
  184. for h in [0.1,0.5,1,5.5,10]:
  185. for n in [33,65,55]:
  186. f = random((n,))
  187. af = sum(f,axis=0)/n
  188. f = f-af
  189. assert_almost_equal(sum(f,axis=0),0.0)
  190. assert_array_almost_equal(itilbert(tilbert(f,h),h),f)
  191. assert_array_almost_equal(tilbert(itilbert(f,h),h),f)
  192. class TestITilbert:
  193. def test_definition(self):
  194. for h in [0.1,0.5,1,5.5,10]:
  195. for n in [16,17,64,127]:
  196. x = arange(n)*2*pi/n
  197. y = itilbert(sin(x),h)
  198. y1 = direct_itilbert(sin(x),h)
  199. assert_array_almost_equal(y,y1)
  200. assert_array_almost_equal(itilbert(sin(x),h),
  201. direct_itilbert(sin(x),h))
  202. assert_array_almost_equal(itilbert(sin(2*x),h),
  203. direct_itilbert(sin(2*x),h))
  204. class TestHilbert:
  205. def test_definition(self):
  206. for n in [16,17,64,127]:
  207. x = arange(n)*2*pi/n
  208. y = hilbert(sin(x))
  209. y1 = direct_hilbert(sin(x))
  210. assert_array_almost_equal(y,y1)
  211. assert_array_almost_equal(hilbert(sin(2*x)),
  212. direct_hilbert(sin(2*x)))
  213. def test_tilbert_relation(self):
  214. for n in [16,17,64,127]:
  215. x = arange(n)*2*pi/n
  216. f = sin(x)+cos(2*x)*sin(x)
  217. y = hilbert(f)
  218. y1 = direct_hilbert(f)
  219. assert_array_almost_equal(y,y1)
  220. y2 = tilbert(f,h=10)
  221. assert_array_almost_equal(y,y2)
  222. def test_random_odd(self):
  223. for n in [33,65,55]:
  224. f = random((n,))
  225. af = sum(f,axis=0)/n
  226. f = f-af
  227. assert_almost_equal(sum(f,axis=0),0.0)
  228. assert_array_almost_equal(ihilbert(hilbert(f)),f)
  229. assert_array_almost_equal(hilbert(ihilbert(f)),f)
  230. def test_random_even(self):
  231. for n in [32,64,56]:
  232. f = random((n,))
  233. af = sum(f,axis=0)/n
  234. f = f-af
  235. # zeroing Nyquist mode:
  236. f = diff(diff(f,1),-1)
  237. assert_almost_equal(sum(f,axis=0),0.0)
  238. assert_array_almost_equal(direct_hilbert(direct_ihilbert(f)),f)
  239. assert_array_almost_equal(hilbert(ihilbert(f)),f)
  240. class TestIHilbert:
  241. def test_definition(self):
  242. for n in [16,17,64,127]:
  243. x = arange(n)*2*pi/n
  244. y = ihilbert(sin(x))
  245. y1 = direct_ihilbert(sin(x))
  246. assert_array_almost_equal(y,y1)
  247. assert_array_almost_equal(ihilbert(sin(2*x)),
  248. direct_ihilbert(sin(2*x)))
  249. def test_itilbert_relation(self):
  250. for n in [16,17,64,127]:
  251. x = arange(n)*2*pi/n
  252. f = sin(x)+cos(2*x)*sin(x)
  253. y = ihilbert(f)
  254. y1 = direct_ihilbert(f)
  255. assert_array_almost_equal(y,y1)
  256. y2 = itilbert(f,h=10)
  257. assert_array_almost_equal(y,y2)
  258. class TestShift:
  259. def test_definition(self):
  260. for n in [18,17,64,127,32,2048,256]:
  261. x = arange(n)*2*pi/n
  262. for a in [0.1,3]:
  263. assert_array_almost_equal(shift(sin(x),a),direct_shift(sin(x),a))
  264. assert_array_almost_equal(shift(sin(x),a),sin(x+a))
  265. assert_array_almost_equal(shift(cos(x),a),cos(x+a))
  266. assert_array_almost_equal(shift(cos(2*x)+sin(x),a),
  267. cos(2*(x+a))+sin(x+a))
  268. assert_array_almost_equal(shift(exp(sin(x)),a),exp(sin(x+a)))
  269. assert_array_almost_equal(shift(sin(x),2*pi),sin(x))
  270. assert_array_almost_equal(shift(sin(x),pi),-sin(x))
  271. assert_array_almost_equal(shift(sin(x),pi/2),cos(x))
  272. class TestOverwrite:
  273. """Check input overwrite behavior """
  274. real_dtypes = (np.float32, np.float64)
  275. dtypes = real_dtypes + (np.complex64, np.complex128)
  276. def _check(self, x, routine, *args, **kwargs):
  277. x2 = x.copy()
  278. routine(x2, *args, **kwargs)
  279. sig = routine.__name__
  280. if args:
  281. sig += repr(args)
  282. if kwargs:
  283. sig += repr(kwargs)
  284. assert_equal(x2, x, err_msg="spurious overwrite in %s" % sig)
  285. def _check_1d(self, routine, dtype, shape, *args, **kwargs):
  286. np.random.seed(1234)
  287. if np.issubdtype(dtype, np.complexfloating):
  288. data = np.random.randn(*shape) + 1j*np.random.randn(*shape)
  289. else:
  290. data = np.random.randn(*shape)
  291. data = data.astype(dtype)
  292. self._check(data, routine, *args, **kwargs)
  293. def test_diff(self):
  294. for dtype in self.dtypes:
  295. self._check_1d(diff, dtype, (16,))
  296. def test_tilbert(self):
  297. for dtype in self.dtypes:
  298. self._check_1d(tilbert, dtype, (16,), 1.6)
  299. def test_itilbert(self):
  300. for dtype in self.dtypes:
  301. self._check_1d(itilbert, dtype, (16,), 1.6)
  302. def test_hilbert(self):
  303. for dtype in self.dtypes:
  304. self._check_1d(hilbert, dtype, (16,))
  305. def test_cs_diff(self):
  306. for dtype in self.dtypes:
  307. self._check_1d(cs_diff, dtype, (16,), 1.0, 4.0)
  308. def test_sc_diff(self):
  309. for dtype in self.dtypes:
  310. self._check_1d(sc_diff, dtype, (16,), 1.0, 4.0)
  311. def test_ss_diff(self):
  312. for dtype in self.dtypes:
  313. self._check_1d(ss_diff, dtype, (16,), 1.0, 4.0)
  314. def test_cc_diff(self):
  315. for dtype in self.dtypes:
  316. self._check_1d(cc_diff, dtype, (16,), 1.0, 4.0)
  317. def test_shift(self):
  318. for dtype in self.dtypes:
  319. self._check_1d(shift, dtype, (16,), 1.0)