test_rgi.py 40 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019
  1. import itertools
  2. import pytest
  3. import numpy as np
  4. from numpy.testing import (assert_allclose, assert_equal, assert_warns,
  5. assert_array_almost_equal, assert_array_equal)
  6. from pytest import raises as assert_raises
  7. from scipy.interpolate import (RegularGridInterpolator, interpn,
  8. RectBivariateSpline,
  9. NearestNDInterpolator, LinearNDInterpolator)
  10. from scipy.sparse._sputils import matrix
  11. parametrize_rgi_interp_methods = pytest.mark.parametrize(
  12. "method", ['linear', 'nearest', 'slinear', 'cubic', 'quintic', 'pchip']
  13. )
  14. class TestRegularGridInterpolator:
  15. def _get_sample_4d(self):
  16. # create a 4-D grid of 3 points in each dimension
  17. points = [(0., .5, 1.)] * 4
  18. values = np.asarray([0., .5, 1.])
  19. values0 = values[:, np.newaxis, np.newaxis, np.newaxis]
  20. values1 = values[np.newaxis, :, np.newaxis, np.newaxis]
  21. values2 = values[np.newaxis, np.newaxis, :, np.newaxis]
  22. values3 = values[np.newaxis, np.newaxis, np.newaxis, :]
  23. values = (values0 + values1 * 10 + values2 * 100 + values3 * 1000)
  24. return points, values
  25. def _get_sample_4d_2(self):
  26. # create another 4-D grid of 3 points in each dimension
  27. points = [(0., .5, 1.)] * 2 + [(0., 5., 10.)] * 2
  28. values = np.asarray([0., .5, 1.])
  29. values0 = values[:, np.newaxis, np.newaxis, np.newaxis]
  30. values1 = values[np.newaxis, :, np.newaxis, np.newaxis]
  31. values2 = values[np.newaxis, np.newaxis, :, np.newaxis]
  32. values3 = values[np.newaxis, np.newaxis, np.newaxis, :]
  33. values = (values0 + values1 * 10 + values2 * 100 + values3 * 1000)
  34. return points, values
  35. def _get_sample_4d_3(self):
  36. # create another 4-D grid of 7 points in each dimension
  37. points = [(0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0)] * 4
  38. values = np.asarray([0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0])
  39. values0 = values[:, np.newaxis, np.newaxis, np.newaxis]
  40. values1 = values[np.newaxis, :, np.newaxis, np.newaxis]
  41. values2 = values[np.newaxis, np.newaxis, :, np.newaxis]
  42. values3 = values[np.newaxis, np.newaxis, np.newaxis, :]
  43. values = (values0 + values1 * 10 + values2 * 100 + values3 * 1000)
  44. return points, values
  45. def _get_sample_4d_4(self):
  46. # create another 4-D grid of 2 points in each dimension
  47. points = [(0.0, 1.0)] * 4
  48. values = np.asarray([0.0, 1.0])
  49. values0 = values[:, np.newaxis, np.newaxis, np.newaxis]
  50. values1 = values[np.newaxis, :, np.newaxis, np.newaxis]
  51. values2 = values[np.newaxis, np.newaxis, :, np.newaxis]
  52. values3 = values[np.newaxis, np.newaxis, np.newaxis, :]
  53. values = (values0 + values1 * 10 + values2 * 100 + values3 * 1000)
  54. return points, values
  55. @parametrize_rgi_interp_methods
  56. def test_list_input(self, method):
  57. points, values = self._get_sample_4d_3()
  58. sample = np.asarray([[0.1, 0.1, 1., .9], [0.2, 0.1, .45, .8],
  59. [0.5, 0.5, .5, .5]])
  60. interp = RegularGridInterpolator(points,
  61. values.tolist(),
  62. method=method)
  63. v1 = interp(sample.tolist())
  64. interp = RegularGridInterpolator(points,
  65. values,
  66. method=method)
  67. v2 = interp(sample)
  68. assert_allclose(v1, v2)
  69. @pytest.mark.parametrize('method', ['cubic', 'quintic', 'pchip'])
  70. def test_spline_dim_error(self, method):
  71. points, values = self._get_sample_4d_4()
  72. match = "points in dimension"
  73. # Check error raise when creating interpolator
  74. with pytest.raises(ValueError, match=match):
  75. RegularGridInterpolator(points, values, method=method)
  76. # Check error raise when creating interpolator
  77. interp = RegularGridInterpolator(points, values)
  78. sample = np.asarray([[0.1, 0.1, 1., .9], [0.2, 0.1, .45, .8],
  79. [0.5, 0.5, .5, .5]])
  80. with pytest.raises(ValueError, match=match):
  81. interp(sample, method=method)
  82. @pytest.mark.parametrize(
  83. "points_values, sample",
  84. [
  85. (
  86. _get_sample_4d,
  87. np.asarray(
  88. [[0.1, 0.1, 1.0, 0.9],
  89. [0.2, 0.1, 0.45, 0.8],
  90. [0.5, 0.5, 0.5, 0.5]]
  91. ),
  92. ),
  93. (_get_sample_4d_2, np.asarray([0.1, 0.1, 10.0, 9.0])),
  94. ],
  95. )
  96. def test_linear_and_slinear_close(self, points_values, sample):
  97. points, values = points_values(self)
  98. interp = RegularGridInterpolator(points, values, method="linear")
  99. v1 = interp(sample)
  100. interp = RegularGridInterpolator(points, values, method="slinear")
  101. v2 = interp(sample)
  102. assert_allclose(v1, v2)
  103. @parametrize_rgi_interp_methods
  104. def test_complex(self, method):
  105. points, values = self._get_sample_4d_3()
  106. values = values - 2j*values
  107. sample = np.asarray([[0.1, 0.1, 1., .9], [0.2, 0.1, .45, .8],
  108. [0.5, 0.5, .5, .5]])
  109. interp = RegularGridInterpolator(points, values, method=method)
  110. rinterp = RegularGridInterpolator(points, values.real, method=method)
  111. iinterp = RegularGridInterpolator(points, values.imag, method=method)
  112. v1 = interp(sample)
  113. v2 = rinterp(sample) + 1j*iinterp(sample)
  114. assert_allclose(v1, v2)
  115. def test_cubic_vs_pchip(self):
  116. x, y = [1, 2, 3, 4], [1, 2, 3, 4]
  117. xg, yg = np.meshgrid(x, y, indexing='ij')
  118. values = (lambda x, y: x**4 * y**4)(xg, yg)
  119. cubic = RegularGridInterpolator((x, y), values, method='cubic')
  120. pchip = RegularGridInterpolator((x, y), values, method='pchip')
  121. vals_cubic = cubic([1.5, 2])
  122. vals_pchip = pchip([1.5, 2])
  123. assert not np.allclose(vals_cubic, vals_pchip, atol=1e-14, rtol=0)
  124. def test_linear_xi1d(self):
  125. points, values = self._get_sample_4d_2()
  126. interp = RegularGridInterpolator(points, values)
  127. sample = np.asarray([0.1, 0.1, 10., 9.])
  128. wanted = 1001.1
  129. assert_array_almost_equal(interp(sample), wanted)
  130. def test_linear_xi3d(self):
  131. points, values = self._get_sample_4d()
  132. interp = RegularGridInterpolator(points, values)
  133. sample = np.asarray([[0.1, 0.1, 1., .9], [0.2, 0.1, .45, .8],
  134. [0.5, 0.5, .5, .5]])
  135. wanted = np.asarray([1001.1, 846.2, 555.5])
  136. assert_array_almost_equal(interp(sample), wanted)
  137. @pytest.mark.parametrize(
  138. "sample, wanted",
  139. [
  140. (np.asarray([0.1, 0.1, 0.9, 0.9]), 1100.0),
  141. (np.asarray([0.1, 0.1, 0.1, 0.1]), 0.0),
  142. (np.asarray([0.0, 0.0, 0.0, 0.0]), 0.0),
  143. (np.asarray([1.0, 1.0, 1.0, 1.0]), 1111.0),
  144. (np.asarray([0.1, 0.4, 0.6, 0.9]), 1055.0),
  145. ],
  146. )
  147. def test_nearest(self, sample, wanted):
  148. points, values = self._get_sample_4d()
  149. interp = RegularGridInterpolator(points, values, method="nearest")
  150. assert_array_almost_equal(interp(sample), wanted)
  151. def test_linear_edges(self):
  152. points, values = self._get_sample_4d()
  153. interp = RegularGridInterpolator(points, values)
  154. sample = np.asarray([[0., 0., 0., 0.], [1., 1., 1., 1.]])
  155. wanted = np.asarray([0., 1111.])
  156. assert_array_almost_equal(interp(sample), wanted)
  157. def test_valid_create(self):
  158. # create a 2-D grid of 3 points in each dimension
  159. points = [(0., .5, 1.), (0., 1., .5)]
  160. values = np.asarray([0., .5, 1.])
  161. values0 = values[:, np.newaxis]
  162. values1 = values[np.newaxis, :]
  163. values = (values0 + values1 * 10)
  164. assert_raises(ValueError, RegularGridInterpolator, points, values)
  165. points = [((0., .5, 1.), ), (0., .5, 1.)]
  166. assert_raises(ValueError, RegularGridInterpolator, points, values)
  167. points = [(0., .5, .75, 1.), (0., .5, 1.)]
  168. assert_raises(ValueError, RegularGridInterpolator, points, values)
  169. points = [(0., .5, 1.), (0., .5, 1.), (0., .5, 1.)]
  170. assert_raises(ValueError, RegularGridInterpolator, points, values)
  171. points = [(0., .5, 1.), (0., .5, 1.)]
  172. assert_raises(ValueError, RegularGridInterpolator, points, values,
  173. method="undefmethod")
  174. def test_valid_call(self):
  175. points, values = self._get_sample_4d()
  176. interp = RegularGridInterpolator(points, values)
  177. sample = np.asarray([[0., 0., 0., 0.], [1., 1., 1., 1.]])
  178. assert_raises(ValueError, interp, sample, "undefmethod")
  179. sample = np.asarray([[0., 0., 0.], [1., 1., 1.]])
  180. assert_raises(ValueError, interp, sample)
  181. sample = np.asarray([[0., 0., 0., 0.], [1., 1., 1., 1.1]])
  182. assert_raises(ValueError, interp, sample)
  183. def test_out_of_bounds_extrap(self):
  184. points, values = self._get_sample_4d()
  185. interp = RegularGridInterpolator(points, values, bounds_error=False,
  186. fill_value=None)
  187. sample = np.asarray([[-.1, -.1, -.1, -.1], [1.1, 1.1, 1.1, 1.1],
  188. [21, 2.1, -1.1, -11], [2.1, 2.1, -1.1, -1.1]])
  189. wanted = np.asarray([0., 1111., 11., 11.])
  190. assert_array_almost_equal(interp(sample, method="nearest"), wanted)
  191. wanted = np.asarray([-111.1, 1222.1, -11068., -1186.9])
  192. assert_array_almost_equal(interp(sample, method="linear"), wanted)
  193. def test_out_of_bounds_extrap2(self):
  194. points, values = self._get_sample_4d_2()
  195. interp = RegularGridInterpolator(points, values, bounds_error=False,
  196. fill_value=None)
  197. sample = np.asarray([[-.1, -.1, -.1, -.1], [1.1, 1.1, 1.1, 1.1],
  198. [21, 2.1, -1.1, -11], [2.1, 2.1, -1.1, -1.1]])
  199. wanted = np.asarray([0., 11., 11., 11.])
  200. assert_array_almost_equal(interp(sample, method="nearest"), wanted)
  201. wanted = np.asarray([-12.1, 133.1, -1069., -97.9])
  202. assert_array_almost_equal(interp(sample, method="linear"), wanted)
  203. def test_out_of_bounds_fill(self):
  204. points, values = self._get_sample_4d()
  205. interp = RegularGridInterpolator(points, values, bounds_error=False,
  206. fill_value=np.nan)
  207. sample = np.asarray([[-.1, -.1, -.1, -.1], [1.1, 1.1, 1.1, 1.1],
  208. [2.1, 2.1, -1.1, -1.1]])
  209. wanted = np.asarray([np.nan, np.nan, np.nan])
  210. assert_array_almost_equal(interp(sample, method="nearest"), wanted)
  211. assert_array_almost_equal(interp(sample, method="linear"), wanted)
  212. sample = np.asarray([[0.1, 0.1, 1., .9], [0.2, 0.1, .45, .8],
  213. [0.5, 0.5, .5, .5]])
  214. wanted = np.asarray([1001.1, 846.2, 555.5])
  215. assert_array_almost_equal(interp(sample), wanted)
  216. def test_nearest_compare_qhull(self):
  217. points, values = self._get_sample_4d()
  218. interp = RegularGridInterpolator(points, values, method="nearest")
  219. points_qhull = itertools.product(*points)
  220. points_qhull = [p for p in points_qhull]
  221. points_qhull = np.asarray(points_qhull)
  222. values_qhull = values.reshape(-1)
  223. interp_qhull = NearestNDInterpolator(points_qhull, values_qhull)
  224. sample = np.asarray([[0.1, 0.1, 1., .9], [0.2, 0.1, .45, .8],
  225. [0.5, 0.5, .5, .5]])
  226. assert_array_almost_equal(interp(sample), interp_qhull(sample))
  227. def test_linear_compare_qhull(self):
  228. points, values = self._get_sample_4d()
  229. interp = RegularGridInterpolator(points, values)
  230. points_qhull = itertools.product(*points)
  231. points_qhull = [p for p in points_qhull]
  232. points_qhull = np.asarray(points_qhull)
  233. values_qhull = values.reshape(-1)
  234. interp_qhull = LinearNDInterpolator(points_qhull, values_qhull)
  235. sample = np.asarray([[0.1, 0.1, 1., .9], [0.2, 0.1, .45, .8],
  236. [0.5, 0.5, .5, .5]])
  237. assert_array_almost_equal(interp(sample), interp_qhull(sample))
  238. @pytest.mark.parametrize("method", ["nearest", "linear"])
  239. def test_duck_typed_values(self, method):
  240. x = np.linspace(0, 2, 5)
  241. y = np.linspace(0, 1, 7)
  242. values = MyValue((5, 7))
  243. interp = RegularGridInterpolator((x, y), values, method=method)
  244. v1 = interp([0.4, 0.7])
  245. interp = RegularGridInterpolator((x, y), values._v, method=method)
  246. v2 = interp([0.4, 0.7])
  247. assert_allclose(v1, v2)
  248. def test_invalid_fill_value(self):
  249. np.random.seed(1234)
  250. x = np.linspace(0, 2, 5)
  251. y = np.linspace(0, 1, 7)
  252. values = np.random.rand(5, 7)
  253. # integers can be cast to floats
  254. RegularGridInterpolator((x, y), values, fill_value=1)
  255. # complex values cannot
  256. assert_raises(ValueError, RegularGridInterpolator,
  257. (x, y), values, fill_value=1+2j)
  258. def test_fillvalue_type(self):
  259. # from #3703; test that interpolator object construction succeeds
  260. values = np.ones((10, 20, 30), dtype='>f4')
  261. points = [np.arange(n) for n in values.shape]
  262. # xi = [(1, 1, 1)]
  263. RegularGridInterpolator(points, values)
  264. RegularGridInterpolator(points, values, fill_value=0.)
  265. def test_length_one_axis(self):
  266. # gh-5890, gh-9524 : length-1 axis is legal for method='linear'.
  267. # Along the axis it's linear interpolation; away from the length-1
  268. # axis, it's an extrapolation, so fill_value should be used.
  269. def f(x, y):
  270. return x + y
  271. x = np.linspace(1, 1, 1)
  272. y = np.linspace(1, 10, 10)
  273. data = f(*np.meshgrid(x, y, indexing="ij", sparse=True))
  274. interp = RegularGridInterpolator((x, y), data, method="linear",
  275. bounds_error=False, fill_value=101)
  276. # check values at the grid
  277. assert_allclose(interp(np.array([[1, 1], [1, 5], [1, 10]])),
  278. [2, 6, 11],
  279. atol=1e-14)
  280. # check off-grid interpolation is indeed linear
  281. assert_allclose(interp(np.array([[1, 1.4], [1, 5.3], [1, 10]])),
  282. [2.4, 6.3, 11],
  283. atol=1e-14)
  284. # check exrapolation w/ fill_value
  285. assert_allclose(interp(np.array([1.1, 2.4])),
  286. interp.fill_value,
  287. atol=1e-14)
  288. # check extrapolation: linear along the `y` axis, const along `x`
  289. interp.fill_value = None
  290. assert_allclose(interp([[1, 0.3], [1, 11.5]]),
  291. [1.3, 12.5], atol=1e-15)
  292. assert_allclose(interp([[1.5, 0.3], [1.9, 11.5]]),
  293. [1.3, 12.5], atol=1e-15)
  294. # extrapolation with method='nearest'
  295. interp = RegularGridInterpolator((x, y), data, method="nearest",
  296. bounds_error=False, fill_value=None)
  297. assert_allclose(interp([[1.5, 1.8], [-4, 5.1]]),
  298. [3, 6],
  299. atol=1e-15)
  300. @pytest.mark.parametrize("fill_value", [None, np.nan, np.pi])
  301. @pytest.mark.parametrize("method", ['linear', 'nearest'])
  302. def test_length_one_axis2(self, fill_value, method):
  303. options = {"fill_value": fill_value, "bounds_error": False,
  304. "method": method}
  305. x = np.linspace(0, 2*np.pi, 20)
  306. z = np.sin(x)
  307. fa = RegularGridInterpolator((x,), z[:], **options)
  308. fb = RegularGridInterpolator((x, [0]), z[:, None], **options)
  309. x1a = np.linspace(-1, 2*np.pi+1, 100)
  310. za = fa(x1a)
  311. # evaluated at provided y-value, fb should behave exactly as fa
  312. y1b = np.zeros(100)
  313. zb = fb(np.vstack([x1a, y1b]).T)
  314. assert_allclose(zb, za)
  315. # evaluated at a different y-value, fb should return fill value
  316. y1b = np.ones(100)
  317. zb = fb(np.vstack([x1a, y1b]).T)
  318. if fill_value is None:
  319. assert_allclose(zb, za)
  320. else:
  321. assert_allclose(zb, fill_value)
  322. @pytest.mark.parametrize("method", ['nearest', 'linear'])
  323. def test_nan_x_1d(self, method):
  324. # gh-6624 : if x is nan, result should be nan
  325. f = RegularGridInterpolator(([1, 2, 3],), [10, 20, 30], fill_value=1,
  326. bounds_error=False, method=method)
  327. assert np.isnan(f([np.nan]))
  328. # test arbitrary nan pattern
  329. rng = np.random.default_rng(8143215468)
  330. x = rng.random(size=100)*4
  331. i = rng.random(size=100) > 0.5
  332. x[i] = np.nan
  333. with np.errstate(invalid='ignore'):
  334. # out-of-bounds comparisons, `out_of_bounds += x < grid[0]`,
  335. # generate numpy warnings if `x` contains nans.
  336. # These warnings should propagate to user (since `x` is user
  337. # input) and we simply filter them out.
  338. res = f(x)
  339. assert_equal(res[i], np.nan)
  340. assert_equal(res[~i], f(x[~i]))
  341. # also test the length-one axis f(nan)
  342. x = [1, 2, 3]
  343. y = [1, ]
  344. data = np.ones((3, 1))
  345. f = RegularGridInterpolator((x, y), data, fill_value=1,
  346. bounds_error=False, method=method)
  347. assert np.isnan(f([np.nan, 1]))
  348. assert np.isnan(f([1, np.nan]))
  349. @pytest.mark.parametrize("method", ['nearest', 'linear'])
  350. def test_nan_x_2d(self, method):
  351. x, y = np.array([0, 1, 2]), np.array([1, 3, 7])
  352. def f(x, y):
  353. return x**2 + y**2
  354. xg, yg = np.meshgrid(x, y, indexing='ij', sparse=True)
  355. data = f(xg, yg)
  356. interp = RegularGridInterpolator((x, y), data,
  357. method=method, bounds_error=False)
  358. with np.errstate(invalid='ignore'):
  359. res = interp([[1.5, np.nan], [1, 1]])
  360. assert_allclose(res[1], 2, atol=1e-14)
  361. assert np.isnan(res[0])
  362. # test arbitrary nan pattern
  363. rng = np.random.default_rng(8143215468)
  364. x = rng.random(size=100)*4-1
  365. y = rng.random(size=100)*8
  366. i1 = rng.random(size=100) > 0.5
  367. i2 = rng.random(size=100) > 0.5
  368. i = i1 | i2
  369. x[i1] = np.nan
  370. y[i2] = np.nan
  371. z = np.array([x, y]).T
  372. with np.errstate(invalid='ignore'):
  373. # out-of-bounds comparisons, `out_of_bounds += x < grid[0]`,
  374. # generate numpy warnings if `x` contains nans.
  375. # These warnings should propagate to user (since `x` is user
  376. # input) and we simply filter them out.
  377. res = interp(z)
  378. assert_equal(res[i], np.nan)
  379. assert_equal(res[~i], interp(z[~i]))
  380. @parametrize_rgi_interp_methods
  381. def test_descending_points(self, method):
  382. def val_func_3d(x, y, z):
  383. return 2 * x ** 3 + 3 * y ** 2 - z
  384. x = np.linspace(1, 4, 11)
  385. y = np.linspace(4, 7, 22)
  386. z = np.linspace(7, 9, 33)
  387. points = (x, y, z)
  388. values = val_func_3d(
  389. *np.meshgrid(*points, indexing='ij', sparse=True))
  390. my_interpolating_function = RegularGridInterpolator(points,
  391. values,
  392. method=method)
  393. pts = np.array([[2.1, 6.2, 8.3], [3.3, 5.2, 7.1]])
  394. correct_result = my_interpolating_function(pts)
  395. # descending data
  396. x_descending = x[::-1]
  397. y_descending = y[::-1]
  398. z_descending = z[::-1]
  399. points_shuffled = (x_descending, y_descending, z_descending)
  400. values_shuffled = val_func_3d(
  401. *np.meshgrid(*points_shuffled, indexing='ij', sparse=True))
  402. my_interpolating_function = RegularGridInterpolator(
  403. points_shuffled, values_shuffled, method=method)
  404. test_result = my_interpolating_function(pts)
  405. assert_array_equal(correct_result, test_result)
  406. def test_invalid_points_order(self):
  407. def val_func_2d(x, y):
  408. return 2 * x ** 3 + 3 * y ** 2
  409. x = np.array([.5, 2., 0., 4., 5.5]) # not ascending or descending
  410. y = np.array([.5, 2., 3., 4., 5.5])
  411. points = (x, y)
  412. values = val_func_2d(*np.meshgrid(*points, indexing='ij',
  413. sparse=True))
  414. match = "must be strictly ascending or descending"
  415. with pytest.raises(ValueError, match=match):
  416. RegularGridInterpolator(points, values)
  417. @parametrize_rgi_interp_methods
  418. def test_fill_value(self, method):
  419. interp = RegularGridInterpolator([np.arange(6)], np.ones(6),
  420. method=method, bounds_error=False)
  421. assert np.isnan(interp([10]))
  422. @parametrize_rgi_interp_methods
  423. def test_nonscalar_values(self, method):
  424. # Verify that non-scalar valued values also works
  425. points = [(0.0, 0.5, 1.0, 1.5, 2.0, 2.5)] * 2 + [
  426. (0.0, 5.0, 10.0, 15.0, 20, 25.0)
  427. ] * 2
  428. rng = np.random.default_rng(1234)
  429. values = rng.random((6, 6, 6, 6, 8))
  430. sample = rng.random((7, 3, 4))
  431. interp = RegularGridInterpolator(points, values, method=method,
  432. bounds_error=False)
  433. v = interp(sample)
  434. assert_equal(v.shape, (7, 3, 8), err_msg=method)
  435. vs = []
  436. for j in range(8):
  437. interp = RegularGridInterpolator(points, values[..., j],
  438. method=method,
  439. bounds_error=False)
  440. vs.append(interp(sample))
  441. v2 = np.array(vs).transpose(1, 2, 0)
  442. assert_allclose(v, v2, atol=1e-14, err_msg=method)
  443. @parametrize_rgi_interp_methods
  444. @pytest.mark.parametrize("flip_points", [False, True])
  445. def test_nonscalar_values_2(self, method, flip_points):
  446. # Verify that non-scalar valued values also work : use different
  447. # lengths of axes to simplify tracing the internals
  448. points = [(0.0, 0.5, 1.0, 1.5, 2.0, 2.5),
  449. (0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0),
  450. (0.0, 5.0, 10.0, 15.0, 20, 25.0, 35.0, 36.0),
  451. (0.0, 5.0, 10.0, 15.0, 20, 25.0, 35.0, 36.0, 47)]
  452. # verify, that strictly decreasing dimensions work
  453. if flip_points:
  454. points = [tuple(reversed(p)) for p in points]
  455. rng = np.random.default_rng(1234)
  456. trailing_points = (3, 2)
  457. # NB: values has a `num_trailing_dims` trailing dimension
  458. values = rng.random((6, 7, 8, 9, *trailing_points))
  459. sample = rng.random(4) # a single sample point !
  460. interp = RegularGridInterpolator(points, values, method=method,
  461. bounds_error=False)
  462. v = interp(sample)
  463. # v has a single sample point *per entry in the trailing dimensions*
  464. assert v.shape == (1, *trailing_points)
  465. # check the values, too : manually loop over the trailing dimensions
  466. vs = np.empty((values.shape[-2:]))
  467. for i in range(values.shape[-2]):
  468. for j in range(values.shape[-1]):
  469. interp = RegularGridInterpolator(points, values[..., i, j],
  470. method=method,
  471. bounds_error=False)
  472. vs[i, j] = interp(sample)
  473. v2 = np.expand_dims(vs, axis=0)
  474. assert_allclose(v, v2, atol=1e-14, err_msg=method)
  475. def test_nonscalar_values_linear_2D(self):
  476. # Verify that non-scalar values work in the 2D fast path
  477. method = 'linear'
  478. points = [(0.0, 0.5, 1.0, 1.5, 2.0, 2.5),
  479. (0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0), ]
  480. rng = np.random.default_rng(1234)
  481. trailing_points = (3, 4)
  482. # NB: values has a `num_trailing_dims` trailing dimension
  483. values = rng.random((6, 7, *trailing_points))
  484. sample = rng.random(2) # a single sample point !
  485. interp = RegularGridInterpolator(points, values, method=method,
  486. bounds_error=False)
  487. v = interp(sample)
  488. # v has a single sample point *per entry in the trailing dimensions*
  489. assert v.shape == (1, *trailing_points)
  490. # check the values, too : manually loop over the trailing dimensions
  491. vs = np.empty((values.shape[-2:]))
  492. for i in range(values.shape[-2]):
  493. for j in range(values.shape[-1]):
  494. interp = RegularGridInterpolator(points, values[..., i, j],
  495. method=method,
  496. bounds_error=False)
  497. vs[i, j] = interp(sample)
  498. v2 = np.expand_dims(vs, axis=0)
  499. assert_allclose(v, v2, atol=1e-14, err_msg=method)
  500. @pytest.mark.parametrize(
  501. "dtype",
  502. [np.float32, np.float64, np.complex64, np.complex128]
  503. )
  504. @pytest.mark.parametrize("xi_dtype", [np.float32, np.float64])
  505. def test_float32_values(self, dtype, xi_dtype):
  506. # regression test for gh-17718: values.dtype=float32 fails
  507. def f(x, y):
  508. return 2 * x**3 + 3 * y**2
  509. x = np.linspace(1, 4, 11)
  510. y = np.linspace(4, 7, 22)
  511. xg, yg = np.meshgrid(x, y, indexing='ij', sparse=True)
  512. data = f(xg, yg)
  513. data = data.astype(dtype)
  514. interp = RegularGridInterpolator((x, y), data)
  515. pts = np.array([[2.1, 6.2],
  516. [3.3, 5.2]], dtype=xi_dtype)
  517. # the values here are just what the call returns; the test checks that
  518. # that the call succeeds at all, instead of failing with cython not
  519. # having a float32 kernel
  520. assert_allclose(interp(pts), [134.10469388, 153.40069388], atol=1e-7)
  521. class MyValue:
  522. """
  523. Minimal indexable object
  524. """
  525. def __init__(self, shape):
  526. self.ndim = 2
  527. self.shape = shape
  528. self._v = np.arange(np.prod(shape)).reshape(shape)
  529. def __getitem__(self, idx):
  530. return self._v[idx]
  531. def __array_interface__(self):
  532. return None
  533. def __array__(self):
  534. raise RuntimeError("No array representation")
  535. class TestInterpN:
  536. def _sample_2d_data(self):
  537. x = np.array([.5, 2., 3., 4., 5.5, 6.])
  538. y = np.array([.5, 2., 3., 4., 5.5, 6.])
  539. z = np.array(
  540. [
  541. [1, 2, 1, 2, 1, 1],
  542. [1, 2, 1, 2, 1, 1],
  543. [1, 2, 3, 2, 1, 1],
  544. [1, 2, 2, 2, 1, 1],
  545. [1, 2, 1, 2, 1, 1],
  546. [1, 2, 2, 2, 1, 1],
  547. ]
  548. )
  549. return x, y, z
  550. def test_spline_2d(self):
  551. x, y, z = self._sample_2d_data()
  552. lut = RectBivariateSpline(x, y, z)
  553. xi = np.array([[1, 2.3, 5.3, 0.5, 3.3, 1.2, 3],
  554. [1, 3.3, 1.2, 4.0, 5.0, 1.0, 3]]).T
  555. assert_array_almost_equal(interpn((x, y), z, xi, method="splinef2d"),
  556. lut.ev(xi[:, 0], xi[:, 1]))
  557. @parametrize_rgi_interp_methods
  558. def test_list_input(self, method):
  559. x, y, z = self._sample_2d_data()
  560. xi = np.array([[1, 2.3, 5.3, 0.5, 3.3, 1.2, 3],
  561. [1, 3.3, 1.2, 4.0, 5.0, 1.0, 3]]).T
  562. v1 = interpn((x, y), z, xi, method=method)
  563. v2 = interpn(
  564. (x.tolist(), y.tolist()), z.tolist(), xi.tolist(), method=method
  565. )
  566. assert_allclose(v1, v2, err_msg=method)
  567. def test_spline_2d_outofbounds(self):
  568. x = np.array([.5, 2., 3., 4., 5.5])
  569. y = np.array([.5, 2., 3., 4., 5.5])
  570. z = np.array([[1, 2, 1, 2, 1], [1, 2, 1, 2, 1], [1, 2, 3, 2, 1],
  571. [1, 2, 2, 2, 1], [1, 2, 1, 2, 1]])
  572. lut = RectBivariateSpline(x, y, z)
  573. xi = np.array([[1, 2.3, 6.3, 0.5, 3.3, 1.2, 3],
  574. [1, 3.3, 1.2, -4.0, 5.0, 1.0, 3]]).T
  575. actual = interpn((x, y), z, xi, method="splinef2d",
  576. bounds_error=False, fill_value=999.99)
  577. expected = lut.ev(xi[:, 0], xi[:, 1])
  578. expected[2:4] = 999.99
  579. assert_array_almost_equal(actual, expected)
  580. # no extrapolation for splinef2d
  581. assert_raises(ValueError, interpn, (x, y), z, xi, method="splinef2d",
  582. bounds_error=False, fill_value=None)
  583. def _sample_4d_data(self):
  584. points = [(0., .5, 1.)] * 2 + [(0., 5., 10.)] * 2
  585. values = np.asarray([0., .5, 1.])
  586. values0 = values[:, np.newaxis, np.newaxis, np.newaxis]
  587. values1 = values[np.newaxis, :, np.newaxis, np.newaxis]
  588. values2 = values[np.newaxis, np.newaxis, :, np.newaxis]
  589. values3 = values[np.newaxis, np.newaxis, np.newaxis, :]
  590. values = (values0 + values1 * 10 + values2 * 100 + values3 * 1000)
  591. return points, values
  592. def test_linear_4d(self):
  593. # create a 4-D grid of 3 points in each dimension
  594. points, values = self._sample_4d_data()
  595. interp_rg = RegularGridInterpolator(points, values)
  596. sample = np.asarray([[0.1, 0.1, 10., 9.]])
  597. wanted = interpn(points, values, sample, method="linear")
  598. assert_array_almost_equal(interp_rg(sample), wanted)
  599. def test_4d_linear_outofbounds(self):
  600. # create a 4-D grid of 3 points in each dimension
  601. points, values = self._sample_4d_data()
  602. sample = np.asarray([[0.1, -0.1, 10.1, 9.]])
  603. wanted = 999.99
  604. actual = interpn(points, values, sample, method="linear",
  605. bounds_error=False, fill_value=999.99)
  606. assert_array_almost_equal(actual, wanted)
  607. def test_nearest_4d(self):
  608. # create a 4-D grid of 3 points in each dimension
  609. points, values = self._sample_4d_data()
  610. interp_rg = RegularGridInterpolator(points, values, method="nearest")
  611. sample = np.asarray([[0.1, 0.1, 10., 9.]])
  612. wanted = interpn(points, values, sample, method="nearest")
  613. assert_array_almost_equal(interp_rg(sample), wanted)
  614. def test_4d_nearest_outofbounds(self):
  615. # create a 4-D grid of 3 points in each dimension
  616. points, values = self._sample_4d_data()
  617. sample = np.asarray([[0.1, -0.1, 10.1, 9.]])
  618. wanted = 999.99
  619. actual = interpn(points, values, sample, method="nearest",
  620. bounds_error=False, fill_value=999.99)
  621. assert_array_almost_equal(actual, wanted)
  622. def test_xi_1d(self):
  623. # verify that 1-D xi works as expected
  624. points, values = self._sample_4d_data()
  625. sample = np.asarray([0.1, 0.1, 10., 9.])
  626. v1 = interpn(points, values, sample, bounds_error=False)
  627. v2 = interpn(points, values, sample[None,:], bounds_error=False)
  628. assert_allclose(v1, v2)
  629. def test_xi_nd(self):
  630. # verify that higher-d xi works as expected
  631. points, values = self._sample_4d_data()
  632. np.random.seed(1234)
  633. sample = np.random.rand(2, 3, 4)
  634. v1 = interpn(points, values, sample, method='nearest',
  635. bounds_error=False)
  636. assert_equal(v1.shape, (2, 3))
  637. v2 = interpn(points, values, sample.reshape(-1, 4),
  638. method='nearest', bounds_error=False)
  639. assert_allclose(v1, v2.reshape(v1.shape))
  640. @parametrize_rgi_interp_methods
  641. def test_xi_broadcast(self, method):
  642. # verify that the interpolators broadcast xi
  643. x, y, values = self._sample_2d_data()
  644. points = (x, y)
  645. xi = np.linspace(0, 1, 2)
  646. yi = np.linspace(0, 3, 3)
  647. sample = (xi[:, None], yi[None, :])
  648. v1 = interpn(points, values, sample, method=method, bounds_error=False)
  649. assert_equal(v1.shape, (2, 3))
  650. xx, yy = np.meshgrid(xi, yi)
  651. sample = np.c_[xx.T.ravel(), yy.T.ravel()]
  652. v2 = interpn(points, values, sample,
  653. method=method, bounds_error=False)
  654. assert_allclose(v1, v2.reshape(v1.shape))
  655. @parametrize_rgi_interp_methods
  656. def test_nonscalar_values(self, method):
  657. # Verify that non-scalar valued values also works
  658. points = [(0.0, 0.5, 1.0, 1.5, 2.0, 2.5)] * 2 + [
  659. (0.0, 5.0, 10.0, 15.0, 20, 25.0)
  660. ] * 2
  661. rng = np.random.default_rng(1234)
  662. values = rng.random((6, 6, 6, 6, 8))
  663. sample = rng.random((7, 3, 4))
  664. v = interpn(points, values, sample, method=method,
  665. bounds_error=False)
  666. assert_equal(v.shape, (7, 3, 8), err_msg=method)
  667. vs = [interpn(points, values[..., j], sample, method=method,
  668. bounds_error=False) for j in range(8)]
  669. v2 = np.array(vs).transpose(1, 2, 0)
  670. assert_allclose(v, v2, atol=1e-14, err_msg=method)
  671. @parametrize_rgi_interp_methods
  672. def test_nonscalar_values_2(self, method):
  673. # Verify that non-scalar valued values also work : use different
  674. # lengths of axes to simplify tracing the internals
  675. points = [(0.0, 0.5, 1.0, 1.5, 2.0, 2.5),
  676. (0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0),
  677. (0.0, 5.0, 10.0, 15.0, 20, 25.0, 35.0, 36.0),
  678. (0.0, 5.0, 10.0, 15.0, 20, 25.0, 35.0, 36.0, 47)]
  679. rng = np.random.default_rng(1234)
  680. trailing_points = (3, 2)
  681. # NB: values has a `num_trailing_dims` trailing dimension
  682. values = rng.random((6, 7, 8, 9, *trailing_points))
  683. sample = rng.random(4) # a single sample point !
  684. v = interpn(points, values, sample, method=method, bounds_error=False)
  685. # v has a single sample point *per entry in the trailing dimensions*
  686. assert v.shape == (1, *trailing_points)
  687. # check the values, too : manually loop over the trailing dimensions
  688. vs = [[
  689. interpn(points, values[..., i, j], sample, method=method,
  690. bounds_error=False) for i in range(values.shape[-2])
  691. ] for j in range(values.shape[-1])]
  692. assert_allclose(v, np.asarray(vs).T, atol=1e-14, err_msg=method)
  693. def test_non_scalar_values_splinef2d(self):
  694. # Vector-valued splines supported with fitpack
  695. points, values = self._sample_4d_data()
  696. np.random.seed(1234)
  697. values = np.random.rand(3, 3, 3, 3, 6)
  698. sample = np.random.rand(7, 11, 4)
  699. assert_raises(ValueError, interpn, points, values, sample,
  700. method='splinef2d')
  701. @parametrize_rgi_interp_methods
  702. def test_complex(self, method):
  703. x, y, values = self._sample_2d_data()
  704. points = (x, y)
  705. values = values - 2j*values
  706. sample = np.array([[1, 2.3, 5.3, 0.5, 3.3, 1.2, 3],
  707. [1, 3.3, 1.2, 4.0, 5.0, 1.0, 3]]).T
  708. v1 = interpn(points, values, sample, method=method)
  709. v2r = interpn(points, values.real, sample, method=method)
  710. v2i = interpn(points, values.imag, sample, method=method)
  711. v2 = v2r + 1j*v2i
  712. assert_allclose(v1, v2)
  713. def test_complex_spline2fd(self):
  714. # Complex-valued data not supported by spline2fd
  715. x, y, values = self._sample_2d_data()
  716. points = (x, y)
  717. values = values - 2j*values
  718. sample = np.array([[1, 2.3, 5.3, 0.5, 3.3, 1.2, 3],
  719. [1, 3.3, 1.2, 4.0, 5.0, 1.0, 3]]).T
  720. with assert_warns(np.ComplexWarning):
  721. interpn(points, values, sample, method='splinef2d')
  722. @pytest.mark.parametrize(
  723. "method",
  724. ["linear", "nearest"]
  725. )
  726. def test_duck_typed_values(self, method):
  727. x = np.linspace(0, 2, 5)
  728. y = np.linspace(0, 1, 7)
  729. values = MyValue((5, 7))
  730. v1 = interpn((x, y), values, [0.4, 0.7], method=method)
  731. v2 = interpn((x, y), values._v, [0.4, 0.7], method=method)
  732. assert_allclose(v1, v2)
  733. @parametrize_rgi_interp_methods
  734. def test_matrix_input(self, method):
  735. x = np.linspace(0, 2, 6)
  736. y = np.linspace(0, 1, 7)
  737. values = matrix(np.random.rand(6, 7))
  738. sample = np.random.rand(3, 7, 2)
  739. v1 = interpn((x, y), values, sample, method=method)
  740. v2 = interpn((x, y), np.asarray(values), sample, method=method)
  741. assert_allclose(v1, v2)
  742. def test_length_one_axis(self):
  743. # gh-5890, gh-9524 : length-1 axis is legal for method='linear'.
  744. # Along the axis it's linear interpolation; away from the length-1
  745. # axis, it's an extrapolation, so fill_value should be used.
  746. values = np.array([[0.1, 1, 10]])
  747. xi = np.array([[1, 2.2], [1, 3.2], [1, 3.8]])
  748. res = interpn(([1], [2, 3, 4]), values, xi)
  749. wanted = [0.9*0.2 + 0.1, # on [2, 3) it's 0.9*(x-2) + 0.1
  750. 9*0.2 + 1, # on [3, 4] it's 9*(x-3) + 1
  751. 9*0.8 + 1]
  752. assert_allclose(res, wanted, atol=1e-15)
  753. # check extrapolation
  754. xi = np.array([[1.1, 2.2], [1.5, 3.2], [-2.3, 3.8]])
  755. res = interpn(([1], [2, 3, 4]), values, xi,
  756. bounds_error=False, fill_value=None)
  757. assert_allclose(res, wanted, atol=1e-15)
  758. def test_descending_points(self):
  759. def value_func_4d(x, y, z, a):
  760. return 2 * x ** 3 + 3 * y ** 2 - z - a
  761. x1 = np.array([0, 1, 2, 3])
  762. x2 = np.array([0, 10, 20, 30])
  763. x3 = np.array([0, 10, 20, 30])
  764. x4 = np.array([0, .1, .2, .30])
  765. points = (x1, x2, x3, x4)
  766. values = value_func_4d(
  767. *np.meshgrid(*points, indexing='ij', sparse=True))
  768. pts = (0.1, 0.3, np.transpose(np.linspace(0, 30, 4)),
  769. np.linspace(0, 0.3, 4))
  770. correct_result = interpn(points, values, pts)
  771. x1_descend = x1[::-1]
  772. x2_descend = x2[::-1]
  773. x3_descend = x3[::-1]
  774. x4_descend = x4[::-1]
  775. points_shuffled = (x1_descend, x2_descend, x3_descend, x4_descend)
  776. values_shuffled = value_func_4d(
  777. *np.meshgrid(*points_shuffled, indexing='ij', sparse=True))
  778. test_result = interpn(points_shuffled, values_shuffled, pts)
  779. assert_array_equal(correct_result, test_result)
  780. def test_invalid_points_order(self):
  781. x = np.array([.5, 2., 0., 4., 5.5]) # not ascending or descending
  782. y = np.array([.5, 2., 3., 4., 5.5])
  783. z = np.array([[1, 2, 1, 2, 1], [1, 2, 1, 2, 1], [1, 2, 3, 2, 1],
  784. [1, 2, 2, 2, 1], [1, 2, 1, 2, 1]])
  785. xi = np.array([[1, 2.3, 6.3, 0.5, 3.3, 1.2, 3],
  786. [1, 3.3, 1.2, -4.0, 5.0, 1.0, 3]]).T
  787. match = "must be strictly ascending or descending"
  788. with pytest.raises(ValueError, match=match):
  789. interpn((x, y), z, xi)
  790. def test_invalid_xi_dimensions(self):
  791. # https://github.com/scipy/scipy/issues/16519
  792. points = [(0, 1)]
  793. values = [0, 1]
  794. xi = np.ones((1, 1, 3))
  795. msg = ("The requested sample points xi have dimension 3, but this "
  796. "RegularGridInterpolator has dimension 1")
  797. with assert_raises(ValueError, match=msg):
  798. interpn(points, values, xi)
  799. def test_readonly_grid(self):
  800. # https://github.com/scipy/scipy/issues/17716
  801. x = np.linspace(0, 4, 5)
  802. y = np.linspace(0, 5, 6)
  803. z = np.linspace(0, 6, 7)
  804. points = (x, y, z)
  805. values = np.ones((5, 6, 7))
  806. point = np.array([2.21, 3.12, 1.15])
  807. for d in points:
  808. d.flags.writeable = False
  809. values.flags.writeable = False
  810. point.flags.writeable = False
  811. interpn(points, values, point)
  812. RegularGridInterpolator(points, values)(point)
  813. def test_2d_readonly_grid(self):
  814. # https://github.com/scipy/scipy/issues/17716
  815. # test special 2d case
  816. x = np.linspace(0, 4, 5)
  817. y = np.linspace(0, 5, 6)
  818. points = (x, y)
  819. values = np.ones((5, 6))
  820. point = np.array([2.21, 3.12])
  821. for d in points:
  822. d.flags.writeable = False
  823. values.flags.writeable = False
  824. point.flags.writeable = False
  825. interpn(points, values, point)
  826. RegularGridInterpolator(points, values)(point)
  827. def test_non_c_contiguous_grid(self):
  828. # https://github.com/scipy/scipy/issues/17716
  829. x = np.linspace(0, 4, 5)
  830. x = np.vstack((x, np.empty_like(x))).T.copy()[:, 0]
  831. assert not x.flags.c_contiguous
  832. y = np.linspace(0, 5, 6)
  833. z = np.linspace(0, 6, 7)
  834. points = (x, y, z)
  835. values = np.ones((5, 6, 7))
  836. point = np.array([2.21, 3.12, 1.15])
  837. interpn(points, values, point)
  838. RegularGridInterpolator(points, values)(point)
  839. @pytest.mark.parametrize("dtype", ['>f8', '<f8'])
  840. def test_endianness(self, dtype):
  841. # https://github.com/scipy/scipy/issues/17716
  842. # test special 2d case
  843. x = np.linspace(0, 4, 5, dtype=dtype)
  844. y = np.linspace(0, 5, 6, dtype=dtype)
  845. points = (x, y)
  846. values = np.ones((5, 6), dtype=dtype)
  847. point = np.array([2.21, 3.12], dtype=dtype)
  848. interpn(points, values, point)
  849. RegularGridInterpolator(points, values)(point)