test_orthogonal.py 31 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791
  1. import numpy as np
  2. from numpy import array, sqrt
  3. from numpy.testing import (assert_array_almost_equal, assert_equal,
  4. assert_almost_equal, assert_allclose)
  5. from pytest import raises as assert_raises
  6. from scipy import integrate
  7. import scipy.special as sc
  8. from scipy.special import gamma
  9. import scipy.special._orthogonal as orth
  10. class TestCheby:
  11. def test_chebyc(self):
  12. C0 = orth.chebyc(0)
  13. C1 = orth.chebyc(1)
  14. with np.errstate(all='ignore'):
  15. C2 = orth.chebyc(2)
  16. C3 = orth.chebyc(3)
  17. C4 = orth.chebyc(4)
  18. C5 = orth.chebyc(5)
  19. assert_array_almost_equal(C0.c,[2],13)
  20. assert_array_almost_equal(C1.c,[1,0],13)
  21. assert_array_almost_equal(C2.c,[1,0,-2],13)
  22. assert_array_almost_equal(C3.c,[1,0,-3,0],13)
  23. assert_array_almost_equal(C4.c,[1,0,-4,0,2],13)
  24. assert_array_almost_equal(C5.c,[1,0,-5,0,5,0],13)
  25. def test_chebys(self):
  26. S0 = orth.chebys(0)
  27. S1 = orth.chebys(1)
  28. S2 = orth.chebys(2)
  29. S3 = orth.chebys(3)
  30. S4 = orth.chebys(4)
  31. S5 = orth.chebys(5)
  32. assert_array_almost_equal(S0.c,[1],13)
  33. assert_array_almost_equal(S1.c,[1,0],13)
  34. assert_array_almost_equal(S2.c,[1,0,-1],13)
  35. assert_array_almost_equal(S3.c,[1,0,-2,0],13)
  36. assert_array_almost_equal(S4.c,[1,0,-3,0,1],13)
  37. assert_array_almost_equal(S5.c,[1,0,-4,0,3,0],13)
  38. def test_chebyt(self):
  39. T0 = orth.chebyt(0)
  40. T1 = orth.chebyt(1)
  41. T2 = orth.chebyt(2)
  42. T3 = orth.chebyt(3)
  43. T4 = orth.chebyt(4)
  44. T5 = orth.chebyt(5)
  45. assert_array_almost_equal(T0.c,[1],13)
  46. assert_array_almost_equal(T1.c,[1,0],13)
  47. assert_array_almost_equal(T2.c,[2,0,-1],13)
  48. assert_array_almost_equal(T3.c,[4,0,-3,0],13)
  49. assert_array_almost_equal(T4.c,[8,0,-8,0,1],13)
  50. assert_array_almost_equal(T5.c,[16,0,-20,0,5,0],13)
  51. def test_chebyu(self):
  52. U0 = orth.chebyu(0)
  53. U1 = orth.chebyu(1)
  54. U2 = orth.chebyu(2)
  55. U3 = orth.chebyu(3)
  56. U4 = orth.chebyu(4)
  57. U5 = orth.chebyu(5)
  58. assert_array_almost_equal(U0.c,[1],13)
  59. assert_array_almost_equal(U1.c,[2,0],13)
  60. assert_array_almost_equal(U2.c,[4,0,-1],13)
  61. assert_array_almost_equal(U3.c,[8,0,-4,0],13)
  62. assert_array_almost_equal(U4.c,[16,0,-12,0,1],13)
  63. assert_array_almost_equal(U5.c,[32,0,-32,0,6,0],13)
  64. class TestGegenbauer:
  65. def test_gegenbauer(self):
  66. a = 5*np.random.random() - 0.5
  67. if np.any(a == 0):
  68. a = -0.2
  69. Ca0 = orth.gegenbauer(0,a)
  70. Ca1 = orth.gegenbauer(1,a)
  71. Ca2 = orth.gegenbauer(2,a)
  72. Ca3 = orth.gegenbauer(3,a)
  73. Ca4 = orth.gegenbauer(4,a)
  74. Ca5 = orth.gegenbauer(5,a)
  75. assert_array_almost_equal(Ca0.c,array([1]),13)
  76. assert_array_almost_equal(Ca1.c,array([2*a,0]),13)
  77. assert_array_almost_equal(Ca2.c,array([2*a*(a+1),0,-a]),13)
  78. assert_array_almost_equal(Ca3.c,array([4*sc.poch(a,3),0,-6*a*(a+1),
  79. 0])/3.0,11)
  80. assert_array_almost_equal(Ca4.c,array([4*sc.poch(a,4),0,-12*sc.poch(a,3),
  81. 0,3*a*(a+1)])/6.0,11)
  82. assert_array_almost_equal(Ca5.c,array([4*sc.poch(a,5),0,-20*sc.poch(a,4),
  83. 0,15*sc.poch(a,3),0])/15.0,11)
  84. class TestHermite:
  85. def test_hermite(self):
  86. H0 = orth.hermite(0)
  87. H1 = orth.hermite(1)
  88. H2 = orth.hermite(2)
  89. H3 = orth.hermite(3)
  90. H4 = orth.hermite(4)
  91. H5 = orth.hermite(5)
  92. assert_array_almost_equal(H0.c,[1],13)
  93. assert_array_almost_equal(H1.c,[2,0],13)
  94. assert_array_almost_equal(H2.c,[4,0,-2],13)
  95. assert_array_almost_equal(H3.c,[8,0,-12,0],13)
  96. assert_array_almost_equal(H4.c,[16,0,-48,0,12],12)
  97. assert_array_almost_equal(H5.c,[32,0,-160,0,120,0],12)
  98. def test_hermitenorm(self):
  99. # He_n(x) = 2**(-n/2) H_n(x/sqrt(2))
  100. psub = np.poly1d([1.0/sqrt(2),0])
  101. H0 = orth.hermitenorm(0)
  102. H1 = orth.hermitenorm(1)
  103. H2 = orth.hermitenorm(2)
  104. H3 = orth.hermitenorm(3)
  105. H4 = orth.hermitenorm(4)
  106. H5 = orth.hermitenorm(5)
  107. he0 = orth.hermite(0)(psub)
  108. he1 = orth.hermite(1)(psub) / sqrt(2)
  109. he2 = orth.hermite(2)(psub) / 2.0
  110. he3 = orth.hermite(3)(psub) / (2*sqrt(2))
  111. he4 = orth.hermite(4)(psub) / 4.0
  112. he5 = orth.hermite(5)(psub) / (4.0*sqrt(2))
  113. assert_array_almost_equal(H0.c,he0.c,13)
  114. assert_array_almost_equal(H1.c,he1.c,13)
  115. assert_array_almost_equal(H2.c,he2.c,13)
  116. assert_array_almost_equal(H3.c,he3.c,13)
  117. assert_array_almost_equal(H4.c,he4.c,13)
  118. assert_array_almost_equal(H5.c,he5.c,13)
  119. class _test_sh_legendre:
  120. def test_sh_legendre(self):
  121. # P*_n(x) = P_n(2x-1)
  122. psub = np.poly1d([2,-1])
  123. Ps0 = orth.sh_legendre(0)
  124. Ps1 = orth.sh_legendre(1)
  125. Ps2 = orth.sh_legendre(2)
  126. Ps3 = orth.sh_legendre(3)
  127. Ps4 = orth.sh_legendre(4)
  128. Ps5 = orth.sh_legendre(5)
  129. pse0 = orth.legendre(0)(psub)
  130. pse1 = orth.legendre(1)(psub)
  131. pse2 = orth.legendre(2)(psub)
  132. pse3 = orth.legendre(3)(psub)
  133. pse4 = orth.legendre(4)(psub)
  134. pse5 = orth.legendre(5)(psub)
  135. assert_array_almost_equal(Ps0.c,pse0.c,13)
  136. assert_array_almost_equal(Ps1.c,pse1.c,13)
  137. assert_array_almost_equal(Ps2.c,pse2.c,13)
  138. assert_array_almost_equal(Ps3.c,pse3.c,13)
  139. assert_array_almost_equal(Ps4.c,pse4.c,12)
  140. assert_array_almost_equal(Ps5.c,pse5.c,12)
  141. class _test_sh_chebyt:
  142. def test_sh_chebyt(self):
  143. # T*_n(x) = T_n(2x-1)
  144. psub = np.poly1d([2,-1])
  145. Ts0 = orth.sh_chebyt(0)
  146. Ts1 = orth.sh_chebyt(1)
  147. Ts2 = orth.sh_chebyt(2)
  148. Ts3 = orth.sh_chebyt(3)
  149. Ts4 = orth.sh_chebyt(4)
  150. Ts5 = orth.sh_chebyt(5)
  151. tse0 = orth.chebyt(0)(psub)
  152. tse1 = orth.chebyt(1)(psub)
  153. tse2 = orth.chebyt(2)(psub)
  154. tse3 = orth.chebyt(3)(psub)
  155. tse4 = orth.chebyt(4)(psub)
  156. tse5 = orth.chebyt(5)(psub)
  157. assert_array_almost_equal(Ts0.c,tse0.c,13)
  158. assert_array_almost_equal(Ts1.c,tse1.c,13)
  159. assert_array_almost_equal(Ts2.c,tse2.c,13)
  160. assert_array_almost_equal(Ts3.c,tse3.c,13)
  161. assert_array_almost_equal(Ts4.c,tse4.c,12)
  162. assert_array_almost_equal(Ts5.c,tse5.c,12)
  163. class _test_sh_chebyu:
  164. def test_sh_chebyu(self):
  165. # U*_n(x) = U_n(2x-1)
  166. psub = np.poly1d([2,-1])
  167. Us0 = orth.sh_chebyu(0)
  168. Us1 = orth.sh_chebyu(1)
  169. Us2 = orth.sh_chebyu(2)
  170. Us3 = orth.sh_chebyu(3)
  171. Us4 = orth.sh_chebyu(4)
  172. Us5 = orth.sh_chebyu(5)
  173. use0 = orth.chebyu(0)(psub)
  174. use1 = orth.chebyu(1)(psub)
  175. use2 = orth.chebyu(2)(psub)
  176. use3 = orth.chebyu(3)(psub)
  177. use4 = orth.chebyu(4)(psub)
  178. use5 = orth.chebyu(5)(psub)
  179. assert_array_almost_equal(Us0.c,use0.c,13)
  180. assert_array_almost_equal(Us1.c,use1.c,13)
  181. assert_array_almost_equal(Us2.c,use2.c,13)
  182. assert_array_almost_equal(Us3.c,use3.c,13)
  183. assert_array_almost_equal(Us4.c,use4.c,12)
  184. assert_array_almost_equal(Us5.c,use5.c,11)
  185. class _test_sh_jacobi:
  186. def test_sh_jacobi(self):
  187. # G^(p,q)_n(x) = n! gamma(n+p)/gamma(2*n+p) * P^(p-q,q-1)_n(2*x-1)
  188. conv = lambda n,p: gamma(n+1)*gamma(n+p)/gamma(2*n+p)
  189. psub = np.poly1d([2,-1])
  190. q = 4 * np.random.random()
  191. p = q-1 + 2*np.random.random()
  192. # print("shifted jacobi p,q = ", p, q)
  193. G0 = orth.sh_jacobi(0,p,q)
  194. G1 = orth.sh_jacobi(1,p,q)
  195. G2 = orth.sh_jacobi(2,p,q)
  196. G3 = orth.sh_jacobi(3,p,q)
  197. G4 = orth.sh_jacobi(4,p,q)
  198. G5 = orth.sh_jacobi(5,p,q)
  199. ge0 = orth.jacobi(0,p-q,q-1)(psub) * conv(0,p)
  200. ge1 = orth.jacobi(1,p-q,q-1)(psub) * conv(1,p)
  201. ge2 = orth.jacobi(2,p-q,q-1)(psub) * conv(2,p)
  202. ge3 = orth.jacobi(3,p-q,q-1)(psub) * conv(3,p)
  203. ge4 = orth.jacobi(4,p-q,q-1)(psub) * conv(4,p)
  204. ge5 = orth.jacobi(5,p-q,q-1)(psub) * conv(5,p)
  205. assert_array_almost_equal(G0.c,ge0.c,13)
  206. assert_array_almost_equal(G1.c,ge1.c,13)
  207. assert_array_almost_equal(G2.c,ge2.c,13)
  208. assert_array_almost_equal(G3.c,ge3.c,13)
  209. assert_array_almost_equal(G4.c,ge4.c,13)
  210. assert_array_almost_equal(G5.c,ge5.c,13)
  211. class TestCall:
  212. def test_call(self):
  213. poly = []
  214. for n in range(5):
  215. poly.extend([x.strip() for x in
  216. ("""
  217. orth.jacobi(%(n)d,0.3,0.9)
  218. orth.sh_jacobi(%(n)d,0.3,0.9)
  219. orth.genlaguerre(%(n)d,0.3)
  220. orth.laguerre(%(n)d)
  221. orth.hermite(%(n)d)
  222. orth.hermitenorm(%(n)d)
  223. orth.gegenbauer(%(n)d,0.3)
  224. orth.chebyt(%(n)d)
  225. orth.chebyu(%(n)d)
  226. orth.chebyc(%(n)d)
  227. orth.chebys(%(n)d)
  228. orth.sh_chebyt(%(n)d)
  229. orth.sh_chebyu(%(n)d)
  230. orth.legendre(%(n)d)
  231. orth.sh_legendre(%(n)d)
  232. """ % dict(n=n)).split()
  233. ])
  234. with np.errstate(all='ignore'):
  235. for pstr in poly:
  236. p = eval(pstr)
  237. assert_almost_equal(p(0.315), np.poly1d(p.coef)(0.315),
  238. err_msg=pstr)
  239. class TestGenlaguerre:
  240. def test_regression(self):
  241. assert_equal(orth.genlaguerre(1, 1, monic=False)(0), 2.)
  242. assert_equal(orth.genlaguerre(1, 1, monic=True)(0), -2.)
  243. assert_equal(orth.genlaguerre(1, 1, monic=False), np.poly1d([-1, 2]))
  244. assert_equal(orth.genlaguerre(1, 1, monic=True), np.poly1d([1, -2]))
  245. def verify_gauss_quad(root_func, eval_func, weight_func, a, b, N,
  246. rtol=1e-15, atol=5e-14):
  247. # this test is copied from numpy's TestGauss in test_hermite.py
  248. x, w, mu = root_func(N, True)
  249. n = np.arange(N)
  250. v = eval_func(n[:,np.newaxis], x)
  251. vv = np.dot(v*w, v.T)
  252. vd = 1 / np.sqrt(vv.diagonal())
  253. vv = vd[:, np.newaxis] * vv * vd
  254. assert_allclose(vv, np.eye(N), rtol, atol)
  255. # check that the integral of 1 is correct
  256. assert_allclose(w.sum(), mu, rtol, atol)
  257. # compare the results of integrating a function with quad.
  258. f = lambda x: x**3 - 3*x**2 + x - 2
  259. resI = integrate.quad(lambda x: f(x)*weight_func(x), a, b)
  260. resG = np.vdot(f(x), w)
  261. rtol = 1e-6 if 1e-6 < resI[1] else resI[1] * 10
  262. assert_allclose(resI[0], resG, rtol=rtol)
  263. def test_roots_jacobi():
  264. rf = lambda a, b: lambda n, mu: sc.roots_jacobi(n, a, b, mu)
  265. ef = lambda a, b: lambda n, x: sc.eval_jacobi(n, a, b, x)
  266. wf = lambda a, b: lambda x: (1 - x)**a * (1 + x)**b
  267. vgq = verify_gauss_quad
  268. vgq(rf(-0.5, -0.75), ef(-0.5, -0.75), wf(-0.5, -0.75), -1., 1., 5)
  269. vgq(rf(-0.5, -0.75), ef(-0.5, -0.75), wf(-0.5, -0.75), -1., 1.,
  270. 25, atol=1e-12)
  271. vgq(rf(-0.5, -0.75), ef(-0.5, -0.75), wf(-0.5, -0.75), -1., 1.,
  272. 100, atol=1e-11)
  273. vgq(rf(0.5, -0.5), ef(0.5, -0.5), wf(0.5, -0.5), -1., 1., 5)
  274. vgq(rf(0.5, -0.5), ef(0.5, -0.5), wf(0.5, -0.5), -1., 1., 25, atol=1.5e-13)
  275. vgq(rf(0.5, -0.5), ef(0.5, -0.5), wf(0.5, -0.5), -1., 1., 100, atol=2e-12)
  276. vgq(rf(1, 0.5), ef(1, 0.5), wf(1, 0.5), -1., 1., 5, atol=2e-13)
  277. vgq(rf(1, 0.5), ef(1, 0.5), wf(1, 0.5), -1., 1., 25, atol=2e-13)
  278. vgq(rf(1, 0.5), ef(1, 0.5), wf(1, 0.5), -1., 1., 100, atol=1e-12)
  279. vgq(rf(0.9, 2), ef(0.9, 2), wf(0.9, 2), -1., 1., 5)
  280. vgq(rf(0.9, 2), ef(0.9, 2), wf(0.9, 2), -1., 1., 25, atol=1e-13)
  281. vgq(rf(0.9, 2), ef(0.9, 2), wf(0.9, 2), -1., 1., 100, atol=3e-13)
  282. vgq(rf(18.24, 27.3), ef(18.24, 27.3), wf(18.24, 27.3), -1., 1., 5)
  283. vgq(rf(18.24, 27.3), ef(18.24, 27.3), wf(18.24, 27.3), -1., 1., 25,
  284. atol=1.1e-14)
  285. vgq(rf(18.24, 27.3), ef(18.24, 27.3), wf(18.24, 27.3), -1., 1.,
  286. 100, atol=1e-13)
  287. vgq(rf(47.1, -0.2), ef(47.1, -0.2), wf(47.1, -0.2), -1., 1., 5, atol=1e-13)
  288. vgq(rf(47.1, -0.2), ef(47.1, -0.2), wf(47.1, -0.2), -1., 1., 25, atol=2e-13)
  289. vgq(rf(47.1, -0.2), ef(47.1, -0.2), wf(47.1, -0.2), -1., 1.,
  290. 100, atol=1e-11)
  291. vgq(rf(1., 658.), ef(1., 658.), wf(1., 658.), -1., 1., 5, atol=2e-13)
  292. vgq(rf(1., 658.), ef(1., 658.), wf(1., 658.), -1., 1., 25, atol=1e-12)
  293. vgq(rf(1., 658.), ef(1., 658.), wf(1., 658.), -1., 1., 100, atol=1e-11)
  294. vgq(rf(1., 658.), ef(1., 658.), wf(1., 658.), -1., 1., 250, atol=1e-11)
  295. vgq(rf(511., 511.), ef(511., 511.), wf(511., 511.), -1., 1., 5,
  296. atol=1e-12)
  297. vgq(rf(511., 511.), ef(511., 511.), wf(511., 511.), -1., 1., 25,
  298. atol=1e-11)
  299. vgq(rf(511., 511.), ef(511., 511.), wf(511., 511.), -1., 1., 100,
  300. atol=1e-10)
  301. vgq(rf(511., 512.), ef(511., 512.), wf(511., 512.), -1., 1., 5,
  302. atol=1e-12)
  303. vgq(rf(511., 512.), ef(511., 512.), wf(511., 512.), -1., 1., 25,
  304. atol=1e-11)
  305. vgq(rf(511., 512.), ef(511., 512.), wf(511., 512.), -1., 1., 100,
  306. atol=1e-10)
  307. vgq(rf(1000., 500.), ef(1000., 500.), wf(1000., 500.), -1., 1., 5,
  308. atol=1e-12)
  309. vgq(rf(1000., 500.), ef(1000., 500.), wf(1000., 500.), -1., 1., 25,
  310. atol=1e-11)
  311. vgq(rf(1000., 500.), ef(1000., 500.), wf(1000., 500.), -1., 1., 100,
  312. atol=1e-10)
  313. vgq(rf(2.25, 68.9), ef(2.25, 68.9), wf(2.25, 68.9), -1., 1., 5)
  314. vgq(rf(2.25, 68.9), ef(2.25, 68.9), wf(2.25, 68.9), -1., 1., 25,
  315. atol=1e-13)
  316. vgq(rf(2.25, 68.9), ef(2.25, 68.9), wf(2.25, 68.9), -1., 1., 100,
  317. atol=1e-13)
  318. # when alpha == beta == 0, P_n^{a,b}(x) == P_n(x)
  319. xj, wj = sc.roots_jacobi(6, 0.0, 0.0)
  320. xl, wl = sc.roots_legendre(6)
  321. assert_allclose(xj, xl, 1e-14, 1e-14)
  322. assert_allclose(wj, wl, 1e-14, 1e-14)
  323. # when alpha == beta != 0, P_n^{a,b}(x) == C_n^{alpha+0.5}(x)
  324. xj, wj = sc.roots_jacobi(6, 4.0, 4.0)
  325. xc, wc = sc.roots_gegenbauer(6, 4.5)
  326. assert_allclose(xj, xc, 1e-14, 1e-14)
  327. assert_allclose(wj, wc, 1e-14, 1e-14)
  328. x, w = sc.roots_jacobi(5, 2, 3, False)
  329. y, v, m = sc.roots_jacobi(5, 2, 3, True)
  330. assert_allclose(x, y, 1e-14, 1e-14)
  331. assert_allclose(w, v, 1e-14, 1e-14)
  332. muI, muI_err = integrate.quad(wf(2,3), -1, 1)
  333. assert_allclose(m, muI, rtol=muI_err)
  334. assert_raises(ValueError, sc.roots_jacobi, 0, 1, 1)
  335. assert_raises(ValueError, sc.roots_jacobi, 3.3, 1, 1)
  336. assert_raises(ValueError, sc.roots_jacobi, 3, -2, 1)
  337. assert_raises(ValueError, sc.roots_jacobi, 3, 1, -2)
  338. assert_raises(ValueError, sc.roots_jacobi, 3, -2, -2)
  339. def test_roots_sh_jacobi():
  340. rf = lambda a, b: lambda n, mu: sc.roots_sh_jacobi(n, a, b, mu)
  341. ef = lambda a, b: lambda n, x: sc.eval_sh_jacobi(n, a, b, x)
  342. wf = lambda a, b: lambda x: (1. - x)**(a - b) * (x)**(b - 1.)
  343. vgq = verify_gauss_quad
  344. vgq(rf(-0.5, 0.25), ef(-0.5, 0.25), wf(-0.5, 0.25), 0., 1., 5)
  345. vgq(rf(-0.5, 0.25), ef(-0.5, 0.25), wf(-0.5, 0.25), 0., 1.,
  346. 25, atol=1e-12)
  347. vgq(rf(-0.5, 0.25), ef(-0.5, 0.25), wf(-0.5, 0.25), 0., 1.,
  348. 100, atol=1e-11)
  349. vgq(rf(0.5, 0.5), ef(0.5, 0.5), wf(0.5, 0.5), 0., 1., 5)
  350. vgq(rf(0.5, 0.5), ef(0.5, 0.5), wf(0.5, 0.5), 0., 1., 25, atol=1e-13)
  351. vgq(rf(0.5, 0.5), ef(0.5, 0.5), wf(0.5, 0.5), 0., 1., 100, atol=1e-12)
  352. vgq(rf(1, 0.5), ef(1, 0.5), wf(1, 0.5), 0., 1., 5)
  353. vgq(rf(1, 0.5), ef(1, 0.5), wf(1, 0.5), 0., 1., 25, atol=1.5e-13)
  354. vgq(rf(1, 0.5), ef(1, 0.5), wf(1, 0.5), 0., 1., 100, atol=2e-12)
  355. vgq(rf(2, 0.9), ef(2, 0.9), wf(2, 0.9), 0., 1., 5)
  356. vgq(rf(2, 0.9), ef(2, 0.9), wf(2, 0.9), 0., 1., 25, atol=1e-13)
  357. vgq(rf(2, 0.9), ef(2, 0.9), wf(2, 0.9), 0., 1., 100, atol=1e-12)
  358. vgq(rf(27.3, 18.24), ef(27.3, 18.24), wf(27.3, 18.24), 0., 1., 5)
  359. vgq(rf(27.3, 18.24), ef(27.3, 18.24), wf(27.3, 18.24), 0., 1., 25)
  360. vgq(rf(27.3, 18.24), ef(27.3, 18.24), wf(27.3, 18.24), 0., 1.,
  361. 100, atol=1e-13)
  362. vgq(rf(47.1, 0.2), ef(47.1, 0.2), wf(47.1, 0.2), 0., 1., 5, atol=1e-12)
  363. vgq(rf(47.1, 0.2), ef(47.1, 0.2), wf(47.1, 0.2), 0., 1., 25, atol=1e-11)
  364. vgq(rf(47.1, 0.2), ef(47.1, 0.2), wf(47.1, 0.2), 0., 1., 100, atol=1e-10)
  365. vgq(rf(68.9, 2.25), ef(68.9, 2.25), wf(68.9, 2.25), 0., 1., 5, atol=3.5e-14)
  366. vgq(rf(68.9, 2.25), ef(68.9, 2.25), wf(68.9, 2.25), 0., 1., 25, atol=2e-13)
  367. vgq(rf(68.9, 2.25), ef(68.9, 2.25), wf(68.9, 2.25), 0., 1.,
  368. 100, atol=1e-12)
  369. x, w = sc.roots_sh_jacobi(5, 3, 2, False)
  370. y, v, m = sc.roots_sh_jacobi(5, 3, 2, True)
  371. assert_allclose(x, y, 1e-14, 1e-14)
  372. assert_allclose(w, v, 1e-14, 1e-14)
  373. muI, muI_err = integrate.quad(wf(3,2), 0, 1)
  374. assert_allclose(m, muI, rtol=muI_err)
  375. assert_raises(ValueError, sc.roots_sh_jacobi, 0, 1, 1)
  376. assert_raises(ValueError, sc.roots_sh_jacobi, 3.3, 1, 1)
  377. assert_raises(ValueError, sc.roots_sh_jacobi, 3, 1, 2) # p - q <= -1
  378. assert_raises(ValueError, sc.roots_sh_jacobi, 3, 2, -1) # q <= 0
  379. assert_raises(ValueError, sc.roots_sh_jacobi, 3, -2, -1) # both
  380. def test_roots_hermite():
  381. rootf = sc.roots_hermite
  382. evalf = sc.eval_hermite
  383. weightf = orth.hermite(5).weight_func
  384. verify_gauss_quad(rootf, evalf, weightf, -np.inf, np.inf, 5)
  385. verify_gauss_quad(rootf, evalf, weightf, -np.inf, np.inf, 25, atol=1e-13)
  386. verify_gauss_quad(rootf, evalf, weightf, -np.inf, np.inf, 100, atol=1e-12)
  387. # Golub-Welsch branch
  388. x, w = sc.roots_hermite(5, False)
  389. y, v, m = sc.roots_hermite(5, True)
  390. assert_allclose(x, y, 1e-14, 1e-14)
  391. assert_allclose(w, v, 1e-14, 1e-14)
  392. muI, muI_err = integrate.quad(weightf, -np.inf, np.inf)
  393. assert_allclose(m, muI, rtol=muI_err)
  394. # Asymptotic branch (switch over at n >= 150)
  395. x, w = sc.roots_hermite(200, False)
  396. y, v, m = sc.roots_hermite(200, True)
  397. assert_allclose(x, y, 1e-14, 1e-14)
  398. assert_allclose(w, v, 1e-14, 1e-14)
  399. assert_allclose(sum(v), m, 1e-14, 1e-14)
  400. assert_raises(ValueError, sc.roots_hermite, 0)
  401. assert_raises(ValueError, sc.roots_hermite, 3.3)
  402. def test_roots_hermite_asy():
  403. # Recursion for Hermite functions
  404. def hermite_recursion(n, nodes):
  405. H = np.zeros((n, nodes.size))
  406. H[0,:] = np.pi**(-0.25) * np.exp(-0.5*nodes**2)
  407. if n > 1:
  408. H[1,:] = sqrt(2.0) * nodes * H[0,:]
  409. for k in range(2, n):
  410. H[k,:] = sqrt(2.0/k) * nodes * H[k-1,:] - sqrt((k-1.0)/k) * H[k-2,:]
  411. return H
  412. # This tests only the nodes
  413. def test(N, rtol=1e-15, atol=1e-14):
  414. x, w = orth._roots_hermite_asy(N)
  415. H = hermite_recursion(N+1, x)
  416. assert_allclose(H[-1,:], np.zeros(N), rtol, atol)
  417. assert_allclose(sum(w), sqrt(np.pi), rtol, atol)
  418. test(150, atol=1e-12)
  419. test(151, atol=1e-12)
  420. test(300, atol=1e-12)
  421. test(301, atol=1e-12)
  422. test(500, atol=1e-12)
  423. test(501, atol=1e-12)
  424. test(999, atol=1e-12)
  425. test(1000, atol=1e-12)
  426. test(2000, atol=1e-12)
  427. test(5000, atol=1e-12)
  428. def test_roots_hermitenorm():
  429. rootf = sc.roots_hermitenorm
  430. evalf = sc.eval_hermitenorm
  431. weightf = orth.hermitenorm(5).weight_func
  432. verify_gauss_quad(rootf, evalf, weightf, -np.inf, np.inf, 5)
  433. verify_gauss_quad(rootf, evalf, weightf, -np.inf, np.inf, 25, atol=1e-13)
  434. verify_gauss_quad(rootf, evalf, weightf, -np.inf, np.inf, 100, atol=1e-12)
  435. x, w = sc.roots_hermitenorm(5, False)
  436. y, v, m = sc.roots_hermitenorm(5, True)
  437. assert_allclose(x, y, 1e-14, 1e-14)
  438. assert_allclose(w, v, 1e-14, 1e-14)
  439. muI, muI_err = integrate.quad(weightf, -np.inf, np.inf)
  440. assert_allclose(m, muI, rtol=muI_err)
  441. assert_raises(ValueError, sc.roots_hermitenorm, 0)
  442. assert_raises(ValueError, sc.roots_hermitenorm, 3.3)
  443. def test_roots_gegenbauer():
  444. rootf = lambda a: lambda n, mu: sc.roots_gegenbauer(n, a, mu)
  445. evalf = lambda a: lambda n, x: sc.eval_gegenbauer(n, a, x)
  446. weightf = lambda a: lambda x: (1 - x**2)**(a - 0.5)
  447. vgq = verify_gauss_quad
  448. vgq(rootf(-0.25), evalf(-0.25), weightf(-0.25), -1., 1., 5)
  449. vgq(rootf(-0.25), evalf(-0.25), weightf(-0.25), -1., 1., 25, atol=1e-12)
  450. vgq(rootf(-0.25), evalf(-0.25), weightf(-0.25), -1., 1., 100, atol=1e-11)
  451. vgq(rootf(0.1), evalf(0.1), weightf(0.1), -1., 1., 5)
  452. vgq(rootf(0.1), evalf(0.1), weightf(0.1), -1., 1., 25, atol=1e-13)
  453. vgq(rootf(0.1), evalf(0.1), weightf(0.1), -1., 1., 100, atol=1e-12)
  454. vgq(rootf(1), evalf(1), weightf(1), -1., 1., 5)
  455. vgq(rootf(1), evalf(1), weightf(1), -1., 1., 25, atol=1e-13)
  456. vgq(rootf(1), evalf(1), weightf(1), -1., 1., 100, atol=1e-12)
  457. vgq(rootf(10), evalf(10), weightf(10), -1., 1., 5)
  458. vgq(rootf(10), evalf(10), weightf(10), -1., 1., 25, atol=1e-13)
  459. vgq(rootf(10), evalf(10), weightf(10), -1., 1., 100, atol=1e-12)
  460. vgq(rootf(50), evalf(50), weightf(50), -1., 1., 5, atol=1e-13)
  461. vgq(rootf(50), evalf(50), weightf(50), -1., 1., 25, atol=1e-12)
  462. vgq(rootf(50), evalf(50), weightf(50), -1., 1., 100, atol=1e-11)
  463. # Alpha=170 is where the approximation used in roots_gegenbauer changes
  464. vgq(rootf(170), evalf(170), weightf(170), -1., 1., 5, atol=1e-13)
  465. vgq(rootf(170), evalf(170), weightf(170), -1., 1., 25, atol=1e-12)
  466. vgq(rootf(170), evalf(170), weightf(170), -1., 1., 100, atol=1e-11)
  467. vgq(rootf(170.5), evalf(170.5), weightf(170.5), -1., 1., 5, atol=1.25e-13)
  468. vgq(rootf(170.5), evalf(170.5), weightf(170.5), -1., 1., 25, atol=1e-12)
  469. vgq(rootf(170.5), evalf(170.5), weightf(170.5), -1., 1., 100, atol=1e-11)
  470. # Test for failures, e.g. overflows, resulting from large alphas
  471. vgq(rootf(238), evalf(238), weightf(238), -1., 1., 5, atol=1e-13)
  472. vgq(rootf(238), evalf(238), weightf(238), -1., 1., 25, atol=1e-12)
  473. vgq(rootf(238), evalf(238), weightf(238), -1., 1., 100, atol=1e-11)
  474. vgq(rootf(512.5), evalf(512.5), weightf(512.5), -1., 1., 5, atol=1e-12)
  475. vgq(rootf(512.5), evalf(512.5), weightf(512.5), -1., 1., 25, atol=1e-11)
  476. vgq(rootf(512.5), evalf(512.5), weightf(512.5), -1., 1., 100, atol=1e-10)
  477. # this is a special case that the old code supported.
  478. # when alpha = 0, the gegenbauer polynomial is uniformly 0. but it goes
  479. # to a scaled down copy of T_n(x) there.
  480. vgq(rootf(0), sc.eval_chebyt, weightf(0), -1., 1., 5)
  481. vgq(rootf(0), sc.eval_chebyt, weightf(0), -1., 1., 25)
  482. vgq(rootf(0), sc.eval_chebyt, weightf(0), -1., 1., 100, atol=1e-12)
  483. x, w = sc.roots_gegenbauer(5, 2, False)
  484. y, v, m = sc.roots_gegenbauer(5, 2, True)
  485. assert_allclose(x, y, 1e-14, 1e-14)
  486. assert_allclose(w, v, 1e-14, 1e-14)
  487. muI, muI_err = integrate.quad(weightf(2), -1, 1)
  488. assert_allclose(m, muI, rtol=muI_err)
  489. assert_raises(ValueError, sc.roots_gegenbauer, 0, 2)
  490. assert_raises(ValueError, sc.roots_gegenbauer, 3.3, 2)
  491. assert_raises(ValueError, sc.roots_gegenbauer, 3, -.75)
  492. def test_roots_chebyt():
  493. weightf = orth.chebyt(5).weight_func
  494. verify_gauss_quad(sc.roots_chebyt, sc.eval_chebyt, weightf, -1., 1., 5)
  495. verify_gauss_quad(sc.roots_chebyt, sc.eval_chebyt, weightf, -1., 1., 25)
  496. verify_gauss_quad(sc.roots_chebyt, sc.eval_chebyt, weightf, -1., 1., 100, atol=1e-12)
  497. x, w = sc.roots_chebyt(5, False)
  498. y, v, m = sc.roots_chebyt(5, True)
  499. assert_allclose(x, y, 1e-14, 1e-14)
  500. assert_allclose(w, v, 1e-14, 1e-14)
  501. muI, muI_err = integrate.quad(weightf, -1, 1)
  502. assert_allclose(m, muI, rtol=muI_err)
  503. assert_raises(ValueError, sc.roots_chebyt, 0)
  504. assert_raises(ValueError, sc.roots_chebyt, 3.3)
  505. def test_chebyt_symmetry():
  506. x, w = sc.roots_chebyt(21)
  507. pos, neg = x[:10], x[11:]
  508. assert_equal(neg, -pos[::-1])
  509. assert_equal(x[10], 0)
  510. def test_roots_chebyu():
  511. weightf = orth.chebyu(5).weight_func
  512. verify_gauss_quad(sc.roots_chebyu, sc.eval_chebyu, weightf, -1., 1., 5)
  513. verify_gauss_quad(sc.roots_chebyu, sc.eval_chebyu, weightf, -1., 1., 25)
  514. verify_gauss_quad(sc.roots_chebyu, sc.eval_chebyu, weightf, -1., 1., 100)
  515. x, w = sc.roots_chebyu(5, False)
  516. y, v, m = sc.roots_chebyu(5, True)
  517. assert_allclose(x, y, 1e-14, 1e-14)
  518. assert_allclose(w, v, 1e-14, 1e-14)
  519. muI, muI_err = integrate.quad(weightf, -1, 1)
  520. assert_allclose(m, muI, rtol=muI_err)
  521. assert_raises(ValueError, sc.roots_chebyu, 0)
  522. assert_raises(ValueError, sc.roots_chebyu, 3.3)
  523. def test_roots_chebyc():
  524. weightf = orth.chebyc(5).weight_func
  525. verify_gauss_quad(sc.roots_chebyc, sc.eval_chebyc, weightf, -2., 2., 5)
  526. verify_gauss_quad(sc.roots_chebyc, sc.eval_chebyc, weightf, -2., 2., 25)
  527. verify_gauss_quad(sc.roots_chebyc, sc.eval_chebyc, weightf, -2., 2., 100, atol=1e-12)
  528. x, w = sc.roots_chebyc(5, False)
  529. y, v, m = sc.roots_chebyc(5, True)
  530. assert_allclose(x, y, 1e-14, 1e-14)
  531. assert_allclose(w, v, 1e-14, 1e-14)
  532. muI, muI_err = integrate.quad(weightf, -2, 2)
  533. assert_allclose(m, muI, rtol=muI_err)
  534. assert_raises(ValueError, sc.roots_chebyc, 0)
  535. assert_raises(ValueError, sc.roots_chebyc, 3.3)
  536. def test_roots_chebys():
  537. weightf = orth.chebys(5).weight_func
  538. verify_gauss_quad(sc.roots_chebys, sc.eval_chebys, weightf, -2., 2., 5)
  539. verify_gauss_quad(sc.roots_chebys, sc.eval_chebys, weightf, -2., 2., 25)
  540. verify_gauss_quad(sc.roots_chebys, sc.eval_chebys, weightf, -2., 2., 100)
  541. x, w = sc.roots_chebys(5, False)
  542. y, v, m = sc.roots_chebys(5, True)
  543. assert_allclose(x, y, 1e-14, 1e-14)
  544. assert_allclose(w, v, 1e-14, 1e-14)
  545. muI, muI_err = integrate.quad(weightf, -2, 2)
  546. assert_allclose(m, muI, rtol=muI_err)
  547. assert_raises(ValueError, sc.roots_chebys, 0)
  548. assert_raises(ValueError, sc.roots_chebys, 3.3)
  549. def test_roots_sh_chebyt():
  550. weightf = orth.sh_chebyt(5).weight_func
  551. verify_gauss_quad(sc.roots_sh_chebyt, sc.eval_sh_chebyt, weightf, 0., 1., 5)
  552. verify_gauss_quad(sc.roots_sh_chebyt, sc.eval_sh_chebyt, weightf, 0., 1., 25)
  553. verify_gauss_quad(sc.roots_sh_chebyt, sc.eval_sh_chebyt, weightf, 0., 1.,
  554. 100, atol=1e-13)
  555. x, w = sc.roots_sh_chebyt(5, False)
  556. y, v, m = sc.roots_sh_chebyt(5, True)
  557. assert_allclose(x, y, 1e-14, 1e-14)
  558. assert_allclose(w, v, 1e-14, 1e-14)
  559. muI, muI_err = integrate.quad(weightf, 0, 1)
  560. assert_allclose(m, muI, rtol=muI_err)
  561. assert_raises(ValueError, sc.roots_sh_chebyt, 0)
  562. assert_raises(ValueError, sc.roots_sh_chebyt, 3.3)
  563. def test_roots_sh_chebyu():
  564. weightf = orth.sh_chebyu(5).weight_func
  565. verify_gauss_quad(sc.roots_sh_chebyu, sc.eval_sh_chebyu, weightf, 0., 1., 5)
  566. verify_gauss_quad(sc.roots_sh_chebyu, sc.eval_sh_chebyu, weightf, 0., 1., 25)
  567. verify_gauss_quad(sc.roots_sh_chebyu, sc.eval_sh_chebyu, weightf, 0., 1.,
  568. 100, atol=1e-13)
  569. x, w = sc.roots_sh_chebyu(5, False)
  570. y, v, m = sc.roots_sh_chebyu(5, True)
  571. assert_allclose(x, y, 1e-14, 1e-14)
  572. assert_allclose(w, v, 1e-14, 1e-14)
  573. muI, muI_err = integrate.quad(weightf, 0, 1)
  574. assert_allclose(m, muI, rtol=muI_err)
  575. assert_raises(ValueError, sc.roots_sh_chebyu, 0)
  576. assert_raises(ValueError, sc.roots_sh_chebyu, 3.3)
  577. def test_roots_legendre():
  578. weightf = orth.legendre(5).weight_func
  579. verify_gauss_quad(sc.roots_legendre, sc.eval_legendre, weightf, -1., 1., 5)
  580. verify_gauss_quad(sc.roots_legendre, sc.eval_legendre, weightf, -1., 1.,
  581. 25, atol=1e-13)
  582. verify_gauss_quad(sc.roots_legendre, sc.eval_legendre, weightf, -1., 1.,
  583. 100, atol=1e-12)
  584. x, w = sc.roots_legendre(5, False)
  585. y, v, m = sc.roots_legendre(5, True)
  586. assert_allclose(x, y, 1e-14, 1e-14)
  587. assert_allclose(w, v, 1e-14, 1e-14)
  588. muI, muI_err = integrate.quad(weightf, -1, 1)
  589. assert_allclose(m, muI, rtol=muI_err)
  590. assert_raises(ValueError, sc.roots_legendre, 0)
  591. assert_raises(ValueError, sc.roots_legendre, 3.3)
  592. def test_roots_sh_legendre():
  593. weightf = orth.sh_legendre(5).weight_func
  594. verify_gauss_quad(sc.roots_sh_legendre, sc.eval_sh_legendre, weightf, 0., 1., 5)
  595. verify_gauss_quad(sc.roots_sh_legendre, sc.eval_sh_legendre, weightf, 0., 1.,
  596. 25, atol=1e-13)
  597. verify_gauss_quad(sc.roots_sh_legendre, sc.eval_sh_legendre, weightf, 0., 1.,
  598. 100, atol=1e-12)
  599. x, w = sc.roots_sh_legendre(5, False)
  600. y, v, m = sc.roots_sh_legendre(5, True)
  601. assert_allclose(x, y, 1e-14, 1e-14)
  602. assert_allclose(w, v, 1e-14, 1e-14)
  603. muI, muI_err = integrate.quad(weightf, 0, 1)
  604. assert_allclose(m, muI, rtol=muI_err)
  605. assert_raises(ValueError, sc.roots_sh_legendre, 0)
  606. assert_raises(ValueError, sc.roots_sh_legendre, 3.3)
  607. def test_roots_laguerre():
  608. weightf = orth.laguerre(5).weight_func
  609. verify_gauss_quad(sc.roots_laguerre, sc.eval_laguerre, weightf, 0., np.inf, 5)
  610. verify_gauss_quad(sc.roots_laguerre, sc.eval_laguerre, weightf, 0., np.inf,
  611. 25, atol=1e-13)
  612. verify_gauss_quad(sc.roots_laguerre, sc.eval_laguerre, weightf, 0., np.inf,
  613. 100, atol=1e-12)
  614. x, w = sc.roots_laguerre(5, False)
  615. y, v, m = sc.roots_laguerre(5, True)
  616. assert_allclose(x, y, 1e-14, 1e-14)
  617. assert_allclose(w, v, 1e-14, 1e-14)
  618. muI, muI_err = integrate.quad(weightf, 0, np.inf)
  619. assert_allclose(m, muI, rtol=muI_err)
  620. assert_raises(ValueError, sc.roots_laguerre, 0)
  621. assert_raises(ValueError, sc.roots_laguerre, 3.3)
  622. def test_roots_genlaguerre():
  623. rootf = lambda a: lambda n, mu: sc.roots_genlaguerre(n, a, mu)
  624. evalf = lambda a: lambda n, x: sc.eval_genlaguerre(n, a, x)
  625. weightf = lambda a: lambda x: x**a * np.exp(-x)
  626. vgq = verify_gauss_quad
  627. vgq(rootf(-0.5), evalf(-0.5), weightf(-0.5), 0., np.inf, 5)
  628. vgq(rootf(-0.5), evalf(-0.5), weightf(-0.5), 0., np.inf, 25, atol=1e-13)
  629. vgq(rootf(-0.5), evalf(-0.5), weightf(-0.5), 0., np.inf, 100, atol=1e-12)
  630. vgq(rootf(0.1), evalf(0.1), weightf(0.1), 0., np.inf, 5)
  631. vgq(rootf(0.1), evalf(0.1), weightf(0.1), 0., np.inf, 25, atol=1e-13)
  632. vgq(rootf(0.1), evalf(0.1), weightf(0.1), 0., np.inf, 100, atol=1.6e-13)
  633. vgq(rootf(1), evalf(1), weightf(1), 0., np.inf, 5)
  634. vgq(rootf(1), evalf(1), weightf(1), 0., np.inf, 25, atol=1e-13)
  635. vgq(rootf(1), evalf(1), weightf(1), 0., np.inf, 100, atol=1.03e-13)
  636. vgq(rootf(10), evalf(10), weightf(10), 0., np.inf, 5)
  637. vgq(rootf(10), evalf(10), weightf(10), 0., np.inf, 25, atol=1e-13)
  638. vgq(rootf(10), evalf(10), weightf(10), 0., np.inf, 100, atol=1e-12)
  639. vgq(rootf(50), evalf(50), weightf(50), 0., np.inf, 5)
  640. vgq(rootf(50), evalf(50), weightf(50), 0., np.inf, 25, atol=1e-13)
  641. vgq(rootf(50), evalf(50), weightf(50), 0., np.inf, 100, rtol=1e-14, atol=2e-13)
  642. x, w = sc.roots_genlaguerre(5, 2, False)
  643. y, v, m = sc.roots_genlaguerre(5, 2, True)
  644. assert_allclose(x, y, 1e-14, 1e-14)
  645. assert_allclose(w, v, 1e-14, 1e-14)
  646. muI, muI_err = integrate.quad(weightf(2.), 0., np.inf)
  647. assert_allclose(m, muI, rtol=muI_err)
  648. assert_raises(ValueError, sc.roots_genlaguerre, 0, 2)
  649. assert_raises(ValueError, sc.roots_genlaguerre, 3.3, 2)
  650. assert_raises(ValueError, sc.roots_genlaguerre, 3, -1.1)
  651. def test_gh_6721():
  652. # Regresssion test for gh_6721. This should not raise.
  653. sc.chebyt(65)(0.2)