scipy.interpolate.PchipInterpolator.from_spline#
- classmethod PchipInterpolator.from_spline(tck, extrapolate=None)[source]#
- Construct a piecewise polynomial from a spline - Parameters:
- tck
- A spline, as returned by - splrepor a BSpline object.
- extrapolatebool or ‘periodic’, optional
- If bool, determines whether to extrapolate to out-of-bounds points based on first and last intervals, or to return NaNs. If ‘periodic’, periodic extrapolation is used. Default is True. 
 
 - Examples - Construct an interpolating spline and convert it to a - PPolyinstance- >>> import numpy as np >>> from scipy.interpolate import splrep, PPoly >>> x = np.linspace(0, 1, 11) >>> y = np.sin(2*np.pi*x) >>> tck = splrep(x, y, s=0) >>> p = PPoly.from_spline(tck) >>> isinstance(p, PPoly) True - Note that this function only supports 1D splines out of the box. - If the - tckobject represents a parametric spline (e.g. constructed by- splprepor a- BSplinewith- c.ndim > 1), you will need to loop over the dimensions manually.- >>> from scipy.interpolate import splprep, splev >>> t = np.linspace(0, 1, 11) >>> x = np.sin(2*np.pi*t) >>> y = np.cos(2*np.pi*t) >>> (t, c, k), u = splprep([x, y], s=0) - Note that - cis a list of two arrays of length 11.- >>> unew = np.arange(0, 1.01, 0.01) >>> out = splev(unew, (t, c, k)) - To convert this spline to the power basis, we convert each component of the list of b-spline coefficients, - c, into the corresponding cubic polynomial.- >>> polys = [PPoly.from_spline((t, cj, k)) for cj in c] >>> polys[0].c.shape (4, 14) - Note that the coefficients of the polynomials polys are in the power basis and their dimensions reflect just that: here 4 is the order (degree+1), and 14 is the number of intervals—which is nothing but the length of the knot array of the original tck minus one. - Optionally, we can stack the components into a single - PPolyalong the third dimension:- >>> cc = np.dstack([p.c for p in polys]) # has shape = (4, 14, 2) >>> poly = PPoly(cc, polys[0].x) >>> np.allclose(poly(unew).T, # note the transpose to match `splev` ... out, atol=1e-15) True