_fitpack_repro.py 51 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339
  1. """ Replicate FITPACK's logic for constructing smoothing spline functions and curves.
  2. Currently provides analogs of splrep and splprep python routines, i.e.
  3. curfit.f and parcur.f routines (the drivers are fpcurf.f and fppara.f, respectively)
  4. The Fortran sources are from
  5. https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/
  6. .. [1] P. Dierckx, "Algorithms for smoothing data with periodic and
  7. parametric splines, Computer Graphics and Image Processing",
  8. 20 (1982) 171-184.
  9. :doi:`10.1016/0146-664X(82)90043-0`.
  10. .. [2] P. Dierckx, "Curve and surface fitting with splines", Monographs on
  11. Numerical Analysis, Oxford University Press, 1993.
  12. .. [3] P. Dierckx, "An algorithm for smoothing, differentiation and integration
  13. of experimental data using spline functions",
  14. Journal of Computational and Applied Mathematics, vol. I, no 3, p. 165 (1975).
  15. https://doi.org/10.1016/0771-050X(75)90034-0
  16. """
  17. import warnings
  18. import operator
  19. import numpy as np
  20. from scipy._lib._array_api import array_namespace, concat_1d, xp_capabilities
  21. from ._bsplines import (
  22. _not_a_knot, make_interp_spline, BSpline, fpcheck, _lsq_solve_qr,
  23. _lsq_solve_qr_for_root_rati_periodic, _periodic_knots
  24. )
  25. from . import _dierckx # type: ignore[attr-defined]
  26. # cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
  27. # c part 1: determination of the number of knots and their position c
  28. # c ************************************************************** c
  29. #
  30. # https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpcurf.f#L31
  31. # Hardcoded in curfit.f
  32. TOL = 0.001
  33. MAXIT = 20
  34. def _get_residuals(x, y, t, k, w, periodic=False):
  35. # inline the relevant part of
  36. # >>> spl = make_lsq_spline(x, y, w=w2, t=t, k=k)
  37. # NB:
  38. # 1. y is assumed to be 2D here. For 1D case (parametric=False),
  39. # the call must have been preceded by y = y[:, None] (cf _validate_inputs)
  40. # 2. We always sum the squares across axis=1:
  41. # * For 1D (parametric=False), the last dimension has size one,
  42. # so the summation is a no-op.
  43. # * For 2D (parametric=True), the summation is actually how the
  44. # 'residuals' are defined, see Eq. (42) in Dierckx1982
  45. # (the reference is in the docstring of `class F`) below.
  46. _, _, _, fp, residuals = _lsq_solve_qr(x, y, t, k, w, periodic=periodic)
  47. if np.isnan(residuals.sum()):
  48. raise ValueError(_iermesg[1])
  49. return residuals, fp
  50. def _validate_bc_type(bc_type):
  51. if bc_type is None:
  52. return "not-a-knot"
  53. if bc_type not in ("not-a-knot", "periodic"):
  54. raise ValueError("Only 'not-a-knot' and 'periodic' "
  55. "boundary conditions are recognised, "
  56. f"found {bc_type}")
  57. return bc_type
  58. def add_knot(x, t, k, residuals):
  59. """Add a new knot.
  60. (Approximately) replicate FITPACK's logic:
  61. 1. split the `x` array into knot intervals, ``t(j+k) <= x(i) <= t(j+k+1)``
  62. 2. find the interval with the maximum sum of residuals
  63. 3. insert a new knot into the middle of that interval.
  64. NB: a new knot is in fact an `x` value at the middle of the interval.
  65. So *the knots are a subset of `x`*.
  66. This routine is an analog of
  67. https://github.com/scipy/scipy/blob/v1.11.4/scipy/interpolate/fitpack/fpcurf.f#L190-L215
  68. (cf _split function)
  69. and https://github.com/scipy/scipy/blob/v1.11.4/scipy/interpolate/fitpack/fpknot.f
  70. """
  71. new_knot = _dierckx.fpknot(x, t, k, residuals)
  72. idx_t = np.searchsorted(t, new_knot)
  73. t_new = np.r_[t[:idx_t], new_knot, t[idx_t:]]
  74. return t_new
  75. def _validate_inputs(x, y, w, k, s, xb, xe, parametric, periodic=False):
  76. """Common input validations for generate_knots and make_splrep.
  77. """
  78. x = np.asarray(x, dtype=float)
  79. y = np.asarray(y, dtype=float)
  80. if not x.flags.c_contiguous:
  81. x = x.copy()
  82. if not y.flags.c_contiguous:
  83. y = y.copy()
  84. if w is None:
  85. w = np.ones_like(x, dtype=float)
  86. else:
  87. w = np.asarray(w, dtype=float)
  88. if not w.flags.c_contiguous:
  89. w = w.copy()
  90. if w.ndim != 1:
  91. raise ValueError(f"{w.ndim = } not implemented yet.")
  92. if (w < 0).any():
  93. raise ValueError("Weights must be non-negative")
  94. if w.sum() == 0:
  95. raise ValueError("All weights are zero.")
  96. if y.ndim == 0 or y.ndim > 2:
  97. raise ValueError(f"{y.ndim = } not supported (must be 1 or 2.)")
  98. parametric = bool(parametric)
  99. if parametric:
  100. if y.ndim != 2:
  101. raise ValueError(f"{y.ndim = } != 2 not supported with {parametric =}.")
  102. else:
  103. if y.ndim != 1:
  104. raise ValueError(f"{y.ndim = } != 1 not supported with {parametric =}.")
  105. # all _impl functions expect y.ndim = 2
  106. y = y[:, None]
  107. if w.shape[0] != x.shape[0]:
  108. raise ValueError(f"Weights is incompatible: {w.shape =} != {x.shape}.")
  109. if x.shape[0] != y.shape[0]:
  110. raise ValueError(f"Data is incompatible: {x.shape = } and {y.shape = }.")
  111. if x.ndim != 1 or (x[1:] < x[:-1]).any():
  112. raise ValueError("Expect `x` to be an ordered 1D sequence.")
  113. k = operator.index(k)
  114. if s < 0:
  115. raise ValueError(f"`s` must be non-negative. Got {s = }")
  116. if xb is None:
  117. xb = min(x)
  118. if xe is None:
  119. xe = max(x)
  120. if periodic and not np.allclose(y[0], y[-1], atol=1e-15):
  121. raise ValueError("First and last points does not match which is required "
  122. "for `bc_type='periodic'`.")
  123. return x, y, w, k, s, xb, xe
  124. @xp_capabilities(cpu_only=True, jax_jit=False, allow_dask_compute=True)
  125. def generate_knots(x, y, *, w=None, xb=None, xe=None,
  126. k=3, s=0, nest=None, bc_type=None):
  127. """Generate knot vectors until the Least SQuares (LSQ) criterion is satified.
  128. Parameters
  129. ----------
  130. x, y : array_like
  131. The data points defining the curve ``y = f(x)``.
  132. w : array_like, optional
  133. Weights.
  134. xb : float, optional
  135. The boundary of the approximation interval. If None (default),
  136. is set to ``x[0]``.
  137. xe : float, optional
  138. The boundary of the approximation interval. If None (default),
  139. is set to ``x[-1]``.
  140. k : int, optional
  141. The spline degree. Default is cubic, ``k = 3``.
  142. s : float, optional
  143. The smoothing factor. Default is ``s = 0``.
  144. nest : int, optional
  145. Stop when at least this many knots are placed.
  146. ends are equivalent.
  147. bc_type : str, optional
  148. Boundary conditions.
  149. Default is `"not-a-knot"`.
  150. The following boundary conditions are recognized:
  151. * ``"not-a-knot"`` (default): The first and second segments are the
  152. same polynomial. This is equivalent to having ``bc_type=None``.
  153. * ``"periodic"``: The values and the first ``k-1`` derivatives at the
  154. ends are equivalent.
  155. Yields
  156. ------
  157. t : ndarray
  158. Knot vectors with an increasing number of knots.
  159. The generator is finite: it stops when the smoothing critetion is
  160. satisfied, or when then number of knots exceeds the maximum value:
  161. the user-provided `nest` or `x.size + k + 1` --- which is the knot vector
  162. for the interpolating spline.
  163. Examples
  164. --------
  165. Generate some noisy data and fit a sequence of LSQ splines:
  166. >>> import numpy as np
  167. >>> import matplotlib.pyplot as plt
  168. >>> from scipy.interpolate import make_lsq_spline, generate_knots
  169. >>> rng = np.random.default_rng(12345)
  170. >>> x = np.linspace(-3, 3, 50)
  171. >>> y = np.exp(-x**2) + 0.1 * rng.standard_normal(size=50)
  172. >>> knots = list(generate_knots(x, y, s=1e-10))
  173. >>> for t in knots[::3]:
  174. ... spl = make_lsq_spline(x, y, t)
  175. ... xs = xs = np.linspace(-3, 3, 201)
  176. ... plt.plot(xs, spl(xs), '-', label=f'n = {len(t)}', lw=3, alpha=0.7)
  177. >>> plt.plot(x, y, 'o', label='data')
  178. >>> plt.plot(xs, np.exp(-xs**2), '--')
  179. >>> plt.legend()
  180. Note that increasing the number of knots make the result follow the data
  181. more and more closely.
  182. Also note that a step of the generator may add multiple knots:
  183. >>> [len(t) for t in knots]
  184. [8, 9, 10, 12, 16, 24, 40, 48, 52, 54]
  185. Notes
  186. -----
  187. The routine generates successive knots vectors of increasing length, starting
  188. from ``2*(k+1)`` to ``len(x) + k + 1``, trying to make knots more dense
  189. in the regions where the deviation of the LSQ spline from data is large.
  190. When the maximum number of knots, ``len(x) + k + 1`` is reached
  191. (this happens when ``s`` is small and ``nest`` is large), the generator
  192. stops, and the last output is the knots for the interpolation with the
  193. not-a-knot boundary condition.
  194. Knots are located at data sites, unless ``k`` is even and the number of knots
  195. is ``len(x) + k + 1``. In that case, the last output of the generator
  196. has internal knots at Greville sites, ``(x[1:] + x[:-1]) / 2``.
  197. .. versionadded:: 1.15.0
  198. """
  199. xp = array_namespace(x, y, w)
  200. bc_type = _validate_bc_type(bc_type)
  201. periodic = bc_type == 'periodic'
  202. if s == 0:
  203. if nest is not None or w is not None:
  204. raise ValueError("s == 0 is interpolation only")
  205. # For special-case k=1 (e.g., Lyche and Morken, Eq.(2.16)),
  206. # _not_a_knot produces desired knot vector
  207. if periodic and k > 1:
  208. t = _periodic_knots(x, k)
  209. else:
  210. t = _not_a_knot(x, k)
  211. yield xp.asarray(t)
  212. return
  213. x, y, w, k, s, xb, xe = _validate_inputs(
  214. x, y, w, k, s, xb, xe, parametric=np.ndim(y) == 2,
  215. periodic=periodic
  216. )
  217. yield from _generate_knots_impl(x, y, w, xb, xe, k, s, nest, periodic, xp=xp)
  218. def _generate_knots_impl(x, y, w, xb, xe, k, s, nest, periodic, xp=np):
  219. acc = s * TOL
  220. m = x.size # the number of data points
  221. if nest is None:
  222. # the max number of knots. This is set in _fitpack_impl.py line 274
  223. # and fitpack.pyf line 198
  224. # Ref: https://github.com/scipy/scipy/blob/596b586e25e34bd842b575bac134b4d6924c6556/scipy/interpolate/_fitpack_impl.py#L260-L263
  225. if periodic:
  226. nest = max(m + 2*k, 2*k + 3)
  227. else:
  228. nest = max(m + k + 1, 2*k + 3)
  229. else:
  230. if nest < 2*(k + 1):
  231. raise ValueError(f"`nest` too small: {nest = } < 2*(k+1) = {2*(k+1)}.")
  232. if not periodic:
  233. nmin = 2*(k + 1) # the number of knots for an LSQ polynomial approximation
  234. nmax = m + k + 1 # the number of knots for the spline interpolation
  235. else:
  236. # Ref: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpperi.f#L54
  237. nmin = 2*(k + 1) # the number of knots for an LSQ polynomial approximation
  238. # Ref: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpperi.f#L61
  239. nmax = m + 2*k # the number of knots for the spline interpolation
  240. per = xe - xb
  241. # Ref: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpperi.f#L107-L123
  242. # Computes fp0 for constant function
  243. if periodic:
  244. t = np.zeros(nmin, dtype=float)
  245. for i in range(0, k + 1):
  246. t[i] = x[0] - (k - i) * per
  247. t[i + k + 1] = x[m - 1] + i * per
  248. _, fp = _get_residuals(x, y, t, k, w, periodic=periodic)
  249. # For periodic splines, check whether constant function
  250. # satisfies accuracy criterion
  251. # Also if maximal number of nodes is equal to the minimal
  252. # then constant function is the direct solution
  253. # Ref: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpperi.f#L600-L610
  254. if fp - s < acc or nmax == nmin:
  255. yield t
  256. return
  257. else:
  258. fp = 0.0
  259. fpold = 0.0
  260. # start from no internal knots
  261. if not periodic:
  262. t = np.asarray([xb]*(k+1) + [xe]*(k+1), dtype=float)
  263. else:
  264. # Ref: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpperi.f#L131
  265. # Initialize knot vector `t` of size (2k + 3) with zeros.
  266. # The central knot `t[k + 1]` is seeded with the midpoint value from `x`.
  267. # Note that, in the `if periodic:` block (in the main loop below),
  268. # the boundary knots `t[k]` and `t[n - k - 1]` are set to the endpoints `xb` and
  269. # `xe`. Then, the surrounding knots on both ends are updated to ensure
  270. # periodicity.
  271. # Specifically:
  272. # - Left-side knots are mirrored from the right end minus the period (`per`).
  273. # - Right-side knots are mirrored from the left end plus the period.
  274. # These updates ensure that the knot vector wraps around correctly for periodic
  275. # B-spline fitting.
  276. t = np.zeros(2*k + 3, dtype=float)
  277. t[k + 1] = x[(m + 1)//2 - 1]
  278. nplus = 1
  279. n = t.shape[0]
  280. # c main loop for the different sets of knots. m is a safe upper bound
  281. # c for the number of trials.
  282. for iter in range(m):
  283. # Ref: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpperi.f#L147-L158
  284. if periodic:
  285. n = t.shape[0]
  286. t[k] = xb
  287. t[n - k - 1] = xe
  288. for j in range(1, k + 1):
  289. t[k - j] = t[n - k - j - 1] - per
  290. t[n - k + j - 1] = t[k + j] + per
  291. yield xp.asarray(t)
  292. # construct the LSQ spline with this set of knots
  293. fpold = fp
  294. residuals, fp = _get_residuals(x, y, t, k, w=w,
  295. periodic=periodic)
  296. fpms = fp - s
  297. # c test whether the approximation sinf(x) is an acceptable solution.
  298. # c if f(p=inf) < s accept the choice of knots.
  299. if (abs(fpms) < acc) or (fpms < 0):
  300. return
  301. # ### c increase the number of knots. ###
  302. # c determine the number of knots nplus we are going to add.
  303. if n == nmin:
  304. # the first iteration
  305. nplus = 1
  306. else:
  307. delta = fpold - fp
  308. npl1 = int(nplus * fpms / delta) if delta > acc else nplus*2
  309. nplus = min(nplus*2, max(npl1, nplus//2, 1))
  310. # actually add knots
  311. for j in range(nplus):
  312. t = add_knot(x, t, k, residuals)
  313. # check if we have enough knots already
  314. n = t.shape[0]
  315. # c if n = nmax, sinf(x) is an interpolating spline.
  316. # c if n=nmax we locate the knots as for interpolation.
  317. if n >= nmax:
  318. if not periodic:
  319. t = _not_a_knot(x, k)
  320. else:
  321. t = _periodic_knots(x, k)
  322. yield xp.asarray(t)
  323. return
  324. # c if n=nest we cannot increase the number of knots because of
  325. # c the storage capacity limitation.
  326. if n >= nest:
  327. yield xp.asarray(t)
  328. return
  329. # recompute if needed
  330. if j < nplus - 1:
  331. residuals, _ = _get_residuals(x, y, t, k, w=w, periodic=periodic)
  332. # this should never be reached
  333. return
  334. # cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
  335. # c part 2: determination of the smoothing spline sp(x). c
  336. # c *************************************************** c
  337. # c we have determined the number of knots and their position. c
  338. # c we now compute the b-spline coefficients of the smoothing spline c
  339. # c sp(x). the observation matrix a is extended by the rows of matrix c
  340. # c b expressing that the kth derivative discontinuities of sp(x) at c
  341. # c the interior knots t(k+2),...t(n-k-1) must be zero. the corres- c
  342. # c ponding weights of these additional rows are set to 1/p. c
  343. # c iteratively we then have to determine the value of p such that c
  344. # c f(p)=sum((w(i)*(y(i)-sp(x(i))))**2) be = s. we already know that c
  345. # c the least-squares kth degree polynomial corresponds to p=0, and c
  346. # c that the least-squares spline corresponds to p=infinity. the c
  347. # c iteration process which is proposed here, makes use of rational c
  348. # c interpolation. since f(p) is a convex and strictly decreasing c
  349. # c function of p, it can be approximated by a rational function c
  350. # c r(p) = (u*p+v)/(p+w). three values of p(p1,p2,p3) with correspond- c
  351. # c ing values of f(p) (f1=f(p1)-s,f2=f(p2)-s,f3=f(p3)-s) are used c
  352. # c to calculate the new value of p such that r(p)=s. convergence is c
  353. # c guaranteed by taking f1>0 and f3<0. c
  354. # cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
  355. def prodd(t, i, j, k):
  356. res = 1.0
  357. for s in range(k+2):
  358. if i + s != j:
  359. res *= (t[j] - t[i+s])
  360. return res
  361. def disc(t, k):
  362. """Discontinuity matrix: jumps of k-th derivatives of b-splines at internal knots.
  363. See Eqs. (9)-(10) of Ref. [1], or, equivalently, Eq. (3.43) of Ref. [2].
  364. This routine assumes internal knots are all simple (have multiplicity =1).
  365. Parameters
  366. ----------
  367. t : ndarray, 1D, shape(n,)
  368. Knots.
  369. k : int
  370. The spline degree
  371. Returns
  372. -------
  373. disc : ndarray, shape(n-2*k-1, k+2)
  374. The jumps of the k-th derivatives of b-splines at internal knots,
  375. ``t[k+1], ...., t[n-k-1]``.
  376. offset : ndarray, shape(2-2*k-1,)
  377. Offsets
  378. nc : int
  379. Notes
  380. -----
  381. The normalization here follows FITPACK:
  382. (https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpdisc.f#L36)
  383. The k-th derivative jumps are multiplied by a factor::
  384. (delta / nrint)**k / k!
  385. where ``delta`` is the length of the interval spanned by internal knots, and
  386. ``nrint`` is one less the number of internal knots (i.e., the number of
  387. subintervals between them).
  388. References
  389. ----------
  390. .. [1] Paul Dierckx, Algorithms for smoothing data with periodic and parametric
  391. splines, Computer Graphics and Image Processing, vol. 20, p. 171 (1982).
  392. :doi:`10.1016/0146-664X(82)90043-0`
  393. .. [2] Tom Lyche and Knut Morken, Spline methods,
  394. http://www.uio.no/studier/emner/matnat/ifi/INF-MAT5340/v05/undervisningsmateriale/
  395. """
  396. n = t.shape[0]
  397. # the length of the base interval spanned by internal knots & the number
  398. # of subintervas between these internal knots
  399. delta = t[n - k - 1] - t[k]
  400. nrint = n - 2*k - 1
  401. matr = np.empty((nrint - 1, k + 2), dtype=float)
  402. for jj in range(nrint - 1):
  403. j = jj + k + 1
  404. for ii in range(k + 2):
  405. i = jj + ii
  406. matr[jj, ii] = (t[i + k + 1] - t[i]) / prodd(t, i, j, k)
  407. # NB: equivalent to
  408. # row = [(t[i + k + 1] - t[i]) / prodd(t, i, j, k) for i in range(j-k-1, j+1)]
  409. # assert (matr[j-k-1, :] == row).all()
  410. # follow FITPACK
  411. matr *= (delta/ nrint)**k
  412. # make it packed
  413. offset = np.array([i for i in range(nrint-1)], dtype=np.int64)
  414. nc = n - k - 1
  415. return matr, offset, nc
  416. class F:
  417. """ The r.h.s. of ``f(p) = s``.
  418. Given scalar `p`, we solve the system of equations in the LSQ sense:
  419. | A | @ | c | = | y |
  420. | B / p | | 0 | | 0 |
  421. where `A` is the matrix of b-splines and `b` is the discontinuity matrix
  422. (the jumps of the k-th derivatives of b-spline basis elements at knots).
  423. Since we do that repeatedly while minimizing over `p`, we QR-factorize
  424. `A` only once and update the QR factorization only of the `B` rows of the
  425. augmented matrix |A, B/p|.
  426. The system of equations is Eq. (15) Ref. [1]_, the strategy and implementation
  427. follows that of FITPACK, see specific links below.
  428. References
  429. ----------
  430. [1] P. Dierckx, Algorithms for Smoothing Data with Periodic and Parametric Splines,
  431. COMPUTER GRAPHICS AND IMAGE PROCESSING vol. 20, pp 171-184 (1982.)
  432. https://doi.org/10.1016/0146-664X(82)90043-0
  433. """
  434. def __init__(self, x, y, t, k, s, w=None, *, R=None, Y=None):
  435. self.x = x
  436. self.y = y
  437. self.t = t
  438. self.k = k
  439. w = np.ones_like(x, dtype=float) if w is None else w
  440. if w.ndim != 1:
  441. raise ValueError(f"{w.ndim = } != 1.")
  442. self.w = w
  443. self.s = s
  444. if y.ndim != 2:
  445. raise ValueError(f"F: expected y.ndim == 2, got {y.ndim = } instead.")
  446. # ### precompute what we can ###
  447. # https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpcurf.f#L250
  448. # c evaluate the discontinuity jump of the kth derivative of the
  449. # c b-splines at the knots t(l),l=k+2,...n-k-1 and store in b.
  450. b, b_offset, b_nc = disc(t, k)
  451. # the QR factorization of the data matrix, if not provided
  452. # NB: otherwise, must be consistent with x,y & s, but this is not checked
  453. if R is None and Y is None:
  454. R, Y, _, _, _ = _lsq_solve_qr(x, y, t, k, w)
  455. # prepare to combine R and the discontinuity matrix (AB); also r.h.s. (YY)
  456. # https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpcurf.f#L269
  457. # c the rows of matrix b with weight 1/p are rotated into the
  458. # c triangularised observation matrix a which is stored in g.
  459. nc = t.shape[0] - k - 1
  460. nz = k + 1
  461. if R.shape[1] != nz:
  462. raise ValueError(f"Internal error: {R.shape[1] =} != {k+1 =}.")
  463. # r.h.s. of the augmented system
  464. z = np.zeros((b.shape[0], Y.shape[1]), dtype=float)
  465. self.YY = np.r_[Y[:nc], z]
  466. # l.h.s. of the augmented system
  467. AA = np.zeros((nc + b.shape[0], self.k+2), dtype=float)
  468. AA[:nc, :nz] = R[:nc, :]
  469. # AA[nc:, :] = b.a / p # done in __call__(self, p)
  470. self.AA = AA
  471. self.offset = np.r_[np.arange(nc, dtype=np.int64), b_offset]
  472. self.nc = nc
  473. self.b = b
  474. def __call__(self, p):
  475. # https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpcurf.f#L279
  476. # c the row of matrix b is rotated into triangle by givens transformation
  477. # copy the precomputed matrices over for in-place work
  478. # R = PackedMatrix(self.AB.a.copy(), self.AB.offset.copy(), nc)
  479. AB = self.AA.copy()
  480. offset = self.offset.copy()
  481. nc = self.nc
  482. AB[nc:, :] = self.b / p
  483. QY = self.YY.copy()
  484. # heavy lifting happens here, in-place
  485. _dierckx.qr_reduce(AB, offset, nc, QY, startrow=nc)
  486. # solve for the coefficients
  487. c, _, fp = _dierckx.fpback(AB, nc, self.x, self.y, self.t, self.k, self.w, QY)
  488. spl = BSpline(self.t, c, self.k)
  489. self.spl = spl # store it
  490. return fp - self.s
  491. class Fperiodic:
  492. """
  493. Fit a smooth periodic B-spline curve to given data points.
  494. This class fits a periodic B-spline curve S(t) of degree k through data points
  495. (x, y) with knots t. The spline is smooth and repeats itself at the start and
  496. end, meaning the function and its derivatives up to order k-1 are equal
  497. at the boundaries.
  498. We want to find spline coefficients c that minimize the difference between
  499. the spline and the data, while also keeping the spline smooth. This is done
  500. by solving:
  501. minimize || W^{1/2} (Y - B c) ||^2 + s * c^T @ R @ c
  502. subject to periodic constraints on c.
  503. where:
  504. - Y is the data values,
  505. - B is the matrix of B-spline basis functions at points x,
  506. - W is a weighting matrix for the data points,
  507. - s is the smoothing parameter (larger s means smoother curve),
  508. - R is a matrix that penalizes wiggliness of the spline,
  509. - c spline coefficients to be solved for
  510. - periodic constraints ensure the spline repeats smoothly.
  511. The solution is obtained by forming augmented matrices and performing
  512. a QR factorization that incorporates these constraints, following the
  513. approach in FITPACK's `fpperi.f`.
  514. Parameters:
  515. -----------
  516. x : array_like, shape (n,)
  517. y : array_like, shape (n, m)
  518. t : array_like, shape (nt,)
  519. Knot vector for the spline
  520. k : int
  521. Degree of the spline.
  522. s : float
  523. Controls smoothness: bigger s means smoother curve, smaller s fits data closer.
  524. w : array_like, shape (n,), optional
  525. Weights for data points. Defaults to all ones.
  526. R, Y, A1, A2, Z : arrays, optional
  527. Precomputed matrices from least squares and QR factorization steps to speed up
  528. repeated fits with the same knots and data.
  529. Attributes:
  530. -----------
  531. G1, G2 : arrays
  532. Augmented matrices combining the original QR factors and constraints related to
  533. the spline basis and data. G1 is roughly the "upper-triangular" part;
  534. G2 contains additional constraint information for periodicity.
  535. H1, H2 : arrays
  536. Matrices associated with the discontinuity jump constraints of the k-th
  537. derivative of B-splines at the knots. These encode the periodicity
  538. conditions and are scaled by the smoothing parameter.
  539. offset : array
  540. Offset indices used for efficient indexing during QR reduction.
  541. Methods:
  542. --------
  543. __call__(p):
  544. Perform QR reduction of augmented matrices scaled by 1/p, solve for spline
  545. coefficients, and return residual difference fp - s.
  546. References:
  547. -----------
  548. - FITPACK's fpperi.f and fpcurf.f Fortran routines for periodic spline fitting.
  549. """
  550. def __init__(self, x, y, t, k, s, w=None, *,
  551. R=None, Y=None, A1=None, A2=None, Z=None):
  552. # Initialize the class with input data points x, y,
  553. # knot vector t, spline degree k, smoothing factor s,
  554. # optional weights w, and optionally precomputed matrices.
  555. self.x = x
  556. self.y = y
  557. self.t = t
  558. self.k = k
  559. # If weights not provided, default to uniform weights (all ones)
  560. w = np.ones_like(x, dtype=float) if w is None else w
  561. if w.ndim != 1:
  562. raise ValueError(f"{w.ndim = } != 1.")
  563. self.w = w
  564. self.s = s
  565. if y.ndim != 2:
  566. raise ValueError(f"F: expected y.ndim == 2, got {y.ndim = } instead.")
  567. # ### precompute matrices and factors needed for spline fitting ###
  568. # Compute the discontinuity jump vector 'b' for the k-th derivative
  569. # of B-splines at internal knots. This is needed for enforcing
  570. # periodicity constraints. The jump discontinuities are stored in b.
  571. # Refer: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpcurf.f#L252
  572. b, b_offset, b_nc = disc(t, k)
  573. # If QR factorization or auxiliary matrices (A1, A2, Z) are not provided,
  574. # compute them via least squares on the data (x,y) with weights w.
  575. # These matrices come from fitting B-spline basis to data.
  576. if ((A1 is None or A2 is None or Z is None) or
  577. (R is None and Y is None)):
  578. # _lsq_solve_qr_for_root_rati_periodic computes QR factorization of
  579. # data matrix and returns related matrices.
  580. # Refer: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpperi.f#L171-L215
  581. R, A1, A2, Z, _, _, _, _ = _lsq_solve_qr_for_root_rati_periodic(
  582. x, y, t, k, w)
  583. # Initialize augmented matrices G1, G2, H1, H2 and offset used
  584. # for periodic B-spline fitting with constraints.
  585. # This calls the C++ function `init_augmented_matrices` to set
  586. # up these matrices based on A1, A2, and discontinuity vector b.
  587. # Refer: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpperi.f#L441-L493
  588. G1, G2, H1, H2, offset = _dierckx.init_augmented_matrices(A1, A2, b, len(t), k)
  589. # Store these matrices as class attributes for use in evaluation
  590. self.G1_ = G1
  591. self.G2_ = G2
  592. self.H1_ = H1
  593. self.H2_ = H2
  594. self.Z_ = Z
  595. self.offset_ = offset
  596. def __call__(self, p):
  597. # Create copies of the augmented matrices to avoid overwriting originals
  598. G1 = self.G1_.copy()
  599. G2 = self.G2_.copy()
  600. H1 = self.H1_.copy()
  601. H2 = self.H2_.copy()
  602. Z = self.Z_.copy()
  603. # Scale H1 and H2 by the inverse of p (1/p), applying the pinv multiplier
  604. # This scaling is required before QR reduction step to control smoothing.
  605. pinv = 1/p
  606. H1 = H1 * pinv
  607. H2 = H2 * pinv
  608. # Initialize vector c for coefficients with shape compatible with matrices.
  609. # The first part of c is set from Z, which is related to least squares fit.
  610. c = np.empty((len(self.t) - self.k - 1, Z.shape[1]), dtype=Z.dtype)
  611. c[:len(self.t) - 2*self.k - 1, :] = Z[:len(self.t) - 2*self.k - 1, :]
  612. # Perform QR factorization reduction on the augmented matrices
  613. # This corresponds to the C++ function qr_reduce_augmented_matrices.
  614. # It applies Givens rotations to eliminate entries and
  615. # reduces the problem dimension by accounting for constraints.
  616. # Refer: https://github.com/scipy/scipy/blob/maintenance/1.16.x/scipy/interpolate/fitpack/fpperi.f#L498-L534
  617. _dierckx.qr_reduce_augmented_matrices(
  618. G1, G2, H1, H2, c, self.offset_, len(self.t), self.k)
  619. # Solve for the B-spline coefficients by backward substitution
  620. # using the reduced matrices. This corresponds to the fpbacp routine in fitpack.
  621. c, _, fp = _dierckx.fpbacp(
  622. G1, G2, c,
  623. self.k, self.k + 1, self.x[:-1], self.y[:-1, :],
  624. self.t, self.w[:-1])
  625. # Construct a BSpline object using knot vector t, coefficients c, and degree k
  626. spl = BSpline(self.t, c, self.k)
  627. self.spl = spl # store it in the object
  628. # Return the difference between the final residual fp and smoothing parameter s
  629. # This quantity is used to assess fit quality in root-finding procedures.
  630. return fp - self.s
  631. def fprati(p1, f1, p2, f2, p3, f3):
  632. """The root of r(p) = (u*p + v) / (p + w) given three points and values,
  633. (p1, f2), (p2, f2) and (p3, f3).
  634. The FITPACK analog adjusts the bounds, and we do not
  635. https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fprati.f
  636. NB: FITPACK uses p < 0 to encode p=infinity. We just use the infinity itself.
  637. Since the bracket is ``p1 <= p2 <= p3``, ``p3`` can be infinite (in fact,
  638. this is what the minimizer starts with, ``p3=inf``).
  639. """
  640. h1 = f1 * (f2 - f3)
  641. h2 = f2 * (f3 - f1)
  642. h3 = f3 * (f1 - f2)
  643. if p3 == np.inf:
  644. return -(p2*h1 + p1*h2) / h3
  645. return -(p1*p2*h3 + p2*p3*h1 + p1*p3*h2) / (p1*h1 + p2*h2 + p3*h3)
  646. class Bunch:
  647. def __init__(self, **kwargs):
  648. self.__dict__.update(**kwargs)
  649. _iermesg1 = """error. a theoretically impossible result was found during
  650. the iteration process for finding a smoothing spline with
  651. fp = s. probably causes : s too small.
  652. """
  653. _iermesg = {
  654. 1: _iermesg1 + """the weighted sum of squared residuals is becoming NaN
  655. """,
  656. 2: _iermesg1 + """there is an approximation returned but the corresponding
  657. weighted sum of squared residuals does not satisfy the
  658. condition abs(fp-s)/s < tol.
  659. """,
  660. 3: """error. the maximal number of iterations maxit (set to 20
  661. by the program) allowed for finding a smoothing spline
  662. with fp=s has been reached. probably causes : s too small
  663. there is an approximation returned but the corresponding
  664. weighted sum of squared residuals does not satisfy the
  665. condition abs(fp-s)/s < tol.
  666. """
  667. }
  668. def root_rati(f, p0, bracket, acc):
  669. """Solve `f(p) = 0` using a rational function approximation.
  670. In a nutshell, since the function f(p) is known to be monotonically decreasing, we
  671. - maintain the bracket (p1, f1), (p2, f2) and (p3, f3)
  672. - at each iteration step, approximate f(p) by a rational function
  673. r(p) = (u*p + v) / (p + w)
  674. and make a step to p_new to the root of f(p): r(p_new) = 0.
  675. The coefficients u, v and w are found from the bracket values p1..3 and f1...3
  676. The algorithm and implementation follows
  677. https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpcurf.f#L229
  678. and
  679. https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fppara.f#L290
  680. Note that the latter is for parametric splines and the former is for 1D spline
  681. functions. The minimization is indentical though [modulo a summation over the
  682. dimensions in the computation of f(p)], so we reuse the minimizer for both
  683. d=1 and d>1.
  684. """
  685. # Magic values from
  686. # https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpcurf.f#L27
  687. con1 = 0.1
  688. con9 = 0.9
  689. con4 = 0.04
  690. # bracketing flags (follow FITPACK)
  691. # https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fppara.f#L365
  692. ich1, ich3 = 0, 0
  693. (p1, f1), (p3, f3) = bracket
  694. p = p0
  695. for it in range(MAXIT):
  696. p2, f2 = p, f(p)
  697. # c test whether the approximation sp(x) is an acceptable solution.
  698. if abs(f2) < acc:
  699. ier, converged = 0, True
  700. break
  701. # c carry out one more step of the iteration process.
  702. if ich3 == 0:
  703. if f2 - f3 <= acc:
  704. # c our initial choice of p is too large.
  705. p3 = p2
  706. f3 = f2
  707. p = p*con4
  708. if p <= p1:
  709. p = p1*con9 + p2*con1
  710. continue
  711. else:
  712. if f2 < 0:
  713. ich3 = 1
  714. if ich1 == 0:
  715. if f1 - f2 <= acc:
  716. # c our initial choice of p is too small
  717. p1 = p2
  718. f1 = f2
  719. p = p/con4
  720. if p3 != np.inf and p <= p3:
  721. p = p2*con1 + p3*con9
  722. continue
  723. else:
  724. if f2 > 0:
  725. ich1 = 1
  726. # c test whether the iteration process proceeds as theoretically expected.
  727. # [f(p) should be monotonically decreasing]
  728. if f1 <= f2 or f2 <= f3:
  729. ier, converged = 2, False
  730. break
  731. # actually make the iteration step
  732. p = fprati(p1, f1, p2, f2, p3, f3)
  733. # c adjust the value of p1,f1,p3 and f3 such that f1 > 0 and f3 < 0.
  734. if f2 < 0:
  735. p3, f3 = p2, f2
  736. else:
  737. p1, f1 = p2, f2
  738. else:
  739. # not converged in MAXIT iterations
  740. ier, converged = 3, False
  741. if ier != 0:
  742. warnings.warn(RuntimeWarning(_iermesg[ier]), stacklevel=2)
  743. return Bunch(converged=converged, root=p, iterations=it, ier=ier)
  744. def _make_splrep_impl(x, y, w, xb, xe, k, s, t, nest, periodic, xp=np):
  745. """Shared infra for make_splrep and make_splprep.
  746. """
  747. acc = s * TOL
  748. m = x.size # the number of data points
  749. if nest is None:
  750. # the max number of knots. This is set in _fitpack_impl.py line 274
  751. # and fitpack.pyf line 198
  752. # Ref: https://github.com/scipy/scipy/blob/596b586e25e34bd842b575bac134b4d6924c6556/scipy/interpolate/_fitpack_impl.py#L260-L263
  753. if periodic:
  754. nest = max(m + 2*k, 2*k + 3)
  755. else:
  756. nest = max(m + k + 1, 2*k + 3)
  757. else:
  758. if nest < 2*(k + 1):
  759. raise ValueError(f"`nest` too small: {nest = } < 2*(k+1) = {2*(k+1)}.")
  760. if t is not None:
  761. raise ValueError("Either supply `t` or `nest`.")
  762. if t is None:
  763. gen = _generate_knots_impl(x, y, w, xb, xe, k, s, nest, periodic)
  764. t = list(gen)[-1]
  765. else:
  766. fpcheck(x, t, k, periodic=periodic)
  767. if t.shape[0] == 2 * (k + 1):
  768. # nothing to optimize
  769. _, _, c, _, _ = _lsq_solve_qr(x, y, t, k, w, periodic=periodic)
  770. t, c = xp.asarray(t), xp.asarray(c)
  771. return BSpline(t, c, k)
  772. ### solve ###
  773. # c initial value for p.
  774. # https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpcurf.f#L253
  775. if periodic:
  776. # N.B. - Check _lsq_solve_qr computation
  777. # of p for periodic splines
  778. R, A1, A2, Z, Y, _, p, _ = _lsq_solve_qr_for_root_rati_periodic(x, y, t, k, w)
  779. else:
  780. R, Y, _, _, _ = _lsq_solve_qr(x, y, t, k, w, periodic=periodic)
  781. nc = t.shape[0] -k -1
  782. if not periodic:
  783. p = nc / R[:, 0].sum()
  784. # ### bespoke solver ####
  785. # initial conditions
  786. # f(p=inf) : LSQ spline with knots t (XXX: reuse R, c)
  787. # N.B. - Check _lsq_solve_qr which is called
  788. # via _get_residuals for logic behind
  789. # computation of fp for periodic splines
  790. _, fp = _get_residuals(x, y, t, k, w, periodic=periodic)
  791. fpinf = fp - s
  792. # f(p=0): LSQ spline without internal knots
  793. if not periodic:
  794. _, fp0 = _get_residuals(x, y, np.array([xb]*(k+1) + [xe]*(k+1)), k, w)
  795. fp0 = fp0 - s
  796. else:
  797. # f(p=0) is fp for constant function
  798. # in case of periodic splines
  799. per = xe - xb
  800. tc = np.zeros(2*(k + 1), dtype=float)
  801. for i in range(0, k + 1):
  802. tc[i] = x[0] - (k - i) * per
  803. tc[i + k + 1] = x[m - 1] + i * per
  804. _, fp0 = _get_residuals(x, y, tc, k, w, periodic=periodic)
  805. fp0 = fp0 - s
  806. # solve
  807. bracket = (0, fp0), (np.inf, fpinf)
  808. if not periodic:
  809. f = F(x, y, t, k=k, s=s, w=w, R=R, Y=Y)
  810. else:
  811. f = Fperiodic(x, y, t, k=k, s=s, w=w, R=R, Y=Y, A1=A1, A2=A2, Z=Z)
  812. _ = root_rati(f, p, bracket, acc)
  813. # solve ALTERNATIVE: is roughly equivalent, gives slightly different results
  814. # starting from scratch, that would have probably been tolerable;
  815. # backwards compatibility dictates that we replicate the FITPACK minimizer though.
  816. # f = F(x, y, t, k=k, s=s, w=w, R=R, Y=Y)
  817. # from scipy.optimize import root_scalar
  818. # res_ = root_scalar(f, x0=p, rtol=acc)
  819. # assert res_.converged
  820. # f.spl is the spline corresponding to the found `p` value
  821. t, c, k = f.spl.tck
  822. axis, extrap = f.spl.axis, f.spl.extrapolate
  823. t, c = xp.asarray(t), xp.asarray(c)
  824. spl = BSpline.construct_fast(t, c, k, axis=axis, extrapolate=extrap)
  825. return spl
  826. @xp_capabilities(cpu_only=True, jax_jit=False, allow_dask_compute=True)
  827. def make_splrep(x, y, *, w=None, xb=None, xe=None,
  828. k=3, s=0, t=None, nest=None, bc_type=None):
  829. r"""Create a smoothing B-spline function with bounded error, minimizing derivative jumps.
  830. Given the set of data points ``(x[i], y[i])``, determine a smooth spline
  831. approximation of degree ``k`` on the interval ``xb <= x <= xe``.
  832. Parameters
  833. ----------
  834. x, y : array_like, shape (m,)
  835. The data points defining a curve ``y = f(x)``.
  836. w : array_like, shape (m,), optional
  837. Strictly positive 1D array of weights, of the same length as `x` and `y`.
  838. The weights are used in computing the weighted least-squares spline
  839. fit. If the errors in the y values have standard-deviation given by the
  840. vector ``d``, then `w` should be ``1/d``.
  841. Default is ``np.ones(m)``.
  842. xb, xe : float, optional
  843. The interval to fit. If None, these default to ``x[0]`` and ``x[-1]``,
  844. respectively.
  845. k : int, optional
  846. The degree of the spline fit. It is recommended to use cubic splines,
  847. ``k=3``, which is the default. Even values of `k` should be avoided,
  848. especially with small `s` values.
  849. s : float, optional
  850. The smoothing condition. The amount of smoothness is determined by
  851. satisfying the LSQ (least-squares) constraint::
  852. sum((w * (g(x) - y))**2 ) <= s
  853. where ``g(x)`` is the smoothed fit to ``(x, y)``. The user can use `s`
  854. to control the tradeoff between closeness to data and smoothness of fit.
  855. Larger `s` means more smoothing while smaller values of `s` indicate less
  856. smoothing.
  857. Recommended values of `s` depend on the weights, `w`. If the weights
  858. represent the inverse of the standard deviation of `y`, then a good `s`
  859. value should be found in the range ``(m-sqrt(2*m), m+sqrt(2*m))`` where
  860. ``m`` is the number of datapoints in `x`, `y`, and `w`.
  861. Default is ``s = 0.0``, i.e. interpolation.
  862. t : array_like, optional
  863. The spline knots. If None (default), the knots will be constructed
  864. automatically.
  865. There must be at least ``2*k + 2`` and at most ``m + k + 1`` knots.
  866. nest : int, optional
  867. The target length of the knot vector. Should be between ``2*(k + 1)``
  868. (the minimum number of knots for a degree-``k`` spline), and
  869. ``m + k + 1`` (the number of knots of the interpolating spline).
  870. The actual number of knots returned by this routine may be slightly
  871. larger than `nest`.
  872. Default is None (no limit, add up to ``m + k + 1`` knots).
  873. periodic : bool, optional
  874. If True, data points are considered periodic with period ``x[m-1]`` -
  875. ``x[0]`` and a smooth periodic spline approximation is returned. Values of
  876. ``y[m-1]`` and ``w[m-1]`` are not used.
  877. The default is False, corresponding to boundary condition 'not-a-knot'.
  878. bc_type : str, optional
  879. Boundary conditions.
  880. Default is `"not-a-knot"`.
  881. The following boundary conditions are recognized:
  882. * ``"not-a-knot"`` (default): The first and second segments are the
  883. same polynomial. This is equivalent to having ``bc_type=None``.
  884. * ``"periodic"``: The values and the first ``k-1`` derivatives at the
  885. ends are equivalent.
  886. Returns
  887. -------
  888. spl : a `BSpline` instance
  889. For `s=0`, ``spl(x) == y``.
  890. For non-zero values of `s` the `spl` represents the smoothed approximation
  891. to `(x, y)`, generally with fewer knots.
  892. See Also
  893. --------
  894. generate_knots : is used under the hood for generating the knots
  895. make_splprep : the analog of this routine for parametric curves
  896. make_interp_spline : construct an interpolating spline (``s = 0``)
  897. make_lsq_spline : construct the least-squares spline given the knot vector
  898. splrep : a FITPACK analog of this routine
  899. References
  900. ----------
  901. .. [1] P. Dierckx, "Algorithms for smoothing data with periodic and
  902. parametric splines, Computer Graphics and Image Processing",
  903. 20 (1982) 171-184.
  904. .. [2] P. Dierckx, "Curve and surface fitting with splines", Monographs on
  905. Numerical Analysis, Oxford University Press, 1993.
  906. Notes
  907. -----
  908. This routine constructs the smoothing spline function, :math:`g(x)`, to
  909. minimize the sum of jumps, :math:`D_j`, of the ``k``-th derivative at the
  910. internal knots (:math:`x_b < t_i < x_e`), where
  911. .. math::
  912. D_i = g^{(k)}(t_i + 0) - g^{(k)}(t_i - 0)
  913. Specifically, the routine constructs the spline function :math:`g(x)` which
  914. minimizes
  915. .. math::
  916. \sum_i | D_i |^2 \to \mathrm{min}
  917. provided that
  918. .. math::
  919. \sum_{j=1}^m (w_j \times (g(x_j) - y_j))^2 \leqslant s ,
  920. where :math:`s > 0` is the input parameter.
  921. In other words, we balance maximizing the smoothness (measured as the jumps
  922. of the derivative, the first criterion), and the deviation of :math:`g(x_j)`
  923. from the data :math:`y_j` (the second criterion).
  924. Note that the summation in the second criterion is over all data points,
  925. and in the first criterion it is over the internal spline knots (i.e.
  926. those with ``xb < t[i] < xe``). The spline knots are in general a subset
  927. of data, see `generate_knots` for details.
  928. Also note the difference of this routine to `make_lsq_spline`: the latter
  929. routine does not consider smoothness and simply solves a least-squares
  930. problem
  931. .. math::
  932. \sum w_j \times (g(x_j) - y_j)^2 \to \mathrm{min}
  933. for a spline function :math:`g(x)` with a _fixed_ knot vector ``t``.
  934. .. versionadded:: 1.15.0
  935. """ # noqa:E501
  936. xp = array_namespace(x, y, w, t)
  937. if t is not None:
  938. t = xp.asarray(t)
  939. bc_type = _validate_bc_type(bc_type)
  940. if s == 0:
  941. if t is not None or w is not None or nest is not None:
  942. raise ValueError("s==0 is for interpolation only")
  943. return make_interp_spline(x, y, k=k, bc_type=bc_type)
  944. periodic = (bc_type == 'periodic')
  945. x, y, w, k, s, xb, xe = _validate_inputs(x, y, w, k, s, xb, xe,
  946. parametric=False, periodic=periodic)
  947. spl = _make_splrep_impl(x, y, w, xb, xe, k, s, t, nest, periodic, xp=xp)
  948. # postprocess: squeeze out the last dimension: was added to simplify the internals.
  949. spl.c = spl.c[:, 0]
  950. return spl
  951. @xp_capabilities(cpu_only=True, jax_jit=False, allow_dask_compute=True)
  952. def make_splprep(x, *, w=None, u=None, ub=None, ue=None,
  953. k=3, s=0, t=None, nest=None, bc_type=None):
  954. r"""
  955. Create a smoothing parametric B-spline curve with bounded error, minimizing derivative jumps.
  956. Given a list of N 1D arrays, `x`, which represent a curve in
  957. N-dimensional space parametrized by `u`, find a smooth approximating
  958. spline curve ``g(u)``.
  959. Parameters
  960. ----------
  961. x : array_like, shape (ndim, m)
  962. Sampled data points representing the curve in ``ndim`` dimensions.
  963. The typical use is a list of 1D arrays, each of length ``m``.
  964. w : array_like, shape(m,), optional
  965. Strictly positive 1D array of weights.
  966. The weights are used in computing the weighted least-squares spline
  967. fit. If the errors in the `x` values have standard deviation given by
  968. the vector d, then `w` should be 1/d. Default is ``np.ones(m)``.
  969. u : array_like, optional
  970. An array of parameter values for the curve in the parametric form.
  971. If not given, these values are calculated automatically, according to::
  972. v[0] = 0
  973. v[i] = v[i-1] + distance(x[i], x[i-1])
  974. u[i] = v[i] / v[-1]
  975. ub, ue : float, optional
  976. The end-points of the parameters interval. Default to ``u[0]`` and ``u[-1]``.
  977. k : int, optional
  978. Degree of the spline. Cubic splines, ``k=3``, are recommended.
  979. Even values of `k` should be avoided especially with a small ``s`` value.
  980. Default is ``k=3``
  981. s : float, optional
  982. A smoothing condition. The amount of smoothness is determined by
  983. satisfying the conditions::
  984. sum((w * (g(u) - x))**2) <= s,
  985. where ``g(u)`` is the smoothed approximation to ``x``. The user can
  986. use `s` to control the trade-off between closeness and smoothness
  987. of fit. Larger ``s`` means more smoothing while smaller values of ``s``
  988. indicate less smoothing.
  989. Recommended values of ``s`` depend on the weights, ``w``. If the weights
  990. represent the inverse of the standard deviation of ``x``, then a good
  991. ``s`` value should be found in the range ``(m - sqrt(2*m), m + sqrt(2*m))``,
  992. where ``m`` is the number of data points in ``x`` and ``w``.
  993. t : array_like, optional
  994. The spline knots. If None (default), the knots will be constructed
  995. automatically.
  996. There must be at least ``2*k + 2`` and at most ``m + k + 1`` knots.
  997. nest : int, optional
  998. The target length of the knot vector. Should be between ``2*(k + 1)``
  999. (the minimum number of knots for a degree-``k`` spline), and
  1000. ``m + k + 1`` (the number of knots of the interpolating spline).
  1001. The actual number of knots returned by this routine may be slightly
  1002. larger than `nest`.
  1003. Default is None (no limit, add up to ``m + k + 1`` knots).
  1004. bc_type : str, optional
  1005. Boundary conditions.
  1006. Default is `"not-a-knot"`.
  1007. The following boundary conditions are recognized:
  1008. * ``"not-a-knot"`` (default): The first and second segments are the
  1009. same polynomial. This is equivalent to having ``bc_type=None``.
  1010. * ``"periodic"``: The values and the first ``k-1`` derivatives at the
  1011. ends are equivalent.
  1012. Returns
  1013. -------
  1014. spl : a `BSpline` instance
  1015. For `s=0`, ``spl(u) == x``.
  1016. For non-zero values of ``s``, `spl` represents the smoothed approximation
  1017. to ``x``, generally with fewer knots.
  1018. u : ndarray
  1019. The values of the parameters
  1020. See Also
  1021. --------
  1022. generate_knots : is used under the hood for generating the knots
  1023. make_splrep : the analog of this routine 1D functions
  1024. make_interp_spline : construct an interpolating spline (``s = 0``)
  1025. make_lsq_spline : construct the least-squares spline given the knot vector
  1026. splprep : a FITPACK analog of this routine
  1027. Notes
  1028. -----
  1029. Given a set of :math:`m` data points in :math:`D` dimensions, :math:`\vec{x}_j`,
  1030. with :math:`j=1, ..., m` and :math:`\vec{x}_j = (x_{j; 1}, ..., x_{j; D})`,
  1031. this routine constructs the parametric spline curve :math:`g_a(u)` with
  1032. :math:`a=1, ..., D`, to minimize the sum of jumps, :math:`D_{i; a}`, of the
  1033. ``k``-th derivative at the internal knots (:math:`u_b < t_i < u_e`), where
  1034. .. math::
  1035. D_{i; a} = g_a^{(k)}(t_i + 0) - g_a^{(k)}(t_i - 0)
  1036. Specifically, the routine constructs the spline function :math:`g(u)` which
  1037. minimizes
  1038. .. math::
  1039. \sum_i \sum_{a=1}^D | D_{i; a} |^2 \to \mathrm{min}
  1040. provided that
  1041. .. math::
  1042. \sum_{j=1}^m \sum_{a=1}^D (w_j \times (g_a(u_j) - x_{j; a}))^2 \leqslant s
  1043. where :math:`u_j` is the value of the parameter corresponding to the data point
  1044. :math:`(x_{j; 1}, ..., x_{j; D})`, and :math:`s > 0` is the input parameter.
  1045. In other words, we balance maximizing the smoothness (measured as the jumps
  1046. of the derivative, the first criterion), and the deviation of :math:`g(u_j)`
  1047. from the data :math:`x_j` (the second criterion).
  1048. Note that the summation in the second criterion is over all data points,
  1049. and in the first criterion it is over the internal spline knots (i.e.
  1050. those with ``ub < t[i] < ue``). The spline knots are in general a subset
  1051. of data, see `generate_knots` for details.
  1052. .. versionadded:: 1.15.0
  1053. References
  1054. ----------
  1055. .. [1] P. Dierckx, "Algorithms for smoothing data with periodic and
  1056. parametric splines, Computer Graphics and Image Processing",
  1057. 20 (1982) 171-184.
  1058. .. [2] P. Dierckx, "Curve and surface fitting with splines", Monographs on
  1059. Numerical Analysis, Oxford University Press, 1993.
  1060. """ # noqa:E501
  1061. bc_type = _validate_bc_type(bc_type)
  1062. periodic = (bc_type == 'periodic')
  1063. # x can be either a 2D array or a list of 1D arrays
  1064. if isinstance(x, list):
  1065. xp = array_namespace(*x, w, u, t)
  1066. else:
  1067. xp = array_namespace(x, w, u, t)
  1068. x = xp.stack(x, axis=1)
  1069. # construct the default parametrization of the curve
  1070. if u is None:
  1071. dp = (x[1:, :] - x[:-1, :])**2
  1072. u = xp.cumulative_sum( xp.sqrt(xp.sum(dp, axis=1)) )
  1073. u = concat_1d(xp, 0., u / u[-1])
  1074. if s == 0:
  1075. if t is not None or w is not None or nest is not None:
  1076. raise ValueError("s==0 is for interpolation only")
  1077. return make_interp_spline(u, x.T, k=k, axis=1), u
  1078. u, x, w, k, s, ub, ue = _validate_inputs(u, x, w, k, s, ub, ue,
  1079. parametric=True, periodic=periodic)
  1080. if t is not None:
  1081. t = xp.asarray(t)
  1082. spl = _make_splrep_impl(u, x, w, ub, ue, k, s, t,
  1083. nest, periodic=periodic, xp=xp)
  1084. # posprocess: `axis=1` so that spl(u).shape == np.shape(x)
  1085. # when `x` is a list of 1D arrays (cf original splPrep)
  1086. cc = spl.c.T
  1087. spl1 = BSpline(spl.t, cc, spl.k, axis=1)
  1088. return spl1, xp.asarray(u)