test_spherical_bessel.py 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407
  1. #
  2. # Tests of spherical Bessel functions.
  3. #
  4. import warnings
  5. import numpy as np
  6. from numpy.testing import assert_allclose
  7. import pytest
  8. from numpy import sin, cos, sinh, cosh, exp, inf, nan, r_, pi
  9. from scipy.special import spherical_jn, spherical_yn, spherical_in, spherical_kn
  10. from scipy.integrate import quad
  11. class TestSphericalJn:
  12. def test_spherical_jn_exact(self):
  13. # https://dlmf.nist.gov/10.49.E3
  14. # Note: exact expression is numerically stable only for small
  15. # n or z >> n.
  16. x = np.array([0.12, 1.23, 12.34, 123.45, 1234.5])
  17. assert_allclose(spherical_jn(2, x),
  18. (-1/x + 3/x**3)*sin(x) - 3/x**2*cos(x))
  19. def test_spherical_jn_recurrence_complex(self):
  20. # https://dlmf.nist.gov/10.51.E1
  21. n = np.array([1, 2, 3, 7, 12])
  22. x = 1.1 + 1.5j
  23. assert_allclose(spherical_jn(n - 1, x) + spherical_jn(n + 1, x),
  24. (2*n + 1)/x*spherical_jn(n, x))
  25. def test_spherical_jn_recurrence_real(self):
  26. # https://dlmf.nist.gov/10.51.E1
  27. n = np.array([1, 2, 3, 7, 12])
  28. x = 0.12
  29. assert_allclose(spherical_jn(n - 1, x) + spherical_jn(n + 1,x),
  30. (2*n + 1)/x*spherical_jn(n, x))
  31. def test_spherical_jn_inf_real(self):
  32. # https://dlmf.nist.gov/10.52.E3
  33. n = 6
  34. x = np.array([-inf, inf])
  35. assert_allclose(spherical_jn(n, x), np.array([0, 0]))
  36. def test_spherical_jn_inf_complex(self):
  37. # https://dlmf.nist.gov/10.52.E3
  38. n = 7
  39. x = np.array([-inf + 0j, inf + 0j, inf*(1+1j)])
  40. with warnings.catch_warnings():
  41. warnings.filterwarnings(
  42. "ignore", "invalid value encountered in multiply", RuntimeWarning)
  43. assert_allclose(spherical_jn(n, x), np.array([0, 0, inf*(1+1j)]))
  44. def test_spherical_jn_large_arg_1(self):
  45. # https://github.com/scipy/scipy/issues/2165
  46. # Reference value computed using mpmath, via
  47. # besselj(n + mpf(1)/2, z)*sqrt(pi/(2*z))
  48. assert_allclose(spherical_jn(2, 3350.507), -0.00029846226538040747)
  49. def test_spherical_jn_large_arg_2(self):
  50. # https://github.com/scipy/scipy/issues/1641
  51. # Reference value computed using mpmath, via
  52. # besselj(n + mpf(1)/2, z)*sqrt(pi/(2*z))
  53. assert_allclose(spherical_jn(2, 10000), 3.0590002633029811e-05)
  54. def test_spherical_jn_at_zero(self):
  55. # https://dlmf.nist.gov/10.52.E1
  56. # But note that n = 0 is a special case: j0 = sin(x)/x -> 1
  57. n = np.array([0, 1, 2, 5, 10, 100])
  58. x = 0
  59. assert_allclose(spherical_jn(n, x), np.array([1, 0, 0, 0, 0, 0]))
  60. class TestSphericalYn:
  61. def test_spherical_yn_exact(self):
  62. # https://dlmf.nist.gov/10.49.E5
  63. # Note: exact expression is numerically stable only for small
  64. # n or z >> n.
  65. x = np.array([0.12, 1.23, 12.34, 123.45, 1234.5])
  66. assert_allclose(spherical_yn(2, x),
  67. (1/x - 3/x**3)*cos(x) - 3/x**2*sin(x))
  68. def test_spherical_yn_recurrence_real(self):
  69. # https://dlmf.nist.gov/10.51.E1
  70. n = np.array([1, 2, 3, 7, 12])
  71. x = 0.12
  72. assert_allclose(spherical_yn(n - 1, x) + spherical_yn(n + 1,x),
  73. (2*n + 1)/x*spherical_yn(n, x))
  74. def test_spherical_yn_recurrence_complex(self):
  75. # https://dlmf.nist.gov/10.51.E1
  76. n = np.array([1, 2, 3, 7, 12])
  77. x = 1.1 + 1.5j
  78. assert_allclose(spherical_yn(n - 1, x) + spherical_yn(n + 1, x),
  79. (2*n + 1)/x*spherical_yn(n, x))
  80. def test_spherical_yn_inf_real(self):
  81. # https://dlmf.nist.gov/10.52.E3
  82. n = 6
  83. x = np.array([-inf, inf])
  84. assert_allclose(spherical_yn(n, x), np.array([0, 0]))
  85. def test_spherical_yn_inf_complex(self):
  86. # https://dlmf.nist.gov/10.52.E3
  87. n = 7
  88. x = np.array([-inf + 0j, inf + 0j, inf*(1+1j)])
  89. with warnings.catch_warnings():
  90. warnings.filterwarnings(
  91. "ignore", "invalid value encountered in multiply", RuntimeWarning)
  92. assert_allclose(spherical_yn(n, x), np.array([0, 0, inf*(1+1j)]))
  93. def test_spherical_yn_at_zero(self):
  94. # https://dlmf.nist.gov/10.52.E2
  95. n = np.array([0, 1, 2, 5, 10, 100])
  96. x = 0
  97. assert_allclose(spherical_yn(n, x), np.full(n.shape, -inf))
  98. def test_spherical_yn_at_zero_complex(self):
  99. # Consistently with numpy:
  100. # >>> -np.cos(0)/0
  101. # -inf
  102. # >>> -np.cos(0+0j)/(0+0j)
  103. # (-inf + nan*j)
  104. n = np.array([0, 1, 2, 5, 10, 100])
  105. x = 0 + 0j
  106. assert_allclose(spherical_yn(n, x), np.full(n.shape, nan))
  107. class TestSphericalJnYnCrossProduct:
  108. def test_spherical_jn_yn_cross_product_1(self):
  109. # https://dlmf.nist.gov/10.50.E3
  110. n = np.array([1, 5, 8])
  111. x = np.array([0.1, 1, 10])
  112. left = (spherical_jn(n + 1, x) * spherical_yn(n, x) -
  113. spherical_jn(n, x) * spherical_yn(n + 1, x))
  114. right = 1/x**2
  115. assert_allclose(left, right)
  116. def test_spherical_jn_yn_cross_product_2(self):
  117. # https://dlmf.nist.gov/10.50.E3
  118. n = np.array([1, 5, 8])
  119. x = np.array([0.1, 1, 10])
  120. left = (spherical_jn(n + 2, x) * spherical_yn(n, x) -
  121. spherical_jn(n, x) * spherical_yn(n + 2, x))
  122. right = (2*n + 3)/x**3
  123. assert_allclose(left, right)
  124. class TestSphericalIn:
  125. def test_spherical_in_exact(self):
  126. # https://dlmf.nist.gov/10.49.E9
  127. x = np.array([0.12, 1.23, 12.34, 123.45])
  128. assert_allclose(spherical_in(2, x),
  129. (1/x + 3/x**3)*sinh(x) - 3/x**2*cosh(x))
  130. def test_spherical_in_recurrence_real(self):
  131. # https://dlmf.nist.gov/10.51.E4
  132. n = np.array([1, 2, 3, 7, 12])
  133. x = 0.12
  134. assert_allclose(spherical_in(n - 1, x) - spherical_in(n + 1,x),
  135. (2*n + 1)/x*spherical_in(n, x))
  136. def test_spherical_in_recurrence_complex(self):
  137. # https://dlmf.nist.gov/10.51.E1
  138. n = np.array([1, 2, 3, 7, 12])
  139. x = 1.1 + 1.5j
  140. assert_allclose(spherical_in(n - 1, x) - spherical_in(n + 1,x),
  141. (2*n + 1)/x*spherical_in(n, x))
  142. def test_spherical_in_inf_real(self):
  143. # https://dlmf.nist.gov/10.52.E3
  144. n = 5
  145. x = np.array([-inf, inf])
  146. assert_allclose(spherical_in(n, x), np.array([-inf, inf]))
  147. def test_spherical_in_inf_complex(self):
  148. # https://dlmf.nist.gov/10.52.E5
  149. # Ideally, i1n(n, 1j*inf) = 0 and i1n(n, (1+1j)*inf) = (1+1j)*inf, but
  150. # this appears impossible to achieve because C99 regards any complex
  151. # value with at least one infinite part as a complex infinity, so
  152. # 1j*inf cannot be distinguished from (1+1j)*inf. Therefore, nan is
  153. # the correct return value.
  154. n = 7
  155. x = np.array([-inf + 0j, inf + 0j, inf*(1+1j)])
  156. assert_allclose(spherical_in(n, x), np.array([-inf, inf, nan]))
  157. def test_spherical_in_at_zero(self):
  158. # https://dlmf.nist.gov/10.52.E1
  159. # But note that n = 0 is a special case: i0 = sinh(x)/x -> 1
  160. n = np.array([0, 1, 2, 5, 10, 100])
  161. x = 0
  162. assert_allclose(spherical_in(n, x), np.array([1, 0, 0, 0, 0, 0]))
  163. class TestSphericalKn:
  164. def test_spherical_kn_exact(self):
  165. # https://dlmf.nist.gov/10.49.E13
  166. x = np.array([0.12, 1.23, 12.34, 123.45])
  167. assert_allclose(spherical_kn(2, x),
  168. pi/2*exp(-x)*(1/x + 3/x**2 + 3/x**3))
  169. def test_spherical_kn_recurrence_real(self):
  170. # https://dlmf.nist.gov/10.51.E4
  171. n = np.array([1, 2, 3, 7, 12])
  172. x = 0.12
  173. assert_allclose(
  174. (-1)**(n - 1)*spherical_kn(n - 1, x) - (-1)**(n + 1)*spherical_kn(n + 1,x),
  175. (-1)**n*(2*n + 1)/x*spherical_kn(n, x)
  176. )
  177. def test_spherical_kn_recurrence_complex(self):
  178. # https://dlmf.nist.gov/10.51.E4
  179. n = np.array([1, 2, 3, 7, 12])
  180. x = 1.1 + 1.5j
  181. assert_allclose(
  182. (-1)**(n - 1)*spherical_kn(n - 1, x) - (-1)**(n + 1)*spherical_kn(n + 1,x),
  183. (-1)**n*(2*n + 1)/x*spherical_kn(n, x)
  184. )
  185. def test_spherical_kn_inf_real(self):
  186. # https://dlmf.nist.gov/10.52.E6
  187. n = 5
  188. x = np.array([-inf, inf])
  189. assert_allclose(spherical_kn(n, x), np.array([-inf, 0]))
  190. def test_spherical_kn_inf_complex(self):
  191. # https://dlmf.nist.gov/10.52.E6
  192. # The behavior at complex infinity depends on the sign of the real
  193. # part: if Re(z) >= 0, then the limit is 0; if Re(z) < 0, then it's
  194. # z*inf. This distinction cannot be captured, so we return nan.
  195. n = 7
  196. x = np.array([-inf + 0j, inf + 0j, inf*(1+1j)])
  197. assert_allclose(spherical_kn(n, x), np.array([-inf, 0, nan]))
  198. def test_spherical_kn_at_zero(self):
  199. # https://dlmf.nist.gov/10.52.E2
  200. n = np.array([0, 1, 2, 5, 10, 100])
  201. x = 0
  202. assert_allclose(spherical_kn(n, x), np.full(n.shape, inf))
  203. def test_spherical_kn_at_zero_complex(self):
  204. # https://dlmf.nist.gov/10.52.E2
  205. n = np.array([0, 1, 2, 5, 10, 100])
  206. x = 0 + 0j
  207. assert_allclose(spherical_kn(n, x), np.full(n.shape, nan))
  208. class SphericalDerivativesTestCase:
  209. def fundamental_theorem(self, n, a, b):
  210. integral, tolerance = quad(lambda z: self.df(n, z), a, b)
  211. assert_allclose(integral,
  212. self.f(n, b) - self.f(n, a),
  213. atol=tolerance)
  214. @pytest.mark.slow
  215. def test_fundamental_theorem_0(self):
  216. self.fundamental_theorem(0, 3.0, 15.0)
  217. @pytest.mark.slow
  218. def test_fundamental_theorem_7(self):
  219. self.fundamental_theorem(7, 0.5, 1.2)
  220. class TestSphericalJnDerivatives(SphericalDerivativesTestCase):
  221. def f(self, n, z):
  222. return spherical_jn(n, z)
  223. def df(self, n, z):
  224. return spherical_jn(n, z, derivative=True)
  225. def test_spherical_jn_d_zero(self):
  226. n = np.array([0, 1, 2, 3, 7, 15])
  227. assert_allclose(spherical_jn(n, 0, derivative=True),
  228. np.array([0, 1/3, 0, 0, 0, 0]))
  229. class TestSphericalYnDerivatives(SphericalDerivativesTestCase):
  230. def f(self, n, z):
  231. return spherical_yn(n, z)
  232. def df(self, n, z):
  233. return spherical_yn(n, z, derivative=True)
  234. class TestSphericalInDerivatives(SphericalDerivativesTestCase):
  235. def f(self, n, z):
  236. return spherical_in(n, z)
  237. def df(self, n, z):
  238. return spherical_in(n, z, derivative=True)
  239. def test_spherical_in_d_zero(self):
  240. n = np.array([0, 1, 2, 3, 7, 15])
  241. spherical_in(n, 0, derivative=False)
  242. assert_allclose(spherical_in(n, 0, derivative=True),
  243. np.array([0, 1/3, 0, 0, 0, 0]))
  244. class TestSphericalKnDerivatives(SphericalDerivativesTestCase):
  245. def f(self, n, z):
  246. return spherical_kn(n, z)
  247. def df(self, n, z):
  248. return spherical_kn(n, z, derivative=True)
  249. class TestSphericalOld:
  250. # These are tests from the TestSpherical class of test_basic.py,
  251. # rewritten to use spherical_* instead of sph_* but otherwise unchanged.
  252. def test_sph_in(self):
  253. # This test reproduces test_basic.TestSpherical.test_sph_in.
  254. i1n = np.empty((2,2))
  255. x = 0.2
  256. i1n[0][0] = spherical_in(0, x)
  257. i1n[0][1] = spherical_in(1, x)
  258. i1n[1][0] = spherical_in(0, x, derivative=True)
  259. i1n[1][1] = spherical_in(1, x, derivative=True)
  260. inp0 = (i1n[0][1])
  261. inp1 = (i1n[0][0] - 2.0/0.2 * i1n[0][1])
  262. assert_allclose(i1n[0], np.array([1.0066800127054699381,
  263. 0.066933714568029540839]),
  264. atol=1.5e-12, rtol=0.0)
  265. assert_allclose(i1n[1], [inp0, inp1], atol=1.5e-12, rtol=0)
  266. def test_sph_in_kn_order0(self):
  267. x = 1.
  268. sph_i0 = np.empty((2,))
  269. sph_i0[0] = spherical_in(0, x)
  270. sph_i0[1] = spherical_in(0, x, derivative=True)
  271. sph_i0_expected = np.array([np.sinh(x)/x,
  272. np.cosh(x)/x-np.sinh(x)/x**2])
  273. assert_allclose(r_[sph_i0], sph_i0_expected, atol=1.5e-7, rtol=0)
  274. sph_k0 = np.empty((2,))
  275. sph_k0[0] = spherical_kn(0, x)
  276. sph_k0[1] = spherical_kn(0, x, derivative=True)
  277. sph_k0_expected = np.array([0.5*pi*exp(-x)/x,
  278. -0.5*pi*exp(-x)*(1/x+1/x**2)])
  279. assert_allclose(r_[sph_k0], sph_k0_expected, atol=1.5e-7, rtol=0)
  280. def test_sph_jn(self):
  281. s1 = np.empty((2,3))
  282. x = 0.2
  283. s1[0][0] = spherical_jn(0, x)
  284. s1[0][1] = spherical_jn(1, x)
  285. s1[0][2] = spherical_jn(2, x)
  286. s1[1][0] = spherical_jn(0, x, derivative=True)
  287. s1[1][1] = spherical_jn(1, x, derivative=True)
  288. s1[1][2] = spherical_jn(2, x, derivative=True)
  289. s10 = -s1[0][1]
  290. s11 = s1[0][0]-2.0/0.2*s1[0][1]
  291. s12 = s1[0][1]-3.0/0.2*s1[0][2]
  292. assert_allclose(s1[0], [0.99334665397530607731,
  293. 0.066400380670322230863,
  294. 0.0026590560795273856680],
  295. atol=1.5e-12, rtol=0)
  296. assert_allclose(s1[1], [s10, s11, s12], atol=1.5e-12, rtol=0)
  297. def test_sph_kn(self):
  298. kn = np.empty((2,3))
  299. x = 0.2
  300. kn[0][0] = spherical_kn(0, x)
  301. kn[0][1] = spherical_kn(1, x)
  302. kn[0][2] = spherical_kn(2, x)
  303. kn[1][0] = spherical_kn(0, x, derivative=True)
  304. kn[1][1] = spherical_kn(1, x, derivative=True)
  305. kn[1][2] = spherical_kn(2, x, derivative=True)
  306. kn0 = -kn[0][1]
  307. kn1 = -kn[0][0]-2.0/0.2*kn[0][1]
  308. kn2 = -kn[0][1]-3.0/0.2*kn[0][2]
  309. assert_allclose(kn[0], [6.4302962978445670140,
  310. 38.581777787067402086,
  311. 585.15696310385559829],
  312. atol=1.5e-12, rtol=0)
  313. assert_allclose(kn[1], [kn0, kn1, kn2], atol=1.5e-9, rtol=0)
  314. def test_sph_yn(self):
  315. sy1 = spherical_yn(2, 0.2)
  316. sy2 = spherical_yn(0, 0.2)
  317. # previous values in the system
  318. assert_allclose(sy1, -377.52483, atol=1.5e-5, rtol=0)
  319. assert_allclose(sy2, -4.9003329, atol=1.5e-5, rtol=0)
  320. sphpy = (spherical_yn(0, 0.2) - 2*spherical_yn(2, 0.2))/3
  321. sy3 = spherical_yn(1, 0.2, derivative=True)
  322. # compare correct derivative val. (correct =-system val).
  323. assert_allclose(sy3, sphpy, atol=1.5e-4, rtol=0)
  324. @pytest.mark.parametrize('derivative', [False, True])
  325. @pytest.mark.parametrize('fun', [spherical_jn, spherical_in,
  326. spherical_yn, spherical_kn])
  327. def test_negative_real_gh14582(derivative, fun):
  328. # gh-14582 reported that the spherical Bessel functions did not work
  329. # with negative real argument `z`. Check that this is resolved.
  330. rng = np.random.default_rng(3598435982345987234)
  331. size = 25
  332. n = rng.integers(0, 10, size=size)
  333. z = rng.standard_normal(size=size)
  334. res = fun(n, z, derivative=derivative)
  335. ref = fun(n, z+0j, derivative=derivative)
  336. assert_allclose(res, ref.real)