make_splprep#
- scipy.interpolate.make_splprep(x, *, w=None, u=None, ub=None, ue=None, k=3, s=0, t=None, nest=None)[source]#
Find a smoothed B-spline representation of a parametric N-D curve.
Given a list of N 1D arrays, x, which represent a curve in N-dimensional space parametrized by u, find a smooth approximating spline curve
g(u)
.- Parameters:
- xarray_like, shape (m, ndim)
Sampled data points representing the curve in
ndim
dimensions. The typical use is a list of 1D arrays, each of lengthm
.- warray_like, shape(m,), optional
Strictly positive 1D array of weights. The weights are used in computing the weighted least-squares spline fit. If the errors in the x values have standard deviation given by the vector d, then w should be 1/d. Default is
np.ones(m)
.- uarray_like, optional
An array of parameter values for the curve in the parametric form. If not given, these values are calculated automatically, according to:
v[0] = 0 v[i] = v[i-1] + distance(x[i], x[i-1]) u[i] = v[i] / v[-1]
- ub, uefloat, optional
The end-points of the parameters interval. Default to
u[0]
andu[-1]
.- kint, optional
Degree of the spline. Cubic splines,
k=3
, are recommended. Even values of k should be avoided especially with a smalls
value. Default isk=3
- sfloat, optional
A smoothing condition. The amount of smoothness is determined by satisfying the conditions:
sum((w * (g(u) - x))**2) <= s,
where
g(u)
is the smoothed approximation tox
. The user can use s to control the trade-off between closeness and smoothness of fit. Largers
means more smoothing while smaller values ofs
indicate less smoothing. Recommended values ofs
depend on the weights,w
. If the weights represent the inverse of the standard deviation ofx
, then a goods
value should be found in the range(m - sqrt(2*m), m + sqrt(2*m))
, wherem
is the number of data points inx
andw
.- tarray_like, optional
The spline knots. If None (default), the knots will be constructed automatically. There must be at least
2*k + 2
and at mostm + k + 1
knots.- nestint, optional
The target length of the knot vector. Should be between
2*(k + 1)
(the minimum number of knots for a degree-k
spline), andm + k + 1
(the number of knots of the interpolating spline). The actual number of knots returned by this routine may be slightly larger than nest. Default is None (no limit, add up tom + k + 1
knots).
- Returns:
- spla
BSpline
instance For s=0,
spl(u) == x
. For non-zero values ofs
, spl represents the smoothed approximation tox
, generally with fewer knots.- undarray
The values of the parameters
- spla
See also
generate_knots
is used under the hood for generating the knots
make_splrep
the analog of this routine 1D functions
make_interp_spline
construct an interpolating spline (
s = 0
)make_lsq_spline
construct the least-squares spline given the knot vector
splprep
a FITPACK analog of this routine
Notes
Given a set of \(m\) data points in \(D\) dimensions, \(\vec{x}_j\), with \(j=1, ..., m\) and \(\vec{x}_j = (x_{j; 1}, ..., x_{j; D})\), this routine constructs the parametric spline curve \(g_a(u)\) with \(a=1, ..., D\), to minimize the sum of jumps, \(D_{i; a}\), of the
k
-th derivative at the internal knots (\(u_b < t_i < u_e\)), where\[D_{i; a} = g_a^{(k)}(t_i + 0) - g_a^{(k)}(t_i - 0)\]Specifically, the routine constructs the spline function \(g(u)\) which minimizes
\[\sum_i \sum_{a=1}^D | D_{i; a} |^2 \to \mathrm{min}\]provided that
\[\sum_{j=1}^m \sum_{a=1}^D (w_j \times (g_a(u_j) - x_{j; a}))^2 \leqslant s\]where \(u_j\) is the value of the parameter corresponding to the data point \((x_{j; 1}, ..., x_{j; D})\), and \(s > 0\) is the input parameter.
In other words, we balance maximizing the smoothness (measured as the jumps of the derivative, the first criterion), and the deviation of \(g(u_j)\) from the data \(x_j\) (the second criterion).
Note that the summation in the second criterion is over all data points, and in the first criterion it is over the internal spline knots (i.e. those with
ub < t[i] < ue
). The spline knots are in general a subset of data, seegenerate_knots
for details.Added in version 1.15.0.
References
[1]P. Dierckx, “Algorithms for smoothing data with periodic and parametric splines, Computer Graphics and Image Processing”, 20 (1982) 171-184.
[2]P. Dierckx, “Curve and surface fitting with splines”, Monographs on Numerical Analysis, Oxford University Press, 1993.