test__quad_vec.py 6.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204
  1. import pytest
  2. import numpy as np
  3. from numpy.testing import assert_allclose
  4. from scipy.integrate import quad_vec
  5. from multiprocessing.dummy import Pool
  6. quadrature_params = pytest.mark.parametrize(
  7. 'quadrature', [None, "gk15", "gk21", "trapezoid"])
  8. @quadrature_params
  9. def test_quad_vec_simple(quadrature):
  10. n = np.arange(10)
  11. f = lambda x: x**n
  12. for epsabs in [0.1, 1e-3, 1e-6]:
  13. if quadrature == 'trapezoid' and epsabs < 1e-4:
  14. # slow: skip
  15. continue
  16. kwargs = dict(epsabs=epsabs, quadrature=quadrature)
  17. exact = 2**(n+1)/(n + 1)
  18. res, err = quad_vec(f, 0, 2, norm='max', **kwargs)
  19. assert_allclose(res, exact, rtol=0, atol=epsabs)
  20. res, err = quad_vec(f, 0, 2, norm='2', **kwargs)
  21. assert np.linalg.norm(res - exact) < epsabs
  22. res, err = quad_vec(f, 0, 2, norm='max', points=(0.5, 1.0), **kwargs)
  23. assert_allclose(res, exact, rtol=0, atol=epsabs)
  24. res, err, *rest = quad_vec(f, 0, 2, norm='max',
  25. epsrel=1e-8,
  26. full_output=True,
  27. limit=10000,
  28. **kwargs)
  29. assert_allclose(res, exact, rtol=0, atol=epsabs)
  30. @quadrature_params
  31. def test_quad_vec_simple_inf(quadrature):
  32. f = lambda x: 1 / (1 + np.float64(x)**2)
  33. for epsabs in [0.1, 1e-3, 1e-6]:
  34. if quadrature == 'trapezoid' and epsabs < 1e-4:
  35. # slow: skip
  36. continue
  37. kwargs = dict(norm='max', epsabs=epsabs, quadrature=quadrature)
  38. res, err = quad_vec(f, 0, np.inf, **kwargs)
  39. assert_allclose(res, np.pi/2, rtol=0, atol=max(epsabs, err))
  40. res, err = quad_vec(f, 0, -np.inf, **kwargs)
  41. assert_allclose(res, -np.pi/2, rtol=0, atol=max(epsabs, err))
  42. res, err = quad_vec(f, -np.inf, 0, **kwargs)
  43. assert_allclose(res, np.pi/2, rtol=0, atol=max(epsabs, err))
  44. res, err = quad_vec(f, np.inf, 0, **kwargs)
  45. assert_allclose(res, -np.pi/2, rtol=0, atol=max(epsabs, err))
  46. res, err = quad_vec(f, -np.inf, np.inf, **kwargs)
  47. assert_allclose(res, np.pi, rtol=0, atol=max(epsabs, err))
  48. res, err = quad_vec(f, np.inf, -np.inf, **kwargs)
  49. assert_allclose(res, -np.pi, rtol=0, atol=max(epsabs, err))
  50. res, err = quad_vec(f, np.inf, np.inf, **kwargs)
  51. assert_allclose(res, 0, rtol=0, atol=max(epsabs, err))
  52. res, err = quad_vec(f, -np.inf, -np.inf, **kwargs)
  53. assert_allclose(res, 0, rtol=0, atol=max(epsabs, err))
  54. res, err = quad_vec(f, 0, np.inf, points=(1.0, 2.0), **kwargs)
  55. assert_allclose(res, np.pi/2, rtol=0, atol=max(epsabs, err))
  56. f = lambda x: np.sin(x + 2) / (1 + x**2)
  57. exact = np.pi / np.e * np.sin(2)
  58. epsabs = 1e-5
  59. res, err, info = quad_vec(f, -np.inf, np.inf, limit=1000, norm='max', epsabs=epsabs,
  60. quadrature=quadrature, full_output=True)
  61. assert info.status == 1
  62. assert_allclose(res, exact, rtol=0, atol=max(epsabs, 1.5 * err))
  63. def test_quad_vec_args():
  64. f = lambda x, a: x * (x + a) * np.arange(3)
  65. a = 2
  66. exact = np.array([0, 4/3, 8/3])
  67. res, err = quad_vec(f, 0, 1, args=(a,))
  68. assert_allclose(res, exact, rtol=0, atol=1e-4)
  69. def _lorenzian(x):
  70. return 1 / (1 + x**2)
  71. def test_quad_vec_pool():
  72. f = _lorenzian
  73. res, err = quad_vec(f, -np.inf, np.inf, norm='max', epsabs=1e-4, workers=4)
  74. assert_allclose(res, np.pi, rtol=0, atol=1e-4)
  75. with Pool(10) as pool:
  76. f = lambda x: 1 / (1 + x**2)
  77. res, err = quad_vec(f, -np.inf, np.inf, norm='max', epsabs=1e-4, workers=pool.map)
  78. assert_allclose(res, np.pi, rtol=0, atol=1e-4)
  79. def _func_with_args(x, a):
  80. return x * (x + a) * np.arange(3)
  81. @pytest.mark.parametrize('extra_args', [2, (2,)])
  82. @pytest.mark.parametrize('workers', [1, 10])
  83. def test_quad_vec_pool_args(extra_args, workers):
  84. f = _func_with_args
  85. exact = np.array([0, 4/3, 8/3])
  86. res, err = quad_vec(f, 0, 1, args=extra_args, workers=workers)
  87. assert_allclose(res, exact, rtol=0, atol=1e-4)
  88. with Pool(workers) as pool:
  89. res, err = quad_vec(f, 0, 1, args=extra_args, workers=pool.map)
  90. assert_allclose(res, exact, rtol=0, atol=1e-4)
  91. @quadrature_params
  92. def test_num_eval(quadrature):
  93. def f(x):
  94. count[0] += 1
  95. return x**5
  96. count = [0]
  97. res = quad_vec(f, 0, 1, norm='max', full_output=True, quadrature=quadrature)
  98. assert res[2].neval == count[0]
  99. def test_info():
  100. def f(x):
  101. return np.ones((3, 2, 1))
  102. res, err, info = quad_vec(f, 0, 1, norm='max', full_output=True)
  103. assert info.success == True
  104. assert info.status == 0
  105. assert info.message == 'Target precision reached.'
  106. assert info.neval > 0
  107. assert info.intervals.shape[1] == 2
  108. assert info.integrals.shape == (info.intervals.shape[0], 3, 2, 1)
  109. assert info.errors.shape == (info.intervals.shape[0],)
  110. def test_nan_inf():
  111. def f_nan(x):
  112. return np.nan
  113. def f_inf(x):
  114. return np.inf if x < 0.1 else 1/x
  115. res, err, info = quad_vec(f_nan, 0, 1, full_output=True)
  116. assert info.status == 3
  117. res, err, info = quad_vec(f_inf, 0, 1, full_output=True)
  118. assert info.status == 3
  119. @pytest.mark.parametrize('a,b', [(0, 1), (0, np.inf), (np.inf, 0),
  120. (-np.inf, np.inf), (np.inf, -np.inf)])
  121. def test_points(a, b):
  122. # Check that initial interval splitting is done according to
  123. # `points`, by checking that consecutive sets of 15 point (for
  124. # gk15) function evaluations lie between `points`
  125. points = (0, 0.25, 0.5, 0.75, 1.0)
  126. points += tuple(-x for x in points)
  127. quadrature_points = 15
  128. interval_sets = []
  129. count = 0
  130. def f(x):
  131. nonlocal count
  132. if count % quadrature_points == 0:
  133. interval_sets.append(set())
  134. count += 1
  135. interval_sets[-1].add(float(x))
  136. return 0.0
  137. quad_vec(f, a, b, points=points, quadrature='gk15', limit=0)
  138. # Check that all point sets lie in a single `points` interval
  139. for p in interval_sets:
  140. j = np.searchsorted(sorted(points), tuple(p))
  141. assert np.all(j == j[0])