test_quadrature.py 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397
  1. import pytest
  2. import numpy as np
  3. from numpy import cos, sin, pi
  4. from numpy.testing import (assert_equal, assert_almost_equal, assert_allclose,
  5. assert_, suppress_warnings)
  6. from scipy.integrate import (quadrature, romberg, romb, newton_cotes,
  7. cumulative_trapezoid, cumtrapz, trapz, trapezoid,
  8. quad, simpson, simps, fixed_quad, AccuracyWarning)
  9. from scipy.integrate._quadrature import _qmc_quad as qmc_quad
  10. from scipy import stats, special as sc
  11. class TestFixedQuad:
  12. def test_scalar(self):
  13. n = 4
  14. expected = 1/(2*n)
  15. got, _ = fixed_quad(lambda x: x**(2*n - 1), 0, 1, n=n)
  16. # quadrature exact for this input
  17. assert_allclose(got, expected, rtol=1e-12)
  18. def test_vector(self):
  19. n = 4
  20. p = np.arange(1, 2*n)
  21. expected = 1/(p + 1)
  22. got, _ = fixed_quad(lambda x: x**p[:, None], 0, 1, n=n)
  23. assert_allclose(got, expected, rtol=1e-12)
  24. class TestQuadrature:
  25. def quad(self, x, a, b, args):
  26. raise NotImplementedError
  27. def test_quadrature(self):
  28. # Typical function with two extra arguments:
  29. def myfunc(x, n, z): # Bessel function integrand
  30. return cos(n*x-z*sin(x))/pi
  31. val, err = quadrature(myfunc, 0, pi, (2, 1.8))
  32. table_val = 0.30614353532540296487
  33. assert_almost_equal(val, table_val, decimal=7)
  34. def test_quadrature_rtol(self):
  35. def myfunc(x, n, z): # Bessel function integrand
  36. return 1e90 * cos(n*x-z*sin(x))/pi
  37. val, err = quadrature(myfunc, 0, pi, (2, 1.8), rtol=1e-10)
  38. table_val = 1e90 * 0.30614353532540296487
  39. assert_allclose(val, table_val, rtol=1e-10)
  40. def test_quadrature_miniter(self):
  41. # Typical function with two extra arguments:
  42. def myfunc(x, n, z): # Bessel function integrand
  43. return cos(n*x-z*sin(x))/pi
  44. table_val = 0.30614353532540296487
  45. for miniter in [5, 52]:
  46. val, err = quadrature(myfunc, 0, pi, (2, 1.8), miniter=miniter)
  47. assert_almost_equal(val, table_val, decimal=7)
  48. assert_(err < 1.0)
  49. def test_quadrature_single_args(self):
  50. def myfunc(x, n):
  51. return 1e90 * cos(n*x-1.8*sin(x))/pi
  52. val, err = quadrature(myfunc, 0, pi, args=2, rtol=1e-10)
  53. table_val = 1e90 * 0.30614353532540296487
  54. assert_allclose(val, table_val, rtol=1e-10)
  55. def test_romberg(self):
  56. # Typical function with two extra arguments:
  57. def myfunc(x, n, z): # Bessel function integrand
  58. return cos(n*x-z*sin(x))/pi
  59. val = romberg(myfunc, 0, pi, args=(2, 1.8))
  60. table_val = 0.30614353532540296487
  61. assert_almost_equal(val, table_val, decimal=7)
  62. def test_romberg_rtol(self):
  63. # Typical function with two extra arguments:
  64. def myfunc(x, n, z): # Bessel function integrand
  65. return 1e19*cos(n*x-z*sin(x))/pi
  66. val = romberg(myfunc, 0, pi, args=(2, 1.8), rtol=1e-10)
  67. table_val = 1e19*0.30614353532540296487
  68. assert_allclose(val, table_val, rtol=1e-10)
  69. def test_romb(self):
  70. assert_equal(romb(np.arange(17)), 128)
  71. def test_romb_gh_3731(self):
  72. # Check that romb makes maximal use of data points
  73. x = np.arange(2**4+1)
  74. y = np.cos(0.2*x)
  75. val = romb(y)
  76. val2, err = quad(lambda x: np.cos(0.2*x), x.min(), x.max())
  77. assert_allclose(val, val2, rtol=1e-8, atol=0)
  78. # should be equal to romb with 2**k+1 samples
  79. with suppress_warnings() as sup:
  80. sup.filter(AccuracyWarning, "divmax .4. exceeded")
  81. val3 = romberg(lambda x: np.cos(0.2*x), x.min(), x.max(), divmax=4)
  82. assert_allclose(val, val3, rtol=1e-12, atol=0)
  83. def test_non_dtype(self):
  84. # Check that we work fine with functions returning float
  85. import math
  86. valmath = romberg(math.sin, 0, 1)
  87. expected_val = 0.45969769413185085
  88. assert_almost_equal(valmath, expected_val, decimal=7)
  89. def test_newton_cotes(self):
  90. """Test the first few degrees, for evenly spaced points."""
  91. n = 1
  92. wts, errcoff = newton_cotes(n, 1)
  93. assert_equal(wts, n*np.array([0.5, 0.5]))
  94. assert_almost_equal(errcoff, -n**3/12.0)
  95. n = 2
  96. wts, errcoff = newton_cotes(n, 1)
  97. assert_almost_equal(wts, n*np.array([1.0, 4.0, 1.0])/6.0)
  98. assert_almost_equal(errcoff, -n**5/2880.0)
  99. n = 3
  100. wts, errcoff = newton_cotes(n, 1)
  101. assert_almost_equal(wts, n*np.array([1.0, 3.0, 3.0, 1.0])/8.0)
  102. assert_almost_equal(errcoff, -n**5/6480.0)
  103. n = 4
  104. wts, errcoff = newton_cotes(n, 1)
  105. assert_almost_equal(wts, n*np.array([7.0, 32.0, 12.0, 32.0, 7.0])/90.0)
  106. assert_almost_equal(errcoff, -n**7/1935360.0)
  107. def test_newton_cotes2(self):
  108. """Test newton_cotes with points that are not evenly spaced."""
  109. x = np.array([0.0, 1.5, 2.0])
  110. y = x**2
  111. wts, errcoff = newton_cotes(x)
  112. exact_integral = 8.0/3
  113. numeric_integral = np.dot(wts, y)
  114. assert_almost_equal(numeric_integral, exact_integral)
  115. x = np.array([0.0, 1.4, 2.1, 3.0])
  116. y = x**2
  117. wts, errcoff = newton_cotes(x)
  118. exact_integral = 9.0
  119. numeric_integral = np.dot(wts, y)
  120. assert_almost_equal(numeric_integral, exact_integral)
  121. def test_simpson(self):
  122. y = np.arange(17)
  123. assert_equal(simpson(y), 128)
  124. assert_equal(simpson(y, dx=0.5), 64)
  125. assert_equal(simpson(y, x=np.linspace(0, 4, 17)), 32)
  126. y = np.arange(4)
  127. x = 2**y
  128. assert_equal(simpson(y, x=x, even='avg'), 13.875)
  129. assert_equal(simpson(y, x=x, even='first'), 13.75)
  130. assert_equal(simpson(y, x=x, even='last'), 14)
  131. # Tests for checking base case
  132. x = np.array([3])
  133. y = np.power(x, 2)
  134. assert_equal(simpson(y, x=x, axis=0), 0.0)
  135. assert_equal(simpson(y, x=x, axis=-1), 0.0)
  136. x = np.array([3, 3, 3, 3])
  137. y = np.power(x, 2)
  138. assert_equal(simpson(y, x=x, axis=0), 0.0)
  139. assert_equal(simpson(y, x=x, axis=-1), 0.0)
  140. x = np.array([[1, 2, 4, 8], [1, 2, 4, 8], [1, 2, 4, 8]])
  141. y = np.power(x, 2)
  142. zero_axis = [0.0, 0.0, 0.0, 0.0]
  143. default_axis = [175.75, 175.75, 175.75]
  144. assert_equal(simpson(y, x=x, axis=0), zero_axis)
  145. assert_equal(simpson(y, x=x, axis=-1), default_axis)
  146. x = np.array([[1, 2, 4, 8], [1, 2, 4, 8], [1, 8, 16, 32]])
  147. y = np.power(x, 2)
  148. zero_axis = [0.0, 136.0, 1088.0, 8704.0]
  149. default_axis = [175.75, 175.75, 11292.25]
  150. assert_equal(simpson(y, x=x, axis=0), zero_axis)
  151. assert_equal(simpson(y, x=x, axis=-1), default_axis)
  152. @pytest.mark.parametrize('droplast', [False, True])
  153. def test_simpson_2d_integer_no_x(self, droplast):
  154. # The inputs are 2d integer arrays. The results should be
  155. # identical to the results when the inputs are floating point.
  156. y = np.array([[2, 2, 4, 4, 8, 8, -4, 5],
  157. [4, 4, 2, -4, 10, 22, -2, 10]])
  158. if droplast:
  159. y = y[:, :-1]
  160. result = simpson(y, axis=-1)
  161. expected = simpson(np.array(y, dtype=np.float64), axis=-1)
  162. assert_equal(result, expected)
  163. def test_simps(self):
  164. # Basic coverage test for the alias
  165. y = np.arange(4)
  166. x = 2**y
  167. assert_equal(simpson(y, x=x, dx=0.5, even='first'),
  168. simps(y, x=x, dx=0.5, even='first'))
  169. class TestCumulative_trapezoid:
  170. def test_1d(self):
  171. x = np.linspace(-2, 2, num=5)
  172. y = x
  173. y_int = cumulative_trapezoid(y, x, initial=0)
  174. y_expected = [0., -1.5, -2., -1.5, 0.]
  175. assert_allclose(y_int, y_expected)
  176. y_int = cumulative_trapezoid(y, x, initial=None)
  177. assert_allclose(y_int, y_expected[1:])
  178. def test_y_nd_x_nd(self):
  179. x = np.arange(3 * 2 * 4).reshape(3, 2, 4)
  180. y = x
  181. y_int = cumulative_trapezoid(y, x, initial=0)
  182. y_expected = np.array([[[0., 0.5, 2., 4.5],
  183. [0., 4.5, 10., 16.5]],
  184. [[0., 8.5, 18., 28.5],
  185. [0., 12.5, 26., 40.5]],
  186. [[0., 16.5, 34., 52.5],
  187. [0., 20.5, 42., 64.5]]])
  188. assert_allclose(y_int, y_expected)
  189. # Try with all axes
  190. shapes = [(2, 2, 4), (3, 1, 4), (3, 2, 3)]
  191. for axis, shape in zip([0, 1, 2], shapes):
  192. y_int = cumulative_trapezoid(y, x, initial=3.45, axis=axis)
  193. assert_equal(y_int.shape, (3, 2, 4))
  194. y_int = cumulative_trapezoid(y, x, initial=None, axis=axis)
  195. assert_equal(y_int.shape, shape)
  196. def test_y_nd_x_1d(self):
  197. y = np.arange(3 * 2 * 4).reshape(3, 2, 4)
  198. x = np.arange(4)**2
  199. # Try with all axes
  200. ys_expected = (
  201. np.array([[[4., 5., 6., 7.],
  202. [8., 9., 10., 11.]],
  203. [[40., 44., 48., 52.],
  204. [56., 60., 64., 68.]]]),
  205. np.array([[[2., 3., 4., 5.]],
  206. [[10., 11., 12., 13.]],
  207. [[18., 19., 20., 21.]]]),
  208. np.array([[[0.5, 5., 17.5],
  209. [4.5, 21., 53.5]],
  210. [[8.5, 37., 89.5],
  211. [12.5, 53., 125.5]],
  212. [[16.5, 69., 161.5],
  213. [20.5, 85., 197.5]]]))
  214. for axis, y_expected in zip([0, 1, 2], ys_expected):
  215. y_int = cumulative_trapezoid(y, x=x[:y.shape[axis]], axis=axis,
  216. initial=None)
  217. assert_allclose(y_int, y_expected)
  218. def test_x_none(self):
  219. y = np.linspace(-2, 2, num=5)
  220. y_int = cumulative_trapezoid(y)
  221. y_expected = [-1.5, -2., -1.5, 0.]
  222. assert_allclose(y_int, y_expected)
  223. y_int = cumulative_trapezoid(y, initial=1.23)
  224. y_expected = [1.23, -1.5, -2., -1.5, 0.]
  225. assert_allclose(y_int, y_expected)
  226. y_int = cumulative_trapezoid(y, dx=3)
  227. y_expected = [-4.5, -6., -4.5, 0.]
  228. assert_allclose(y_int, y_expected)
  229. y_int = cumulative_trapezoid(y, dx=3, initial=1.23)
  230. y_expected = [1.23, -4.5, -6., -4.5, 0.]
  231. assert_allclose(y_int, y_expected)
  232. def test_cumtrapz(self):
  233. # Basic coverage test for the alias
  234. x = np.arange(3 * 2 * 4).reshape(3, 2, 4)
  235. y = x
  236. assert_allclose(cumulative_trapezoid(y, x, dx=0.5, axis=0, initial=0),
  237. cumtrapz(y, x, dx=0.5, axis=0, initial=0),
  238. rtol=1e-14)
  239. class TestTrapezoid():
  240. """This function is tested in NumPy more extensive, just do some
  241. basic due diligence here."""
  242. def test_trapezoid(self):
  243. y = np.arange(17)
  244. assert_equal(trapezoid(y), 128)
  245. assert_equal(trapezoid(y, dx=0.5), 64)
  246. assert_equal(trapezoid(y, x=np.linspace(0, 4, 17)), 32)
  247. y = np.arange(4)
  248. x = 2**y
  249. assert_equal(trapezoid(y, x=x, dx=0.1), 13.5)
  250. def test_trapz(self):
  251. # Basic coverage test for the alias
  252. y = np.arange(4)
  253. x = 2**y
  254. assert_equal(trapezoid(y, x=x, dx=0.5, axis=0),
  255. trapz(y, x=x, dx=0.5, axis=0))
  256. class TestQMCQuad():
  257. def test_input_validation(self):
  258. message = "`func` must be callable."
  259. with pytest.raises(TypeError, match=message):
  260. qmc_quad("a duck", [0, 0], [1, 1])
  261. message = "`func` must evaluate the integrand at points..."
  262. with pytest.raises(ValueError, match=message):
  263. qmc_quad(lambda: 1, [0, 0], [1, 1])
  264. def func(x):
  265. assert x.ndim == 1
  266. return np.sum(x)
  267. message = "Exception encountered when attempting vectorized call..."
  268. with pytest.warns(UserWarning, match=message):
  269. qmc_quad(func, [0, 0], [1, 1])
  270. message = "`n_points` must be an integer."
  271. with pytest.raises(TypeError, match=message):
  272. qmc_quad(lambda x: 1, [0, 0], [1, 1], n_points=1024.5)
  273. message = "`n_estimates` must be an integer."
  274. with pytest.raises(TypeError, match=message):
  275. qmc_quad(lambda x: 1, [0, 0], [1, 1], n_estimates=8.5)
  276. message = "`qrng` must be an instance of scipy.stats.qmc.QMCEngine."
  277. with pytest.raises(TypeError, match=message):
  278. qmc_quad(lambda x: 1, [0, 0], [1, 1], qrng="a duck")
  279. message = "`qrng` must be initialized with dimensionality equal to "
  280. with pytest.raises(ValueError, match=message):
  281. qmc_quad(lambda x: 1, [0, 0], [1, 1], qrng=stats.qmc.Sobol(1))
  282. message = r"`log` must be boolean \(`True` or `False`\)."
  283. with pytest.raises(TypeError, match=message):
  284. qmc_quad(lambda x: 1, [0, 0], [1, 1], log=10)
  285. def basic_test(self, n_points=2**8, n_estimates=8, signs=np.ones(2)):
  286. ndim = 2
  287. mean = np.zeros(ndim)
  288. cov = np.eye(ndim)
  289. def func(x):
  290. return stats.multivariate_normal.pdf(x, mean, cov)
  291. rng = np.random.default_rng(2879434385674690281)
  292. qrng = stats.qmc.Sobol(ndim, seed=rng)
  293. a = np.zeros(ndim)
  294. b = np.ones(ndim) * signs
  295. res = qmc_quad(func, a, b, n_points=n_points,
  296. n_estimates=n_estimates, args=(mean, cov), qrng=qrng)
  297. ref = stats.multivariate_normal.cdf(b, mean, cov, lower_limit=a)
  298. atol = sc.stdtrit(n_estimates-1, 0.995) * res.standard_error # 99% CI
  299. assert_allclose(res.integral, ref, atol=atol)
  300. assert np.prod(signs)*res.integral > 0
  301. rng = np.random.default_rng(2879434385674690281)
  302. qrng = stats.qmc.Sobol(ndim, seed=rng)
  303. logres = qmc_quad(lambda *args: np.log(func(*args)), a, b,
  304. n_points=n_points, n_estimates=n_estimates,
  305. args=(mean, cov), log=True, qrng=qrng)
  306. assert_allclose(np.exp(logres.integral), res.integral)
  307. assert np.imag(logres.integral) == (np.pi if np.prod(signs) < 0 else 0)
  308. @pytest.mark.parametrize("n_points", [2**8, 2**12])
  309. @pytest.mark.parametrize("n_estimates", [8, 16])
  310. def test_basic(self, n_points, n_estimates):
  311. self.basic_test(n_points, n_estimates)
  312. @pytest.mark.parametrize("signs", [[1, 1], [-1, -1], [-1, 1], [1, -1]])
  313. def test_sign(self, signs):
  314. self.basic_test(signs=signs)
  315. @pytest.mark.parametrize("log", [False, True])
  316. def test_zero(self, log):
  317. message = "A lower limit was equal to an upper limit, so"
  318. with pytest.warns(UserWarning, match=message):
  319. res = qmc_quad(lambda x: 1, [0, 0], [0, 1], log=log)
  320. assert res.integral == (-np.inf if log else 0)
  321. assert res.standard_error == 0
  322. def test_flexible_input(self):
  323. # check that qrng is not required
  324. # also checks that for 1d problems, a and b can be scalars
  325. def func(x):
  326. return stats.norm.pdf(x, scale=2)
  327. res = qmc_quad(func, 0, 1)
  328. ref = stats.norm.cdf(1, scale=2) - stats.norm.cdf(0, scale=2)
  329. assert_allclose(res.integral, ref, 1e-2)