_bsplines.py 68 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871187218731874187518761877187818791880188118821883188418851886188718881889189018911892189318941895189618971898189919001901190219031904190519061907190819091910191119121913191419151916191719181919192019211922192319241925192619271928192919301931193219331934193519361937193819391940194119421943194419451946194719481949195019511952195319541955195619571958195919601961196219631964196519661967196819691970197119721973197419751976197719781979198019811982198319841985198619871988198919901991199219931994199519961997199819992000200120022003200420052006200720082009201020112012201320142015201620172018201920202021202220232024202520262027202820292030
  1. import operator
  2. import numpy as np
  3. from numpy.core.multiarray import normalize_axis_index
  4. from scipy.linalg import (get_lapack_funcs, LinAlgError,
  5. cholesky_banded, cho_solve_banded,
  6. solve, solve_banded)
  7. from scipy.optimize import minimize_scalar
  8. from . import _bspl
  9. from . import _fitpack_impl
  10. from scipy._lib._util import prod
  11. from scipy.sparse import csr_array
  12. from scipy.special import poch
  13. from itertools import combinations
  14. __all__ = ["BSpline", "make_interp_spline", "make_lsq_spline",
  15. "make_smoothing_spline"]
  16. def _get_dtype(dtype):
  17. """Return np.complex128 for complex dtypes, np.float64 otherwise."""
  18. if np.issubdtype(dtype, np.complexfloating):
  19. return np.complex_
  20. else:
  21. return np.float_
  22. def _as_float_array(x, check_finite=False):
  23. """Convert the input into a C contiguous float array.
  24. NB: Upcasts half- and single-precision floats to double precision.
  25. """
  26. x = np.ascontiguousarray(x)
  27. dtyp = _get_dtype(x.dtype)
  28. x = x.astype(dtyp, copy=False)
  29. if check_finite and not np.isfinite(x).all():
  30. raise ValueError("Array must not contain infs or nans.")
  31. return x
  32. def _dual_poly(j, k, t, y):
  33. """
  34. Dual polynomial of the B-spline B_{j,k,t} -
  35. polynomial which is associated with B_{j,k,t}:
  36. $p_{j,k}(y) = (y - t_{j+1})(y - t_{j+2})...(y - t_{j+k})$
  37. """
  38. if k == 0:
  39. return 1
  40. return np.prod([(y - t[j + i]) for i in range(1, k + 1)])
  41. def _diff_dual_poly(j, k, y, d, t):
  42. """
  43. d-th derivative of the dual polynomial $p_{j,k}(y)$
  44. """
  45. if d == 0:
  46. return _dual_poly(j, k, t, y)
  47. if d == k:
  48. return poch(1, k)
  49. comb = list(combinations(range(j + 1, j + k + 1), d))
  50. res = 0
  51. for i in range(len(comb) * len(comb[0])):
  52. res += np.prod([(y - t[j + p]) for p in range(1, k + 1)
  53. if (j + p) not in comb[i//d]])
  54. return res
  55. class BSpline:
  56. r"""Univariate spline in the B-spline basis.
  57. .. math::
  58. S(x) = \sum_{j=0}^{n-1} c_j B_{j, k; t}(x)
  59. where :math:`B_{j, k; t}` are B-spline basis functions of degree `k`
  60. and knots `t`.
  61. Parameters
  62. ----------
  63. t : ndarray, shape (n+k+1,)
  64. knots
  65. c : ndarray, shape (>=n, ...)
  66. spline coefficients
  67. k : int
  68. B-spline degree
  69. extrapolate : bool or 'periodic', optional
  70. whether to extrapolate beyond the base interval, ``t[k] .. t[n]``,
  71. or to return nans.
  72. If True, extrapolates the first and last polynomial pieces of b-spline
  73. functions active on the base interval.
  74. If 'periodic', periodic extrapolation is used.
  75. Default is True.
  76. axis : int, optional
  77. Interpolation axis. Default is zero.
  78. Attributes
  79. ----------
  80. t : ndarray
  81. knot vector
  82. c : ndarray
  83. spline coefficients
  84. k : int
  85. spline degree
  86. extrapolate : bool
  87. If True, extrapolates the first and last polynomial pieces of b-spline
  88. functions active on the base interval.
  89. axis : int
  90. Interpolation axis.
  91. tck : tuple
  92. A read-only equivalent of ``(self.t, self.c, self.k)``
  93. Methods
  94. -------
  95. __call__
  96. basis_element
  97. derivative
  98. antiderivative
  99. integrate
  100. construct_fast
  101. design_matrix
  102. from_power_basis
  103. Notes
  104. -----
  105. B-spline basis elements are defined via
  106. .. math::
  107. B_{i, 0}(x) = 1, \textrm{if $t_i \le x < t_{i+1}$, otherwise $0$,}
  108. B_{i, k}(x) = \frac{x - t_i}{t_{i+k} - t_i} B_{i, k-1}(x)
  109. + \frac{t_{i+k+1} - x}{t_{i+k+1} - t_{i+1}} B_{i+1, k-1}(x)
  110. **Implementation details**
  111. - At least ``k+1`` coefficients are required for a spline of degree `k`,
  112. so that ``n >= k+1``. Additional coefficients, ``c[j]`` with
  113. ``j > n``, are ignored.
  114. - B-spline basis elements of degree `k` form a partition of unity on the
  115. *base interval*, ``t[k] <= x <= t[n]``.
  116. Examples
  117. --------
  118. Translating the recursive definition of B-splines into Python code, we have:
  119. >>> def B(x, k, i, t):
  120. ... if k == 0:
  121. ... return 1.0 if t[i] <= x < t[i+1] else 0.0
  122. ... if t[i+k] == t[i]:
  123. ... c1 = 0.0
  124. ... else:
  125. ... c1 = (x - t[i])/(t[i+k] - t[i]) * B(x, k-1, i, t)
  126. ... if t[i+k+1] == t[i+1]:
  127. ... c2 = 0.0
  128. ... else:
  129. ... c2 = (t[i+k+1] - x)/(t[i+k+1] - t[i+1]) * B(x, k-1, i+1, t)
  130. ... return c1 + c2
  131. >>> def bspline(x, t, c, k):
  132. ... n = len(t) - k - 1
  133. ... assert (n >= k+1) and (len(c) >= n)
  134. ... return sum(c[i] * B(x, k, i, t) for i in range(n))
  135. Note that this is an inefficient (if straightforward) way to
  136. evaluate B-splines --- this spline class does it in an equivalent,
  137. but much more efficient way.
  138. Here we construct a quadratic spline function on the base interval
  139. ``2 <= x <= 4`` and compare with the naive way of evaluating the spline:
  140. >>> from scipy.interpolate import BSpline
  141. >>> k = 2
  142. >>> t = [0, 1, 2, 3, 4, 5, 6]
  143. >>> c = [-1, 2, 0, -1]
  144. >>> spl = BSpline(t, c, k)
  145. >>> spl(2.5)
  146. array(1.375)
  147. >>> bspline(2.5, t, c, k)
  148. 1.375
  149. Note that outside of the base interval results differ. This is because
  150. `BSpline` extrapolates the first and last polynomial pieces of B-spline
  151. functions active on the base interval.
  152. >>> import matplotlib.pyplot as plt
  153. >>> import numpy as np
  154. >>> fig, ax = plt.subplots()
  155. >>> xx = np.linspace(1.5, 4.5, 50)
  156. >>> ax.plot(xx, [bspline(x, t, c ,k) for x in xx], 'r-', lw=3, label='naive')
  157. >>> ax.plot(xx, spl(xx), 'b-', lw=4, alpha=0.7, label='BSpline')
  158. >>> ax.grid(True)
  159. >>> ax.legend(loc='best')
  160. >>> plt.show()
  161. References
  162. ----------
  163. .. [1] Tom Lyche and Knut Morken, Spline methods,
  164. http://www.uio.no/studier/emner/matnat/ifi/INF-MAT5340/v05/undervisningsmateriale/
  165. .. [2] Carl de Boor, A practical guide to splines, Springer, 2001.
  166. """
  167. def __init__(self, t, c, k, extrapolate=True, axis=0):
  168. super().__init__()
  169. self.k = operator.index(k)
  170. self.c = np.asarray(c)
  171. self.t = np.ascontiguousarray(t, dtype=np.float64)
  172. if extrapolate == 'periodic':
  173. self.extrapolate = extrapolate
  174. else:
  175. self.extrapolate = bool(extrapolate)
  176. n = self.t.shape[0] - self.k - 1
  177. axis = normalize_axis_index(axis, self.c.ndim)
  178. # Note that the normalized axis is stored in the object.
  179. self.axis = axis
  180. if axis != 0:
  181. # roll the interpolation axis to be the first one in self.c
  182. # More specifically, the target shape for self.c is (n, ...),
  183. # and axis !=0 means that we have c.shape (..., n, ...)
  184. # ^
  185. # axis
  186. self.c = np.moveaxis(self.c, axis, 0)
  187. if k < 0:
  188. raise ValueError("Spline order cannot be negative.")
  189. if self.t.ndim != 1:
  190. raise ValueError("Knot vector must be one-dimensional.")
  191. if n < self.k + 1:
  192. raise ValueError("Need at least %d knots for degree %d" %
  193. (2*k + 2, k))
  194. if (np.diff(self.t) < 0).any():
  195. raise ValueError("Knots must be in a non-decreasing order.")
  196. if len(np.unique(self.t[k:n+1])) < 2:
  197. raise ValueError("Need at least two internal knots.")
  198. if not np.isfinite(self.t).all():
  199. raise ValueError("Knots should not have nans or infs.")
  200. if self.c.ndim < 1:
  201. raise ValueError("Coefficients must be at least 1-dimensional.")
  202. if self.c.shape[0] < n:
  203. raise ValueError("Knots, coefficients and degree are inconsistent.")
  204. dt = _get_dtype(self.c.dtype)
  205. self.c = np.ascontiguousarray(self.c, dtype=dt)
  206. @classmethod
  207. def construct_fast(cls, t, c, k, extrapolate=True, axis=0):
  208. """Construct a spline without making checks.
  209. Accepts same parameters as the regular constructor. Input arrays
  210. `t` and `c` must of correct shape and dtype.
  211. """
  212. self = object.__new__(cls)
  213. self.t, self.c, self.k = t, c, k
  214. self.extrapolate = extrapolate
  215. self.axis = axis
  216. return self
  217. @property
  218. def tck(self):
  219. """Equivalent to ``(self.t, self.c, self.k)`` (read-only).
  220. """
  221. return self.t, self.c, self.k
  222. @classmethod
  223. def basis_element(cls, t, extrapolate=True):
  224. """Return a B-spline basis element ``B(x | t[0], ..., t[k+1])``.
  225. Parameters
  226. ----------
  227. t : ndarray, shape (k+2,)
  228. internal knots
  229. extrapolate : bool or 'periodic', optional
  230. whether to extrapolate beyond the base interval, ``t[0] .. t[k+1]``,
  231. or to return nans.
  232. If 'periodic', periodic extrapolation is used.
  233. Default is True.
  234. Returns
  235. -------
  236. basis_element : callable
  237. A callable representing a B-spline basis element for the knot
  238. vector `t`.
  239. Notes
  240. -----
  241. The degree of the B-spline, `k`, is inferred from the length of `t` as
  242. ``len(t)-2``. The knot vector is constructed by appending and prepending
  243. ``k+1`` elements to internal knots `t`.
  244. Examples
  245. --------
  246. Construct a cubic B-spline:
  247. >>> import numpy as np
  248. >>> from scipy.interpolate import BSpline
  249. >>> b = BSpline.basis_element([0, 1, 2, 3, 4])
  250. >>> k = b.k
  251. >>> b.t[k:-k]
  252. array([ 0., 1., 2., 3., 4.])
  253. >>> k
  254. 3
  255. Construct a quadratic B-spline on ``[0, 1, 1, 2]``, and compare
  256. to its explicit form:
  257. >>> t = [0, 1, 1, 2]
  258. >>> b = BSpline.basis_element(t)
  259. >>> def f(x):
  260. ... return np.where(x < 1, x*x, (2. - x)**2)
  261. >>> import matplotlib.pyplot as plt
  262. >>> fig, ax = plt.subplots()
  263. >>> x = np.linspace(0, 2, 51)
  264. >>> ax.plot(x, b(x), 'g', lw=3)
  265. >>> ax.plot(x, f(x), 'r', lw=8, alpha=0.4)
  266. >>> ax.grid(True)
  267. >>> plt.show()
  268. """
  269. k = len(t) - 2
  270. t = _as_float_array(t)
  271. t = np.r_[(t[0]-1,) * k, t, (t[-1]+1,) * k]
  272. c = np.zeros_like(t)
  273. c[k] = 1.
  274. return cls.construct_fast(t, c, k, extrapolate)
  275. @classmethod
  276. def design_matrix(cls, x, t, k, extrapolate=False):
  277. """
  278. Returns a design matrix as a CSR format sparse array.
  279. Parameters
  280. ----------
  281. x : array_like, shape (n,)
  282. Points to evaluate the spline at.
  283. t : array_like, shape (nt,)
  284. Sorted 1D array of knots.
  285. k : int
  286. B-spline degree.
  287. extrapolate : bool or 'periodic', optional
  288. Whether to extrapolate based on the first and last intervals
  289. or raise an error. If 'periodic', periodic extrapolation is used.
  290. Default is False.
  291. .. versionadded:: 1.10.0
  292. Returns
  293. -------
  294. design_matrix : `csr_array` object
  295. Sparse matrix in CSR format where each row contains all the basis
  296. elements of the input row (first row = basis elements of x[0],
  297. ..., last row = basis elements x[-1]).
  298. Examples
  299. --------
  300. Construct a design matrix for a B-spline
  301. >>> from scipy.interpolate import make_interp_spline, BSpline
  302. >>> import numpy as np
  303. >>> x = np.linspace(0, np.pi * 2, 4)
  304. >>> y = np.sin(x)
  305. >>> k = 3
  306. >>> bspl = make_interp_spline(x, y, k=k)
  307. >>> design_matrix = bspl.design_matrix(x, bspl.t, k)
  308. >>> design_matrix.toarray()
  309. [[1. , 0. , 0. , 0. ],
  310. [0.2962963 , 0.44444444, 0.22222222, 0.03703704],
  311. [0.03703704, 0.22222222, 0.44444444, 0.2962963 ],
  312. [0. , 0. , 0. , 1. ]]
  313. Construct a design matrix for some vector of knots
  314. >>> k = 2
  315. >>> t = [-1, 0, 1, 2, 3, 4, 5, 6]
  316. >>> x = [1, 2, 3, 4]
  317. >>> design_matrix = BSpline.design_matrix(x, t, k).toarray()
  318. >>> design_matrix
  319. [[0.5, 0.5, 0. , 0. , 0. ],
  320. [0. , 0.5, 0.5, 0. , 0. ],
  321. [0. , 0. , 0.5, 0.5, 0. ],
  322. [0. , 0. , 0. , 0.5, 0.5]]
  323. This result is equivalent to the one created in the sparse format
  324. >>> c = np.eye(len(t) - k - 1)
  325. >>> design_matrix_gh = BSpline(t, c, k)(x)
  326. >>> np.allclose(design_matrix, design_matrix_gh, atol=1e-14)
  327. True
  328. Notes
  329. -----
  330. .. versionadded:: 1.8.0
  331. In each row of the design matrix all the basis elements are evaluated
  332. at the certain point (first row - x[0], ..., last row - x[-1]).
  333. `nt` is a length of the vector of knots: as far as there are
  334. `nt - k - 1` basis elements, `nt` should be not less than `2 * k + 2`
  335. to have at least `k + 1` basis element.
  336. Out of bounds `x` raises a ValueError.
  337. """
  338. x = _as_float_array(x, True)
  339. t = _as_float_array(t, True)
  340. if extrapolate != 'periodic':
  341. extrapolate = bool(extrapolate)
  342. if k < 0:
  343. raise ValueError("Spline order cannot be negative.")
  344. if t.ndim != 1 or np.any(t[1:] < t[:-1]):
  345. raise ValueError(f"Expect t to be a 1-D sorted array_like, but "
  346. f"got t={t}.")
  347. # There are `nt - k - 1` basis elements in a BSpline built on the
  348. # vector of knots with length `nt`, so to have at least `k + 1` basis
  349. # elements we need to have at least `2 * k + 2` elements in the vector
  350. # of knots.
  351. if len(t) < 2 * k + 2:
  352. raise ValueError(f"Length t is not enough for k={k}.")
  353. if extrapolate == 'periodic':
  354. # With periodic extrapolation we map x to the segment
  355. # [t[k], t[n]].
  356. n = t.size - k - 1
  357. x = t[k] + (x - t[k]) % (t[n] - t[k])
  358. extrapolate = False
  359. elif not extrapolate and (
  360. (min(x) < t[k]) or (max(x) > t[t.shape[0] - k - 1])
  361. ):
  362. # Checks from `find_interval` function
  363. raise ValueError(f'Out of bounds w/ x = {x}.')
  364. # Compute number of non-zeros of final CSR array in order to determine
  365. # the dtype of indices and indptr of the CSR array.
  366. n = x.shape[0]
  367. nnz = n * (k + 1)
  368. if nnz < np.iinfo(np.int32).max:
  369. int_dtype = np.int32
  370. else:
  371. int_dtype = np.int64
  372. # Preallocate indptr and indices
  373. indices = np.empty(n * (k + 1), dtype=int_dtype)
  374. indptr = np.arange(0, (n + 1) * (k + 1), k + 1, dtype=int_dtype)
  375. # indptr is not passed to Cython as it is already fully computed
  376. data, indices = _bspl._make_design_matrix(
  377. x, t, k, extrapolate, indices
  378. )
  379. return csr_array(
  380. (data, indices, indptr),
  381. shape=(x.shape[0], t.shape[0] - k - 1)
  382. )
  383. def __call__(self, x, nu=0, extrapolate=None):
  384. """
  385. Evaluate a spline function.
  386. Parameters
  387. ----------
  388. x : array_like
  389. points to evaluate the spline at.
  390. nu : int, optional
  391. derivative to evaluate (default is 0).
  392. extrapolate : bool or 'periodic', optional
  393. whether to extrapolate based on the first and last intervals
  394. or return nans. If 'periodic', periodic extrapolation is used.
  395. Default is `self.extrapolate`.
  396. Returns
  397. -------
  398. y : array_like
  399. Shape is determined by replacing the interpolation axis
  400. in the coefficient array with the shape of `x`.
  401. """
  402. if extrapolate is None:
  403. extrapolate = self.extrapolate
  404. x = np.asarray(x)
  405. x_shape, x_ndim = x.shape, x.ndim
  406. x = np.ascontiguousarray(x.ravel(), dtype=np.float_)
  407. # With periodic extrapolation we map x to the segment
  408. # [self.t[k], self.t[n]].
  409. if extrapolate == 'periodic':
  410. n = self.t.size - self.k - 1
  411. x = self.t[self.k] + (x - self.t[self.k]) % (self.t[n] -
  412. self.t[self.k])
  413. extrapolate = False
  414. out = np.empty((len(x), prod(self.c.shape[1:])), dtype=self.c.dtype)
  415. self._ensure_c_contiguous()
  416. self._evaluate(x, nu, extrapolate, out)
  417. out = out.reshape(x_shape + self.c.shape[1:])
  418. if self.axis != 0:
  419. # transpose to move the calculated values to the interpolation axis
  420. l = list(range(out.ndim))
  421. l = l[x_ndim:x_ndim+self.axis] + l[:x_ndim] + l[x_ndim+self.axis:]
  422. out = out.transpose(l)
  423. return out
  424. def _evaluate(self, xp, nu, extrapolate, out):
  425. _bspl.evaluate_spline(self.t, self.c.reshape(self.c.shape[0], -1),
  426. self.k, xp, nu, extrapolate, out)
  427. def _ensure_c_contiguous(self):
  428. """
  429. c and t may be modified by the user. The Cython code expects
  430. that they are C contiguous.
  431. """
  432. if not self.t.flags.c_contiguous:
  433. self.t = self.t.copy()
  434. if not self.c.flags.c_contiguous:
  435. self.c = self.c.copy()
  436. def derivative(self, nu=1):
  437. """Return a B-spline representing the derivative.
  438. Parameters
  439. ----------
  440. nu : int, optional
  441. Derivative order.
  442. Default is 1.
  443. Returns
  444. -------
  445. b : BSpline object
  446. A new instance representing the derivative.
  447. See Also
  448. --------
  449. splder, splantider
  450. """
  451. c = self.c
  452. # pad the c array if needed
  453. ct = len(self.t) - len(c)
  454. if ct > 0:
  455. c = np.r_[c, np.zeros((ct,) + c.shape[1:])]
  456. tck = _fitpack_impl.splder((self.t, c, self.k), nu)
  457. return self.construct_fast(*tck, extrapolate=self.extrapolate,
  458. axis=self.axis)
  459. def antiderivative(self, nu=1):
  460. """Return a B-spline representing the antiderivative.
  461. Parameters
  462. ----------
  463. nu : int, optional
  464. Antiderivative order. Default is 1.
  465. Returns
  466. -------
  467. b : BSpline object
  468. A new instance representing the antiderivative.
  469. Notes
  470. -----
  471. If antiderivative is computed and ``self.extrapolate='periodic'``,
  472. it will be set to False for the returned instance. This is done because
  473. the antiderivative is no longer periodic and its correct evaluation
  474. outside of the initially given x interval is difficult.
  475. See Also
  476. --------
  477. splder, splantider
  478. """
  479. c = self.c
  480. # pad the c array if needed
  481. ct = len(self.t) - len(c)
  482. if ct > 0:
  483. c = np.r_[c, np.zeros((ct,) + c.shape[1:])]
  484. tck = _fitpack_impl.splantider((self.t, c, self.k), nu)
  485. if self.extrapolate == 'periodic':
  486. extrapolate = False
  487. else:
  488. extrapolate = self.extrapolate
  489. return self.construct_fast(*tck, extrapolate=extrapolate,
  490. axis=self.axis)
  491. def integrate(self, a, b, extrapolate=None):
  492. """Compute a definite integral of the spline.
  493. Parameters
  494. ----------
  495. a : float
  496. Lower limit of integration.
  497. b : float
  498. Upper limit of integration.
  499. extrapolate : bool or 'periodic', optional
  500. whether to extrapolate beyond the base interval,
  501. ``t[k] .. t[-k-1]``, or take the spline to be zero outside of the
  502. base interval. If 'periodic', periodic extrapolation is used.
  503. If None (default), use `self.extrapolate`.
  504. Returns
  505. -------
  506. I : array_like
  507. Definite integral of the spline over the interval ``[a, b]``.
  508. Examples
  509. --------
  510. Construct the linear spline ``x if x < 1 else 2 - x`` on the base
  511. interval :math:`[0, 2]`, and integrate it
  512. >>> from scipy.interpolate import BSpline
  513. >>> b = BSpline.basis_element([0, 1, 2])
  514. >>> b.integrate(0, 1)
  515. array(0.5)
  516. If the integration limits are outside of the base interval, the result
  517. is controlled by the `extrapolate` parameter
  518. >>> b.integrate(-1, 1)
  519. array(0.0)
  520. >>> b.integrate(-1, 1, extrapolate=False)
  521. array(0.5)
  522. >>> import matplotlib.pyplot as plt
  523. >>> fig, ax = plt.subplots()
  524. >>> ax.grid(True)
  525. >>> ax.axvline(0, c='r', lw=5, alpha=0.5) # base interval
  526. >>> ax.axvline(2, c='r', lw=5, alpha=0.5)
  527. >>> xx = [-1, 1, 2]
  528. >>> ax.plot(xx, b(xx))
  529. >>> plt.show()
  530. """
  531. if extrapolate is None:
  532. extrapolate = self.extrapolate
  533. # Prepare self.t and self.c.
  534. self._ensure_c_contiguous()
  535. # Swap integration bounds if needed.
  536. sign = 1
  537. if b < a:
  538. a, b = b, a
  539. sign = -1
  540. n = self.t.size - self.k - 1
  541. if extrapolate != "periodic" and not extrapolate:
  542. # Shrink the integration interval, if needed.
  543. a = max(a, self.t[self.k])
  544. b = min(b, self.t[n])
  545. if self.c.ndim == 1:
  546. # Fast path: use FITPACK's routine
  547. # (cf _fitpack_impl.splint).
  548. integral = _fitpack_impl.splint(a, b, self.tck)
  549. return integral * sign
  550. out = np.empty((2, prod(self.c.shape[1:])), dtype=self.c.dtype)
  551. # Compute the antiderivative.
  552. c = self.c
  553. ct = len(self.t) - len(c)
  554. if ct > 0:
  555. c = np.r_[c, np.zeros((ct,) + c.shape[1:])]
  556. ta, ca, ka = _fitpack_impl.splantider((self.t, c, self.k), 1)
  557. if extrapolate == 'periodic':
  558. # Split the integral into the part over period (can be several
  559. # of them) and the remaining part.
  560. ts, te = self.t[self.k], self.t[n]
  561. period = te - ts
  562. interval = b - a
  563. n_periods, left = divmod(interval, period)
  564. if n_periods > 0:
  565. # Evaluate the difference of antiderivatives.
  566. x = np.asarray([ts, te], dtype=np.float_)
  567. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  568. ka, x, 0, False, out)
  569. integral = out[1] - out[0]
  570. integral *= n_periods
  571. else:
  572. integral = np.zeros((1, prod(self.c.shape[1:])),
  573. dtype=self.c.dtype)
  574. # Map a to [ts, te], b is always a + left.
  575. a = ts + (a - ts) % period
  576. b = a + left
  577. # If b <= te then we need to integrate over [a, b], otherwise
  578. # over [a, te] and from xs to what is remained.
  579. if b <= te:
  580. x = np.asarray([a, b], dtype=np.float_)
  581. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  582. ka, x, 0, False, out)
  583. integral += out[1] - out[0]
  584. else:
  585. x = np.asarray([a, te], dtype=np.float_)
  586. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  587. ka, x, 0, False, out)
  588. integral += out[1] - out[0]
  589. x = np.asarray([ts, ts + b - te], dtype=np.float_)
  590. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  591. ka, x, 0, False, out)
  592. integral += out[1] - out[0]
  593. else:
  594. # Evaluate the difference of antiderivatives.
  595. x = np.asarray([a, b], dtype=np.float_)
  596. _bspl.evaluate_spline(ta, ca.reshape(ca.shape[0], -1),
  597. ka, x, 0, extrapolate, out)
  598. integral = out[1] - out[0]
  599. integral *= sign
  600. return integral.reshape(ca.shape[1:])
  601. @classmethod
  602. def from_power_basis(cls, pp, bc_type='not-a-knot'):
  603. r"""
  604. Construct a polynomial in the B-spline basis
  605. from a piecewise polynomial in the power basis.
  606. For now, accepts ``CubicSpline`` instances only.
  607. Parameters
  608. ----------
  609. pp : CubicSpline
  610. A piecewise polynomial in the power basis, as created
  611. by ``CubicSpline``
  612. bc_type : string, optional
  613. Boundary condition type as in ``CubicSpline``: one of the
  614. ``not-a-knot``, ``natural``, ``clamped``, or ``periodic``.
  615. Necessary for construction an instance of ``BSpline`` class.
  616. Default is ``not-a-knot``.
  617. Returns
  618. -------
  619. b : BSpline object
  620. A new instance representing the initial polynomial
  621. in the B-spline basis.
  622. Notes
  623. -----
  624. .. versionadded:: 1.8.0
  625. Accepts only ``CubicSpline`` instances for now.
  626. The algorithm follows from differentiation
  627. the Marsden's identity [1]: each of coefficients of spline
  628. interpolation function in the B-spline basis is computed as follows:
  629. .. math::
  630. c_j = \sum_{m=0}^{k} \frac{(k-m)!}{k!}
  631. c_{m,i} (-1)^{k-m} D^m p_{j,k}(x_i)
  632. :math:`c_{m, i}` - a coefficient of CubicSpline,
  633. :math:`D^m p_{j, k}(x_i)` - an m-th defivative of a dual polynomial
  634. in :math:`x_i`.
  635. ``k`` always equals 3 for now.
  636. First ``n - 2`` coefficients are computed in :math:`x_i = x_j`, e.g.
  637. .. math::
  638. c_1 = \sum_{m=0}^{k} \frac{(k-1)!}{k!} c_{m,1} D^m p_{j,3}(x_1)
  639. Last ``nod + 2`` coefficients are computed in ``x[-2]``,
  640. ``nod`` - number of derivatives at the ends.
  641. For example, consider :math:`x = [0, 1, 2, 3, 4]`,
  642. :math:`y = [1, 1, 1, 1, 1]` and bc_type = ``natural``
  643. The coefficients of CubicSpline in the power basis:
  644. :math:`[[0, 0, 0, 0, 0], [0, 0, 0, 0, 0],
  645. [0, 0, 0, 0, 0], [1, 1, 1, 1, 1]]`
  646. The knot vector: :math:`t = [0, 0, 0, 0, 1, 2, 3, 4, 4, 4, 4]`
  647. In this case
  648. .. math::
  649. c_j = \frac{0!}{k!} c_{3, i} k! = c_{3, i} = 1,~j = 0, ..., 6
  650. References
  651. ----------
  652. .. [1] Tom Lyche and Knut Morken, Spline Methods, 2005, Section 3.1.2
  653. """
  654. from ._cubic import CubicSpline
  655. if not isinstance(pp, CubicSpline):
  656. raise NotImplementedError("Only CubicSpline objects are accepted"
  657. "for now. Got %s instead." % type(pp))
  658. x = pp.x
  659. coef = pp.c
  660. k = pp.c.shape[0] - 1
  661. n = x.shape[0]
  662. if bc_type == 'not-a-knot':
  663. t = _not_a_knot(x, k)
  664. elif bc_type == 'natural' or bc_type == 'clamped':
  665. t = _augknt(x, k)
  666. elif bc_type == 'periodic':
  667. t = _periodic_knots(x, k)
  668. else:
  669. raise TypeError('Unknown boundary condition: %s' % bc_type)
  670. nod = t.shape[0] - (n + k + 1) # number of derivatives at the ends
  671. c = np.zeros(n + nod, dtype=pp.c.dtype)
  672. for m in range(k + 1):
  673. for i in range(n - 2):
  674. c[i] += poch(k + 1, -m) * coef[m, i]\
  675. * np.power(-1, k - m)\
  676. * _diff_dual_poly(i, k, x[i], m, t)
  677. for j in range(n - 2, n + nod):
  678. c[j] += poch(k + 1, -m) * coef[m, n - 2]\
  679. * np.power(-1, k - m)\
  680. * _diff_dual_poly(j, k, x[n - 2], m, t)
  681. return cls.construct_fast(t, c, k, pp.extrapolate, pp.axis)
  682. #################################
  683. # Interpolating spline helpers #
  684. #################################
  685. def _not_a_knot(x, k):
  686. """Given data x, construct the knot vector w/ not-a-knot BC.
  687. cf de Boor, XIII(12)."""
  688. x = np.asarray(x)
  689. if k % 2 != 1:
  690. raise ValueError("Odd degree for now only. Got %s." % k)
  691. m = (k - 1) // 2
  692. t = x[m+1:-m-1]
  693. t = np.r_[(x[0],)*(k+1), t, (x[-1],)*(k+1)]
  694. return t
  695. def _augknt(x, k):
  696. """Construct a knot vector appropriate for the order-k interpolation."""
  697. return np.r_[(x[0],)*k, x, (x[-1],)*k]
  698. def _convert_string_aliases(deriv, target_shape):
  699. if isinstance(deriv, str):
  700. if deriv == "clamped":
  701. deriv = [(1, np.zeros(target_shape))]
  702. elif deriv == "natural":
  703. deriv = [(2, np.zeros(target_shape))]
  704. else:
  705. raise ValueError("Unknown boundary condition : %s" % deriv)
  706. return deriv
  707. def _process_deriv_spec(deriv):
  708. if deriv is not None:
  709. try:
  710. ords, vals = zip(*deriv)
  711. except TypeError as e:
  712. msg = ("Derivatives, `bc_type`, should be specified as a pair of "
  713. "iterables of pairs of (order, value).")
  714. raise ValueError(msg) from e
  715. else:
  716. ords, vals = [], []
  717. return np.atleast_1d(ords, vals)
  718. def _woodbury_algorithm(A, ur, ll, b, k):
  719. '''
  720. Solve a cyclic banded linear system with upper right
  721. and lower blocks of size ``(k-1) / 2`` using
  722. the Woodbury formula
  723. Parameters
  724. ----------
  725. A : 2-D array, shape(k, n)
  726. Matrix of diagonals of original matrix(see
  727. ``solve_banded`` documentation).
  728. ur : 2-D array, shape(bs, bs)
  729. Upper right block matrix.
  730. ll : 2-D array, shape(bs, bs)
  731. Lower left block matrix.
  732. b : 1-D array, shape(n,)
  733. Vector of constant terms of the system of linear equations.
  734. k : int
  735. B-spline degree.
  736. Returns
  737. -------
  738. c : 1-D array, shape(n,)
  739. Solution of the original system of linear equations.
  740. Notes
  741. -----
  742. This algorithm works only for systems with banded matrix A plus
  743. a correction term U @ V.T, where the matrix U @ V.T gives upper right
  744. and lower left block of A
  745. The system is solved with the following steps:
  746. 1. New systems of linear equations are constructed:
  747. A @ z_i = u_i,
  748. u_i - columnn vector of U,
  749. i = 1, ..., k - 1
  750. 2. Matrix Z is formed from vectors z_i:
  751. Z = [ z_1 | z_2 | ... | z_{k - 1} ]
  752. 3. Matrix H = (1 + V.T @ Z)^{-1}
  753. 4. The system A' @ y = b is solved
  754. 5. x = y - Z @ (H @ V.T @ y)
  755. Also, ``n`` should be greater than ``k``, otherwise corner block
  756. elements will intersect with diagonals.
  757. Examples
  758. --------
  759. Consider the case of n = 8, k = 5 (size of blocks - 2 x 2).
  760. The matrix of a system: U: V:
  761. x x x * * a b a b 0 0 0 0 1 0
  762. x x x x * * c 0 c 0 0 0 0 0 1
  763. x x x x x * * 0 0 0 0 0 0 0 0
  764. * x x x x x * 0 0 0 0 0 0 0 0
  765. * * x x x x x 0 0 0 0 0 0 0 0
  766. d * * x x x x 0 0 d 0 1 0 0 0
  767. e f * * x x x 0 0 e f 0 1 0 0
  768. References
  769. ----------
  770. .. [1] William H. Press, Saul A. Teukolsky, William T. Vetterling
  771. and Brian P. Flannery, Numerical Recipes, 2007, Section 2.7.3
  772. '''
  773. k_mod = k - k % 2
  774. bs = int((k - 1) / 2) + (k + 1) % 2
  775. n = A.shape[1] + 1
  776. U = np.zeros((n - 1, k_mod))
  777. VT = np.zeros((k_mod, n - 1)) # V transpose
  778. # upper right block
  779. U[:bs, :bs] = ur
  780. VT[np.arange(bs), np.arange(bs) - bs] = 1
  781. # lower left block
  782. U[-bs:, -bs:] = ll
  783. VT[np.arange(bs) - bs, np.arange(bs)] = 1
  784. Z = solve_banded((bs, bs), A, U)
  785. H = solve(np.identity(k_mod) + VT @ Z, np.identity(k_mod))
  786. y = solve_banded((bs, bs), A, b)
  787. c = y - Z @ (H @ (VT @ y))
  788. return c
  789. def _periodic_knots(x, k):
  790. '''
  791. returns vector of nodes on circle
  792. '''
  793. xc = np.copy(x)
  794. n = len(xc)
  795. if k % 2 == 0:
  796. dx = np.diff(xc)
  797. xc[1: -1] -= dx[:-1] / 2
  798. dx = np.diff(xc)
  799. t = np.zeros(n + 2 * k)
  800. t[k: -k] = xc
  801. for i in range(0, k):
  802. # filling first `k` elements in descending order
  803. t[k - i - 1] = t[k - i] - dx[-(i % (n - 1)) - 1]
  804. # filling last `k` elements in ascending order
  805. t[-k + i] = t[-k + i - 1] + dx[i % (n - 1)]
  806. return t
  807. def _make_interp_per_full_matr(x, y, t, k):
  808. '''
  809. Returns a solution of a system for B-spline interpolation with periodic
  810. boundary conditions. First ``k - 1`` rows of matrix are condtions of
  811. periodicity (continuity of ``k - 1`` derivatives at the boundary points).
  812. Last ``n`` rows are interpolation conditions.
  813. RHS is ``k - 1`` zeros and ``n`` ordinates in this case.
  814. Parameters
  815. ----------
  816. x : 1-D array, shape (n,)
  817. Values of x - coordinate of a given set of points.
  818. y : 1-D array, shape (n,)
  819. Values of y - coordinate of a given set of points.
  820. t : 1-D array, shape(n+2*k,)
  821. Vector of knots.
  822. k : int
  823. The maximum degree of spline
  824. Returns
  825. -------
  826. c : 1-D array, shape (n+k-1,)
  827. B-spline coefficients
  828. Notes
  829. -----
  830. ``t`` is supposed to be taken on circle.
  831. '''
  832. x, y, t = map(np.asarray, (x, y, t))
  833. n = x.size
  834. # LHS: the collocation matrix + derivatives at edges
  835. matr = np.zeros((n + k - 1, n + k - 1))
  836. # derivatives at x[0] and x[-1]:
  837. for i in range(k - 1):
  838. bb = _bspl.evaluate_all_bspl(t, k, x[0], k, nu=i + 1)
  839. matr[i, : k + 1] += bb
  840. bb = _bspl.evaluate_all_bspl(t, k, x[-1], n + k - 1, nu=i + 1)[:-1]
  841. matr[i, -k:] -= bb
  842. # collocation matrix
  843. for i in range(n):
  844. xval = x[i]
  845. # find interval
  846. if xval == t[k]:
  847. left = k
  848. else:
  849. left = np.searchsorted(t, xval) - 1
  850. # fill a row
  851. bb = _bspl.evaluate_all_bspl(t, k, xval, left)
  852. matr[i + k - 1, left-k:left+1] = bb
  853. # RHS
  854. b = np.r_[[0] * (k - 1), y]
  855. c = solve(matr, b)
  856. return c
  857. def _make_periodic_spline(x, y, t, k, axis):
  858. '''
  859. Compute the (coefficients of) interpolating B-spline with periodic
  860. boundary conditions.
  861. Parameters
  862. ----------
  863. x : array_like, shape (n,)
  864. Abscissas.
  865. y : array_like, shape (n,)
  866. Ordinates.
  867. k : int
  868. B-spline degree.
  869. t : array_like, shape (n + 2 * k,).
  870. Knots taken on a circle, ``k`` on the left and ``k`` on the right
  871. of the vector ``x``.
  872. Returns
  873. -------
  874. b : a BSpline object of the degree ``k`` and with knots ``t``.
  875. Notes
  876. -----
  877. The original system is formed by ``n + k - 1`` equations where the first
  878. ``k - 1`` of them stand for the ``k - 1`` derivatives continuity on the
  879. edges while the other equations correspond to an interpolating case
  880. (matching all the input points). Due to a special form of knot vector, it
  881. can be proved that in the original system the first and last ``k``
  882. coefficients of a spline function are the same, respectively. It follows
  883. from the fact that all ``k - 1`` derivatives are equal term by term at ends
  884. and that the matrix of the original system of linear equations is
  885. non-degenerate. So, we can reduce the number of equations to ``n - 1``
  886. (first ``k - 1`` equations could be reduced). Another trick of this
  887. implementation is cyclic shift of values of B-splines due to equality of
  888. ``k`` unknown coefficients. With this we can receive matrix of the system
  889. with upper right and lower left blocks, and ``k`` diagonals. It allows
  890. to use Woodbury formula to optimize the computations.
  891. '''
  892. n = y.shape[0]
  893. extradim = prod(y.shape[1:])
  894. y_new = y.reshape(n, extradim)
  895. c = np.zeros((n + k - 1, extradim))
  896. # n <= k case is solved with full matrix
  897. if n <= k:
  898. for i in range(extradim):
  899. c[:, i] = _make_interp_per_full_matr(x, y_new[:, i], t, k)
  900. c = np.ascontiguousarray(c.reshape((n + k - 1,) + y.shape[1:]))
  901. return BSpline.construct_fast(t, c, k, extrapolate='periodic', axis=axis)
  902. nt = len(t) - k - 1
  903. # size of block elements
  904. kul = int(k / 2)
  905. # kl = ku = k
  906. ab = np.zeros((3 * k + 1, nt), dtype=np.float_, order='F')
  907. # upper right and lower left blocks
  908. ur = np.zeros((kul, kul))
  909. ll = np.zeros_like(ur)
  910. # `offset` is made to shift all the non-zero elements to the end of the
  911. # matrix
  912. _bspl._colloc(x, t, k, ab, offset=k)
  913. # remove zeros before the matrix
  914. ab = ab[-k - (k + 1) % 2:, :]
  915. # The least elements in rows (except repetitions) are diagonals
  916. # of block matrices. Upper right matrix is an upper triangular
  917. # matrix while lower left is a lower triangular one.
  918. for i in range(kul):
  919. ur += np.diag(ab[-i - 1, i: kul], k=i)
  920. ll += np.diag(ab[i, -kul - (k % 2): n - 1 + 2 * kul - i], k=-i)
  921. # remove elements that occur in the last point
  922. # (first and last points are equivalent)
  923. A = ab[:, kul: -k + kul]
  924. for i in range(extradim):
  925. cc = _woodbury_algorithm(A, ur, ll, y_new[:, i][:-1], k)
  926. c[:, i] = np.concatenate((cc[-kul:], cc, cc[:kul + k % 2]))
  927. c = np.ascontiguousarray(c.reshape((n + k - 1,) + y.shape[1:]))
  928. return BSpline.construct_fast(t, c, k, extrapolate='periodic', axis=axis)
  929. def make_interp_spline(x, y, k=3, t=None, bc_type=None, axis=0,
  930. check_finite=True):
  931. """Compute the (coefficients of) interpolating B-spline.
  932. Parameters
  933. ----------
  934. x : array_like, shape (n,)
  935. Abscissas.
  936. y : array_like, shape (n, ...)
  937. Ordinates.
  938. k : int, optional
  939. B-spline degree. Default is cubic, ``k = 3``.
  940. t : array_like, shape (nt + k + 1,), optional.
  941. Knots.
  942. The number of knots needs to agree with the number of data points and
  943. the number of derivatives at the edges. Specifically, ``nt - n`` must
  944. equal ``len(deriv_l) + len(deriv_r)``.
  945. bc_type : 2-tuple or None
  946. Boundary conditions.
  947. Default is None, which means choosing the boundary conditions
  948. automatically. Otherwise, it must be a length-two tuple where the first
  949. element (``deriv_l``) sets the boundary conditions at ``x[0]`` and
  950. the second element (``deriv_r``) sets the boundary conditions at
  951. ``x[-1]``. Each of these must be an iterable of pairs
  952. ``(order, value)`` which gives the values of derivatives of specified
  953. orders at the given edge of the interpolation interval.
  954. Alternatively, the following string aliases are recognized:
  955. * ``"clamped"``: The first derivatives at the ends are zero. This is
  956. equivalent to ``bc_type=([(1, 0.0)], [(1, 0.0)])``.
  957. * ``"natural"``: The second derivatives at ends are zero. This is
  958. equivalent to ``bc_type=([(2, 0.0)], [(2, 0.0)])``.
  959. * ``"not-a-knot"`` (default): The first and second segments are the
  960. same polynomial. This is equivalent to having ``bc_type=None``.
  961. * ``"periodic"``: The values and the first ``k-1`` derivatives at the
  962. ends are equivalent.
  963. axis : int, optional
  964. Interpolation axis. Default is 0.
  965. check_finite : bool, optional
  966. Whether to check that the input arrays contain only finite numbers.
  967. Disabling may give a performance gain, but may result in problems
  968. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  969. Default is True.
  970. Returns
  971. -------
  972. b : a BSpline object of the degree ``k`` and with knots ``t``.
  973. Examples
  974. --------
  975. Use cubic interpolation on Chebyshev nodes:
  976. >>> import numpy as np
  977. >>> import matplotlib.pyplot as plt
  978. >>> def cheb_nodes(N):
  979. ... jj = 2.*np.arange(N) + 1
  980. ... x = np.cos(np.pi * jj / 2 / N)[::-1]
  981. ... return x
  982. >>> x = cheb_nodes(20)
  983. >>> y = np.sqrt(1 - x**2)
  984. >>> from scipy.interpolate import BSpline, make_interp_spline
  985. >>> b = make_interp_spline(x, y)
  986. >>> np.allclose(b(x), y)
  987. True
  988. Note that the default is a cubic spline with a not-a-knot boundary condition
  989. >>> b.k
  990. 3
  991. Here we use a 'natural' spline, with zero 2nd derivatives at edges:
  992. >>> l, r = [(2, 0.0)], [(2, 0.0)]
  993. >>> b_n = make_interp_spline(x, y, bc_type=(l, r)) # or, bc_type="natural"
  994. >>> np.allclose(b_n(x), y)
  995. True
  996. >>> x0, x1 = x[0], x[-1]
  997. >>> np.allclose([b_n(x0, 2), b_n(x1, 2)], [0, 0])
  998. True
  999. Interpolation of parametric curves is also supported. As an example, we
  1000. compute a discretization of a snail curve in polar coordinates
  1001. >>> phi = np.linspace(0, 2.*np.pi, 40)
  1002. >>> r = 0.3 + np.cos(phi)
  1003. >>> x, y = r*np.cos(phi), r*np.sin(phi) # convert to Cartesian coordinates
  1004. Build an interpolating curve, parameterizing it by the angle
  1005. >>> spl = make_interp_spline(phi, np.c_[x, y])
  1006. Evaluate the interpolant on a finer grid (note that we transpose the result
  1007. to unpack it into a pair of x- and y-arrays)
  1008. >>> phi_new = np.linspace(0, 2.*np.pi, 100)
  1009. >>> x_new, y_new = spl(phi_new).T
  1010. Plot the result
  1011. >>> plt.plot(x, y, 'o')
  1012. >>> plt.plot(x_new, y_new, '-')
  1013. >>> plt.show()
  1014. Build a B-spline curve with 2 dimensional y
  1015. >>> x = np.linspace(0, 2*np.pi, 10)
  1016. >>> y = np.array([np.sin(x), np.cos(x)])
  1017. Periodic condition is satisfied because y coordinates of points on the ends
  1018. are equivalent
  1019. >>> ax = plt.axes(projection='3d')
  1020. >>> xx = np.linspace(0, 2*np.pi, 100)
  1021. >>> bspl = make_interp_spline(x, y, k=5, bc_type='periodic', axis=1)
  1022. >>> ax.plot3D(xx, *bspl(xx))
  1023. >>> ax.scatter3D(x, *y, color='red')
  1024. >>> plt.show()
  1025. See Also
  1026. --------
  1027. BSpline : base class representing the B-spline objects
  1028. CubicSpline : a cubic spline in the polynomial basis
  1029. make_lsq_spline : a similar factory function for spline fitting
  1030. UnivariateSpline : a wrapper over FITPACK spline fitting routines
  1031. splrep : a wrapper over FITPACK spline fitting routines
  1032. """
  1033. # convert string aliases for the boundary conditions
  1034. if bc_type is None or bc_type == 'not-a-knot' or bc_type == 'periodic':
  1035. deriv_l, deriv_r = None, None
  1036. elif isinstance(bc_type, str):
  1037. deriv_l, deriv_r = bc_type, bc_type
  1038. else:
  1039. try:
  1040. deriv_l, deriv_r = bc_type
  1041. except TypeError as e:
  1042. raise ValueError("Unknown boundary condition: %s" % bc_type) from e
  1043. y = np.asarray(y)
  1044. axis = normalize_axis_index(axis, y.ndim)
  1045. x = _as_float_array(x, check_finite)
  1046. y = _as_float_array(y, check_finite)
  1047. y = np.moveaxis(y, axis, 0) # now internally interp axis is zero
  1048. # sanity check the input
  1049. if bc_type == 'periodic' and not np.allclose(y[0], y[-1], atol=1e-15):
  1050. raise ValueError("First and last points does not match while "
  1051. "periodic case expected")
  1052. if x.size != y.shape[0]:
  1053. raise ValueError('Shapes of x {} and y {} are incompatible'
  1054. .format(x.shape, y.shape))
  1055. if np.any(x[1:] == x[:-1]):
  1056. raise ValueError("Expect x to not have duplicates")
  1057. if x.ndim != 1 or np.any(x[1:] < x[:-1]):
  1058. raise ValueError("Expect x to be a 1D strictly increasing sequence.")
  1059. # special-case k=0 right away
  1060. if k == 0:
  1061. if any(_ is not None for _ in (t, deriv_l, deriv_r)):
  1062. raise ValueError("Too much info for k=0: t and bc_type can only "
  1063. "be None.")
  1064. t = np.r_[x, x[-1]]
  1065. c = np.asarray(y)
  1066. c = np.ascontiguousarray(c, dtype=_get_dtype(c.dtype))
  1067. return BSpline.construct_fast(t, c, k, axis=axis)
  1068. # special-case k=1 (e.g., Lyche and Morken, Eq.(2.16))
  1069. if k == 1 and t is None:
  1070. if not (deriv_l is None and deriv_r is None):
  1071. raise ValueError("Too much info for k=1: bc_type can only be None.")
  1072. t = np.r_[x[0], x, x[-1]]
  1073. c = np.asarray(y)
  1074. c = np.ascontiguousarray(c, dtype=_get_dtype(c.dtype))
  1075. return BSpline.construct_fast(t, c, k, axis=axis)
  1076. k = operator.index(k)
  1077. if bc_type == 'periodic' and t is not None:
  1078. raise NotImplementedError("For periodic case t is constructed "
  1079. "automatically and can not be passed manually")
  1080. # come up with a sensible knot vector, if needed
  1081. if t is None:
  1082. if deriv_l is None and deriv_r is None:
  1083. if bc_type == 'periodic':
  1084. t = _periodic_knots(x, k)
  1085. elif k == 2:
  1086. # OK, it's a bit ad hoc: Greville sites + omit
  1087. # 2nd and 2nd-to-last points, a la not-a-knot
  1088. t = (x[1:] + x[:-1]) / 2.
  1089. t = np.r_[(x[0],)*(k+1),
  1090. t[1:-1],
  1091. (x[-1],)*(k+1)]
  1092. else:
  1093. t = _not_a_knot(x, k)
  1094. else:
  1095. t = _augknt(x, k)
  1096. t = _as_float_array(t, check_finite)
  1097. if k < 0:
  1098. raise ValueError("Expect non-negative k.")
  1099. if t.ndim != 1 or np.any(t[1:] < t[:-1]):
  1100. raise ValueError("Expect t to be a 1-D sorted array_like.")
  1101. if t.size < x.size + k + 1:
  1102. raise ValueError('Got %d knots, need at least %d.' %
  1103. (t.size, x.size + k + 1))
  1104. if (x[0] < t[k]) or (x[-1] > t[-k]):
  1105. raise ValueError('Out of bounds w/ x = %s.' % x)
  1106. if bc_type == 'periodic':
  1107. return _make_periodic_spline(x, y, t, k, axis)
  1108. # Here : deriv_l, r = [(nu, value), ...]
  1109. deriv_l = _convert_string_aliases(deriv_l, y.shape[1:])
  1110. deriv_l_ords, deriv_l_vals = _process_deriv_spec(deriv_l)
  1111. nleft = deriv_l_ords.shape[0]
  1112. deriv_r = _convert_string_aliases(deriv_r, y.shape[1:])
  1113. deriv_r_ords, deriv_r_vals = _process_deriv_spec(deriv_r)
  1114. nright = deriv_r_ords.shape[0]
  1115. # have `n` conditions for `nt` coefficients; need nt-n derivatives
  1116. n = x.size
  1117. nt = t.size - k - 1
  1118. if nt - n != nleft + nright:
  1119. raise ValueError("The number of derivatives at boundaries does not "
  1120. "match: expected %s, got %s+%s" % (nt-n, nleft, nright))
  1121. # bail out if the `y` array is zero-sized
  1122. if y.size == 0:
  1123. c = np.zeros((nt,) + y.shape[1:], dtype=float)
  1124. return BSpline.construct_fast(t, c, k, axis=axis)
  1125. # set up the LHS: the collocation matrix + derivatives at boundaries
  1126. kl = ku = k
  1127. ab = np.zeros((2*kl + ku + 1, nt), dtype=np.float_, order='F')
  1128. _bspl._colloc(x, t, k, ab, offset=nleft)
  1129. if nleft > 0:
  1130. _bspl._handle_lhs_derivatives(t, k, x[0], ab, kl, ku, deriv_l_ords)
  1131. if nright > 0:
  1132. _bspl._handle_lhs_derivatives(t, k, x[-1], ab, kl, ku, deriv_r_ords,
  1133. offset=nt-nright)
  1134. # set up the RHS: values to interpolate (+ derivative values, if any)
  1135. extradim = prod(y.shape[1:])
  1136. rhs = np.empty((nt, extradim), dtype=y.dtype)
  1137. if nleft > 0:
  1138. rhs[:nleft] = deriv_l_vals.reshape(-1, extradim)
  1139. rhs[nleft:nt - nright] = y.reshape(-1, extradim)
  1140. if nright > 0:
  1141. rhs[nt - nright:] = deriv_r_vals.reshape(-1, extradim)
  1142. # solve Ab @ x = rhs; this is the relevant part of linalg.solve_banded
  1143. if check_finite:
  1144. ab, rhs = map(np.asarray_chkfinite, (ab, rhs))
  1145. gbsv, = get_lapack_funcs(('gbsv',), (ab, rhs))
  1146. lu, piv, c, info = gbsv(kl, ku, ab, rhs,
  1147. overwrite_ab=True, overwrite_b=True)
  1148. if info > 0:
  1149. raise LinAlgError("Collocation matrix is singular.")
  1150. elif info < 0:
  1151. raise ValueError('illegal value in %d-th argument of internal gbsv' % -info)
  1152. c = np.ascontiguousarray(c.reshape((nt,) + y.shape[1:]))
  1153. return BSpline.construct_fast(t, c, k, axis=axis)
  1154. def make_lsq_spline(x, y, t, k=3, w=None, axis=0, check_finite=True):
  1155. r"""Compute the (coefficients of) an LSQ (Least SQuared) based
  1156. fitting B-spline.
  1157. The result is a linear combination
  1158. .. math::
  1159. S(x) = \sum_j c_j B_j(x; t)
  1160. of the B-spline basis elements, :math:`B_j(x; t)`, which minimizes
  1161. .. math::
  1162. \sum_{j} \left( w_j \times (S(x_j) - y_j) \right)^2
  1163. Parameters
  1164. ----------
  1165. x : array_like, shape (m,)
  1166. Abscissas.
  1167. y : array_like, shape (m, ...)
  1168. Ordinates.
  1169. t : array_like, shape (n + k + 1,).
  1170. Knots.
  1171. Knots and data points must satisfy Schoenberg-Whitney conditions.
  1172. k : int, optional
  1173. B-spline degree. Default is cubic, ``k = 3``.
  1174. w : array_like, shape (m,), optional
  1175. Weights for spline fitting. Must be positive. If ``None``,
  1176. then weights are all equal.
  1177. Default is ``None``.
  1178. axis : int, optional
  1179. Interpolation axis. Default is zero.
  1180. check_finite : bool, optional
  1181. Whether to check that the input arrays contain only finite numbers.
  1182. Disabling may give a performance gain, but may result in problems
  1183. (crashes, non-termination) if the inputs do contain infinities or NaNs.
  1184. Default is True.
  1185. Returns
  1186. -------
  1187. b : a BSpline object of the degree ``k`` with knots ``t``.
  1188. Notes
  1189. -----
  1190. The number of data points must be larger than the spline degree ``k``.
  1191. Knots ``t`` must satisfy the Schoenberg-Whitney conditions,
  1192. i.e., there must be a subset of data points ``x[j]`` such that
  1193. ``t[j] < x[j] < t[j+k+1]``, for ``j=0, 1,...,n-k-2``.
  1194. Examples
  1195. --------
  1196. Generate some noisy data:
  1197. >>> import numpy as np
  1198. >>> import matplotlib.pyplot as plt
  1199. >>> rng = np.random.default_rng()
  1200. >>> x = np.linspace(-3, 3, 50)
  1201. >>> y = np.exp(-x**2) + 0.1 * rng.standard_normal(50)
  1202. Now fit a smoothing cubic spline with a pre-defined internal knots.
  1203. Here we make the knot vector (k+1)-regular by adding boundary knots:
  1204. >>> from scipy.interpolate import make_lsq_spline, BSpline
  1205. >>> t = [-1, 0, 1]
  1206. >>> k = 3
  1207. >>> t = np.r_[(x[0],)*(k+1),
  1208. ... t,
  1209. ... (x[-1],)*(k+1)]
  1210. >>> spl = make_lsq_spline(x, y, t, k)
  1211. For comparison, we also construct an interpolating spline for the same
  1212. set of data:
  1213. >>> from scipy.interpolate import make_interp_spline
  1214. >>> spl_i = make_interp_spline(x, y)
  1215. Plot both:
  1216. >>> xs = np.linspace(-3, 3, 100)
  1217. >>> plt.plot(x, y, 'ro', ms=5)
  1218. >>> plt.plot(xs, spl(xs), 'g-', lw=3, label='LSQ spline')
  1219. >>> plt.plot(xs, spl_i(xs), 'b-', lw=3, alpha=0.7, label='interp spline')
  1220. >>> plt.legend(loc='best')
  1221. >>> plt.show()
  1222. **NaN handling**: If the input arrays contain ``nan`` values, the result is
  1223. not useful since the underlying spline fitting routines cannot deal with
  1224. ``nan``. A workaround is to use zero weights for not-a-number data points:
  1225. >>> y[8] = np.nan
  1226. >>> w = np.isnan(y)
  1227. >>> y[w] = 0.
  1228. >>> tck = make_lsq_spline(x, y, t, w=~w)
  1229. Notice the need to replace a ``nan`` by a numerical value (precise value
  1230. does not matter as long as the corresponding weight is zero.)
  1231. See Also
  1232. --------
  1233. BSpline : base class representing the B-spline objects
  1234. make_interp_spline : a similar factory function for interpolating splines
  1235. LSQUnivariateSpline : a FITPACK-based spline fitting routine
  1236. splrep : a FITPACK-based fitting routine
  1237. """
  1238. x = _as_float_array(x, check_finite)
  1239. y = _as_float_array(y, check_finite)
  1240. t = _as_float_array(t, check_finite)
  1241. if w is not None:
  1242. w = _as_float_array(w, check_finite)
  1243. else:
  1244. w = np.ones_like(x)
  1245. k = operator.index(k)
  1246. axis = normalize_axis_index(axis, y.ndim)
  1247. y = np.moveaxis(y, axis, 0) # now internally interp axis is zero
  1248. if x.ndim != 1 or np.any(x[1:] - x[:-1] <= 0):
  1249. raise ValueError("Expect x to be a 1-D sorted array_like.")
  1250. if x.shape[0] < k+1:
  1251. raise ValueError("Need more x points.")
  1252. if k < 0:
  1253. raise ValueError("Expect non-negative k.")
  1254. if t.ndim != 1 or np.any(t[1:] - t[:-1] < 0):
  1255. raise ValueError("Expect t to be a 1-D sorted array_like.")
  1256. if x.size != y.shape[0]:
  1257. raise ValueError('Shapes of x {} and y {} are incompatible'
  1258. .format(x.shape, y.shape))
  1259. if k > 0 and np.any((x < t[k]) | (x > t[-k])):
  1260. raise ValueError('Out of bounds w/ x = %s.' % x)
  1261. if x.size != w.size:
  1262. raise ValueError('Shapes of x {} and w {} are incompatible'
  1263. .format(x.shape, w.shape))
  1264. # number of coefficients
  1265. n = t.size - k - 1
  1266. # construct A.T @ A and rhs with A the collocation matrix, and
  1267. # rhs = A.T @ y for solving the LSQ problem ``A.T @ A @ c = A.T @ y``
  1268. lower = True
  1269. extradim = prod(y.shape[1:])
  1270. ab = np.zeros((k+1, n), dtype=np.float_, order='F')
  1271. rhs = np.zeros((n, extradim), dtype=y.dtype, order='F')
  1272. _bspl._norm_eq_lsq(x, t, k,
  1273. y.reshape(-1, extradim),
  1274. w,
  1275. ab, rhs)
  1276. rhs = rhs.reshape((n,) + y.shape[1:])
  1277. # have observation matrix & rhs, can solve the LSQ problem
  1278. cho_decomp = cholesky_banded(ab, overwrite_ab=True, lower=lower,
  1279. check_finite=check_finite)
  1280. c = cho_solve_banded((cho_decomp, lower), rhs, overwrite_b=True,
  1281. check_finite=check_finite)
  1282. c = np.ascontiguousarray(c)
  1283. return BSpline.construct_fast(t, c, k, axis=axis)
  1284. #############################
  1285. # Smoothing spline helpers #
  1286. #############################
  1287. def _compute_optimal_gcv_parameter(X, wE, y, w):
  1288. """
  1289. Returns an optimal regularization parameter from the GCV criteria [1].
  1290. Parameters
  1291. ----------
  1292. X : array, shape (5, n)
  1293. 5 bands of the design matrix ``X`` stored in LAPACK banded storage.
  1294. wE : array, shape (5, n)
  1295. 5 bands of the penalty matrix :math:`W^{-1} E` stored in LAPACK banded
  1296. storage.
  1297. y : array, shape (n,)
  1298. Ordinates.
  1299. w : array, shape (n,)
  1300. Vector of weights.
  1301. Returns
  1302. -------
  1303. lam : float
  1304. An optimal from the GCV criteria point of view regularization
  1305. parameter.
  1306. Notes
  1307. -----
  1308. No checks are performed.
  1309. References
  1310. ----------
  1311. .. [1] G. Wahba, "Estimating the smoothing parameter" in Spline models
  1312. for observational data, Philadelphia, Pennsylvania: Society for
  1313. Industrial and Applied Mathematics, 1990, pp. 45-65.
  1314. :doi:`10.1137/1.9781611970128`
  1315. """
  1316. def compute_banded_symmetric_XT_W_Y(X, w, Y):
  1317. """
  1318. Assuming that the product :math:`X^T W Y` is symmetric and both ``X``
  1319. and ``Y`` are 5-banded, compute the unique bands of the product.
  1320. Parameters
  1321. ----------
  1322. X : array, shape (5, n)
  1323. 5 bands of the matrix ``X`` stored in LAPACK banded storage.
  1324. w : array, shape (n,)
  1325. Array of weights
  1326. Y : array, shape (5, n)
  1327. 5 bands of the matrix ``Y`` stored in LAPACK banded storage.
  1328. Returns
  1329. -------
  1330. res : array, shape (4, n)
  1331. The result of the product :math:`X^T Y` stored in the banded way.
  1332. Notes
  1333. -----
  1334. As far as the matrices ``X`` and ``Y`` are 5-banded, their product
  1335. :math:`X^T W Y` is 7-banded. It is also symmetric, so we can store only
  1336. unique diagonals.
  1337. """
  1338. # compute W Y
  1339. W_Y = np.copy(Y)
  1340. W_Y[2] *= w
  1341. for i in range(2):
  1342. W_Y[i, 2 - i:] *= w[:-2 + i]
  1343. W_Y[3 + i, :-1 - i] *= w[1 + i:]
  1344. n = X.shape[1]
  1345. res = np.zeros((4, n))
  1346. for i in range(n):
  1347. for j in range(min(n-i, 4)):
  1348. res[-j-1, i + j] = sum(X[j:, i] * W_Y[:5-j, i + j])
  1349. return res
  1350. def compute_b_inv(A):
  1351. """
  1352. Inverse 3 central bands of matrix :math:`A=U^T D^{-1} U` assuming that
  1353. ``U`` is a unit upper triangular banded matrix using an algorithm
  1354. proposed in [1].
  1355. Parameters
  1356. ----------
  1357. A : array, shape (4, n)
  1358. Matrix to inverse, stored in LAPACK banded storage.
  1359. Returns
  1360. -------
  1361. B : array, shape (4, n)
  1362. 3 unique bands of the symmetric matrix that is an inverse to ``A``.
  1363. The first row is filled with zeros.
  1364. Notes
  1365. -----
  1366. The algorithm is based on the cholesky decomposition and, therefore,
  1367. in case matrix ``A`` is close to not positive defined, the function
  1368. raises LinalgError.
  1369. Both matrices ``A`` and ``B`` are stored in LAPACK banded storage.
  1370. References
  1371. ----------
  1372. .. [1] M. F. Hutchinson and F. R. de Hoog, "Smoothing noisy data with
  1373. spline functions," Numerische Mathematik, vol. 47, no. 1,
  1374. pp. 99-106, 1985.
  1375. :doi:`10.1007/BF01389878`
  1376. """
  1377. def find_b_inv_elem(i, j, U, D, B):
  1378. rng = min(3, n - i - 1)
  1379. rng_sum = 0.
  1380. if j == 0:
  1381. # use 2-nd formula from [1]
  1382. for k in range(1, rng + 1):
  1383. rng_sum -= U[-k - 1, i + k] * B[-k - 1, i + k]
  1384. rng_sum += D[i]
  1385. B[-1, i] = rng_sum
  1386. else:
  1387. # use 1-st formula from [1]
  1388. for k in range(1, rng + 1):
  1389. diag = abs(k - j)
  1390. ind = i + min(k, j)
  1391. rng_sum -= U[-k - 1, i + k] * B[-diag - 1, ind + diag]
  1392. B[-j - 1, i + j] = rng_sum
  1393. U = cholesky_banded(A)
  1394. for i in range(2, 5):
  1395. U[-i, i-1:] /= U[-1, :-i+1]
  1396. D = 1. / (U[-1])**2
  1397. U[-1] /= U[-1]
  1398. n = U.shape[1]
  1399. B = np.zeros(shape=(4, n))
  1400. for i in range(n - 1, -1, -1):
  1401. for j in range(min(3, n - i - 1), -1, -1):
  1402. find_b_inv_elem(i, j, U, D, B)
  1403. # the first row contains garbage and should be removed
  1404. B[0] = [0.] * n
  1405. return B
  1406. def _gcv(lam, X, XtWX, wE, XtE):
  1407. r"""
  1408. Computes the generalized cross-validation criteria [1].
  1409. Parameters
  1410. ----------
  1411. lam : float, (:math:`\lambda \geq 0`)
  1412. Regularization parameter.
  1413. X : array, shape (5, n)
  1414. Matrix is stored in LAPACK banded storage.
  1415. XtWX : array, shape (4, n)
  1416. Product :math:`X^T W X` stored in LAPACK banded storage.
  1417. wE : array, shape (5, n)
  1418. Matrix :math:`W^{-1} E` stored in LAPACK banded storage.
  1419. XtE : array, shape (4, n)
  1420. Product :math:`X^T E` stored in LAPACK banded storage.
  1421. Returns
  1422. -------
  1423. res : float
  1424. Value of the GCV criteria with the regularization parameter
  1425. :math:`\lambda`.
  1426. Notes
  1427. -----
  1428. Criteria is computed from the formula (1.3.2) [3]:
  1429. .. math:
  1430. GCV(\lambda) = \dfrac{1}{n} \sum\limits_{k = 1}^{n} \dfrac{ \left(
  1431. y_k - f_{\lambda}(x_k) \right)^2}{\left( 1 - \Tr{A}/n\right)^2}$.
  1432. The criteria is discussed in section 1.3 [3].
  1433. The numerator is computed using (2.2.4) [3] and the denominator is
  1434. computed using an algorithm from [2] (see in the ``compute_b_inv``
  1435. function).
  1436. References
  1437. ----------
  1438. .. [1] G. Wahba, "Estimating the smoothing parameter" in Spline models
  1439. for observational data, Philadelphia, Pennsylvania: Society for
  1440. Industrial and Applied Mathematics, 1990, pp. 45-65.
  1441. :doi:`10.1137/1.9781611970128`
  1442. .. [2] M. F. Hutchinson and F. R. de Hoog, "Smoothing noisy data with
  1443. spline functions," Numerische Mathematik, vol. 47, no. 1,
  1444. pp. 99-106, 1985.
  1445. :doi:`10.1007/BF01389878`
  1446. .. [3] E. Zemlyanoy, "Generalized cross-validation smoothing splines",
  1447. BSc thesis, 2022. Might be available (in Russian)
  1448. `here <https://www.hse.ru/ba/am/students/diplomas/620910604>`_
  1449. """
  1450. # Compute the numerator from (2.2.4) [3]
  1451. n = X.shape[1]
  1452. c = solve_banded((2, 2), X + lam * wE, y)
  1453. res = np.zeros(n)
  1454. # compute ``W^{-1} E c`` with respect to banded-storage of ``E``
  1455. tmp = wE * c
  1456. for i in range(n):
  1457. for j in range(max(0, i - n + 3), min(5, i + 3)):
  1458. res[i] += tmp[j, i + 2 - j]
  1459. numer = np.linalg.norm(lam * res)**2 / n
  1460. # compute the denominator
  1461. lhs = XtWX + lam * XtE
  1462. try:
  1463. b_banded = compute_b_inv(lhs)
  1464. # compute the trace of the product b_banded @ XtX
  1465. tr = b_banded * XtWX
  1466. tr[:-1] *= 2
  1467. # find the denominator
  1468. denom = (1 - sum(sum(tr)) / n)**2
  1469. except LinAlgError:
  1470. # cholesky decomposition cannot be performed
  1471. raise ValueError('Seems like the problem is ill-posed')
  1472. res = numer / denom
  1473. return res
  1474. n = X.shape[1]
  1475. XtWX = compute_banded_symmetric_XT_W_Y(X, w, X)
  1476. XtE = compute_banded_symmetric_XT_W_Y(X, w, wE)
  1477. def fun(lam):
  1478. return _gcv(lam, X, XtWX, wE, XtE)
  1479. gcv_est = minimize_scalar(fun, bounds=(0, n), method='Bounded')
  1480. if gcv_est.success:
  1481. return gcv_est.x
  1482. raise ValueError(f"Unable to find minimum of the GCV "
  1483. f"function: {gcv_est.message}")
  1484. def _coeff_of_divided_diff(x):
  1485. """
  1486. Returns the coefficients of the divided difference.
  1487. Parameters
  1488. ----------
  1489. x : array, shape (n,)
  1490. Array which is used for the computation of divided difference.
  1491. Returns
  1492. -------
  1493. res : array_like, shape (n,)
  1494. Coefficients of the divided difference.
  1495. Notes
  1496. -----
  1497. Vector ``x`` should have unique elements, otherwise an error division by
  1498. zero might be raised.
  1499. No checks are performed.
  1500. """
  1501. n = x.shape[0]
  1502. res = np.zeros(n)
  1503. for i in range(n):
  1504. pp = 1.
  1505. for k in range(n):
  1506. if k != i:
  1507. pp *= (x[i] - x[k])
  1508. res[i] = 1. / pp
  1509. return res
  1510. def make_smoothing_spline(x, y, w=None, lam=None):
  1511. r"""
  1512. Compute the (coefficients of) smoothing cubic spline function using
  1513. ``lam`` to control the tradeoff between the amount of smoothness of the
  1514. curve and its proximity to the data. In case ``lam`` is None, using the
  1515. GCV criteria [1] to find it.
  1516. A smoothing spline is found as a solution to the regularized weighted
  1517. linear regression problem:
  1518. .. math::
  1519. \sum\limits_{i=1}^n w_i\lvert y_i - f(x_i) \rvert^2 +
  1520. \lambda\int\limits_{x_1}^{x_n} (f^{(2)}(u))^2 d u
  1521. where :math:`f` is a spline function, :math:`w` is a vector of weights and
  1522. :math:`\lambda` is a regularization parameter.
  1523. If ``lam`` is None, we use the GCV criteria to find an optimal
  1524. regularization parameter, otherwise we solve the regularized weighted
  1525. linear regression problem with given parameter. The parameter controls
  1526. the tradeoff in the following way: the larger the parameter becomes, the
  1527. smoother the function gets.
  1528. Parameters
  1529. ----------
  1530. x : array_like, shape (n,)
  1531. Abscissas.
  1532. y : array_like, shape (n,)
  1533. Ordinates.
  1534. w : array_like, shape (n,), optional
  1535. Vector of weights. Default is ``np.ones_like(x)``.
  1536. lam : float, (:math:`\lambda \geq 0`), optional
  1537. Regularization parameter. If ``lam`` is None, then it is found from
  1538. the GCV criteria. Default is None.
  1539. Returns
  1540. -------
  1541. func : a BSpline object.
  1542. A callable representing a spline in the B-spline basis
  1543. as a solution of the problem of smoothing splines using
  1544. the GCV criteria [1] in case ``lam`` is None, otherwise using the
  1545. given parameter ``lam``.
  1546. Notes
  1547. -----
  1548. This algorithm is a clean room reimplementation of the algorithm
  1549. introduced by Woltring in FORTRAN [2]. The original version cannot be used
  1550. in SciPy source code because of the license issues. The details of the
  1551. reimplementation are discussed here (available only in Russian) [4].
  1552. If the vector of weights ``w`` is None, we assume that all the points are
  1553. equal in terms of weights, and vector of weights is vector of ones.
  1554. Note that in weighted residual sum of squares, weights are not squared:
  1555. :math:`\sum\limits_{i=1}^n w_i\lvert y_i - f(x_i) \rvert^2` while in
  1556. ``splrep`` the sum is built from the squared weights.
  1557. In cases when the initial problem is ill-posed (for example, the product
  1558. :math:`X^T W X` where :math:`X` is a design matrix is not a positive
  1559. defined matrix) a ValueError is raised.
  1560. References
  1561. ----------
  1562. .. [1] G. Wahba, "Estimating the smoothing parameter" in Spline models for
  1563. observational data, Philadelphia, Pennsylvania: Society for Industrial
  1564. and Applied Mathematics, 1990, pp. 45-65.
  1565. :doi:`10.1137/1.9781611970128`
  1566. .. [2] H. J. Woltring, A Fortran package for generalized, cross-validatory
  1567. spline smoothing and differentiation, Advances in Engineering
  1568. Software, vol. 8, no. 2, pp. 104-113, 1986.
  1569. :doi:`10.1016/0141-1195(86)90098-7`
  1570. .. [3] T. Hastie, J. Friedman, and R. Tisbshirani, "Smoothing Splines" in
  1571. The elements of Statistical Learning: Data Mining, Inference, and
  1572. prediction, New York: Springer, 2017, pp. 241-249.
  1573. :doi:`10.1007/978-0-387-84858-7`
  1574. .. [4] E. Zemlyanoy, "Generalized cross-validation smoothing splines",
  1575. BSc thesis, 2022.
  1576. `<https://www.hse.ru/ba/am/students/diplomas/620910604>`_ (in
  1577. Russian)
  1578. Examples
  1579. --------
  1580. Generate some noisy data
  1581. >>> import numpy as np
  1582. >>> np.random.seed(1234)
  1583. >>> n = 200
  1584. >>> def func(x):
  1585. ... return x**3 + x**2 * np.sin(4 * x)
  1586. >>> x = np.sort(np.random.random_sample(n) * 4 - 2)
  1587. >>> y = func(x) + np.random.normal(scale=1.5, size=n)
  1588. Make a smoothing spline function
  1589. >>> from scipy.interpolate import make_smoothing_spline
  1590. >>> spl = make_smoothing_spline(x, y)
  1591. Plot both
  1592. >>> import matplotlib.pyplot as plt
  1593. >>> grid = np.linspace(x[0], x[-1], 400)
  1594. >>> plt.plot(grid, spl(grid), label='Spline')
  1595. >>> plt.plot(grid, func(grid), label='Original function')
  1596. >>> plt.scatter(x, y, marker='.')
  1597. >>> plt.legend(loc='best')
  1598. >>> plt.show()
  1599. """
  1600. x = np.ascontiguousarray(x, dtype=float)
  1601. y = np.ascontiguousarray(y, dtype=float)
  1602. if any(x[1:] - x[:-1] <= 0):
  1603. raise ValueError('``x`` should be an ascending array')
  1604. if x.ndim != 1 or y.ndim != 1 or x.shape[0] != y.shape[0]:
  1605. raise ValueError('``x`` and ``y`` should be one dimensional and the'
  1606. ' same size')
  1607. if w is None:
  1608. w = np.ones(len(x))
  1609. else:
  1610. w = np.ascontiguousarray(w)
  1611. if any(w <= 0):
  1612. raise ValueError('Invalid vector of weights')
  1613. t = np.r_[[x[0]] * 3, x, [x[-1]] * 3]
  1614. n = x.shape[0]
  1615. # It is known that the solution to the stated minimization problem exists
  1616. # and is a natural cubic spline with vector of knots equal to the unique
  1617. # elements of ``x`` [3], so we will solve the problem in the basis of
  1618. # natural splines.
  1619. # create design matrix in the B-spline basis
  1620. X_bspl = BSpline.design_matrix(x, t, 3)
  1621. # move from B-spline basis to the basis of natural splines using equations
  1622. # (2.1.7) [4]
  1623. # central elements
  1624. X = np.zeros((5, n))
  1625. for i in range(1, 4):
  1626. X[i, 2: -2] = X_bspl[i: i - 4, 3: -3][np.diag_indices(n - 4)]
  1627. # first elements
  1628. X[1, 1] = X_bspl[0, 0]
  1629. X[2, :2] = ((x[2] + x[1] - 2 * x[0]) * X_bspl[0, 0],
  1630. X_bspl[1, 1] + X_bspl[1, 2])
  1631. X[3, :2] = ((x[2] - x[0]) * X_bspl[1, 1], X_bspl[2, 2])
  1632. # last elements
  1633. X[1, -2:] = (X_bspl[-3, -3], (x[-1] - x[-3]) * X_bspl[-2, -2])
  1634. X[2, -2:] = (X_bspl[-2, -3] + X_bspl[-2, -2],
  1635. (2 * x[-1] - x[-2] - x[-3]) * X_bspl[-1, -1])
  1636. X[3, -2] = X_bspl[-1, -1]
  1637. # create penalty matrix and divide it by vector of weights: W^{-1} E
  1638. wE = np.zeros((5, n))
  1639. wE[2:, 0] = _coeff_of_divided_diff(x[:3]) / w[:3]
  1640. wE[1:, 1] = _coeff_of_divided_diff(x[:4]) / w[:4]
  1641. for j in range(2, n - 2):
  1642. wE[:, j] = (x[j+2] - x[j-2]) * _coeff_of_divided_diff(x[j-2:j+3])\
  1643. / w[j-2: j+3]
  1644. wE[:-1, -2] = -_coeff_of_divided_diff(x[-4:]) / w[-4:]
  1645. wE[:-2, -1] = _coeff_of_divided_diff(x[-3:]) / w[-3:]
  1646. wE *= 6
  1647. if lam is None:
  1648. lam = _compute_optimal_gcv_parameter(X, wE, y, w)
  1649. elif lam < 0.:
  1650. raise ValueError('Regularization parameter should be non-negative')
  1651. # solve the initial problem in the basis of natural splines
  1652. c = solve_banded((2, 2), X + lam * wE, y)
  1653. # move back to B-spline basis using equations (2.2.10) [4]
  1654. c_ = np.r_[c[0] * (t[5] + t[4] - 2 * t[3]) + c[1],
  1655. c[0] * (t[5] - t[3]) + c[1],
  1656. c[1: -1],
  1657. c[-1] * (t[-4] - t[-6]) + c[-2],
  1658. c[-1] * (2 * t[-4] - t[-5] - t[-6]) + c[-2]]
  1659. return BSpline.construct_fast(t, c_, 3)