_cubic.py 33 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864
  1. """Interpolation algorithms using piecewise cubic polynomials."""
  2. import numpy as np
  3. from . import PPoly
  4. from ._polyint import _isscalar
  5. from scipy.linalg import solve_banded, solve
  6. __all__ = ["CubicHermiteSpline", "PchipInterpolator", "pchip_interpolate",
  7. "Akima1DInterpolator", "CubicSpline"]
  8. def prepare_input(x, y, axis, dydx=None):
  9. """Prepare input for cubic spline interpolators.
  10. All data are converted to numpy arrays and checked for correctness.
  11. Axes equal to `axis` of arrays `y` and `dydx` are moved to be the 0th
  12. axis. The value of `axis` is converted to lie in
  13. [0, number of dimensions of `y`).
  14. """
  15. x, y = map(np.asarray, (x, y))
  16. if np.issubdtype(x.dtype, np.complexfloating):
  17. raise ValueError("`x` must contain real values.")
  18. x = x.astype(float)
  19. if np.issubdtype(y.dtype, np.complexfloating):
  20. dtype = complex
  21. else:
  22. dtype = float
  23. if dydx is not None:
  24. dydx = np.asarray(dydx)
  25. if y.shape != dydx.shape:
  26. raise ValueError("The shapes of `y` and `dydx` must be identical.")
  27. if np.issubdtype(dydx.dtype, np.complexfloating):
  28. dtype = complex
  29. dydx = dydx.astype(dtype, copy=False)
  30. y = y.astype(dtype, copy=False)
  31. axis = axis % y.ndim
  32. if x.ndim != 1:
  33. raise ValueError("`x` must be 1-dimensional.")
  34. if x.shape[0] < 2:
  35. raise ValueError("`x` must contain at least 2 elements.")
  36. if x.shape[0] != y.shape[axis]:
  37. raise ValueError("The length of `y` along `axis`={0} doesn't "
  38. "match the length of `x`".format(axis))
  39. if not np.all(np.isfinite(x)):
  40. raise ValueError("`x` must contain only finite values.")
  41. if not np.all(np.isfinite(y)):
  42. raise ValueError("`y` must contain only finite values.")
  43. if dydx is not None and not np.all(np.isfinite(dydx)):
  44. raise ValueError("`dydx` must contain only finite values.")
  45. dx = np.diff(x)
  46. if np.any(dx <= 0):
  47. raise ValueError("`x` must be strictly increasing sequence.")
  48. y = np.moveaxis(y, axis, 0)
  49. if dydx is not None:
  50. dydx = np.moveaxis(dydx, axis, 0)
  51. return x, dx, y, axis, dydx
  52. class CubicHermiteSpline(PPoly):
  53. """Piecewise-cubic interpolator matching values and first derivatives.
  54. The result is represented as a `PPoly` instance.
  55. Parameters
  56. ----------
  57. x : array_like, shape (n,)
  58. 1-D array containing values of the independent variable.
  59. Values must be real, finite and in strictly increasing order.
  60. y : array_like
  61. Array containing values of the dependent variable. It can have
  62. arbitrary number of dimensions, but the length along ``axis``
  63. (see below) must match the length of ``x``. Values must be finite.
  64. dydx : array_like
  65. Array containing derivatives of the dependent variable. It can have
  66. arbitrary number of dimensions, but the length along ``axis``
  67. (see below) must match the length of ``x``. Values must be finite.
  68. axis : int, optional
  69. Axis along which `y` is assumed to be varying. Meaning that for
  70. ``x[i]`` the corresponding values are ``np.take(y, i, axis=axis)``.
  71. Default is 0.
  72. extrapolate : {bool, 'periodic', None}, optional
  73. If bool, determines whether to extrapolate to out-of-bounds points
  74. based on first and last intervals, or to return NaNs. If 'periodic',
  75. periodic extrapolation is used. If None (default), it is set to True.
  76. Attributes
  77. ----------
  78. x : ndarray, shape (n,)
  79. Breakpoints. The same ``x`` which was passed to the constructor.
  80. c : ndarray, shape (4, n-1, ...)
  81. Coefficients of the polynomials on each segment. The trailing
  82. dimensions match the dimensions of `y`, excluding ``axis``.
  83. For example, if `y` is 1-D, then ``c[k, i]`` is a coefficient for
  84. ``(x-x[i])**(3-k)`` on the segment between ``x[i]`` and ``x[i+1]``.
  85. axis : int
  86. Interpolation axis. The same axis which was passed to the
  87. constructor.
  88. Methods
  89. -------
  90. __call__
  91. derivative
  92. antiderivative
  93. integrate
  94. roots
  95. See Also
  96. --------
  97. Akima1DInterpolator : Akima 1D interpolator.
  98. PchipInterpolator : PCHIP 1-D monotonic cubic interpolator.
  99. CubicSpline : Cubic spline data interpolator.
  100. PPoly : Piecewise polynomial in terms of coefficients and breakpoints
  101. Notes
  102. -----
  103. If you want to create a higher-order spline matching higher-order
  104. derivatives, use `BPoly.from_derivatives`.
  105. References
  106. ----------
  107. .. [1] `Cubic Hermite spline
  108. <https://en.wikipedia.org/wiki/Cubic_Hermite_spline>`_
  109. on Wikipedia.
  110. """
  111. def __init__(self, x, y, dydx, axis=0, extrapolate=None):
  112. if extrapolate is None:
  113. extrapolate = True
  114. x, dx, y, axis, dydx = prepare_input(x, y, axis, dydx)
  115. dxr = dx.reshape([dx.shape[0]] + [1] * (y.ndim - 1))
  116. slope = np.diff(y, axis=0) / dxr
  117. t = (dydx[:-1] + dydx[1:] - 2 * slope) / dxr
  118. c = np.empty((4, len(x) - 1) + y.shape[1:], dtype=t.dtype)
  119. c[0] = t / dxr
  120. c[1] = (slope - dydx[:-1]) / dxr - t
  121. c[2] = dydx[:-1]
  122. c[3] = y[:-1]
  123. super().__init__(c, x, extrapolate=extrapolate)
  124. self.axis = axis
  125. class PchipInterpolator(CubicHermiteSpline):
  126. r"""PCHIP 1-D monotonic cubic interpolation.
  127. ``x`` and ``y`` are arrays of values used to approximate some function f,
  128. with ``y = f(x)``. The interpolant uses monotonic cubic splines
  129. to find the value of new points. (PCHIP stands for Piecewise Cubic
  130. Hermite Interpolating Polynomial).
  131. Parameters
  132. ----------
  133. x : ndarray
  134. A 1-D array of monotonically increasing real values. ``x`` cannot
  135. include duplicate values (otherwise f is overspecified)
  136. y : ndarray
  137. A 1-D array of real values. ``y``'s length along the interpolation
  138. axis must be equal to the length of ``x``. If N-D array, use ``axis``
  139. parameter to select correct axis.
  140. axis : int, optional
  141. Axis in the y array corresponding to the x-coordinate values.
  142. extrapolate : bool, optional
  143. Whether to extrapolate to out-of-bounds points based on first
  144. and last intervals, or to return NaNs.
  145. Methods
  146. -------
  147. __call__
  148. derivative
  149. antiderivative
  150. roots
  151. See Also
  152. --------
  153. CubicHermiteSpline : Piecewise-cubic interpolator.
  154. Akima1DInterpolator : Akima 1D interpolator.
  155. CubicSpline : Cubic spline data interpolator.
  156. PPoly : Piecewise polynomial in terms of coefficients and breakpoints.
  157. Notes
  158. -----
  159. The interpolator preserves monotonicity in the interpolation data and does
  160. not overshoot if the data is not smooth.
  161. The first derivatives are guaranteed to be continuous, but the second
  162. derivatives may jump at :math:`x_k`.
  163. Determines the derivatives at the points :math:`x_k`, :math:`f'_k`,
  164. by using PCHIP algorithm [1]_.
  165. Let :math:`h_k = x_{k+1} - x_k`, and :math:`d_k = (y_{k+1} - y_k) / h_k`
  166. are the slopes at internal points :math:`x_k`.
  167. If the signs of :math:`d_k` and :math:`d_{k-1}` are different or either of
  168. them equals zero, then :math:`f'_k = 0`. Otherwise, it is given by the
  169. weighted harmonic mean
  170. .. math::
  171. \frac{w_1 + w_2}{f'_k} = \frac{w_1}{d_{k-1}} + \frac{w_2}{d_k}
  172. where :math:`w_1 = 2 h_k + h_{k-1}` and :math:`w_2 = h_k + 2 h_{k-1}`.
  173. The end slopes are set using a one-sided scheme [2]_.
  174. References
  175. ----------
  176. .. [1] F. N. Fritsch and J. Butland,
  177. A method for constructing local
  178. monotone piecewise cubic interpolants,
  179. SIAM J. Sci. Comput., 5(2), 300-304 (1984).
  180. :doi:`10.1137/0905021`.
  181. .. [2] see, e.g., C. Moler, Numerical Computing with Matlab, 2004.
  182. :doi:`10.1137/1.9780898717952`
  183. """
  184. def __init__(self, x, y, axis=0, extrapolate=None):
  185. x, _, y, axis, _ = prepare_input(x, y, axis)
  186. xp = x.reshape((x.shape[0],) + (1,)*(y.ndim-1))
  187. dk = self._find_derivatives(xp, y)
  188. super().__init__(x, y, dk, axis=0, extrapolate=extrapolate)
  189. self.axis = axis
  190. @staticmethod
  191. def _edge_case(h0, h1, m0, m1):
  192. # one-sided three-point estimate for the derivative
  193. d = ((2*h0 + h1)*m0 - h0*m1) / (h0 + h1)
  194. # try to preserve shape
  195. mask = np.sign(d) != np.sign(m0)
  196. mask2 = (np.sign(m0) != np.sign(m1)) & (np.abs(d) > 3.*np.abs(m0))
  197. mmm = (~mask) & mask2
  198. d[mask] = 0.
  199. d[mmm] = 3.*m0[mmm]
  200. return d
  201. @staticmethod
  202. def _find_derivatives(x, y):
  203. # Determine the derivatives at the points y_k, d_k, by using
  204. # PCHIP algorithm is:
  205. # We choose the derivatives at the point x_k by
  206. # Let m_k be the slope of the kth segment (between k and k+1)
  207. # If m_k=0 or m_{k-1}=0 or sgn(m_k) != sgn(m_{k-1}) then d_k == 0
  208. # else use weighted harmonic mean:
  209. # w_1 = 2h_k + h_{k-1}, w_2 = h_k + 2h_{k-1}
  210. # 1/d_k = 1/(w_1 + w_2)*(w_1 / m_k + w_2 / m_{k-1})
  211. # where h_k is the spacing between x_k and x_{k+1}
  212. y_shape = y.shape
  213. if y.ndim == 1:
  214. # So that _edge_case doesn't end up assigning to scalars
  215. x = x[:, None]
  216. y = y[:, None]
  217. hk = x[1:] - x[:-1]
  218. mk = (y[1:] - y[:-1]) / hk
  219. if y.shape[0] == 2:
  220. # edge case: only have two points, use linear interpolation
  221. dk = np.zeros_like(y)
  222. dk[0] = mk
  223. dk[1] = mk
  224. return dk.reshape(y_shape)
  225. smk = np.sign(mk)
  226. condition = (smk[1:] != smk[:-1]) | (mk[1:] == 0) | (mk[:-1] == 0)
  227. w1 = 2*hk[1:] + hk[:-1]
  228. w2 = hk[1:] + 2*hk[:-1]
  229. # values where division by zero occurs will be excluded
  230. # by 'condition' afterwards
  231. with np.errstate(divide='ignore', invalid='ignore'):
  232. whmean = (w1/mk[:-1] + w2/mk[1:]) / (w1 + w2)
  233. dk = np.zeros_like(y)
  234. dk[1:-1][condition] = 0.0
  235. dk[1:-1][~condition] = 1.0 / whmean[~condition]
  236. # special case endpoints, as suggested in
  237. # Cleve Moler, Numerical Computing with MATLAB, Chap 3.6 (pchiptx.m)
  238. dk[0] = PchipInterpolator._edge_case(hk[0], hk[1], mk[0], mk[1])
  239. dk[-1] = PchipInterpolator._edge_case(hk[-1], hk[-2], mk[-1], mk[-2])
  240. return dk.reshape(y_shape)
  241. def pchip_interpolate(xi, yi, x, der=0, axis=0):
  242. """
  243. Convenience function for pchip interpolation.
  244. xi and yi are arrays of values used to approximate some function f,
  245. with ``yi = f(xi)``. The interpolant uses monotonic cubic splines
  246. to find the value of new points x and the derivatives there.
  247. See `scipy.interpolate.PchipInterpolator` for details.
  248. Parameters
  249. ----------
  250. xi : array_like
  251. A sorted list of x-coordinates, of length N.
  252. yi : array_like
  253. A 1-D array of real values. `yi`'s length along the interpolation
  254. axis must be equal to the length of `xi`. If N-D array, use axis
  255. parameter to select correct axis.
  256. x : scalar or array_like
  257. Of length M.
  258. der : int or list, optional
  259. Derivatives to extract. The 0th derivative can be included to
  260. return the function value.
  261. axis : int, optional
  262. Axis in the yi array corresponding to the x-coordinate values.
  263. See Also
  264. --------
  265. PchipInterpolator : PCHIP 1-D monotonic cubic interpolator.
  266. Returns
  267. -------
  268. y : scalar or array_like
  269. The result, of length R or length M or M by R,
  270. Examples
  271. --------
  272. We can interpolate 2D observed data using pchip interpolation:
  273. >>> import numpy as np
  274. >>> import matplotlib.pyplot as plt
  275. >>> from scipy.interpolate import pchip_interpolate
  276. >>> x_observed = np.linspace(0.0, 10.0, 11)
  277. >>> y_observed = np.sin(x_observed)
  278. >>> x = np.linspace(min(x_observed), max(x_observed), num=100)
  279. >>> y = pchip_interpolate(x_observed, y_observed, x)
  280. >>> plt.plot(x_observed, y_observed, "o", label="observation")
  281. >>> plt.plot(x, y, label="pchip interpolation")
  282. >>> plt.legend()
  283. >>> plt.show()
  284. """
  285. P = PchipInterpolator(xi, yi, axis=axis)
  286. if der == 0:
  287. return P(x)
  288. elif _isscalar(der):
  289. return P.derivative(der)(x)
  290. else:
  291. return [P.derivative(nu)(x) for nu in der]
  292. class Akima1DInterpolator(CubicHermiteSpline):
  293. """
  294. Akima interpolator
  295. Fit piecewise cubic polynomials, given vectors x and y. The interpolation
  296. method by Akima uses a continuously differentiable sub-spline built from
  297. piecewise cubic polynomials. The resultant curve passes through the given
  298. data points and will appear smooth and natural.
  299. Parameters
  300. ----------
  301. x : ndarray, shape (m, )
  302. 1-D array of monotonically increasing real values.
  303. y : ndarray, shape (m, ...)
  304. N-D array of real values. The length of ``y`` along the first axis
  305. must be equal to the length of ``x``.
  306. axis : int, optional
  307. Specifies the axis of ``y`` along which to interpolate. Interpolation
  308. defaults to the first axis of ``y``.
  309. Methods
  310. -------
  311. __call__
  312. derivative
  313. antiderivative
  314. roots
  315. See Also
  316. --------
  317. PchipInterpolator : PCHIP 1-D monotonic cubic interpolator.
  318. CubicSpline : Cubic spline data interpolator.
  319. PPoly : Piecewise polynomial in terms of coefficients and breakpoints
  320. Notes
  321. -----
  322. .. versionadded:: 0.14
  323. Use only for precise data, as the fitted curve passes through the given
  324. points exactly. This routine is useful for plotting a pleasingly smooth
  325. curve through a few given points for purposes of plotting.
  326. References
  327. ----------
  328. [1] A new method of interpolation and smooth curve fitting based
  329. on local procedures. Hiroshi Akima, J. ACM, October 1970, 17(4),
  330. 589-602.
  331. """
  332. def __init__(self, x, y, axis=0):
  333. # Original implementation in MATLAB by N. Shamsundar (BSD licensed), see
  334. # https://www.mathworks.com/matlabcentral/fileexchange/1814-akima-interpolation
  335. x, dx, y, axis, _ = prepare_input(x, y, axis)
  336. # determine slopes between breakpoints
  337. m = np.empty((x.size + 3, ) + y.shape[1:])
  338. dx = dx[(slice(None), ) + (None, ) * (y.ndim - 1)]
  339. m[2:-2] = np.diff(y, axis=0) / dx
  340. # add two additional points on the left ...
  341. m[1] = 2. * m[2] - m[3]
  342. m[0] = 2. * m[1] - m[2]
  343. # ... and on the right
  344. m[-2] = 2. * m[-3] - m[-4]
  345. m[-1] = 2. * m[-2] - m[-3]
  346. # if m1 == m2 != m3 == m4, the slope at the breakpoint is not
  347. # defined. This is the fill value:
  348. t = .5 * (m[3:] + m[:-3])
  349. # get the denominator of the slope t
  350. dm = np.abs(np.diff(m, axis=0))
  351. f1 = dm[2:]
  352. f2 = dm[:-2]
  353. f12 = f1 + f2
  354. # These are the mask of where the slope at breakpoint is defined:
  355. ind = np.nonzero(f12 > 1e-9 * np.max(f12, initial=-np.inf))
  356. x_ind, y_ind = ind[0], ind[1:]
  357. # Set the slope at breakpoint
  358. t[ind] = (f1[ind] * m[(x_ind + 1,) + y_ind] +
  359. f2[ind] * m[(x_ind + 2,) + y_ind]) / f12[ind]
  360. super().__init__(x, y, t, axis=0, extrapolate=False)
  361. self.axis = axis
  362. def extend(self, c, x, right=True):
  363. raise NotImplementedError("Extending a 1-D Akima interpolator is not "
  364. "yet implemented")
  365. # These are inherited from PPoly, but they do not produce an Akima
  366. # interpolator. Hence stub them out.
  367. @classmethod
  368. def from_spline(cls, tck, extrapolate=None):
  369. raise NotImplementedError("This method does not make sense for "
  370. "an Akima interpolator.")
  371. @classmethod
  372. def from_bernstein_basis(cls, bp, extrapolate=None):
  373. raise NotImplementedError("This method does not make sense for "
  374. "an Akima interpolator.")
  375. class CubicSpline(CubicHermiteSpline):
  376. """Cubic spline data interpolator.
  377. Interpolate data with a piecewise cubic polynomial which is twice
  378. continuously differentiable [1]_. The result is represented as a `PPoly`
  379. instance with breakpoints matching the given data.
  380. Parameters
  381. ----------
  382. x : array_like, shape (n,)
  383. 1-D array containing values of the independent variable.
  384. Values must be real, finite and in strictly increasing order.
  385. y : array_like
  386. Array containing values of the dependent variable. It can have
  387. arbitrary number of dimensions, but the length along ``axis``
  388. (see below) must match the length of ``x``. Values must be finite.
  389. axis : int, optional
  390. Axis along which `y` is assumed to be varying. Meaning that for
  391. ``x[i]`` the corresponding values are ``np.take(y, i, axis=axis)``.
  392. Default is 0.
  393. bc_type : string or 2-tuple, optional
  394. Boundary condition type. Two additional equations, given by the
  395. boundary conditions, are required to determine all coefficients of
  396. polynomials on each segment [2]_.
  397. If `bc_type` is a string, then the specified condition will be applied
  398. at both ends of a spline. Available conditions are:
  399. * 'not-a-knot' (default): The first and second segment at a curve end
  400. are the same polynomial. It is a good default when there is no
  401. information on boundary conditions.
  402. * 'periodic': The interpolated functions is assumed to be periodic
  403. of period ``x[-1] - x[0]``. The first and last value of `y` must be
  404. identical: ``y[0] == y[-1]``. This boundary condition will result in
  405. ``y'[0] == y'[-1]`` and ``y''[0] == y''[-1]``.
  406. * 'clamped': The first derivative at curves ends are zero. Assuming
  407. a 1D `y`, ``bc_type=((1, 0.0), (1, 0.0))`` is the same condition.
  408. * 'natural': The second derivative at curve ends are zero. Assuming
  409. a 1D `y`, ``bc_type=((2, 0.0), (2, 0.0))`` is the same condition.
  410. If `bc_type` is a 2-tuple, the first and the second value will be
  411. applied at the curve start and end respectively. The tuple values can
  412. be one of the previously mentioned strings (except 'periodic') or a
  413. tuple `(order, deriv_values)` allowing to specify arbitrary
  414. derivatives at curve ends:
  415. * `order`: the derivative order, 1 or 2.
  416. * `deriv_value`: array_like containing derivative values, shape must
  417. be the same as `y`, excluding ``axis`` dimension. For example, if
  418. `y` is 1-D, then `deriv_value` must be a scalar. If `y` is 3-D with
  419. the shape (n0, n1, n2) and axis=2, then `deriv_value` must be 2-D
  420. and have the shape (n0, n1).
  421. extrapolate : {bool, 'periodic', None}, optional
  422. If bool, determines whether to extrapolate to out-of-bounds points
  423. based on first and last intervals, or to return NaNs. If 'periodic',
  424. periodic extrapolation is used. If None (default), ``extrapolate`` is
  425. set to 'periodic' for ``bc_type='periodic'`` and to True otherwise.
  426. Attributes
  427. ----------
  428. x : ndarray, shape (n,)
  429. Breakpoints. The same ``x`` which was passed to the constructor.
  430. c : ndarray, shape (4, n-1, ...)
  431. Coefficients of the polynomials on each segment. The trailing
  432. dimensions match the dimensions of `y`, excluding ``axis``.
  433. For example, if `y` is 1-d, then ``c[k, i]`` is a coefficient for
  434. ``(x-x[i])**(3-k)`` on the segment between ``x[i]`` and ``x[i+1]``.
  435. axis : int
  436. Interpolation axis. The same axis which was passed to the
  437. constructor.
  438. Methods
  439. -------
  440. __call__
  441. derivative
  442. antiderivative
  443. integrate
  444. roots
  445. See Also
  446. --------
  447. Akima1DInterpolator : Akima 1D interpolator.
  448. PchipInterpolator : PCHIP 1-D monotonic cubic interpolator.
  449. PPoly : Piecewise polynomial in terms of coefficients and breakpoints.
  450. Notes
  451. -----
  452. Parameters `bc_type` and ``extrapolate`` work independently, i.e. the
  453. former controls only construction of a spline, and the latter only
  454. evaluation.
  455. When a boundary condition is 'not-a-knot' and n = 2, it is replaced by
  456. a condition that the first derivative is equal to the linear interpolant
  457. slope. When both boundary conditions are 'not-a-knot' and n = 3, the
  458. solution is sought as a parabola passing through given points.
  459. When 'not-a-knot' boundary conditions is applied to both ends, the
  460. resulting spline will be the same as returned by `splrep` (with ``s=0``)
  461. and `InterpolatedUnivariateSpline`, but these two methods use a
  462. representation in B-spline basis.
  463. .. versionadded:: 0.18.0
  464. Examples
  465. --------
  466. In this example the cubic spline is used to interpolate a sampled sinusoid.
  467. You can see that the spline continuity property holds for the first and
  468. second derivatives and violates only for the third derivative.
  469. >>> import numpy as np
  470. >>> from scipy.interpolate import CubicSpline
  471. >>> import matplotlib.pyplot as plt
  472. >>> x = np.arange(10)
  473. >>> y = np.sin(x)
  474. >>> cs = CubicSpline(x, y)
  475. >>> xs = np.arange(-0.5, 9.6, 0.1)
  476. >>> fig, ax = plt.subplots(figsize=(6.5, 4))
  477. >>> ax.plot(x, y, 'o', label='data')
  478. >>> ax.plot(xs, np.sin(xs), label='true')
  479. >>> ax.plot(xs, cs(xs), label="S")
  480. >>> ax.plot(xs, cs(xs, 1), label="S'")
  481. >>> ax.plot(xs, cs(xs, 2), label="S''")
  482. >>> ax.plot(xs, cs(xs, 3), label="S'''")
  483. >>> ax.set_xlim(-0.5, 9.5)
  484. >>> ax.legend(loc='lower left', ncol=2)
  485. >>> plt.show()
  486. In the second example, the unit circle is interpolated with a spline. A
  487. periodic boundary condition is used. You can see that the first derivative
  488. values, ds/dx=0, ds/dy=1 at the periodic point (1, 0) are correctly
  489. computed. Note that a circle cannot be exactly represented by a cubic
  490. spline. To increase precision, more breakpoints would be required.
  491. >>> theta = 2 * np.pi * np.linspace(0, 1, 5)
  492. >>> y = np.c_[np.cos(theta), np.sin(theta)]
  493. >>> cs = CubicSpline(theta, y, bc_type='periodic')
  494. >>> print("ds/dx={:.1f} ds/dy={:.1f}".format(cs(0, 1)[0], cs(0, 1)[1]))
  495. ds/dx=0.0 ds/dy=1.0
  496. >>> xs = 2 * np.pi * np.linspace(0, 1, 100)
  497. >>> fig, ax = plt.subplots(figsize=(6.5, 4))
  498. >>> ax.plot(y[:, 0], y[:, 1], 'o', label='data')
  499. >>> ax.plot(np.cos(xs), np.sin(xs), label='true')
  500. >>> ax.plot(cs(xs)[:, 0], cs(xs)[:, 1], label='spline')
  501. >>> ax.axes.set_aspect('equal')
  502. >>> ax.legend(loc='center')
  503. >>> plt.show()
  504. The third example is the interpolation of a polynomial y = x**3 on the
  505. interval 0 <= x<= 1. A cubic spline can represent this function exactly.
  506. To achieve that we need to specify values and first derivatives at
  507. endpoints of the interval. Note that y' = 3 * x**2 and thus y'(0) = 0 and
  508. y'(1) = 3.
  509. >>> cs = CubicSpline([0, 1], [0, 1], bc_type=((1, 0), (1, 3)))
  510. >>> x = np.linspace(0, 1)
  511. >>> np.allclose(x**3, cs(x))
  512. True
  513. References
  514. ----------
  515. .. [1] `Cubic Spline Interpolation
  516. <https://en.wikiversity.org/wiki/Cubic_Spline_Interpolation>`_
  517. on Wikiversity.
  518. .. [2] Carl de Boor, "A Practical Guide to Splines", Springer-Verlag, 1978.
  519. """
  520. def __init__(self, x, y, axis=0, bc_type='not-a-knot', extrapolate=None):
  521. x, dx, y, axis, _ = prepare_input(x, y, axis)
  522. n = len(x)
  523. bc, y = self._validate_bc(bc_type, y, y.shape[1:], axis)
  524. if extrapolate is None:
  525. if bc[0] == 'periodic':
  526. extrapolate = 'periodic'
  527. else:
  528. extrapolate = True
  529. if y.size == 0:
  530. # bail out early for zero-sized arrays
  531. s = np.zeros_like(y)
  532. else:
  533. dxr = dx.reshape([dx.shape[0]] + [1] * (y.ndim - 1))
  534. slope = np.diff(y, axis=0) / dxr
  535. # If bc is 'not-a-knot' this change is just a convention.
  536. # If bc is 'periodic' then we already checked that y[0] == y[-1],
  537. # and the spline is just a constant, we handle this case in the
  538. # same way by setting the first derivatives to slope, which is 0.
  539. if n == 2:
  540. if bc[0] in ['not-a-knot', 'periodic']:
  541. bc[0] = (1, slope[0])
  542. if bc[1] in ['not-a-knot', 'periodic']:
  543. bc[1] = (1, slope[0])
  544. # This is a special case, when both conditions are 'not-a-knot'
  545. # and n == 3. In this case 'not-a-knot' can't be handled regularly
  546. # as the both conditions are identical. We handle this case by
  547. # constructing a parabola passing through given points.
  548. if n == 3 and bc[0] == 'not-a-knot' and bc[1] == 'not-a-knot':
  549. A = np.zeros((3, 3)) # This is a standard matrix.
  550. b = np.empty((3,) + y.shape[1:], dtype=y.dtype)
  551. A[0, 0] = 1
  552. A[0, 1] = 1
  553. A[1, 0] = dx[1]
  554. A[1, 1] = 2 * (dx[0] + dx[1])
  555. A[1, 2] = dx[0]
  556. A[2, 1] = 1
  557. A[2, 2] = 1
  558. b[0] = 2 * slope[0]
  559. b[1] = 3 * (dxr[0] * slope[1] + dxr[1] * slope[0])
  560. b[2] = 2 * slope[1]
  561. s = solve(A, b, overwrite_a=True, overwrite_b=True,
  562. check_finite=False)
  563. elif n == 3 and bc[0] == 'periodic':
  564. # In case when number of points is 3 we compute the derivatives
  565. # manually
  566. s = np.empty((n,) + y.shape[1:], dtype=y.dtype)
  567. t = (slope / dxr).sum() / (1. / dxr).sum()
  568. s.fill(t)
  569. else:
  570. # Find derivative values at each x[i] by solving a tridiagonal
  571. # system.
  572. A = np.zeros((3, n)) # This is a banded matrix representation.
  573. b = np.empty((n,) + y.shape[1:], dtype=y.dtype)
  574. # Filling the system for i=1..n-2
  575. # (x[i-1] - x[i]) * s[i-1] +\
  576. # 2 * ((x[i] - x[i-1]) + (x[i+1] - x[i])) * s[i] +\
  577. # (x[i] - x[i-1]) * s[i+1] =\
  578. # 3 * ((x[i+1] - x[i])*(y[i] - y[i-1])/(x[i] - x[i-1]) +\
  579. # (x[i] - x[i-1])*(y[i+1] - y[i])/(x[i+1] - x[i]))
  580. A[1, 1:-1] = 2 * (dx[:-1] + dx[1:]) # The diagonal
  581. A[0, 2:] = dx[:-1] # The upper diagonal
  582. A[-1, :-2] = dx[1:] # The lower diagonal
  583. b[1:-1] = 3 * (dxr[1:] * slope[:-1] + dxr[:-1] * slope[1:])
  584. bc_start, bc_end = bc
  585. if bc_start == 'periodic':
  586. # Due to the periodicity, and because y[-1] = y[0], the
  587. # linear system has (n-1) unknowns/equations instead of n:
  588. A = A[:, 0:-1]
  589. A[1, 0] = 2 * (dx[-1] + dx[0])
  590. A[0, 1] = dx[-1]
  591. b = b[:-1]
  592. # Also, due to the periodicity, the system is not tri-diagonal.
  593. # We need to compute a "condensed" matrix of shape (n-2, n-2).
  594. # See https://web.archive.org/web/20151220180652/http://www.cfm.brown.edu/people/gk/chap6/node14.html
  595. # for more explanations.
  596. # The condensed matrix is obtained by removing the last column
  597. # and last row of the (n-1, n-1) system matrix. The removed
  598. # values are saved in scalar variables with the (n-1, n-1)
  599. # system matrix indices forming their names:
  600. a_m1_0 = dx[-2] # lower left corner value: A[-1, 0]
  601. a_m1_m2 = dx[-1]
  602. a_m1_m1 = 2 * (dx[-1] + dx[-2])
  603. a_m2_m1 = dx[-3]
  604. a_0_m1 = dx[0]
  605. b[0] = 3 * (dxr[0] * slope[-1] + dxr[-1] * slope[0])
  606. b[-1] = 3 * (dxr[-1] * slope[-2] + dxr[-2] * slope[-1])
  607. Ac = A[:, :-1]
  608. b1 = b[:-1]
  609. b2 = np.zeros_like(b1)
  610. b2[0] = -a_0_m1
  611. b2[-1] = -a_m2_m1
  612. # s1 and s2 are the solutions of (n-2, n-2) system
  613. s1 = solve_banded((1, 1), Ac, b1, overwrite_ab=False,
  614. overwrite_b=False, check_finite=False)
  615. s2 = solve_banded((1, 1), Ac, b2, overwrite_ab=False,
  616. overwrite_b=False, check_finite=False)
  617. # computing the s[n-2] solution:
  618. s_m1 = ((b[-1] - a_m1_0 * s1[0] - a_m1_m2 * s1[-1]) /
  619. (a_m1_m1 + a_m1_0 * s2[0] + a_m1_m2 * s2[-1]))
  620. # s is the solution of the (n, n) system:
  621. s = np.empty((n,) + y.shape[1:], dtype=y.dtype)
  622. s[:-2] = s1 + s_m1 * s2
  623. s[-2] = s_m1
  624. s[-1] = s[0]
  625. else:
  626. if bc_start == 'not-a-knot':
  627. A[1, 0] = dx[1]
  628. A[0, 1] = x[2] - x[0]
  629. d = x[2] - x[0]
  630. b[0] = ((dxr[0] + 2*d) * dxr[1] * slope[0] +
  631. dxr[0]**2 * slope[1]) / d
  632. elif bc_start[0] == 1:
  633. A[1, 0] = 1
  634. A[0, 1] = 0
  635. b[0] = bc_start[1]
  636. elif bc_start[0] == 2:
  637. A[1, 0] = 2 * dx[0]
  638. A[0, 1] = dx[0]
  639. b[0] = -0.5 * bc_start[1] * dx[0]**2 + 3 * (y[1] - y[0])
  640. if bc_end == 'not-a-knot':
  641. A[1, -1] = dx[-2]
  642. A[-1, -2] = x[-1] - x[-3]
  643. d = x[-1] - x[-3]
  644. b[-1] = ((dxr[-1]**2*slope[-2] +
  645. (2*d + dxr[-1])*dxr[-2]*slope[-1]) / d)
  646. elif bc_end[0] == 1:
  647. A[1, -1] = 1
  648. A[-1, -2] = 0
  649. b[-1] = bc_end[1]
  650. elif bc_end[0] == 2:
  651. A[1, -1] = 2 * dx[-1]
  652. A[-1, -2] = dx[-1]
  653. b[-1] = 0.5 * bc_end[1] * dx[-1]**2 + 3 * (y[-1] - y[-2])
  654. s = solve_banded((1, 1), A, b, overwrite_ab=True,
  655. overwrite_b=True, check_finite=False)
  656. super().__init__(x, y, s, axis=0, extrapolate=extrapolate)
  657. self.axis = axis
  658. @staticmethod
  659. def _validate_bc(bc_type, y, expected_deriv_shape, axis):
  660. """Validate and prepare boundary conditions.
  661. Returns
  662. -------
  663. validated_bc : 2-tuple
  664. Boundary conditions for a curve start and end.
  665. y : ndarray
  666. y casted to complex dtype if one of the boundary conditions has
  667. complex dtype.
  668. """
  669. if isinstance(bc_type, str):
  670. if bc_type == 'periodic':
  671. if not np.allclose(y[0], y[-1], rtol=1e-15, atol=1e-15):
  672. raise ValueError(
  673. "The first and last `y` point along axis {} must "
  674. "be identical (within machine precision) when "
  675. "bc_type='periodic'.".format(axis))
  676. bc_type = (bc_type, bc_type)
  677. else:
  678. if len(bc_type) != 2:
  679. raise ValueError("`bc_type` must contain 2 elements to "
  680. "specify start and end conditions.")
  681. if 'periodic' in bc_type:
  682. raise ValueError("'periodic' `bc_type` is defined for both "
  683. "curve ends and cannot be used with other "
  684. "boundary conditions.")
  685. validated_bc = []
  686. for bc in bc_type:
  687. if isinstance(bc, str):
  688. if bc == 'clamped':
  689. validated_bc.append((1, np.zeros(expected_deriv_shape)))
  690. elif bc == 'natural':
  691. validated_bc.append((2, np.zeros(expected_deriv_shape)))
  692. elif bc in ['not-a-knot', 'periodic']:
  693. validated_bc.append(bc)
  694. else:
  695. raise ValueError("bc_type={} is not allowed.".format(bc))
  696. else:
  697. try:
  698. deriv_order, deriv_value = bc
  699. except Exception as e:
  700. raise ValueError(
  701. "A specified derivative value must be "
  702. "given in the form (order, value)."
  703. ) from e
  704. if deriv_order not in [1, 2]:
  705. raise ValueError("The specified derivative order must "
  706. "be 1 or 2.")
  707. deriv_value = np.asarray(deriv_value)
  708. if deriv_value.shape != expected_deriv_shape:
  709. raise ValueError(
  710. "`deriv_value` shape {} is not the expected one {}."
  711. .format(deriv_value.shape, expected_deriv_shape))
  712. if np.issubdtype(deriv_value.dtype, np.complexfloating):
  713. y = y.astype(complex, copy=False)
  714. validated_bc.append((deriv_order, deriv_value))
  715. return validated_bc, y