_rgi.py 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675
  1. __all__ = ['RegularGridInterpolator', 'interpn']
  2. import itertools
  3. import numpy as np
  4. from .interpnd import _ndim_coords_from_arrays
  5. from ._cubic import PchipInterpolator
  6. from ._rgi_cython import evaluate_linear_2d, find_indices
  7. from ._bsplines import make_interp_spline
  8. from ._fitpack2 import RectBivariateSpline
  9. def _check_points(points):
  10. descending_dimensions = []
  11. grid = []
  12. for i, p in enumerate(points):
  13. # early make points float
  14. # see https://github.com/scipy/scipy/pull/17230
  15. p = np.asarray(p, dtype=float)
  16. if not np.all(p[1:] > p[:-1]):
  17. if np.all(p[1:] < p[:-1]):
  18. # input is descending, so make it ascending
  19. descending_dimensions.append(i)
  20. p = np.flip(p)
  21. else:
  22. raise ValueError(
  23. "The points in dimension %d must be strictly "
  24. "ascending or descending" % i)
  25. # see https://github.com/scipy/scipy/issues/17716
  26. p = np.ascontiguousarray(p)
  27. grid.append(p)
  28. return tuple(grid), tuple(descending_dimensions)
  29. def _check_dimensionality(points, values):
  30. if len(points) > values.ndim:
  31. raise ValueError("There are %d point arrays, but values has %d "
  32. "dimensions" % (len(points), values.ndim))
  33. for i, p in enumerate(points):
  34. if not np.asarray(p).ndim == 1:
  35. raise ValueError("The points in dimension %d must be "
  36. "1-dimensional" % i)
  37. if not values.shape[i] == len(p):
  38. raise ValueError("There are %d points and %d values in "
  39. "dimension %d" % (len(p), values.shape[i], i))
  40. class RegularGridInterpolator:
  41. """
  42. Interpolation on a regular or rectilinear grid in arbitrary dimensions.
  43. The data must be defined on a rectilinear grid; that is, a rectangular
  44. grid with even or uneven spacing. Linear, nearest-neighbor, spline
  45. interpolations are supported. After setting up the interpolator object,
  46. the interpolation method may be chosen at each evaluation.
  47. Parameters
  48. ----------
  49. points : tuple of ndarray of float, with shapes (m1, ), ..., (mn, )
  50. The points defining the regular grid in n dimensions. The points in
  51. each dimension (i.e. every elements of the points tuple) must be
  52. strictly ascending or descending.
  53. values : array_like, shape (m1, ..., mn, ...)
  54. The data on the regular grid in n dimensions. Complex data can be
  55. acceptable.
  56. method : str, optional
  57. The method of interpolation to perform. Supported are "linear",
  58. "nearest", "slinear", "cubic", "quintic" and "pchip". This
  59. parameter will become the default for the object's ``__call__``
  60. method. Default is "linear".
  61. bounds_error : bool, optional
  62. If True, when interpolated values are requested outside of the
  63. domain of the input data, a ValueError is raised.
  64. If False, then `fill_value` is used.
  65. Default is True.
  66. fill_value : float or None, optional
  67. The value to use for points outside of the interpolation domain.
  68. If None, values outside the domain are extrapolated.
  69. Default is ``np.nan``.
  70. Methods
  71. -------
  72. __call__
  73. Attributes
  74. ----------
  75. grid : tuple of ndarrays
  76. The points defining the regular grid in n dimensions.
  77. This tuple defines the full grid via
  78. ``np.meshgrid(*grid, indexing='ij')``
  79. values : ndarray
  80. Data values at the grid.
  81. method : str
  82. Interpolation method.
  83. fill_value : float or ``None``
  84. Use this value for out-of-bounds arguments to `__call__`.
  85. bounds_error : bool
  86. If ``True``, out-of-bounds argument raise a ``ValueError``.
  87. Notes
  88. -----
  89. Contrary to `LinearNDInterpolator` and `NearestNDInterpolator`, this class
  90. avoids expensive triangulation of the input data by taking advantage of the
  91. regular grid structure.
  92. In other words, this class assumes that the data is defined on a
  93. *rectilinear* grid.
  94. .. versionadded:: 0.14
  95. The 'slinear'(k=1), 'cubic'(k=3), and 'quintic'(k=5) methods are
  96. tensor-product spline interpolators, where `k` is the spline degree,
  97. If any dimension has fewer points than `k` + 1, an error will be raised.
  98. .. versionadded:: 1.9
  99. If the input data is such that dimensions have incommensurate
  100. units and differ by many orders of magnitude, the interpolant may have
  101. numerical artifacts. Consider rescaling the data before interpolating.
  102. Examples
  103. --------
  104. **Evaluate a function on the points of a 3-D grid**
  105. As a first example, we evaluate a simple example function on the points of
  106. a 3-D grid:
  107. >>> from scipy.interpolate import RegularGridInterpolator
  108. >>> import numpy as np
  109. >>> def f(x, y, z):
  110. ... return 2 * x**3 + 3 * y**2 - z
  111. >>> x = np.linspace(1, 4, 11)
  112. >>> y = np.linspace(4, 7, 22)
  113. >>> z = np.linspace(7, 9, 33)
  114. >>> xg, yg ,zg = np.meshgrid(x, y, z, indexing='ij', sparse=True)
  115. >>> data = f(xg, yg, zg)
  116. ``data`` is now a 3-D array with ``data[i, j, k] = f(x[i], y[j], z[k])``.
  117. Next, define an interpolating function from this data:
  118. >>> interp = RegularGridInterpolator((x, y, z), data)
  119. Evaluate the interpolating function at the two points
  120. ``(x,y,z) = (2.1, 6.2, 8.3)`` and ``(3.3, 5.2, 7.1)``:
  121. >>> pts = np.array([[2.1, 6.2, 8.3],
  122. ... [3.3, 5.2, 7.1]])
  123. >>> interp(pts)
  124. array([ 125.80469388, 146.30069388])
  125. which is indeed a close approximation to
  126. >>> f(2.1, 6.2, 8.3), f(3.3, 5.2, 7.1)
  127. (125.54200000000002, 145.894)
  128. **Interpolate and extrapolate a 2D dataset**
  129. As a second example, we interpolate and extrapolate a 2D data set:
  130. >>> x, y = np.array([-2, 0, 4]), np.array([-2, 0, 2, 5])
  131. >>> def ff(x, y):
  132. ... return x**2 + y**2
  133. >>> xg, yg = np.meshgrid(x, y, indexing='ij')
  134. >>> data = ff(xg, yg)
  135. >>> interp = RegularGridInterpolator((x, y), data,
  136. ... bounds_error=False, fill_value=None)
  137. >>> import matplotlib.pyplot as plt
  138. >>> fig = plt.figure()
  139. >>> ax = fig.add_subplot(projection='3d')
  140. >>> ax.scatter(xg.ravel(), yg.ravel(), data.ravel(),
  141. ... s=60, c='k', label='data')
  142. Evaluate and plot the interpolator on a finer grid
  143. >>> xx = np.linspace(-4, 9, 31)
  144. >>> yy = np.linspace(-4, 9, 31)
  145. >>> X, Y = np.meshgrid(xx, yy, indexing='ij')
  146. >>> # interpolator
  147. >>> ax.plot_wireframe(X, Y, interp((X, Y)), rstride=3, cstride=3,
  148. ... alpha=0.4, color='m', label='linear interp')
  149. >>> # ground truth
  150. >>> ax.plot_wireframe(X, Y, ff(X, Y), rstride=3, cstride=3,
  151. ... alpha=0.4, label='ground truth')
  152. >>> plt.legend()
  153. >>> plt.show()
  154. Other examples are given
  155. :ref:`in the tutorial <tutorial-interpolate_regular_grid_interpolator>`.
  156. See Also
  157. --------
  158. NearestNDInterpolator : Nearest neighbor interpolation on *unstructured*
  159. data in N dimensions
  160. LinearNDInterpolator : Piecewise linear interpolant on *unstructured* data
  161. in N dimensions
  162. interpn : a convenience function which wraps `RegularGridInterpolator`
  163. scipy.ndimage.map_coordinates : interpolation on grids with equal spacing
  164. (suitable for e.g., N-D image resampling)
  165. References
  166. ----------
  167. .. [1] Python package *regulargrid* by Johannes Buchner, see
  168. https://pypi.python.org/pypi/regulargrid/
  169. .. [2] Wikipedia, "Trilinear interpolation",
  170. https://en.wikipedia.org/wiki/Trilinear_interpolation
  171. .. [3] Weiser, Alan, and Sergio E. Zarantonello. "A note on piecewise linear
  172. and multilinear table interpolation in many dimensions." MATH.
  173. COMPUT. 50.181 (1988): 189-196.
  174. https://www.ams.org/journals/mcom/1988-50-181/S0025-5718-1988-0917826-0/S0025-5718-1988-0917826-0.pdf
  175. :doi:`10.1090/S0025-5718-1988-0917826-0`
  176. """
  177. # this class is based on code originally programmed by Johannes Buchner,
  178. # see https://github.com/JohannesBuchner/regulargrid
  179. _SPLINE_DEGREE_MAP = {"slinear": 1, "cubic": 3, "quintic": 5, 'pchip': 3}
  180. _SPLINE_METHODS = list(_SPLINE_DEGREE_MAP.keys())
  181. _ALL_METHODS = ["linear", "nearest"] + _SPLINE_METHODS
  182. def __init__(self, points, values, method="linear", bounds_error=True,
  183. fill_value=np.nan):
  184. if method not in self._ALL_METHODS:
  185. raise ValueError("Method '%s' is not defined" % method)
  186. elif method in self._SPLINE_METHODS:
  187. self._validate_grid_dimensions(points, method)
  188. self.method = method
  189. self.bounds_error = bounds_error
  190. self.grid, self._descending_dimensions = _check_points(points)
  191. self.values = self._check_values(values)
  192. self._check_dimensionality(self.grid, self.values)
  193. self.fill_value = self._check_fill_value(self.values, fill_value)
  194. if self._descending_dimensions:
  195. self.values = np.flip(values, axis=self._descending_dimensions)
  196. def _check_dimensionality(self, grid, values):
  197. _check_dimensionality(grid, values)
  198. def _check_points(self, points):
  199. return _check_points(points)
  200. def _check_values(self, values):
  201. if not hasattr(values, 'ndim'):
  202. # allow reasonable duck-typed values
  203. values = np.asarray(values)
  204. if hasattr(values, 'dtype') and hasattr(values, 'astype'):
  205. if not np.issubdtype(values.dtype, np.inexact):
  206. values = values.astype(float)
  207. return values
  208. def _check_fill_value(self, values, fill_value):
  209. if fill_value is not None:
  210. fill_value_dtype = np.asarray(fill_value).dtype
  211. if (hasattr(values, 'dtype') and not
  212. np.can_cast(fill_value_dtype, values.dtype,
  213. casting='same_kind')):
  214. raise ValueError("fill_value must be either 'None' or "
  215. "of a type compatible with values")
  216. return fill_value
  217. def __call__(self, xi, method=None):
  218. """
  219. Interpolation at coordinates.
  220. Parameters
  221. ----------
  222. xi : ndarray of shape (..., ndim)
  223. The coordinates to evaluate the interpolator at.
  224. method : str, optional
  225. The method of interpolation to perform. Supported are "linear",
  226. "nearest", "slinear", "cubic", "quintic" and "pchip". Default is
  227. the method chosen when the interpolator was created.
  228. Returns
  229. -------
  230. values_x : ndarray, shape xi.shape[:-1] + values.shape[ndim:]
  231. Interpolated values at `xi`. See notes for behaviour when
  232. ``xi.ndim == 1``.
  233. Notes
  234. -----
  235. In the case that ``xi.ndim == 1`` a new axis is inserted into
  236. the 0 position of the returned array, values_x, so its shape is
  237. instead ``(1,) + values.shape[ndim:]``.
  238. Examples
  239. --------
  240. Here we define a nearest-neighbor interpolator of a simple function
  241. >>> import numpy as np
  242. >>> x, y = np.array([0, 1, 2]), np.array([1, 3, 7])
  243. >>> def f(x, y):
  244. ... return x**2 + y**2
  245. >>> data = f(*np.meshgrid(x, y, indexing='ij', sparse=True))
  246. >>> from scipy.interpolate import RegularGridInterpolator
  247. >>> interp = RegularGridInterpolator((x, y), data, method='nearest')
  248. By construction, the interpolator uses the nearest-neighbor
  249. interpolation
  250. >>> interp([[1.5, 1.3], [0.3, 4.5]])
  251. array([2., 9.])
  252. We can however evaluate the linear interpolant by overriding the
  253. `method` parameter
  254. >>> interp([[1.5, 1.3], [0.3, 4.5]], method='linear')
  255. array([ 4.7, 24.3])
  256. """
  257. is_method_changed = self.method != method
  258. method = self.method if method is None else method
  259. if method not in self._ALL_METHODS:
  260. raise ValueError("Method '%s' is not defined" % method)
  261. xi, xi_shape, ndim, nans, out_of_bounds = self._prepare_xi(xi)
  262. if method == "linear":
  263. indices, norm_distances = self._find_indices(xi.T)
  264. if (ndim == 2 and hasattr(self.values, 'dtype') and
  265. self.values.ndim == 2 and self.values.flags.writeable and
  266. self.values.dtype in (np.float64, np.complex128) and
  267. self.values.dtype.byteorder == '='):
  268. # until cython supports const fused types, the fast path
  269. # cannot support non-writeable values
  270. # a fast path
  271. out = np.empty(indices.shape[1], dtype=self.values.dtype)
  272. result = evaluate_linear_2d(self.values,
  273. indices,
  274. norm_distances,
  275. self.grid,
  276. out)
  277. else:
  278. result = self._evaluate_linear(indices, norm_distances)
  279. elif method == "nearest":
  280. indices, norm_distances = self._find_indices(xi.T)
  281. result = self._evaluate_nearest(indices, norm_distances)
  282. elif method in self._SPLINE_METHODS:
  283. if is_method_changed:
  284. self._validate_grid_dimensions(self.grid, method)
  285. result = self._evaluate_spline(xi, method)
  286. if not self.bounds_error and self.fill_value is not None:
  287. result[out_of_bounds] = self.fill_value
  288. # f(nan) = nan, if any
  289. if np.any(nans):
  290. result[nans] = np.nan
  291. return result.reshape(xi_shape[:-1] + self.values.shape[ndim:])
  292. def _prepare_xi(self, xi):
  293. ndim = len(self.grid)
  294. xi = _ndim_coords_from_arrays(xi, ndim=ndim)
  295. if xi.shape[-1] != len(self.grid):
  296. raise ValueError("The requested sample points xi have dimension "
  297. f"{xi.shape[-1]} but this "
  298. f"RegularGridInterpolator has dimension {ndim}")
  299. xi_shape = xi.shape
  300. xi = xi.reshape(-1, xi_shape[-1])
  301. xi = np.asarray(xi, dtype=float)
  302. # find nans in input
  303. nans = np.any(np.isnan(xi), axis=-1)
  304. if self.bounds_error:
  305. for i, p in enumerate(xi.T):
  306. if not np.logical_and(np.all(self.grid[i][0] <= p),
  307. np.all(p <= self.grid[i][-1])):
  308. raise ValueError("One of the requested xi is out of bounds "
  309. "in dimension %d" % i)
  310. out_of_bounds = None
  311. else:
  312. out_of_bounds = self._find_out_of_bounds(xi.T)
  313. return xi, xi_shape, ndim, nans, out_of_bounds
  314. def _evaluate_linear(self, indices, norm_distances):
  315. # slice for broadcasting over trailing dimensions in self.values
  316. vslice = (slice(None),) + (None,)*(self.values.ndim - len(indices))
  317. # Compute shifting up front before zipping everything together
  318. shift_norm_distances = [1 - yi for yi in norm_distances]
  319. shift_indices = [i + 1 for i in indices]
  320. # The formula for linear interpolation in 2d takes the form:
  321. # values = self.values[(i0, i1)] * (1 - y0) * (1 - y1) + \
  322. # self.values[(i0, i1 + 1)] * (1 - y0) * y1 + \
  323. # self.values[(i0 + 1, i1)] * y0 * (1 - y1) + \
  324. # self.values[(i0 + 1, i1 + 1)] * y0 * y1
  325. # We pair i with 1 - yi (zipped1) and i + 1 with yi (zipped2)
  326. zipped1 = zip(indices, shift_norm_distances)
  327. zipped2 = zip(shift_indices, norm_distances)
  328. # Take all products of zipped1 and zipped2 and iterate over them
  329. # to get the terms in the above formula. This corresponds to iterating
  330. # over the vertices of a hypercube.
  331. hypercube = itertools.product(*zip(zipped1, zipped2))
  332. value = np.array([0.])
  333. for h in hypercube:
  334. edge_indices, weights = zip(*h)
  335. weight = np.array([1.])
  336. for w in weights:
  337. weight = weight * w
  338. term = np.asarray(self.values[edge_indices]) * weight[vslice]
  339. value = value + term # cannot use += because broadcasting
  340. return value
  341. def _evaluate_nearest(self, indices, norm_distances):
  342. idx_res = [np.where(yi <= .5, i, i + 1)
  343. for i, yi in zip(indices, norm_distances)]
  344. return self.values[tuple(idx_res)]
  345. def _validate_grid_dimensions(self, points, method):
  346. k = self._SPLINE_DEGREE_MAP[method]
  347. for i, point in enumerate(points):
  348. ndim = len(np.atleast_1d(point))
  349. if ndim <= k:
  350. raise ValueError(f"There are {ndim} points in dimension {i},"
  351. f" but method {method} requires at least "
  352. f" {k+1} points per dimension.")
  353. def _evaluate_spline(self, xi, method):
  354. # ensure xi is 2D list of points to evaluate (`m` is the number of
  355. # points and `n` is the number of interpolation dimensions,
  356. # ``n == len(self.grid)``.)
  357. if xi.ndim == 1:
  358. xi = xi.reshape((1, xi.size))
  359. m, n = xi.shape
  360. # Reorder the axes: n-dimensional process iterates over the
  361. # interpolation axes from the last axis downwards: E.g. for a 4D grid
  362. # the order of axes is 3, 2, 1, 0. Each 1D interpolation works along
  363. # the 0th axis of its argument array (for 1D routine it's its ``y``
  364. # array). Thus permute the interpolation axes of `values` *and keep
  365. # trailing dimensions trailing*.
  366. axes = tuple(range(self.values.ndim))
  367. axx = axes[:n][::-1] + axes[n:]
  368. values = self.values.transpose(axx)
  369. if method == 'pchip':
  370. _eval_func = self._do_pchip
  371. else:
  372. _eval_func = self._do_spline_fit
  373. k = self._SPLINE_DEGREE_MAP[method]
  374. # Non-stationary procedure: difficult to vectorize this part entirely
  375. # into numpy-level operations. Unfortunately this requires explicit
  376. # looping over each point in xi.
  377. # can at least vectorize the first pass across all points in the
  378. # last variable of xi.
  379. last_dim = n - 1
  380. first_values = _eval_func(self.grid[last_dim],
  381. values,
  382. xi[:, last_dim],
  383. k)
  384. # the rest of the dimensions have to be on a per point-in-xi basis
  385. shape = (m, *self.values.shape[n:])
  386. result = np.empty(shape, dtype=self.values.dtype)
  387. for j in range(m):
  388. # Main process: Apply 1D interpolate in each dimension
  389. # sequentially, starting with the last dimension.
  390. # These are then "folded" into the next dimension in-place.
  391. folded_values = first_values[j, ...]
  392. for i in range(last_dim-1, -1, -1):
  393. # Interpolate for each 1D from the last dimensions.
  394. # This collapses each 1D sequence into a scalar.
  395. folded_values = _eval_func(self.grid[i],
  396. folded_values,
  397. xi[j, i],
  398. k)
  399. result[j, ...] = folded_values
  400. return result
  401. @staticmethod
  402. def _do_spline_fit(x, y, pt, k):
  403. local_interp = make_interp_spline(x, y, k=k, axis=0)
  404. values = local_interp(pt)
  405. return values
  406. @staticmethod
  407. def _do_pchip(x, y, pt, k):
  408. local_interp = PchipInterpolator(x, y, axis=0)
  409. values = local_interp(pt)
  410. return values
  411. def _find_indices(self, xi):
  412. return find_indices(self.grid, xi)
  413. def _find_out_of_bounds(self, xi):
  414. # check for out of bounds xi
  415. out_of_bounds = np.zeros((xi.shape[1]), dtype=bool)
  416. # iterate through dimensions
  417. for x, grid in zip(xi, self.grid):
  418. out_of_bounds += x < grid[0]
  419. out_of_bounds += x > grid[-1]
  420. return out_of_bounds
  421. def interpn(points, values, xi, method="linear", bounds_error=True,
  422. fill_value=np.nan):
  423. """
  424. Multidimensional interpolation on regular or rectilinear grids.
  425. Strictly speaking, not all regular grids are supported - this function
  426. works on *rectilinear* grids, that is, a rectangular grid with even or
  427. uneven spacing.
  428. Parameters
  429. ----------
  430. points : tuple of ndarray of float, with shapes (m1, ), ..., (mn, )
  431. The points defining the regular grid in n dimensions. The points in
  432. each dimension (i.e. every elements of the points tuple) must be
  433. strictly ascending or descending.
  434. values : array_like, shape (m1, ..., mn, ...)
  435. The data on the regular grid in n dimensions. Complex data can be
  436. acceptable.
  437. xi : ndarray of shape (..., ndim)
  438. The coordinates to sample the gridded data at
  439. method : str, optional
  440. The method of interpolation to perform. Supported are "linear",
  441. "nearest", "slinear", "cubic", "quintic", "pchip", and "splinef2d".
  442. "splinef2d" is only supported for 2-dimensional data.
  443. bounds_error : bool, optional
  444. If True, when interpolated values are requested outside of the
  445. domain of the input data, a ValueError is raised.
  446. If False, then `fill_value` is used.
  447. fill_value : number, optional
  448. If provided, the value to use for points outside of the
  449. interpolation domain. If None, values outside
  450. the domain are extrapolated. Extrapolation is not supported by method
  451. "splinef2d".
  452. Returns
  453. -------
  454. values_x : ndarray, shape xi.shape[:-1] + values.shape[ndim:]
  455. Interpolated values at `xi`. See notes for behaviour when
  456. ``xi.ndim == 1``.
  457. Notes
  458. -----
  459. .. versionadded:: 0.14
  460. In the case that ``xi.ndim == 1`` a new axis is inserted into
  461. the 0 position of the returned array, values_x, so its shape is
  462. instead ``(1,) + values.shape[ndim:]``.
  463. If the input data is such that input dimensions have incommensurate
  464. units and differ by many orders of magnitude, the interpolant may have
  465. numerical artifacts. Consider rescaling the data before interpolation.
  466. Examples
  467. --------
  468. Evaluate a simple example function on the points of a regular 3-D grid:
  469. >>> import numpy as np
  470. >>> from scipy.interpolate import interpn
  471. >>> def value_func_3d(x, y, z):
  472. ... return 2 * x + 3 * y - z
  473. >>> x = np.linspace(0, 4, 5)
  474. >>> y = np.linspace(0, 5, 6)
  475. >>> z = np.linspace(0, 6, 7)
  476. >>> points = (x, y, z)
  477. >>> values = value_func_3d(*np.meshgrid(*points, indexing='ij'))
  478. Evaluate the interpolating function at a point
  479. >>> point = np.array([2.21, 3.12, 1.15])
  480. >>> print(interpn(points, values, point))
  481. [12.63]
  482. See Also
  483. --------
  484. NearestNDInterpolator : Nearest neighbor interpolation on unstructured
  485. data in N dimensions
  486. LinearNDInterpolator : Piecewise linear interpolant on unstructured data
  487. in N dimensions
  488. RegularGridInterpolator : interpolation on a regular or rectilinear grid
  489. in arbitrary dimensions (`interpn` wraps this
  490. class).
  491. RectBivariateSpline : Bivariate spline approximation over a rectangular mesh
  492. scipy.ndimage.map_coordinates : interpolation on grids with equal spacing
  493. (suitable for e.g., N-D image resampling)
  494. """
  495. # sanity check 'method' kwarg
  496. if method not in ["linear", "nearest", "cubic", "quintic", "pchip",
  497. "splinef2d", "slinear"]:
  498. raise ValueError("interpn only understands the methods 'linear', "
  499. "'nearest', 'slinear', 'cubic', 'quintic', 'pchip', "
  500. f"and 'splinef2d'. You provided {method}.")
  501. if not hasattr(values, 'ndim'):
  502. values = np.asarray(values)
  503. ndim = values.ndim
  504. if ndim > 2 and method == "splinef2d":
  505. raise ValueError("The method splinef2d can only be used for "
  506. "2-dimensional input data")
  507. if not bounds_error and fill_value is None and method == "splinef2d":
  508. raise ValueError("The method splinef2d does not support extrapolation.")
  509. # sanity check consistency of input dimensions
  510. if len(points) > ndim:
  511. raise ValueError("There are %d point arrays, but values has %d "
  512. "dimensions" % (len(points), ndim))
  513. if len(points) != ndim and method == 'splinef2d':
  514. raise ValueError("The method splinef2d can only be used for "
  515. "scalar data with one point per coordinate")
  516. grid, descending_dimensions = _check_points(points)
  517. _check_dimensionality(grid, values)
  518. # sanity check requested xi
  519. xi = _ndim_coords_from_arrays(xi, ndim=len(grid))
  520. if xi.shape[-1] != len(grid):
  521. raise ValueError("The requested sample points xi have dimension "
  522. "%d, but this RegularGridInterpolator has "
  523. "dimension %d" % (xi.shape[-1], len(grid)))
  524. if bounds_error:
  525. for i, p in enumerate(xi.T):
  526. if not np.logical_and(np.all(grid[i][0] <= p),
  527. np.all(p <= grid[i][-1])):
  528. raise ValueError("One of the requested xi is out of bounds "
  529. "in dimension %d" % i)
  530. # perform interpolation
  531. if method in ["linear", "nearest", "slinear", "cubic", "quintic", "pchip"]:
  532. interp = RegularGridInterpolator(points, values, method=method,
  533. bounds_error=bounds_error,
  534. fill_value=fill_value)
  535. return interp(xi)
  536. elif method == "splinef2d":
  537. xi_shape = xi.shape
  538. xi = xi.reshape(-1, xi.shape[-1])
  539. # RectBivariateSpline doesn't support fill_value; we need to wrap here
  540. idx_valid = np.all((grid[0][0] <= xi[:, 0], xi[:, 0] <= grid[0][-1],
  541. grid[1][0] <= xi[:, 1], xi[:, 1] <= grid[1][-1]),
  542. axis=0)
  543. result = np.empty_like(xi[:, 0])
  544. # make a copy of values for RectBivariateSpline
  545. interp = RectBivariateSpline(points[0], points[1], values[:])
  546. result[idx_valid] = interp.ev(xi[idx_valid, 0], xi[idx_valid, 1])
  547. result[np.logical_not(idx_valid)] = fill_value
  548. return result.reshape(xi_shape[:-1])