1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457
|
.. _tutorial-interpolate_1Dsection:
.. currentmodule:: scipy.interpolate
=================
1-D interpolation
=================
Piecewise linear interpolation
==============================
If all you need is a linear (a.k.a. broken line) interpolation, you can use
the `numpy.interp` routine. It takes two arrays of data to interpolate, ``x``,
and ``y``, and a third array, ``xnew``, of points to evaluate the interpolation on:
.. plot::
>>> import numpy as np
>>> x = np.linspace(0, 10, num=11)
>>> y = np.cos(-x**2 / 9.0)
Construct the interpolation
>>> xnew = np.linspace(0, 10, num=1001)
>>> ynew = np.interp(xnew, x, y)
And plot it
>>> import matplotlib.pyplot as plt
>>> plt.plot(xnew, ynew, '-', label='linear interp')
>>> plt.plot(x, y, 'o', label='data')
>>> plt.legend(loc='best')
>>> plt.show()
.. :caption: One-dimensional interpolation using `numpy.interp`
One limitation of `numpy.interp` is that it does not allow controlling the
extrapolation. See the :ref:`interpolation with B-Splines section <tutorial-interpolate_bsplines>`
section for alternative routines which provide this kind of functionality.
Cubic splines
=============
Of course, piecewise linear interpolation produces corners at data points,
where linear pieces join. To produce a smoother curve, you can use cubic
splines, where the interpolating curve is made of cubic pieces with matching
first and second derivatives. In code, these objects are represented via the
``CubicSpline`` class instances. An instance is constructed with the ``x`` and
``y`` arrays of data, and then it can be evaluated using the target ``xnew``
values:
>>> from scipy.interpolate import CubicSpline
>>> spl = CubicSpline([1, 2, 3, 4, 5, 6], [1, 4, 8, 16, 25, 36])
>>> spl(2.5)
5.57
A `CubicSpline` object's ``__call__`` method accepts both scalar values and
arrays. It also accepts a second argument, ``nu``, to evaluate the
derivative of order ``nu``. As an example, we plot the derivatives of a spline:
.. plot::
>>> from scipy.interpolate import CubicSpline
>>> x = np.linspace(0, 10, num=11)
>>> y = np.cos(-x**2 / 9.)
>>> spl = CubicSpline(x, y)
>>> import matplotlib.pyplot as plt
>>> fig, ax = plt.subplots(4, 1, figsize=(5, 7))
>>> xnew = np.linspace(0, 10, num=1001)
>>> ax[0].plot(xnew, spl(xnew))
>>> ax[0].plot(x, y, 'o', label='data')
>>> ax[1].plot(xnew, spl(xnew, nu=1), '--', label='1st derivative')
>>> ax[2].plot(xnew, spl(xnew, nu=2), '--', label='2nd derivative')
>>> ax[3].plot(xnew, spl(xnew, nu=3), '--', label='3rd derivative')
>>> for j in range(4):
... ax[j].legend(loc='best')
>>> plt.tight_layout()
>>> plt.show()
Note that the first and second derivatives are continuous by construction, and
the third derivative jumps at data points.
Monotone interpolants
=====================
Cubic splines are by construction twice continuously differentiable. This may
lead to the spline function oscillating and ''overshooting'' in between the
data points. In these situations, an alternative is to use the so-called
*monotone* cubic interpolants: these are constructed to be only once
continuously differentiable, and attempt to preserve the local shape implied
by the data. `scipy.interpolate` provides two objects of this kind:
`PchipInterpolator` and `Akima1DInterpolator` . To illustrate, let's consider
data with an outlier:
.. plot::
>>> from scipy.interpolate import CubicSpline, PchipInterpolator, Akima1DInterpolator
>>> x = np.array([1., 2., 3., 4., 4.5, 5., 6., 7., 8])
>>> y = x**2
>>> y[4] += 101
>>> import matplotlib.pyplot as plt
>>> xx = np.linspace(1, 8, 51)
>>> plt.plot(xx, CubicSpline(x, y)(xx), '--', label='spline')
>>> plt.plot(xx, Akima1DInterpolator(x, y)(xx), '-', label='Akima1D')
>>> plt.plot(xx, PchipInterpolator(x, y)(xx), '-', label='pchip')
>>> plt.plot(x, y, 'o')
>>> plt.legend()
>>> plt.show()
.. _tutorial-interpolate_bsplines:
Interpolation with B-splines
============================
B-splines form an alternative (if formally equivalent) representation of piecewise polynomials. This basis is generally more computationally stable than the power basis and is useful for a variety of applications which include interpolation, regression and curve representation. Details are given in the :ref:`piecewise polynomials section <tutorial-interpolate_ppoly>`, and here we illustrate their usage by constructing the interpolation of a sine function:
.. plot::
>>> x = np.linspace(0, 3/2, 7)
>>> y = np.sin(np.pi*x)
To construct the interpolating objects given data arrays, ``x`` and ``y``,
we use the `make_interp_spline` function:
>>> from scipy.interpolate import make_interp_spline
>>> bspl = make_interp_spline(x, y, k=3)
This function returns an object which has an interface similar to that
of the `CubicSpline` objects. In particular, it can be evaluated at a data
point and differentiated:
>>> der = bspl.derivative() # a BSpline representing the derivative
>>> import matplotlib.pyplot as plt
>>> xx = np.linspace(0, 3/2, 51)
>>> plt.plot(xx, bspl(xx), '--', label=r'$\sin(\pi x)$ approx')
>>> plt.plot(x, y, 'o', label='data')
>>> plt.plot(xx, der(xx)/np.pi, '--', label=r'$d \sin(\pi x)/dx / \pi$ approx')
>>> plt.legend()
>>> plt.show()
Note that by specifying ``k=3`` in the `make_interp_spline` call, we requested
a cubic spline (this is the default, so ``k=3`` could have been omitted); the
derivative of a cubic is a quadratic:
>>> bspl.k, der.k
(3, 2)
By default, the result of ``make_interp_spline(x, y)`` is equivalent to
``CubicSpline(x, y)``. The difference is that the former allows several optional
capabilities: it can construct splines of various degrees (via the optional
argument ``k``) and predefined knots (via the optional argument ``t``).
Boundary conditions for the spline interpolation can be controlled by
the ``bc_type`` argument to `make_interp_spline` function and `CubicSpline`
constructor. By default, both use the 'not-a-knot' boundary
condition.
.. _tutorial-interpolate_linear_spline:
Non-cubic splines
-----------------
One use of ``make_interp_spline`` is constructing a linear interpolant with
linear extrapolation since ``make_interp_spline`` extrapolates by default. Consider
>>> from scipy.interpolate import make_interp_spline
>>> x = np.linspace(0, 5, 11)
>>> y = 2*x
>>> spl = make_interp_spline(x, y, k=1) # k=1: linear
>>> spl([-1, 6])
[-2., 12.]
>>> np.interp([-1, 6], x, y)
[0., 10.]
See :ref:`the extrapolation section <tutorial-extrapolation>` for more details
and discussion.
.. _tutorial-interpolate_batching:
Batches of ``y``
================
Univariate interpolators accept not only one-dimensional ``y`` arrays, but also
``y.ndim > 1``. The interpretation is that ``y`` is a *batch* of 1D data arrays: by
default, the zeroth dimension of ``y`` is the interpolation axis, and the trailing
dimensions are batch dimensions.
Consider a collection (a *batch*) of functions :math:`f_j` sampled at the points
:math:`x_i`. We can instantiate a single interpolator for all of these functions by
providing a two-dimensional array ``y`` such that ``y[i, j]`` records :math:`f_j(x_i)`.
To illustrate:
>>> import numpy as np
>>> import matplotlib.pyplot as plt
>>> from scipy.interpolate import make_interp_spline
>>> n = 11
>>> x = 2 * np.pi * np.arange(n) / n
>>> x.shape
(11,)
>>> y = np.stack((np.sin(x)**2, np.cos(x)), axis=1)
>>> y.shape
(11, 2)
>>> spl = make_interp_spline(x, y)
>>> xv = np.linspace(0, 2*np.pi, 51)
>>> plt.plot(x, y, 'o')
>>> plt.plot(xv, spl(xv), '-')
>>> plt.show()
Several notes are in order. First and foremost, the behavior here looks similar to
NumPy's broadcasting, but differs in two respects:
1. The ``x`` array is expected to be 1D even if the ``y`` array is not: ``x.ndim == 1``
while ``y.ndim >= 1``. There is no broadcasting of ``x`` vs ``y``.
2. By default, the *trailing* dimensions are used as batch dimensions, in contrast
to the NumPy convention of using the *leading* dimensions as batch dimensions.
Second, the interpolation axis can be controlled by an optional ``axis`` argument. The
example above uses the default value of ``axis=0``. For a non-default values, the
following is true:
- ``y.shape[axis] == x.size`` (otherwise en error is raised)
- the shape of ``spl(xv)`` is ``y.shape[axis:] + xv.shape + y.shape[:axis]``
While we demonstrated the batching behavior with `make_interp_spline`, in fact the
majority of univariate interpolators support this functionality: `PchipInterpolator` and
`Akima1DInterpolator`, `CubicSpline`; low-level polynomial representation classes,
`PPoly`, `BPoly` and `BSpline`; as well as :ref:`least-squares fit and spline smoothing
functions<tutorial-interpolate_fitpack>`, `make_lsq_spline` and `make_smoothing_spline`.
.. _tutorial-interpolate_parametric:
Parametric spline curves
========================
So far we considered spline *functions*, where the data, ``y``, is expected to
depend explicitly on the independent variable ``x``---so that the interpolating
function satisfies :math:`f(x_j) = y_j`. Spline *curves* treat
the ``x`` and ``y`` arrays as coordinates of points, :math:`\mathbf{p}_j` on a
plane, and an interpolating curve which passes through these points is
parameterized by some additional parameter (typically called ``u``). Note that
this construction readily generalizes to higher dimensions where
:math:`\mathbf{p}_j` are points in an N-dimensional space.
Spline curves can be easily constructed using the fact that interpolation
functions handle multidimensional data arrays, as discussed in
:ref:`the previous section <tutorial-interpolate_batching>`.
The values of the parameter, ``u``, corresponding to the data points, need to be
separately supplied by the user.
The choice of parametrization is problem-dependent and different parametrizations
may produce vastly different curves. As an example, we consider three
parametrizations of (a somewhat difficult) dataset, which we take from
Chapter 6 of Ref [1] listed in the `BSpline` docstring:
.. plot::
>>> x = [0, 1, 2, 3, 4, 5, 6]
>>> y = [0, 0, 0, 9, 0, 0, 0]
>>> p = np.stack((x, y))
>>> p
array([[0, 1, 2, 3, 4, 5, 6],
[0, 0, 0, 9, 0, 0, 0]])
We take elements of the ``p`` array as coordinates of seven points on the
plane, where ``p[:, j]`` gives the coordinates of the point
:math:`\mathbf{p}_j`.
First, consider the *uniform* parametrization, :math:`u_j = j`:
>>> u_unif = x
Second, we consider the so-called *cord length* parametrization, which is
nothing but a cumulative length of straight line segments connecting the
data points:
.. math::
u_j = u_{j-1} + |\mathbf{p}_j - \mathbf{p}_{j-1}|
for :math:`j=1, 2, \dots` and :math:`u_0 = 0`. Here :math:`| \cdots |` is the
length between the consecutive points :math:`p_j` on the plane.
>>> dp = p[:, 1:] - p[:, :-1] # 2-vector distances between points
>>> l = (dp**2).sum(axis=0) # squares of lengths of 2-vectors between points
>>> u_cord = np.sqrt(l).cumsum() # cumulative sums of 2-norms
>>> u_cord = np.r_[0, u_cord] # the first point is parameterized at zero
Finally, we consider what is sometimes called the *centripetal*
parametrization: :math:`u_j = u_{j-1} + |\mathbf{p}_j - \mathbf{p}_{j-1}|^{1/2}`.
Due to the extra square root, the difference between consecutive values
:math:`u_j - u_{j-1}` will be smaller than for the cord length parametrization:
>>> u_c = np.r_[0, np.cumsum((dp**2).sum(axis=0)**0.25)]
Now plot the resulting curves:
>>> from scipy.interpolate import make_interp_spline
>>> import matplotlib.pyplot as plt
>>> fig, ax = plt.subplots(1, 3, figsize=(8, 3))
>>> parametrizations = ['uniform', 'cord length', 'centripetal']
>>>
>>> for j, u in enumerate([u_unif, u_cord, u_c]):
... spl = make_interp_spline(u, p, axis=1) # note p is a 2D array
...
... uu = np.linspace(u[0], u[-1], 51)
... xx, yy = spl(uu)
...
... ax[j].plot(xx, yy, '--')
... ax[j].plot(p[0, :], p[1, :], 'o')
... ax[j].set_title(parametrizations[j])
>>> plt.show()
Missing data
============
We note that `scipy.interpolate` does *not* support interpolation with missing
data. Two popular ways of representing missing data are using masked arrays of
the `numpy.ma` library, and encoding missing values as not-a-number, ``NaN``.
Neither of these two approaches is directly supported in `scipy.interpolate`.
Individual routines may offer partial support, and/or workarounds, but in
general, the library firmly adheres to the IEEE 754 semantics where a ``NaN``
means *not-a-number*, i.e. a result of an illegal mathematical operation
(e.g., division by zero), not *missing*.
.. _tutorial-interpolate_interp1d:
Legacy interface for 1-D interpolation (:class:`interp1d`)
==========================================================
.. note::
`interp1d` is considered legacy API and is not recommended for use in new
code. Consider using :ref:`more specific interpolators instead <tutorial-interpolate_interp1d_replacements>`.
The `interp1d` class in `scipy.interpolate` is a convenient method to
create a function based on fixed data points, which can be evaluated
anywhere within the domain defined by the given data using linear
interpolation. An instance of this class is created by passing the 1-D
vectors comprising the data. The instance of this class defines a
``__call__`` method and can therefore be treated like a function which
interpolates between known data values to obtain unknown values.
Behavior at the boundary can be
specified at instantiation time. The following example demonstrates
its use, for linear and cubic spline interpolation:
.. plot::
:alt: "This code generates an X-Y plot of a time-series with amplitude on the Y axis and time on the X axis. The original time-series is shown as a series of blue markers roughly defining some kind of oscillation. An orange trace showing the linear interpolation is drawn atop the data forming a jagged representation of the original signal. A dotted green cubic interpolation is also drawn that appears to smoothly represent the source data."
>>> from scipy.interpolate import interp1d
>>> x = np.linspace(0, 10, num=11, endpoint=True)
>>> y = np.cos(-x**2/9.0)
>>> f = interp1d(x, y)
>>> f2 = interp1d(x, y, kind='cubic')
>>> xnew = np.linspace(0, 10, num=41, endpoint=True)
>>> import matplotlib.pyplot as plt
>>> plt.plot(x, y, 'o', xnew, f(xnew), '-', xnew, f2(xnew), '--')
>>> plt.legend(['data', 'linear', 'cubic'], loc='best')
>>> plt.show()
.. :caption: One-dimensional interpolation using the
.. class :obj:`interpolate.interp1d` with
.. kind equals `linear` and `cubic`.
The 'cubic' kind of `interp1d` is equivalent to `make_interp_spline`, and
the 'linear' kind is equivalent to `numpy.interp` while also allowing
N-dimensional ``y`` arrays.
Another set of interpolations in `interp1d` is `nearest`, `previous`, and
`next`, where they return the nearest, previous, or next point along the
x-axis. Nearest and next can be thought of as a special case of a causal
interpolating filter. The following example demonstrates their use, using the
same data as in the previous example:
.. plot::
:alt: "This code generates an X-Y plot of a time-series with amplitude on the Y axis and time on the X axis. The original time-series is shown as a series of blue markers roughly defining some kind of oscillation. An orange trace showing the nearest neighbor interpolation is drawn atop the original with a stair-like appearance where the original data is right in the middle of each stair step. A green trace showing the previous neighbor interpolation looks similar to the orange trace but the original data is at the back of each stair step. Similarly a dotted red trace showing the next neighbor interpolation goes through each of the previous points, but it is centered at the front edge of each stair."
>>> from scipy.interpolate import interp1d
>>> x = np.linspace(0, 10, num=11, endpoint=True)
>>> y = np.cos(-x**2/9.0)
>>> f1 = interp1d(x, y, kind='nearest')
>>> f2 = interp1d(x, y, kind='previous')
>>> f3 = interp1d(x, y, kind='next')
>>> xnew = np.linspace(0, 10, num=1001, endpoint=True)
>>> import matplotlib.pyplot as plt
>>> plt.plot(x, y, 'o')
>>> plt.plot(xnew, f1(xnew), '-', xnew, f2(xnew), '--', xnew, f3(xnew), ':')
>>> plt.legend(['data', 'nearest', 'previous', 'next'], loc='best')
>>> plt.show()
.. :caption: One-dimensional interpolation using the
.. class :obj:`interpolate.interp1d` with
.. kind equals `nearest`, `previous`, and
.. `next`.
.. _tutorial-interpolate_interp1d_replacements:
Recommended replacements for `interp1d` modes
---------------------------------------------
As mentioned, `interp1d` class is *legacy*: we have no plans to remove it; we
are going to keep supporting its existing usages; however
we believe there are better alternatives which we recommend using in new code.
Here we list specific recommendations, depending on the interpolation ``kind``.
**Linear interpolation,** ``kind="linear"``
The default recommendation is to use `numpy.interp` function. Alternatively,
you can use linear splines, ``make_interp_spline(x, y, k=1)``, see
:ref:`this section for a discussion<tutorial-interpolate_linear_spline>`.
**Spline interpolators,** ``kind="quadratic"`` or ``"cubic"``
Under the hood, `interp1d` delegates to `make_interp_spline`, so we recommend
using the latter directly.
**Piecewise constant modes,** ``kind="nearest", "previous", "next"``
First, we note that ``interp1d(x, y, kind='previous')`` is equivalent to
``make_interp_spline(x, y, k=0)``.
More generally however, all these piecewise constant interpolation modes are
based on `numpy.searchsorted`. For example, the ``"nearest"`` mode is nothing but
>>> x = np.arange(8)
>>> y = x**2
>>> x_new = np.linspace(0, 7, 101) # input points
>>> x_bds = x[:-1] / 2.0 + x[1:] / 2.0 # halfway points
>>> idx = np.searchsorted(x_bds, x_new, side='left')
>>> idx = np.clip(idx, 0, len(x) - 1) # clip the indices so that they are within the range of x indices.
>>> import matplotlib.pyplot as plt
>>> plt.plot(x, y, 'o')
>>> plt.plot(x_new, y[idx], '--')
>>> plt.show()
Other variants are similar, see the ``interp1d`` `source code <https://github.com/scipy/scipy/blob/v1.14.1/scipy/interpolate/_interpolate.py#L486>`_
for details.
|