_bsplines.py 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683
  1. from numpy import (logical_and, asarray, pi, zeros_like,
  2. piecewise, array, arctan2, tan, zeros, arange, floor)
  3. from numpy.core.umath import (sqrt, exp, greater, less, cos, add, sin,
  4. less_equal, greater_equal)
  5. # From splinemodule.c
  6. from ._spline import cspline2d, sepfir2d
  7. from scipy.special import comb
  8. from scipy._lib._util import float_factorial
  9. __all__ = ['spline_filter', 'bspline', 'gauss_spline', 'cubic', 'quadratic',
  10. 'cspline1d', 'qspline1d', 'cspline1d_eval', 'qspline1d_eval']
  11. def spline_filter(Iin, lmbda=5.0):
  12. """Smoothing spline (cubic) filtering of a rank-2 array.
  13. Filter an input data set, `Iin`, using a (cubic) smoothing spline of
  14. fall-off `lmbda`.
  15. Parameters
  16. ----------
  17. Iin : array_like
  18. input data set
  19. lmbda : float, optional
  20. spline smooghing fall-off value, default is `5.0`.
  21. Returns
  22. -------
  23. res : ndarray
  24. filterd input data
  25. Examples
  26. --------
  27. We can filter an multi dimentional signal (ex: 2D image) using cubic
  28. B-spline filter:
  29. >>> import numpy as np
  30. >>> from scipy.signal import spline_filter
  31. >>> import matplotlib.pyplot as plt
  32. >>> orig_img = np.eye(20) # create an image
  33. >>> orig_img[10, :] = 1.0
  34. >>> sp_filter = spline_filter(orig_img, lmbda=0.1)
  35. >>> f, ax = plt.subplots(1, 2, sharex=True)
  36. >>> for ind, data in enumerate([[orig_img, "original image"],
  37. ... [sp_filter, "spline filter"]]):
  38. ... ax[ind].imshow(data[0], cmap='gray_r')
  39. ... ax[ind].set_title(data[1])
  40. >>> plt.tight_layout()
  41. >>> plt.show()
  42. """
  43. intype = Iin.dtype.char
  44. hcol = array([1.0, 4.0, 1.0], 'f') / 6.0
  45. if intype in ['F', 'D']:
  46. Iin = Iin.astype('F')
  47. ckr = cspline2d(Iin.real, lmbda)
  48. cki = cspline2d(Iin.imag, lmbda)
  49. outr = sepfir2d(ckr, hcol, hcol)
  50. outi = sepfir2d(cki, hcol, hcol)
  51. out = (outr + 1j * outi).astype(intype)
  52. elif intype in ['f', 'd']:
  53. ckr = cspline2d(Iin, lmbda)
  54. out = sepfir2d(ckr, hcol, hcol)
  55. out = out.astype(intype)
  56. else:
  57. raise TypeError("Invalid data type for Iin")
  58. return out
  59. _splinefunc_cache = {}
  60. def _bspline_piecefunctions(order):
  61. """Returns the function defined over the left-side pieces for a bspline of
  62. a given order.
  63. The 0th piece is the first one less than 0. The last piece is a function
  64. identical to 0 (returned as the constant 0). (There are order//2 + 2 total
  65. pieces).
  66. Also returns the condition functions that when evaluated return boolean
  67. arrays for use with `numpy.piecewise`.
  68. """
  69. try:
  70. return _splinefunc_cache[order]
  71. except KeyError:
  72. pass
  73. def condfuncgen(num, val1, val2):
  74. if num == 0:
  75. return lambda x: logical_and(less_equal(x, val1),
  76. greater_equal(x, val2))
  77. elif num == 2:
  78. return lambda x: less_equal(x, val2)
  79. else:
  80. return lambda x: logical_and(less(x, val1),
  81. greater_equal(x, val2))
  82. last = order // 2 + 2
  83. if order % 2:
  84. startbound = -1.0
  85. else:
  86. startbound = -0.5
  87. condfuncs = [condfuncgen(0, 0, startbound)]
  88. bound = startbound
  89. for num in range(1, last - 1):
  90. condfuncs.append(condfuncgen(1, bound, bound - 1))
  91. bound = bound - 1
  92. condfuncs.append(condfuncgen(2, 0, -(order + 1) / 2.0))
  93. # final value of bound is used in piecefuncgen below
  94. # the functions to evaluate are taken from the left-hand side
  95. # in the general expression derived from the central difference
  96. # operator (because they involve fewer terms).
  97. fval = float_factorial(order)
  98. def piecefuncgen(num):
  99. Mk = order // 2 - num
  100. if (Mk < 0):
  101. return 0 # final function is 0
  102. coeffs = [(1 - 2 * (k % 2)) * float(comb(order + 1, k, exact=1)) / fval
  103. for k in range(Mk + 1)]
  104. shifts = [-bound - k for k in range(Mk + 1)]
  105. def thefunc(x):
  106. res = 0.0
  107. for k in range(Mk + 1):
  108. res += coeffs[k] * (x + shifts[k]) ** order
  109. return res
  110. return thefunc
  111. funclist = [piecefuncgen(k) for k in range(last)]
  112. _splinefunc_cache[order] = (funclist, condfuncs)
  113. return funclist, condfuncs
  114. def bspline(x, n):
  115. """B-spline basis function of order n.
  116. Parameters
  117. ----------
  118. x : array_like
  119. a knot vector
  120. n : int
  121. The order of the spline. Must be non-negative, i.e., n >= 0
  122. Returns
  123. -------
  124. res : ndarray
  125. B-spline basis function values
  126. See Also
  127. --------
  128. cubic : A cubic B-spline.
  129. quadratic : A quadratic B-spline.
  130. Notes
  131. -----
  132. Uses numpy.piecewise and automatic function-generator.
  133. Examples
  134. --------
  135. We can calculate B-Spline basis function of several orders:
  136. >>> import numpy as np
  137. >>> from scipy.signal import bspline, cubic, quadratic
  138. >>> bspline(0.0, 1)
  139. 1
  140. >>> knots = [-1.0, 0.0, -1.0]
  141. >>> bspline(knots, 2)
  142. array([0.125, 0.75, 0.125])
  143. >>> np.array_equal(bspline(knots, 2), quadratic(knots))
  144. True
  145. >>> np.array_equal(bspline(knots, 3), cubic(knots))
  146. True
  147. """
  148. ax = -abs(asarray(x))
  149. # number of pieces on the left-side is (n+1)/2
  150. funclist, condfuncs = _bspline_piecefunctions(n)
  151. condlist = [func(ax) for func in condfuncs]
  152. return piecewise(ax, condlist, funclist)
  153. def gauss_spline(x, n):
  154. r"""Gaussian approximation to B-spline basis function of order n.
  155. Parameters
  156. ----------
  157. x : array_like
  158. a knot vector
  159. n : int
  160. The order of the spline. Must be non-negative, i.e., n >= 0
  161. Returns
  162. -------
  163. res : ndarray
  164. B-spline basis function values approximated by a zero-mean Gaussian
  165. function.
  166. Notes
  167. -----
  168. The B-spline basis function can be approximated well by a zero-mean
  169. Gaussian function with standard-deviation equal to :math:`\sigma=(n+1)/12`
  170. for large `n` :
  171. .. math:: \frac{1}{\sqrt {2\pi\sigma^2}}exp(-\frac{x^2}{2\sigma})
  172. References
  173. ----------
  174. .. [1] Bouma H., Vilanova A., Bescos J.O., ter Haar Romeny B.M., Gerritsen
  175. F.A. (2007) Fast and Accurate Gaussian Derivatives Based on B-Splines. In:
  176. Sgallari F., Murli A., Paragios N. (eds) Scale Space and Variational
  177. Methods in Computer Vision. SSVM 2007. Lecture Notes in Computer
  178. Science, vol 4485. Springer, Berlin, Heidelberg
  179. .. [2] http://folk.uio.no/inf3330/scripting/doc/python/SciPy/tutorial/old/node24.html
  180. Examples
  181. --------
  182. We can calculate B-Spline basis functions approximated by a gaussian
  183. distribution:
  184. >>> import numpy as np
  185. >>> from scipy.signal import gauss_spline, bspline
  186. >>> knots = np.array([-1.0, 0.0, -1.0])
  187. >>> gauss_spline(knots, 3)
  188. array([0.15418033, 0.6909883, 0.15418033]) # may vary
  189. >>> bspline(knots, 3)
  190. array([0.16666667, 0.66666667, 0.16666667]) # may vary
  191. """
  192. x = asarray(x)
  193. signsq = (n + 1) / 12.0
  194. return 1 / sqrt(2 * pi * signsq) * exp(-x ** 2 / 2 / signsq)
  195. def cubic(x):
  196. """A cubic B-spline.
  197. This is a special case of `bspline`, and equivalent to ``bspline(x, 3)``.
  198. Parameters
  199. ----------
  200. x : array_like
  201. a knot vector
  202. Returns
  203. -------
  204. res : ndarray
  205. Cubic B-spline basis function values
  206. See Also
  207. --------
  208. bspline : B-spline basis function of order n
  209. quadratic : A quadratic B-spline.
  210. Examples
  211. --------
  212. We can calculate B-Spline basis function of several orders:
  213. >>> import numpy as np
  214. >>> from scipy.signal import bspline, cubic, quadratic
  215. >>> bspline(0.0, 1)
  216. 1
  217. >>> knots = [-1.0, 0.0, -1.0]
  218. >>> bspline(knots, 2)
  219. array([0.125, 0.75, 0.125])
  220. >>> np.array_equal(bspline(knots, 2), quadratic(knots))
  221. True
  222. >>> np.array_equal(bspline(knots, 3), cubic(knots))
  223. True
  224. """
  225. ax = abs(asarray(x))
  226. res = zeros_like(ax)
  227. cond1 = less(ax, 1)
  228. if cond1.any():
  229. ax1 = ax[cond1]
  230. res[cond1] = 2.0 / 3 - 1.0 / 2 * ax1 ** 2 * (2 - ax1)
  231. cond2 = ~cond1 & less(ax, 2)
  232. if cond2.any():
  233. ax2 = ax[cond2]
  234. res[cond2] = 1.0 / 6 * (2 - ax2) ** 3
  235. return res
  236. def quadratic(x):
  237. """A quadratic B-spline.
  238. This is a special case of `bspline`, and equivalent to ``bspline(x, 2)``.
  239. Parameters
  240. ----------
  241. x : array_like
  242. a knot vector
  243. Returns
  244. -------
  245. res : ndarray
  246. Quadratic B-spline basis function values
  247. See Also
  248. --------
  249. bspline : B-spline basis function of order n
  250. cubic : A cubic B-spline.
  251. Examples
  252. --------
  253. We can calculate B-Spline basis function of several orders:
  254. >>> import numpy as np
  255. >>> from scipy.signal import bspline, cubic, quadratic
  256. >>> bspline(0.0, 1)
  257. 1
  258. >>> knots = [-1.0, 0.0, -1.0]
  259. >>> bspline(knots, 2)
  260. array([0.125, 0.75, 0.125])
  261. >>> np.array_equal(bspline(knots, 2), quadratic(knots))
  262. True
  263. >>> np.array_equal(bspline(knots, 3), cubic(knots))
  264. True
  265. """
  266. ax = abs(asarray(x))
  267. res = zeros_like(ax)
  268. cond1 = less(ax, 0.5)
  269. if cond1.any():
  270. ax1 = ax[cond1]
  271. res[cond1] = 0.75 - ax1 ** 2
  272. cond2 = ~cond1 & less(ax, 1.5)
  273. if cond2.any():
  274. ax2 = ax[cond2]
  275. res[cond2] = (ax2 - 1.5) ** 2 / 2.0
  276. return res
  277. def _coeff_smooth(lam):
  278. xi = 1 - 96 * lam + 24 * lam * sqrt(3 + 144 * lam)
  279. omeg = arctan2(sqrt(144 * lam - 1), sqrt(xi))
  280. rho = (24 * lam - 1 - sqrt(xi)) / (24 * lam)
  281. rho = rho * sqrt((48 * lam + 24 * lam * sqrt(3 + 144 * lam)) / xi)
  282. return rho, omeg
  283. def _hc(k, cs, rho, omega):
  284. return (cs / sin(omega) * (rho ** k) * sin(omega * (k + 1)) *
  285. greater(k, -1))
  286. def _hs(k, cs, rho, omega):
  287. c0 = (cs * cs * (1 + rho * rho) / (1 - rho * rho) /
  288. (1 - 2 * rho * rho * cos(2 * omega) + rho ** 4))
  289. gamma = (1 - rho * rho) / (1 + rho * rho) / tan(omega)
  290. ak = abs(k)
  291. return c0 * rho ** ak * (cos(omega * ak) + gamma * sin(omega * ak))
  292. def _cubic_smooth_coeff(signal, lamb):
  293. rho, omega = _coeff_smooth(lamb)
  294. cs = 1 - 2 * rho * cos(omega) + rho * rho
  295. K = len(signal)
  296. yp = zeros((K,), signal.dtype.char)
  297. k = arange(K)
  298. yp[0] = (_hc(0, cs, rho, omega) * signal[0] +
  299. add.reduce(_hc(k + 1, cs, rho, omega) * signal))
  300. yp[1] = (_hc(0, cs, rho, omega) * signal[0] +
  301. _hc(1, cs, rho, omega) * signal[1] +
  302. add.reduce(_hc(k + 2, cs, rho, omega) * signal))
  303. for n in range(2, K):
  304. yp[n] = (cs * signal[n] + 2 * rho * cos(omega) * yp[n - 1] -
  305. rho * rho * yp[n - 2])
  306. y = zeros((K,), signal.dtype.char)
  307. y[K - 1] = add.reduce((_hs(k, cs, rho, omega) +
  308. _hs(k + 1, cs, rho, omega)) * signal[::-1])
  309. y[K - 2] = add.reduce((_hs(k - 1, cs, rho, omega) +
  310. _hs(k + 2, cs, rho, omega)) * signal[::-1])
  311. for n in range(K - 3, -1, -1):
  312. y[n] = (cs * yp[n] + 2 * rho * cos(omega) * y[n + 1] -
  313. rho * rho * y[n + 2])
  314. return y
  315. def _cubic_coeff(signal):
  316. zi = -2 + sqrt(3)
  317. K = len(signal)
  318. yplus = zeros((K,), signal.dtype.char)
  319. powers = zi ** arange(K)
  320. yplus[0] = signal[0] + zi * add.reduce(powers * signal)
  321. for k in range(1, K):
  322. yplus[k] = signal[k] + zi * yplus[k - 1]
  323. output = zeros((K,), signal.dtype)
  324. output[K - 1] = zi / (zi - 1) * yplus[K - 1]
  325. for k in range(K - 2, -1, -1):
  326. output[k] = zi * (output[k + 1] - yplus[k])
  327. return output * 6.0
  328. def _quadratic_coeff(signal):
  329. zi = -3 + 2 * sqrt(2.0)
  330. K = len(signal)
  331. yplus = zeros((K,), signal.dtype.char)
  332. powers = zi ** arange(K)
  333. yplus[0] = signal[0] + zi * add.reduce(powers * signal)
  334. for k in range(1, K):
  335. yplus[k] = signal[k] + zi * yplus[k - 1]
  336. output = zeros((K,), signal.dtype.char)
  337. output[K - 1] = zi / (zi - 1) * yplus[K - 1]
  338. for k in range(K - 2, -1, -1):
  339. output[k] = zi * (output[k + 1] - yplus[k])
  340. return output * 8.0
  341. def cspline1d(signal, lamb=0.0):
  342. """
  343. Compute cubic spline coefficients for rank-1 array.
  344. Find the cubic spline coefficients for a 1-D signal assuming
  345. mirror-symmetric boundary conditions. To obtain the signal back from the
  346. spline representation mirror-symmetric-convolve these coefficients with a
  347. length 3 FIR window [1.0, 4.0, 1.0]/ 6.0 .
  348. Parameters
  349. ----------
  350. signal : ndarray
  351. A rank-1 array representing samples of a signal.
  352. lamb : float, optional
  353. Smoothing coefficient, default is 0.0.
  354. Returns
  355. -------
  356. c : ndarray
  357. Cubic spline coefficients.
  358. See Also
  359. --------
  360. cspline1d_eval : Evaluate a cubic spline at the new set of points.
  361. Examples
  362. --------
  363. We can filter a signal to reduce and smooth out high-frequency noise with
  364. a cubic spline:
  365. >>> import numpy as np
  366. >>> import matplotlib.pyplot as plt
  367. >>> from scipy.signal import cspline1d, cspline1d_eval
  368. >>> rng = np.random.default_rng()
  369. >>> sig = np.repeat([0., 1., 0.], 100)
  370. >>> sig += rng.standard_normal(len(sig))*0.05 # add noise
  371. >>> time = np.linspace(0, len(sig))
  372. >>> filtered = cspline1d_eval(cspline1d(sig), time)
  373. >>> plt.plot(sig, label="signal")
  374. >>> plt.plot(time, filtered, label="filtered")
  375. >>> plt.legend()
  376. >>> plt.show()
  377. """
  378. if lamb != 0.0:
  379. return _cubic_smooth_coeff(signal, lamb)
  380. else:
  381. return _cubic_coeff(signal)
  382. def qspline1d(signal, lamb=0.0):
  383. """Compute quadratic spline coefficients for rank-1 array.
  384. Parameters
  385. ----------
  386. signal : ndarray
  387. A rank-1 array representing samples of a signal.
  388. lamb : float, optional
  389. Smoothing coefficient (must be zero for now).
  390. Returns
  391. -------
  392. c : ndarray
  393. Quadratic spline coefficients.
  394. See Also
  395. --------
  396. qspline1d_eval : Evaluate a quadratic spline at the new set of points.
  397. Notes
  398. -----
  399. Find the quadratic spline coefficients for a 1-D signal assuming
  400. mirror-symmetric boundary conditions. To obtain the signal back from the
  401. spline representation mirror-symmetric-convolve these coefficients with a
  402. length 3 FIR window [1.0, 6.0, 1.0]/ 8.0 .
  403. Examples
  404. --------
  405. We can filter a signal to reduce and smooth out high-frequency noise with
  406. a quadratic spline:
  407. >>> import numpy as np
  408. >>> import matplotlib.pyplot as plt
  409. >>> from scipy.signal import qspline1d, qspline1d_eval
  410. >>> rng = np.random.default_rng()
  411. >>> sig = np.repeat([0., 1., 0.], 100)
  412. >>> sig += rng.standard_normal(len(sig))*0.05 # add noise
  413. >>> time = np.linspace(0, len(sig))
  414. >>> filtered = qspline1d_eval(qspline1d(sig), time)
  415. >>> plt.plot(sig, label="signal")
  416. >>> plt.plot(time, filtered, label="filtered")
  417. >>> plt.legend()
  418. >>> plt.show()
  419. """
  420. if lamb != 0.0:
  421. raise ValueError("Smoothing quadratic splines not supported yet.")
  422. else:
  423. return _quadratic_coeff(signal)
  424. def cspline1d_eval(cj, newx, dx=1.0, x0=0):
  425. """Evaluate a cubic spline at the new set of points.
  426. `dx` is the old sample-spacing while `x0` was the old origin. In
  427. other-words the old-sample points (knot-points) for which the `cj`
  428. represent spline coefficients were at equally-spaced points of:
  429. oldx = x0 + j*dx j=0...N-1, with N=len(cj)
  430. Edges are handled using mirror-symmetric boundary conditions.
  431. Parameters
  432. ----------
  433. cj : ndarray
  434. cublic spline coefficients
  435. newx : ndarray
  436. New set of points.
  437. dx : float, optional
  438. Old sample-spacing, the default value is 1.0.
  439. x0 : int, optional
  440. Old origin, the default value is 0.
  441. Returns
  442. -------
  443. res : ndarray
  444. Evaluated a cubic spline points.
  445. See Also
  446. --------
  447. cspline1d : Compute cubic spline coefficients for rank-1 array.
  448. Examples
  449. --------
  450. We can filter a signal to reduce and smooth out high-frequency noise with
  451. a cubic spline:
  452. >>> import numpy as np
  453. >>> import matplotlib.pyplot as plt
  454. >>> from scipy.signal import cspline1d, cspline1d_eval
  455. >>> rng = np.random.default_rng()
  456. >>> sig = np.repeat([0., 1., 0.], 100)
  457. >>> sig += rng.standard_normal(len(sig))*0.05 # add noise
  458. >>> time = np.linspace(0, len(sig))
  459. >>> filtered = cspline1d_eval(cspline1d(sig), time)
  460. >>> plt.plot(sig, label="signal")
  461. >>> plt.plot(time, filtered, label="filtered")
  462. >>> plt.legend()
  463. >>> plt.show()
  464. """
  465. newx = (asarray(newx) - x0) / float(dx)
  466. res = zeros_like(newx, dtype=cj.dtype)
  467. if res.size == 0:
  468. return res
  469. N = len(cj)
  470. cond1 = newx < 0
  471. cond2 = newx > (N - 1)
  472. cond3 = ~(cond1 | cond2)
  473. # handle general mirror-symmetry
  474. res[cond1] = cspline1d_eval(cj, -newx[cond1])
  475. res[cond2] = cspline1d_eval(cj, 2 * (N - 1) - newx[cond2])
  476. newx = newx[cond3]
  477. if newx.size == 0:
  478. return res
  479. result = zeros_like(newx, dtype=cj.dtype)
  480. jlower = floor(newx - 2).astype(int) + 1
  481. for i in range(4):
  482. thisj = jlower + i
  483. indj = thisj.clip(0, N - 1) # handle edge cases
  484. result += cj[indj] * cubic(newx - thisj)
  485. res[cond3] = result
  486. return res
  487. def qspline1d_eval(cj, newx, dx=1.0, x0=0):
  488. """Evaluate a quadratic spline at the new set of points.
  489. Parameters
  490. ----------
  491. cj : ndarray
  492. Quadratic spline coefficients
  493. newx : ndarray
  494. New set of points.
  495. dx : float, optional
  496. Old sample-spacing, the default value is 1.0.
  497. x0 : int, optional
  498. Old origin, the default value is 0.
  499. Returns
  500. -------
  501. res : ndarray
  502. Evaluated a quadratic spline points.
  503. See Also
  504. --------
  505. qspline1d : Compute quadratic spline coefficients for rank-1 array.
  506. Notes
  507. -----
  508. `dx` is the old sample-spacing while `x0` was the old origin. In
  509. other-words the old-sample points (knot-points) for which the `cj`
  510. represent spline coefficients were at equally-spaced points of::
  511. oldx = x0 + j*dx j=0...N-1, with N=len(cj)
  512. Edges are handled using mirror-symmetric boundary conditions.
  513. Examples
  514. --------
  515. We can filter a signal to reduce and smooth out high-frequency noise with
  516. a quadratic spline:
  517. >>> import numpy as np
  518. >>> import matplotlib.pyplot as plt
  519. >>> from scipy.signal import qspline1d, qspline1d_eval
  520. >>> rng = np.random.default_rng()
  521. >>> sig = np.repeat([0., 1., 0.], 100)
  522. >>> sig += rng.standard_normal(len(sig))*0.05 # add noise
  523. >>> time = np.linspace(0, len(sig))
  524. >>> filtered = qspline1d_eval(qspline1d(sig), time)
  525. >>> plt.plot(sig, label="signal")
  526. >>> plt.plot(time, filtered, label="filtered")
  527. >>> plt.legend()
  528. >>> plt.show()
  529. """
  530. newx = (asarray(newx) - x0) / dx
  531. res = zeros_like(newx)
  532. if res.size == 0:
  533. return res
  534. N = len(cj)
  535. cond1 = newx < 0
  536. cond2 = newx > (N - 1)
  537. cond3 = ~(cond1 | cond2)
  538. # handle general mirror-symmetry
  539. res[cond1] = qspline1d_eval(cj, -newx[cond1])
  540. res[cond2] = qspline1d_eval(cj, 2 * (N - 1) - newx[cond2])
  541. newx = newx[cond3]
  542. if newx.size == 0:
  543. return res
  544. result = zeros_like(newx)
  545. jlower = floor(newx - 1.5).astype(int) + 1
  546. for i in range(3):
  547. thisj = jlower + i
  548. indj = thisj.clip(0, N - 1) # handle edge cases
  549. result += cj[indj] * quadratic(newx - thisj)
  550. res[cond3] = result
  551. return res