_bsplines.py 34 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023
  1. from __future__ import division, print_function, absolute_import
  2. import functools
  3. import operator
  4. import numpy as np
  5. from scipy._lib.six import string_types
  6. from scipy.linalg import (get_lapack_funcs, LinAlgError,
  7. cholesky_banded, cho_solve_banded)
  8. from . import _bspl
  9. from . import _fitpack_impl
  10. from . import _fitpack as _dierckx
  11. __all__ = ["BSpline", "make_interp_spline", "make_lsq_spline"]
  12. # copy-paste from interpolate.py
  13. def prod(x):
  14. """Product of a list of numbers; ~40x faster vs np.prod for Python tuples"""
  15. if len(x) == 0:
  16. return 1
  17. return functools.reduce(operator.mul, x)
  18. def _get_dtype(dtype):
  19. """Return np.complex128 for complex dtypes, np.float64 otherwise."""
  20. if np.issubdtype(dtype, np.complexfloating):
  21. return np.complex_
  22. else:
  23. return np.float_
  24. def _as_float_array(x, check_finite=False):
  25. """Convert the input into a C contiguous float array.
  26. NB: Upcasts half- and single-precision floats to double precision.
  27. """
  28. x = np.ascontiguousarray(x)
  29. dtyp = _get_dtype(x.dtype)
  30. x = x.astype(dtyp, copy=False)
  31. if check_finite and not np.isfinite(x).all():
  32. raise ValueError("Array must not contain infs or nans.")
  33. return x
  34. class BSpline(object):
  35. r"""Univariate spline in the B-spline basis.
  36. .. math::
  37. S(x) = \sum_{j=0}^{n-1} c_j B_{j, k; t}(x)
  38. where :math:`B_{j, k; t}` are B-spline basis functions of degree `k`
  39. and knots `t`.
  40. Parameters
  41. ----------
  42. t : ndarray, shape (n+k+1,)
  43. knots
  44. c : ndarray, shape (>=n, ...)
  45. spline coefficients
  46. k : int
  47. B-spline order
  48. extrapolate : bool or 'periodic', optional
  49. whether to extrapolate beyond the base interval, ``t[k] .. t[n]``,
  50. or to return nans.
  51. If True, extrapolates the first and last polynomial pieces of b-spline
  52. functions active on the base interval.
  53. If 'periodic', periodic extrapolation is used.
  54. Default is True.
  55. axis : int, optional
  56. Interpolation axis. Default is zero.
  57. Attributes
  58. ----------
  59. t : ndarray
  60. knot vector
  61. c : ndarray
  62. spline coefficients
  63. k : int
  64. spline degree
  65. extrapolate : bool
  66. If True, extrapolates the first and last polynomial pieces of b-spline
  67. functions active on the base interval.
  68. axis : int
  69. Interpolation axis.
  70. tck : tuple
  71. A read-only equivalent of ``(self.t, self.c, self.k)``
  72. Methods
  73. -------
  74. __call__
  75. basis_element
  76. derivative
  77. antiderivative
  78. integrate
  79. construct_fast
  80. Notes
  81. -----
  82. B-spline basis elements are defined via
  83. .. math::
  84. B_{i, 0}(x) = 1, \textrm{if $t_i \le x < t_{i+1}$, otherwise $0$,}
  85. B_{i, k}(x) = \frac{x - t_i}{t_{i+k} - t_i} B_{i, k-1}(x)
  86. + \frac{t_{i+k+1} - x}{t_{i+k+1} - t_{i+1}} B_{i+1, k-1}(x)
  87. **Implementation details**
  88. - At least ``k+1`` coefficients are required for a spline of degree `k`,
  89. so that ``n >= k+1``. Additional coefficients, ``c[j]`` with
  90. ``j > n``, are ignored.
  91. - B-spline basis elements of degree `k` form a partition of unity on the
  92. *base interval*, ``t[k] <= x <= t[n]``.
  93. Examples
  94. --------
  95. Translating the recursive definition of B-splines into Python code, we have:
  96. >>> def B(x, k, i, t):
  97. ... if k == 0:
  98. ... return 1.0 if t[i] <= x < t[i+1] else 0.0
  99. ... if t[i+k] == t[i]:
  100. ... c1 = 0.0
  101. ... else:
  102. ... c1 = (x - t[i])/(t[i+k] - t[i]) * B(x, k-1, i, t)
  103. ... if t[i+k+1] == t[i+1]:
  104. ... c2 = 0.0
  105. ... else:
  106. ... c2 = (t[i+k+1] - x)/(t[i+k+1] - t[i+1]) * B(x, k-1, i+1, t)
  107. ... return c1 + c2
  108. >>> def bspline(x, t, c, k):
  109. ... n = len(t) - k - 1
  110. ... assert (n >= k+1) and (len(c) >= n)
  111. ... return sum(c[i] * B(x, k, i, t) for i in range(n))
  112. Note that this is an inefficient (if straightforward) way to
  113. evaluate B-splines --- this spline class does it in an equivalent,
  114. but much more efficient way.
  115. Here we construct a quadratic spline function on the base interval
  116. ``2 <= x <= 4`` and compare with the naive way of evaluating the spline:
  117. >>> from scipy.interpolate import BSpline
  118. >>> k = 2
  119. >>> t = [0, 1, 2, 3, 4, 5, 6]
  120. >>> c = [-1, 2, 0, -1]
  121. >>> spl = BSpline(t, c, k)
  122. >>> spl(2.5)
  123. array(1.375)
  124. >>> bspline(2.5, t, c, k)
  125. 1.375
  126. Note that outside of the base interval results differ. This is because
  127. `BSpline` extrapolates the first and last polynomial pieces of b-spline
  128. functions active on the base interval.
  129. >>> import matplotlib.pyplot as plt
  130. >>> fig, ax = plt.subplots()
  131. >>> xx = np.linspace(1.5, 4.5, 50)
  132. >>> ax.plot(xx, [bspline(x, t, c ,k) for x in xx], 'r-', lw=3, label='naive')
  133. >>> ax.plot(xx, spl(xx), 'b-', lw=4, alpha=0.7, label='BSpline')
  134. >>> ax.grid(True)
  135. >>> ax.legend(loc='best')
  136. >>> plt.show()
  137. References
  138. ----------
  139. .. [1] Tom Lyche and Knut Morken, Spline methods,
  140. http://www.uio.no/studier/emner/matnat/ifi/INF-MAT5340/v05/undervisningsmateriale/
  141. .. [2] Carl de Boor, A practical guide to splines, Springer, 2001.
  142. """
  143. def __init__(self, t, c, k, extrapolate=True, axis=0):
  144. super(BSpline, self).__init__()
  145. self.k = operator.index(k)
  146. self.c = np.asarray(c)
  147. self.t = np.ascontiguousarray(t, dtype=np.float64)
  148. if extrapolate == 'periodic':
  149. self.extrapolate = extrapolate
  150. else:
  151. self.extrapolate = bool(extrapolate)
  152. n = self.t.shape[0] - self.k - 1
  153. if not (0 <= axis < self.c.ndim):
  154. raise ValueError("%s must be between 0 and %s" % (axis, c.ndim))
  155. self.axis = axis
  156. if axis != 0:
  157. # roll the interpolation axis to be the first one in self.c
  158. # More specifically, the target shape for self.c is (n, ...),
  159. # and axis !=0 means that we have c.shape (..., n, ...)
  160. # ^
  161. # axis
  162. self.c = np.rollaxis(self.c, axis)
  163. if k < 0:
  164. raise ValueError("Spline order cannot be negative.")
  165. if self.t.ndim != 1:
  166. raise ValueError("Knot vector must be one-dimensional.")
  167. if n < self.k + 1:
  168. raise ValueError("Need at least %d knots for degree %d" %
  169. (2*k + 2, k))
  170. if (np.diff(self.t) < 0).any():
  171. raise ValueError("Knots must be in a non-decreasing order.")
  172. if len(np.unique(self.t[k:n+1])) < 2:
  173. raise ValueError("Need at least two internal knots.")
  174. if not np.isfinite(self.t).all():
  175. raise ValueError("Knots should not have nans or infs.")
  176. if self.c.ndim < 1:
  177. raise ValueError("Coefficients must be at least 1-dimensional.")
  178. if self.c.shape[0] < n:
  179. raise ValueError("Knots, coefficients and degree are inconsistent.")
  180. dt = _get_dtype(self.c.dtype)
  181. self.c = np.ascontiguousarray(self.c, dtype=dt)
  182. @classmethod
  183. def construct_fast(cls, t, c, k, extrapolate=True, axis=0):
  184. """Construct a spline without making checks.
  185. Accepts same parameters as the regular constructor. Input arrays
  186. `t` and `c` must of correct shape and dtype.
  187. """
  188. self = object.__new__(cls)
  189. self.t, self.c, self.k = t, c, k
  190. self.extrapolate = extrapolate
  191. self.axis = axis
  192. return self
  193. @property
  194. def tck(self):
  195. """Equivalent to ``(self.t, self.c, self.k)`` (read-only).
  196. """
  197. return self.t, self.c, self.k
  198. @classmethod
  199. def basis_element(cls, t, extrapolate=True):
  200. """Return a B-spline basis element ``B(x | t[0], ..., t[k+1])``.
  201. Parameters
  202. ----------
  203. t : ndarray, shape (k+1,)
  204. internal knots
  205. extrapolate : bool or 'periodic', optional
  206. whether to extrapolate beyond the base interval, ``t[0] .. t[k+1]``,
  207. or to return nans.
  208. If 'periodic', periodic extrapolation is used.
  209. Default is True.
  210. Returns
  211. -------
  212. basis_element : callable
  213. A callable representing a B-spline basis element for the knot
  214. vector `t`.
  215. Notes
  216. -----
  217. The order of the b-spline, `k`, is inferred from the length of `t` as
  218. ``len(t)-2``. The knot vector is constructed by appending and prepending
  219. ``k+1`` elements to internal knots `t`.
  220. Examples
  221. --------
  222. Construct a cubic b-spline:
  223. >>> from scipy.interpolate import BSpline
  224. >>> b = BSpline.basis_element([0, 1, 2, 3, 4])
  225. >>> k = b.k
  226. >>> b.t[k:-k]
  227. array([ 0., 1., 2., 3., 4.])
  228. >>> k
  229. 3
  230. Construct a second order b-spline on ``[0, 1, 1, 2]``, and compare
  231. to its explicit form:
  232. >>> t = [-1, 0, 1, 1, 2]
  233. >>> b = BSpline.basis_element(t[1:])
  234. >>> def f(x):
  235. ... return np.where(x < 1, x*x, (2. - x)**2)
  236. >>> import matplotlib.pyplot as plt
  237. >>> fig, ax = plt.subplots()
  238. >>> x = np.linspace(0, 2, 51)
  239. >>> ax.plot(x, b(x), 'g', lw=3)
  240. >>> ax.plot(x, f(x), 'r', lw=8, alpha=0.4)
  241. >>> ax.grid(True)
  242. >>> plt.show()
  243. """
  244. k = len(t) - 2
  245. t = _as_float_array(t)
  246. t = np.r_[(t[0]-1,) * k, t, (t[-1]+1,) * k]
  247. c = np.zeros_like(t)
  248. c[k] = 1.
  249. return cls.construct_fast(t, c, k, extrapolate)
  250. def __call__(self, x, nu=0, extrapolate=None):
  251. """
  252. Evaluate a spline function.
  253. Parameters
  254. ----------
  255. x : array_like
  256. points to evaluate the spline at.
  257. nu: int, optional
  258. derivative to evaluate (default is 0).
  259. extrapolate : bool or 'periodic', optional
  260. whether to extrapolate based on the first and last intervals
  261. or return nans. If 'periodic', periodic extrapolation is used.
  262. Default is `self.extrapolate`.
  263. Returns
  264. -------
  265. y : array_like
  266. Shape is determined by replacing the interpolation axis
  267. in the coefficient array with the shape of `x`.
  268. """
  269. if extrapolate is None:
  270. extrapolate = self.extrapolate
  271. x = np.asarray(x)
  272. x_shape, x_ndim = x.shape, x.ndim
  273. x = np.ascontiguousarray(x.ravel(), dtype=np.float_)
  274. # With periodic extrapolation we map x to the segment
  275. # [self.t[k], self.t[n]].
  276. if extrapolate == 'periodic':
  277. n = self.t.size - self.k - 1
  278. x = self.t[self.k] + (x - self.t[self.k]) % (self.t[n] -
  279. self.t[self.k])
  280. extrapolate = False
  281. out = np.empty((len(x), prod(self.c.shape[1:])), dtype=self.c.dtype)
  282. self._ensure_c_contiguous()
  283. self._evaluate(x, nu, extrapolate, out)
  284. out = out.reshape(x_shape + self.c.shape[1:])
  285. if self.axis != 0:
  286. # transpose to move the calculated values to the interpolation axis
  287. l = list(range(out.ndim))
  288. l = l[x_ndim:x_ndim+self.axis] + l[:x_ndim] + l[x_ndim+self.axis:]
  289. out = out.transpose(l)
  290. return out
  291. def _evaluate(self, xp, nu, extrapolate, out):
  292. _bspl.evaluate_spline(self.t, self.c.reshape(self.c.shape[0], -1),
  293. self.k, xp, nu, extrapolate, out)
  294. def _ensure_c_contiguous(self):
  295. """
  296. c and t may be modified by the user. The Cython code expects
  297. that they are C contiguous.
  298. """
  299. if not self.t.flags.c_contiguous:
  300. self.t = self.t.copy()
  301. if not self.c.flags.c_contiguous:
  302. self.c = self.c.copy()
  303. def derivative(self, nu=1):
  304. """Return a b-spline representing the derivative.
  305. Parameters
  306. ----------
  307. nu : int, optional
  308. Derivative order.
  309. Default is 1.
  310. Returns
  311. -------
  312. b : BSpline object
  313. A new instance representing the derivative.
  314. See Also
  315. --------
  316. splder, splantider
  317. """
  318. c = self.c
  319. # pad the c array if needed
  320. ct = len(self.t) - len(c)
  321. if ct > 0:
  322. c = np.r_[c, np.zeros((ct,) + c.shape[1:])]
  323. tck = _fitpack_impl.splder((self.t, c, self.k), nu)
  324. return self.construct_fast(*tck, extrapolate=self.extrapolate,
  325. axis=self.axis)
  326. def antiderivative(self, nu=1):
  327. """Return a b-spline representing the antiderivative.
  328. Parameters
  329. ----------
  330. nu : int, optional
  331. Antiderivative order. Default is 1.
  332. Returns
  333. -------
  334. b : BSpline object
  335. A new instance representing the antiderivative.
  336. Notes
  337. -----
  338. If antiderivative is computed and ``self.extrapolate='periodic'``,
  339. it will be set to False for the returned instance. This is done because
  340. the antiderivative is no longer periodic and its correct evaluation
  341. outside of the initially given x interval is difficult.
  342. See Also
  343. --------
  344. splder, splantider
  345. """
  346. c = self.c
  347. # pad the c array if needed
  348. ct = len(self.t) - len(c)
  349. if ct > 0:
  350. c = np.r_[c, np.zeros((ct,) + c.shape[1:])]
  351. tck = _fitpack_impl.splantider((self.t, c, self.k), nu)
  352. if self.extrapolate == 'periodic':
  353. extrapolate = False
  354. else:
  355. extrapolate = self.extrapolate
  356. return self.construct_fast(*tck, extrapolate=extrapolate,
  357. axis=self.axis)
  358. def integrate(self, a, b, extrapolate=None):
  359. """Compute a definite integral of the spline.
  360. Parameters
  361. ----------
  362. a : float
  363. Lower limit of integration.
  364. b : float
  365. Upper limit of integration.
  366. extrapolate : bool or 'periodic', optional
  367. whether to extrapolate beyond the base interval,
  368. ``t[k] .. t[-k-1]``, or take the spline to be zero outside of the
  369. base interval. If 'periodic', periodic extrapolation is used.
  370. If None (default), use `self.extrapolate`.
  371. Returns
  372. -------
  373. I : array_like
  374. Definite integral of the spline over the interval ``[a, b]``.
  375. Examples
  376. --------
  377. Construct the linear spline ``x if x < 1 else 2 - x`` on the base
  378. interval :math:`[0, 2]`, and integrate it
  379. >>> from scipy.interpolate import BSpline
  380. >>> b = BSpline.basis_element([0, 1, 2])
  381. >>> b.integrate(0, 1)
  382. array(0.5)
  383. If the integration limits are outside of the base interval, the result
  384. is controlled by the `extrapolate` parameter
  385. >>> b.integrate(-1, 1)
  386. array(0.0)
  387. >>> b.integrate(-1, 1, extrapolate=False)
  388. array(0.5)
  389. >>> import matplotlib.pyplot as plt
  390. >>> fig, ax = plt.subplots()
  391. >>> ax.grid(True)
  392. >>> ax.axvline(0, c='r', lw=5, alpha=0.5) # base interval
  393. >>> ax.axvline(2, c='r', lw=5, alpha=0.5)
  394. >>> xx = [-1, 1, 2]
  395. >>> ax.plot(xx, b(xx))
  396. >>> plt.show()
  397. """
  398. if extrapolate is None:
  399. extrapolate = self.extrapolate
  400. # Prepare self.t and self.c.
  401. self._ensure_c_contiguous()
  402. # Swap integration bounds if needed.
  403. sign = 1
  404. if b < a:
  405. a, b = b, a
  406. sign = -1
  407. n = self.t.size - self.k - 1
  408. if extrapolate != "periodic" and not extrapolate:
  409. # Shrink the integration interval, if needed.
  410. a = max(a, self.t[self.k])
  411. b = min(b, self.t[n])
  412. if self.c.ndim == 1:
  413. # Fast path: use FITPACK's routine
  414. # (cf _fitpack_impl.splint).
  415. t, c, k = self.tck
  416. integral, wrk = _dierckx._splint(t, c, k, a, b)
  417. return integral * sign
  418. out = np.empty((2, prod(self.c.shape[1:])), dtype=self.c.dtype)
  419. # Compute the antiderivative.
  420. c = self.c
  421. ct = len(self.t) - len(c)
  422. if ct > 0:
  423. c = np.r_[c, np.zeros((ct,) + c.shape[1:])]
  424. ta, ca, ka = _fitpack_impl.splantider((self.t, c, self.k), 1)
  425. if extrapolate == 'periodic':
  426. # Split the integral into the part over period (can be several
  427. # of them) and the remaining part.
  428. ts, te = self.t[self.k], self.t[n]
  429. period = te - ts
  430. interval = b - a
  431. n_periods, left = divmod(interval, period)
  432. if n_periods > 0:
  433. # Evaluate the difference of antiderivatives.
  434. x = np.asarray([ts, te], dtype=np.float_)
  435. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  436. ka, x, 0, False, out)
  437. integral = out[1] - out[0]
  438. integral *= n_periods
  439. else:
  440. integral = np.zeros((1, prod(self.c.shape[1:])),
  441. dtype=self.c.dtype)
  442. # Map a to [ts, te], b is always a + left.
  443. a = ts + (a - ts) % period
  444. b = a + left
  445. # If b <= te then we need to integrate over [a, b], otherwise
  446. # over [a, te] and from xs to what is remained.
  447. if b <= te:
  448. x = np.asarray([a, b], dtype=np.float_)
  449. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  450. ka, x, 0, False, out)
  451. integral += out[1] - out[0]
  452. else:
  453. x = np.asarray([a, te], dtype=np.float_)
  454. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  455. ka, x, 0, False, out)
  456. integral += out[1] - out[0]
  457. x = np.asarray([ts, ts + b - te], dtype=np.float_)
  458. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  459. ka, x, 0, False, out)
  460. integral += out[1] - out[0]
  461. else:
  462. # Evaluate the difference of antiderivatives.
  463. x = np.asarray([a, b], dtype=np.float_)
  464. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  465. ka, x, 0, extrapolate, out)
  466. integral = out[1] - out[0]
  467. integral *= sign
  468. return integral.reshape(ca.shape[1:])
  469. #################################
  470. # Interpolating spline helpers #
  471. #################################
  472. def _not_a_knot(x, k):
  473. """Given data x, construct the knot vector w/ not-a-knot BC.
  474. cf de Boor, XIII(12)."""
  475. x = np.asarray(x)
  476. if k % 2 != 1:
  477. raise ValueError("Odd degree for now only. Got %s." % k)
  478. m = (k - 1) // 2
  479. t = x[m+1:-m-1]
  480. t = np.r_[(x[0],)*(k+1), t, (x[-1],)*(k+1)]
  481. return t
  482. def _augknt(x, k):
  483. """Construct a knot vector appropriate for the order-k interpolation."""
  484. return np.r_[(x[0],)*k, x, (x[-1],)*k]
  485. def _convert_string_aliases(deriv, target_shape):
  486. if isinstance(deriv, string_types):
  487. if deriv == "clamped":
  488. deriv = [(1, np.zeros(target_shape))]
  489. elif deriv == "natural":
  490. deriv = [(2, np.zeros(target_shape))]
  491. else:
  492. raise ValueError("Unknown boundary condition : %s" % deriv)
  493. return deriv
  494. def _process_deriv_spec(deriv):
  495. if deriv is not None:
  496. try:
  497. ords, vals = zip(*deriv)
  498. except TypeError:
  499. msg = ("Derivatives, `bc_type`, should be specified as a pair of "
  500. "iterables of pairs of (order, value).")
  501. raise ValueError(msg)
  502. else:
  503. ords, vals = [], []
  504. return np.atleast_1d(ords, vals)
  505. def make_interp_spline(x, y, k=3, t=None, bc_type=None, axis=0,
  506. check_finite=True):
  507. """Compute the (coefficients of) interpolating B-spline.
  508. Parameters
  509. ----------
  510. x : array_like, shape (n,)
  511. Abscissas.
  512. y : array_like, shape (n, ...)
  513. Ordinates.
  514. k : int, optional
  515. B-spline degree. Default is cubic, k=3.
  516. t : array_like, shape (nt + k + 1,), optional.
  517. Knots.
  518. The number of knots needs to agree with the number of datapoints and
  519. the number of derivatives at the edges. Specifically, ``nt - n`` must
  520. equal ``len(deriv_l) + len(deriv_r)``.
  521. bc_type : 2-tuple or None
  522. Boundary conditions.
  523. Default is None, which means choosing the boundary conditions
  524. automatically. Otherwise, it must be a length-two tuple where the first
  525. element sets the boundary conditions at ``x[0]`` and the second
  526. element sets the boundary conditions at ``x[-1]``. Each of these must
  527. be an iterable of pairs ``(order, value)`` which gives the values of
  528. derivatives of specified orders at the given edge of the interpolation
  529. interval.
  530. Alternatively, the following string aliases are recognized:
  531. * ``"clamped"``: The first derivatives at the ends are zero. This is
  532. equivalent to ``bc_type=([(1, 0.0)], [(1, 0.0)])``.
  533. * ``"natural"``: The second derivatives at ends are zero. This is
  534. equivalent to ``bc_type=([(2, 0.0)], [(2, 0.0)])``.
  535. * ``"not-a-knot"`` (default): The first and second segments are the same
  536. polynomial. This is equivalent to having ``bc_type=None``.
  537. axis : int, optional
  538. Interpolation axis. Default is 0.
  539. check_finite : bool, optional
  540. Whether to check that the input arrays contain only finite numbers.
  541. Disabling may give a performance gain, but may result in problems
  542. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  543. Default is True.
  544. Returns
  545. -------
  546. b : a BSpline object of the degree ``k`` and with knots ``t``.
  547. Examples
  548. --------
  549. Use cubic interpolation on Chebyshev nodes:
  550. >>> def cheb_nodes(N):
  551. ... jj = 2.*np.arange(N) + 1
  552. ... x = np.cos(np.pi * jj / 2 / N)[::-1]
  553. ... return x
  554. >>> x = cheb_nodes(20)
  555. >>> y = np.sqrt(1 - x**2)
  556. >>> from scipy.interpolate import BSpline, make_interp_spline
  557. >>> b = make_interp_spline(x, y)
  558. >>> np.allclose(b(x), y)
  559. True
  560. Note that the default is a cubic spline with a not-a-knot boundary condition
  561. >>> b.k
  562. 3
  563. Here we use a 'natural' spline, with zero 2nd derivatives at edges:
  564. >>> l, r = [(2, 0.0)], [(2, 0.0)]
  565. >>> b_n = make_interp_spline(x, y, bc_type=(l, r)) # or, bc_type="natural"
  566. >>> np.allclose(b_n(x), y)
  567. True
  568. >>> x0, x1 = x[0], x[-1]
  569. >>> np.allclose([b_n(x0, 2), b_n(x1, 2)], [0, 0])
  570. True
  571. Interpolation of parametric curves is also supported. As an example, we
  572. compute a discretization of a snail curve in polar coordinates
  573. >>> phi = np.linspace(0, 2.*np.pi, 40)
  574. >>> r = 0.3 + np.cos(phi)
  575. >>> x, y = r*np.cos(phi), r*np.sin(phi) # convert to Cartesian coordinates
  576. Build an interpolating curve, parameterizing it by the angle
  577. >>> from scipy.interpolate import make_interp_spline
  578. >>> spl = make_interp_spline(phi, np.c_[x, y])
  579. Evaluate the interpolant on a finer grid (note that we transpose the result
  580. to unpack it into a pair of x- and y-arrays)
  581. >>> phi_new = np.linspace(0, 2.*np.pi, 100)
  582. >>> x_new, y_new = spl(phi_new).T
  583. Plot the result
  584. >>> import matplotlib.pyplot as plt
  585. >>> plt.plot(x, y, 'o')
  586. >>> plt.plot(x_new, y_new, '-')
  587. >>> plt.show()
  588. See Also
  589. --------
  590. BSpline : base class representing the B-spline objects
  591. CubicSpline : a cubic spline in the polynomial basis
  592. make_lsq_spline : a similar factory function for spline fitting
  593. UnivariateSpline : a wrapper over FITPACK spline fitting routines
  594. splrep : a wrapper over FITPACK spline fitting routines
  595. """
  596. # convert string aliases for the boundary conditions
  597. if bc_type is None or bc_type == 'not-a-knot':
  598. deriv_l, deriv_r = None, None
  599. elif isinstance(bc_type, string_types):
  600. deriv_l, deriv_r = bc_type, bc_type
  601. else:
  602. try:
  603. deriv_l, deriv_r = bc_type
  604. except TypeError:
  605. raise ValueError("Unknown boundary condition: %s" % bc_type)
  606. y = np.asarray(y)
  607. if not -y.ndim <= axis < y.ndim:
  608. raise ValueError("axis {} is out of bounds".format(axis))
  609. if axis < 0:
  610. axis += y.ndim
  611. # special-case k=0 right away
  612. if k == 0:
  613. if any(_ is not None for _ in (t, deriv_l, deriv_r)):
  614. raise ValueError("Too much info for k=0: t and bc_type can only "
  615. "be None.")
  616. x = _as_float_array(x, check_finite)
  617. t = np.r_[x, x[-1]]
  618. c = np.asarray(y)
  619. c = np.rollaxis(c, axis)
  620. c = np.ascontiguousarray(c, dtype=_get_dtype(c.dtype))
  621. return BSpline.construct_fast(t, c, k, axis=axis)
  622. # special-case k=1 (e.g., Lyche and Morken, Eq.(2.16))
  623. if k == 1 and t is None:
  624. if not (deriv_l is None and deriv_r is None):
  625. raise ValueError("Too much info for k=1: bc_type can only be None.")
  626. x = _as_float_array(x, check_finite)
  627. t = np.r_[x[0], x, x[-1]]
  628. c = np.asarray(y)
  629. c = np.rollaxis(c, axis)
  630. c = np.ascontiguousarray(c, dtype=_get_dtype(c.dtype))
  631. return BSpline.construct_fast(t, c, k, axis=axis)
  632. x = _as_float_array(x, check_finite)
  633. y = _as_float_array(y, check_finite)
  634. k = operator.index(k)
  635. # come up with a sensible knot vector, if needed
  636. if t is None:
  637. if deriv_l is None and deriv_r is None:
  638. if k == 2:
  639. # OK, it's a bit ad hoc: Greville sites + omit
  640. # 2nd and 2nd-to-last points, a la not-a-knot
  641. t = (x[1:] + x[:-1]) / 2.
  642. t = np.r_[(x[0],)*(k+1),
  643. t[1:-1],
  644. (x[-1],)*(k+1)]
  645. else:
  646. t = _not_a_knot(x, k)
  647. else:
  648. t = _augknt(x, k)
  649. t = _as_float_array(t, check_finite)
  650. y = np.rollaxis(y, axis) # now internally interp axis is zero
  651. if x.ndim != 1 or np.any(x[1:] <= x[:-1]):
  652. raise ValueError("Expect x to be a 1-D sorted array_like.")
  653. if k < 0:
  654. raise ValueError("Expect non-negative k.")
  655. if t.ndim != 1 or np.any(t[1:] < t[:-1]):
  656. raise ValueError("Expect t to be a 1-D sorted array_like.")
  657. if x.size != y.shape[0]:
  658. raise ValueError('x and y are incompatible.')
  659. if t.size < x.size + k + 1:
  660. raise ValueError('Got %d knots, need at least %d.' %
  661. (t.size, x.size + k + 1))
  662. if (x[0] < t[k]) or (x[-1] > t[-k]):
  663. raise ValueError('Out of bounds w/ x = %s.' % x)
  664. # Here : deriv_l, r = [(nu, value), ...]
  665. deriv_l = _convert_string_aliases(deriv_l, y.shape[1:])
  666. deriv_l_ords, deriv_l_vals = _process_deriv_spec(deriv_l)
  667. nleft = deriv_l_ords.shape[0]
  668. deriv_r = _convert_string_aliases(deriv_r, y.shape[1:])
  669. deriv_r_ords, deriv_r_vals = _process_deriv_spec(deriv_r)
  670. nright = deriv_r_ords.shape[0]
  671. # have `n` conditions for `nt` coefficients; need nt-n derivatives
  672. n = x.size
  673. nt = t.size - k - 1
  674. if nt - n != nleft + nright:
  675. raise ValueError("The number of derivatives at boundaries does not "
  676. "match: expected %s, got %s+%s" % (nt-n, nleft, nright))
  677. # set up the LHS: the collocation matrix + derivatives at boundaries
  678. kl = ku = k
  679. ab = np.zeros((2*kl + ku + 1, nt), dtype=np.float_, order='F')
  680. _bspl._colloc(x, t, k, ab, offset=nleft)
  681. if nleft > 0:
  682. _bspl._handle_lhs_derivatives(t, k, x[0], ab, kl, ku, deriv_l_ords)
  683. if nright > 0:
  684. _bspl._handle_lhs_derivatives(t, k, x[-1], ab, kl, ku, deriv_r_ords,
  685. offset=nt-nright)
  686. # set up the RHS: values to interpolate (+ derivative values, if any)
  687. extradim = prod(y.shape[1:])
  688. rhs = np.empty((nt, extradim), dtype=y.dtype)
  689. if nleft > 0:
  690. rhs[:nleft] = deriv_l_vals.reshape(-1, extradim)
  691. rhs[nleft:nt - nright] = y.reshape(-1, extradim)
  692. if nright > 0:
  693. rhs[nt - nright:] = deriv_r_vals.reshape(-1, extradim)
  694. # solve Ab @ x = rhs; this is the relevant part of linalg.solve_banded
  695. if check_finite:
  696. ab, rhs = map(np.asarray_chkfinite, (ab, rhs))
  697. gbsv, = get_lapack_funcs(('gbsv',), (ab, rhs))
  698. lu, piv, c, info = gbsv(kl, ku, ab, rhs,
  699. overwrite_ab=True, overwrite_b=True)
  700. if info > 0:
  701. raise LinAlgError("Collocation matix is singular.")
  702. elif info < 0:
  703. raise ValueError('illegal value in %d-th argument of internal gbsv' % -info)
  704. c = np.ascontiguousarray(c.reshape((nt,) + y.shape[1:]))
  705. return BSpline.construct_fast(t, c, k, axis=axis)
  706. def make_lsq_spline(x, y, t, k=3, w=None, axis=0, check_finite=True):
  707. r"""Compute the (coefficients of) an LSQ B-spline.
  708. The result is a linear combination
  709. .. math::
  710. S(x) = \sum_j c_j B_j(x; t)
  711. of the B-spline basis elements, :math:`B_j(x; t)`, which minimizes
  712. .. math::
  713. \sum_{j} \left( w_j \times (S(x_j) - y_j) \right)^2
  714. Parameters
  715. ----------
  716. x : array_like, shape (m,)
  717. Abscissas.
  718. y : array_like, shape (m, ...)
  719. Ordinates.
  720. t : array_like, shape (n + k + 1,).
  721. Knots.
  722. Knots and data points must satisfy Schoenberg-Whitney conditions.
  723. k : int, optional
  724. B-spline degree. Default is cubic, k=3.
  725. w : array_like, shape (n,), optional
  726. Weights for spline fitting. Must be positive. If ``None``,
  727. then weights are all equal.
  728. Default is ``None``.
  729. axis : int, optional
  730. Interpolation axis. Default is zero.
  731. check_finite : bool, optional
  732. Whether to check that the input arrays contain only finite numbers.
  733. Disabling may give a performance gain, but may result in problems
  734. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  735. Default is True.
  736. Returns
  737. -------
  738. b : a BSpline object of the degree `k` with knots `t`.
  739. Notes
  740. -----
  741. The number of data points must be larger than the spline degree `k`.
  742. Knots `t` must satisfy the Schoenberg-Whitney conditions,
  743. i.e., there must be a subset of data points ``x[j]`` such that
  744. ``t[j] < x[j] < t[j+k+1]``, for ``j=0, 1,...,n-k-2``.
  745. Examples
  746. --------
  747. Generate some noisy data:
  748. >>> x = np.linspace(-3, 3, 50)
  749. >>> y = np.exp(-x**2) + 0.1 * np.random.randn(50)
  750. Now fit a smoothing cubic spline with a pre-defined internal knots.
  751. Here we make the knot vector (k+1)-regular by adding boundary knots:
  752. >>> from scipy.interpolate import make_lsq_spline, BSpline
  753. >>> t = [-1, 0, 1]
  754. >>> k = 3
  755. >>> t = np.r_[(x[0],)*(k+1),
  756. ... t,
  757. ... (x[-1],)*(k+1)]
  758. >>> spl = make_lsq_spline(x, y, t, k)
  759. For comparison, we also construct an interpolating spline for the same
  760. set of data:
  761. >>> from scipy.interpolate import make_interp_spline
  762. >>> spl_i = make_interp_spline(x, y)
  763. Plot both:
  764. >>> import matplotlib.pyplot as plt
  765. >>> xs = np.linspace(-3, 3, 100)
  766. >>> plt.plot(x, y, 'ro', ms=5)
  767. >>> plt.plot(xs, spl(xs), 'g-', lw=3, label='LSQ spline')
  768. >>> plt.plot(xs, spl_i(xs), 'b-', lw=3, alpha=0.7, label='interp spline')
  769. >>> plt.legend(loc='best')
  770. >>> plt.show()
  771. **NaN handling**: If the input arrays contain ``nan`` values, the result is
  772. not useful since the underlying spline fitting routines cannot deal with
  773. ``nan``. A workaround is to use zero weights for not-a-number data points:
  774. >>> y[8] = np.nan
  775. >>> w = np.isnan(y)
  776. >>> y[w] = 0.
  777. >>> tck = make_lsq_spline(x, y, t, w=~w)
  778. Notice the need to replace a ``nan`` by a numerical value (precise value
  779. does not matter as long as the corresponding weight is zero.)
  780. See Also
  781. --------
  782. BSpline : base class representing the B-spline objects
  783. make_interp_spline : a similar factory function for interpolating splines
  784. LSQUnivariateSpline : a FITPACK-based spline fitting routine
  785. splrep : a FITPACK-based fitting routine
  786. """
  787. x = _as_float_array(x, check_finite)
  788. y = _as_float_array(y, check_finite)
  789. t = _as_float_array(t, check_finite)
  790. if w is not None:
  791. w = _as_float_array(w, check_finite)
  792. else:
  793. w = np.ones_like(x)
  794. k = operator.index(k)
  795. if not -y.ndim <= axis < y.ndim:
  796. raise ValueError("axis {} is out of bounds".format(axis))
  797. if axis < 0:
  798. axis += y.ndim
  799. y = np.rollaxis(y, axis) # now internally interp axis is zero
  800. if x.ndim != 1 or np.any(x[1:] - x[:-1] <= 0):
  801. raise ValueError("Expect x to be a 1-D sorted array_like.")
  802. if x.shape[0] < k+1:
  803. raise ValueError("Need more x points.")
  804. if k < 0:
  805. raise ValueError("Expect non-negative k.")
  806. if t.ndim != 1 or np.any(t[1:] - t[:-1] < 0):
  807. raise ValueError("Expect t to be a 1-D sorted array_like.")
  808. if x.size != y.shape[0]:
  809. raise ValueError('x & y are incompatible.')
  810. if k > 0 and np.any((x < t[k]) | (x > t[-k])):
  811. raise ValueError('Out of bounds w/ x = %s.' % x)
  812. if x.size != w.size:
  813. raise ValueError('Incompatible weights.')
  814. # number of coefficients
  815. n = t.size - k - 1
  816. # construct A.T @ A and rhs with A the collocation matrix, and
  817. # rhs = A.T @ y for solving the LSQ problem ``A.T @ A @ c = A.T @ y``
  818. lower = True
  819. extradim = prod(y.shape[1:])
  820. ab = np.zeros((k+1, n), dtype=np.float_, order='F')
  821. rhs = np.zeros((n, extradim), dtype=y.dtype, order='F')
  822. _bspl._norm_eq_lsq(x, t, k,
  823. y.reshape(-1, extradim),
  824. w,
  825. ab, rhs)
  826. rhs = rhs.reshape((n,) + y.shape[1:])
  827. # have observation matrix & rhs, can solve the LSQ problem
  828. cho_decomp = cholesky_banded(ab, overwrite_ab=True, lower=lower,
  829. check_finite=check_finite)
  830. c = cho_solve_banded((cho_decomp, lower), rhs, overwrite_b=True,
  831. check_finite=check_finite)
  832. c = np.ascontiguousarray(c)
  833. return BSpline.construct_fast(t, c, k, axis=axis)