test_quadpack.py 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418
  1. from __future__ import division, print_function, absolute_import
  2. import sys
  3. import math
  4. import numpy as np
  5. from numpy import sqrt, cos, sin, arctan, exp, log, pi, Inf
  6. from numpy.testing import (assert_,
  7. assert_allclose, assert_array_less, assert_almost_equal)
  8. import pytest
  9. from pytest import raises as assert_raises
  10. from scipy.integrate import quad, dblquad, tplquad, nquad
  11. from scipy._lib.six import xrange
  12. from scipy._lib._ccallback import LowLevelCallable
  13. import ctypes
  14. import ctypes.util
  15. from scipy._lib._ccallback_c import sine_ctypes
  16. import scipy.integrate._test_multivariate as clib_test
  17. def assert_quad(value_and_err, tabled_value, errTol=1.5e-8):
  18. value, err = value_and_err
  19. assert_allclose(value, tabled_value, atol=err, rtol=0)
  20. if errTol is not None:
  21. assert_array_less(err, errTol)
  22. def get_clib_test_routine(name, restype, *argtypes):
  23. ptr = getattr(clib_test, name)
  24. return ctypes.cast(ptr, ctypes.CFUNCTYPE(restype, *argtypes))
  25. class TestCtypesQuad(object):
  26. def setup_method(self):
  27. if sys.platform == 'win32':
  28. if sys.version_info < (3, 5):
  29. files = [ctypes.util.find_msvcrt()]
  30. else:
  31. files = ['api-ms-win-crt-math-l1-1-0.dll']
  32. elif sys.platform == 'darwin':
  33. files = ['libm.dylib']
  34. else:
  35. files = ['libm.so', 'libm.so.6']
  36. for file in files:
  37. try:
  38. self.lib = ctypes.CDLL(file)
  39. break
  40. except OSError:
  41. pass
  42. else:
  43. # This test doesn't work on some Linux platforms (Fedora for
  44. # example) that put an ld script in libm.so - see gh-5370
  45. self.skipTest("Ctypes can't import libm.so")
  46. restype = ctypes.c_double
  47. argtypes = (ctypes.c_double,)
  48. for name in ['sin', 'cos', 'tan']:
  49. func = getattr(self.lib, name)
  50. func.restype = restype
  51. func.argtypes = argtypes
  52. def test_typical(self):
  53. assert_quad(quad(self.lib.sin, 0, 5), quad(math.sin, 0, 5)[0])
  54. assert_quad(quad(self.lib.cos, 0, 5), quad(math.cos, 0, 5)[0])
  55. assert_quad(quad(self.lib.tan, 0, 1), quad(math.tan, 0, 1)[0])
  56. def test_ctypes_sine(self):
  57. quad(LowLevelCallable(sine_ctypes), 0, 1)
  58. def test_ctypes_variants(self):
  59. sin_0 = get_clib_test_routine('_sin_0', ctypes.c_double,
  60. ctypes.c_double, ctypes.c_void_p)
  61. sin_1 = get_clib_test_routine('_sin_1', ctypes.c_double,
  62. ctypes.c_int, ctypes.POINTER(ctypes.c_double),
  63. ctypes.c_void_p)
  64. sin_2 = get_clib_test_routine('_sin_2', ctypes.c_double,
  65. ctypes.c_double)
  66. sin_3 = get_clib_test_routine('_sin_3', ctypes.c_double,
  67. ctypes.c_int, ctypes.POINTER(ctypes.c_double))
  68. sin_4 = get_clib_test_routine('_sin_3', ctypes.c_double,
  69. ctypes.c_int, ctypes.c_double)
  70. all_sigs = [sin_0, sin_1, sin_2, sin_3, sin_4]
  71. legacy_sigs = [sin_2, sin_4]
  72. legacy_only_sigs = [sin_4]
  73. # LowLevelCallables work for new signatures
  74. for j, func in enumerate(all_sigs):
  75. callback = LowLevelCallable(func)
  76. if func in legacy_only_sigs:
  77. assert_raises(ValueError, quad, callback, 0, pi)
  78. else:
  79. assert_allclose(quad(callback, 0, pi)[0], 2.0)
  80. # Plain ctypes items work only for legacy signatures
  81. for j, func in enumerate(legacy_sigs):
  82. if func in legacy_sigs:
  83. assert_allclose(quad(func, 0, pi)[0], 2.0)
  84. else:
  85. assert_raises(ValueError, quad, func, 0, pi)
  86. class TestMultivariateCtypesQuad(object):
  87. def setup_method(self):
  88. restype = ctypes.c_double
  89. argtypes = (ctypes.c_int, ctypes.c_double)
  90. for name in ['_multivariate_typical', '_multivariate_indefinite',
  91. '_multivariate_sin']:
  92. func = get_clib_test_routine(name, restype, *argtypes)
  93. setattr(self, name, func)
  94. def test_typical(self):
  95. # 1) Typical function with two extra arguments:
  96. assert_quad(quad(self._multivariate_typical, 0, pi, (2, 1.8)),
  97. 0.30614353532540296487)
  98. def test_indefinite(self):
  99. # 2) Infinite integration limits --- Euler's constant
  100. assert_quad(quad(self._multivariate_indefinite, 0, Inf),
  101. 0.577215664901532860606512)
  102. def test_threadsafety(self):
  103. # Ensure multivariate ctypes are threadsafe
  104. def threadsafety(y):
  105. return y + quad(self._multivariate_sin, 0, 1)[0]
  106. assert_quad(quad(threadsafety, 0, 1), 0.9596976941318602)
  107. class TestQuad(object):
  108. def test_typical(self):
  109. # 1) Typical function with two extra arguments:
  110. def myfunc(x, n, z): # Bessel function integrand
  111. return cos(n*x-z*sin(x))/pi
  112. assert_quad(quad(myfunc, 0, pi, (2, 1.8)), 0.30614353532540296487)
  113. def test_indefinite(self):
  114. # 2) Infinite integration limits --- Euler's constant
  115. def myfunc(x): # Euler's constant integrand
  116. return -exp(-x)*log(x)
  117. assert_quad(quad(myfunc, 0, Inf), 0.577215664901532860606512)
  118. def test_singular(self):
  119. # 3) Singular points in region of integration.
  120. def myfunc(x):
  121. if 0 < x < 2.5:
  122. return sin(x)
  123. elif 2.5 <= x <= 5.0:
  124. return exp(-x)
  125. else:
  126. return 0.0
  127. assert_quad(quad(myfunc, 0, 10, points=[2.5, 5.0]),
  128. 1 - cos(2.5) + exp(-2.5) - exp(-5.0))
  129. def test_sine_weighted_finite(self):
  130. # 4) Sine weighted integral (finite limits)
  131. def myfunc(x, a):
  132. return exp(a*(x-1))
  133. ome = 2.0**3.4
  134. assert_quad(quad(myfunc, 0, 1, args=20, weight='sin', wvar=ome),
  135. (20*sin(ome)-ome*cos(ome)+ome*exp(-20))/(20**2 + ome**2))
  136. def test_sine_weighted_infinite(self):
  137. # 5) Sine weighted integral (infinite limits)
  138. def myfunc(x, a):
  139. return exp(-x*a)
  140. a = 4.0
  141. ome = 3.0
  142. assert_quad(quad(myfunc, 0, Inf, args=a, weight='sin', wvar=ome),
  143. ome/(a**2 + ome**2))
  144. def test_cosine_weighted_infinite(self):
  145. # 6) Cosine weighted integral (negative infinite limits)
  146. def myfunc(x, a):
  147. return exp(x*a)
  148. a = 2.5
  149. ome = 2.3
  150. assert_quad(quad(myfunc, -Inf, 0, args=a, weight='cos', wvar=ome),
  151. a/(a**2 + ome**2))
  152. def test_algebraic_log_weight(self):
  153. # 6) Algebraic-logarithmic weight.
  154. def myfunc(x, a):
  155. return 1/(1+x+2**(-a))
  156. a = 1.5
  157. assert_quad(quad(myfunc, -1, 1, args=a, weight='alg',
  158. wvar=(-0.5, -0.5)),
  159. pi/sqrt((1+2**(-a))**2 - 1))
  160. def test_cauchypv_weight(self):
  161. # 7) Cauchy prinicpal value weighting w(x) = 1/(x-c)
  162. def myfunc(x, a):
  163. return 2.0**(-a)/((x-1)**2+4.0**(-a))
  164. a = 0.4
  165. tabledValue = ((2.0**(-0.4)*log(1.5) -
  166. 2.0**(-1.4)*log((4.0**(-a)+16) / (4.0**(-a)+1)) -
  167. arctan(2.0**(a+2)) -
  168. arctan(2.0**a)) /
  169. (4.0**(-a) + 1))
  170. assert_quad(quad(myfunc, 0, 5, args=0.4, weight='cauchy', wvar=2.0),
  171. tabledValue, errTol=1.9e-8)
  172. def test_b_less_than_a(self):
  173. def f(x, p, q):
  174. return p * np.exp(-q*x)
  175. val_1, err_1 = quad(f, 0, np.inf, args=(2, 3))
  176. val_2, err_2 = quad(f, np.inf, 0, args=(2, 3))
  177. assert_allclose(val_1, -val_2, atol=max(err_1, err_2))
  178. def test_b_less_than_a_2(self):
  179. def f(x, s):
  180. return np.exp(-x**2 / 2 / s) / np.sqrt(2.*s)
  181. val_1, err_1 = quad(f, -np.inf, np.inf, args=(2,))
  182. val_2, err_2 = quad(f, np.inf, -np.inf, args=(2,))
  183. assert_allclose(val_1, -val_2, atol=max(err_1, err_2))
  184. def test_b_less_than_a_3(self):
  185. def f(x):
  186. return 1.0
  187. val_1, err_1 = quad(f, 0, 1, weight='alg', wvar=(0, 0))
  188. val_2, err_2 = quad(f, 1, 0, weight='alg', wvar=(0, 0))
  189. assert_allclose(val_1, -val_2, atol=max(err_1, err_2))
  190. def test_b_less_than_a_full_output(self):
  191. def f(x):
  192. return 1.0
  193. res_1 = quad(f, 0, 1, weight='alg', wvar=(0, 0), full_output=True)
  194. res_2 = quad(f, 1, 0, weight='alg', wvar=(0, 0), full_output=True)
  195. err = max(res_1[1], res_2[1])
  196. assert_allclose(res_1[0], -res_2[0], atol=err)
  197. def test_double_integral(self):
  198. # 8) Double Integral test
  199. def simpfunc(y, x): # Note order of arguments.
  200. return x+y
  201. a, b = 1.0, 2.0
  202. assert_quad(dblquad(simpfunc, a, b, lambda x: x, lambda x: 2*x),
  203. 5/6.0 * (b**3.0-a**3.0))
  204. def test_double_integral2(self):
  205. def func(x0, x1, t0, t1):
  206. return x0 + x1 + t0 + t1
  207. g = lambda x: x
  208. h = lambda x: 2 * x
  209. args = 1, 2
  210. assert_quad(dblquad(func, 1, 2, g, h, args=args),35./6 + 9*.5)
  211. def test_double_integral3(self):
  212. def func(x0, x1):
  213. return x0 + x1 + 1 + 2
  214. assert_quad(dblquad(func, 1, 2, 1, 2),6.)
  215. def test_triple_integral(self):
  216. # 9) Triple Integral test
  217. def simpfunc(z, y, x, t): # Note order of arguments.
  218. return (x+y+z)*t
  219. a, b = 1.0, 2.0
  220. assert_quad(tplquad(simpfunc, a, b,
  221. lambda x: x, lambda x: 2*x,
  222. lambda x, y: x - y, lambda x, y: x + y,
  223. (2.,)),
  224. 2*8/3.0 * (b**4.0 - a**4.0))
  225. class TestNQuad(object):
  226. def test_fixed_limits(self):
  227. def func1(x0, x1, x2, x3):
  228. val = (x0**2 + x1*x2 - x3**3 + np.sin(x0) +
  229. (1 if (x0 - 0.2*x3 - 0.5 - 0.25*x1 > 0) else 0))
  230. return val
  231. def opts_basic(*args):
  232. return {'points': [0.2*args[2] + 0.5 + 0.25*args[0]]}
  233. res = nquad(func1, [[0, 1], [-1, 1], [.13, .8], [-.15, 1]],
  234. opts=[opts_basic, {}, {}, {}], full_output=True)
  235. assert_quad(res[:-1], 1.5267454070738635)
  236. assert_(res[-1]['neval'] > 0 and res[-1]['neval'] < 4e5)
  237. def test_variable_limits(self):
  238. scale = .1
  239. def func2(x0, x1, x2, x3, t0, t1):
  240. val = (x0*x1*x3**2 + np.sin(x2) + 1 +
  241. (1 if x0 + t1*x1 - t0 > 0 else 0))
  242. return val
  243. def lim0(x1, x2, x3, t0, t1):
  244. return [scale * (x1**2 + x2 + np.cos(x3)*t0*t1 + 1) - 1,
  245. scale * (x1**2 + x2 + np.cos(x3)*t0*t1 + 1) + 1]
  246. def lim1(x2, x3, t0, t1):
  247. return [scale * (t0*x2 + t1*x3) - 1,
  248. scale * (t0*x2 + t1*x3) + 1]
  249. def lim2(x3, t0, t1):
  250. return [scale * (x3 + t0**2*t1**3) - 1,
  251. scale * (x3 + t0**2*t1**3) + 1]
  252. def lim3(t0, t1):
  253. return [scale * (t0 + t1) - 1, scale * (t0 + t1) + 1]
  254. def opts0(x1, x2, x3, t0, t1):
  255. return {'points': [t0 - t1*x1]}
  256. def opts1(x2, x3, t0, t1):
  257. return {}
  258. def opts2(x3, t0, t1):
  259. return {}
  260. def opts3(t0, t1):
  261. return {}
  262. res = nquad(func2, [lim0, lim1, lim2, lim3], args=(0, 0),
  263. opts=[opts0, opts1, opts2, opts3])
  264. assert_quad(res, 25.066666666666663)
  265. def test_square_separate_ranges_and_opts(self):
  266. def f(y, x):
  267. return 1.0
  268. assert_quad(nquad(f, [[-1, 1], [-1, 1]], opts=[{}, {}]), 4.0)
  269. def test_square_aliased_ranges_and_opts(self):
  270. def f(y, x):
  271. return 1.0
  272. r = [-1, 1]
  273. opt = {}
  274. assert_quad(nquad(f, [r, r], opts=[opt, opt]), 4.0)
  275. def test_square_separate_fn_ranges_and_opts(self):
  276. def f(y, x):
  277. return 1.0
  278. def fn_range0(*args):
  279. return (-1, 1)
  280. def fn_range1(*args):
  281. return (-1, 1)
  282. def fn_opt0(*args):
  283. return {}
  284. def fn_opt1(*args):
  285. return {}
  286. ranges = [fn_range0, fn_range1]
  287. opts = [fn_opt0, fn_opt1]
  288. assert_quad(nquad(f, ranges, opts=opts), 4.0)
  289. def test_square_aliased_fn_ranges_and_opts(self):
  290. def f(y, x):
  291. return 1.0
  292. def fn_range(*args):
  293. return (-1, 1)
  294. def fn_opt(*args):
  295. return {}
  296. ranges = [fn_range, fn_range]
  297. opts = [fn_opt, fn_opt]
  298. assert_quad(nquad(f, ranges, opts=opts), 4.0)
  299. def test_matching_quad(self):
  300. def func(x):
  301. return x**2 + 1
  302. res, reserr = quad(func, 0, 4)
  303. res2, reserr2 = nquad(func, ranges=[[0, 4]])
  304. assert_almost_equal(res, res2)
  305. assert_almost_equal(reserr, reserr2)
  306. def test_matching_dblquad(self):
  307. def func2d(x0, x1):
  308. return x0**2 + x1**3 - x0 * x1 + 1
  309. res, reserr = dblquad(func2d, -2, 2, lambda x: -3, lambda x: 3)
  310. res2, reserr2 = nquad(func2d, [[-3, 3], (-2, 2)])
  311. assert_almost_equal(res, res2)
  312. assert_almost_equal(reserr, reserr2)
  313. def test_matching_tplquad(self):
  314. def func3d(x0, x1, x2, c0, c1):
  315. return x0**2 + c0 * x1**3 - x0 * x1 + 1 + c1 * np.sin(x2)
  316. res = tplquad(func3d, -1, 2, lambda x: -2, lambda x: 2,
  317. lambda x, y: -np.pi, lambda x, y: np.pi,
  318. args=(2, 3))
  319. res2 = nquad(func3d, [[-np.pi, np.pi], [-2, 2], (-1, 2)], args=(2, 3))
  320. assert_almost_equal(res, res2)
  321. def test_dict_as_opts(self):
  322. try:
  323. out = nquad(lambda x, y: x * y, [[0, 1], [0, 1]], opts={'epsrel': 0.0001})
  324. except(TypeError):
  325. assert False