cupyx.scipy.interpolate.UnivariateSpline#

class cupyx.scipy.interpolate.UnivariateSpline(x, y, w=None, bbox=[None, None], k=3, s=None, ext=0)[source]#

1-D smoothing spline fit to a given set of data points.

Fits a spline y = spl(x) of degree k to the provided x, y data. s specifies the number of knots by specifying a smoothing condition.

Parameters:
  • x ((N,) array_like) – 1-D array of independent input data. Must be increasing; must be strictly increasing if s is 0.

  • y ((N,) array_like) – 1-D array of dependent input data, of the same length as x.

  • w ((N,) array_like, optional) – Weights for spline fitting. Must be positive. If w is None, weights are all 1. Default is None.

  • bbox ((2,) array_like, optional) – 2-sequence specifying the boundary of the approximation interval. If bbox is None, bbox=[x[0], x[-1]]. Default is None.

  • k (int, optional) – Degree of the smoothing spline. k = 3 is a cubic spline. Default is 3.

  • s (float or None, optional) –

    Positive smoothing factor used to choose the number of knots. Number of knots will be increased until the smoothing condition is satisfied:

    sum((w[i] * (y[i]-spl(x[i])))**2, axis=0) <= s
    

    However, because of numerical issues, the actual condition is:

    abs(sum((w[i] * (y[i]-spl(x[i])))**2, axis=0) - s) < 0.001 * s
    

    If s is None, s will be set as len(w) for a smoothing spline that uses all data points. If 0, spline will interpolate through all data points. This is equivalent to InterpolatedUnivariateSpline. Default is None. The user can use the s to control the tradeoff between closeness and smoothness of fit. Larger s means more smoothing while smaller values of s indicate less smoothing. Recommended values of s depend on the weights, w. If the weights represent the inverse of the standard-deviation of y, then a good s value should be found in the range (m-sqrt(2*m),m+sqrt(2*m)) where m is the number of datapoints in x, y, and w. This means s = len(w) should be a good value if 1/w[i] is an estimate of the standard deviation of y[i].

  • ext (int or str, optional) –

    Controls the extrapolation mode for elements not in the interval defined by the knot sequence.

    • if ext=0 or ‘extrapolate’, return the extrapolated value.

    • if ext=1 or ‘zeros’, return 0

    • if ext=2 or ‘raise’, raise a ValueError

    • if ext=3 or ‘const’, return the boundary value.

    Default is 0.

Methods

__call__(x, nu=0, ext=None)[source]#

Evaluate spline (or its nu-th derivative) at positions x.

Parameters:
  • x (ndarray) – A 1-D array of points at which to return the value of the smoothed spline or its derivatives. Note: x can be unordered but the evaluation is more efficient if x is (partially) ordered.

  • nu (int) – The order of derivative of the spline to compute.

  • ext (int) –

    Controls the value returned for elements of x not in the interval defined by the knot sequence.

    • if ext=0 or ‘extrapolate’, return the extrapolated value.

    • if ext=1 or ‘zeros’, return 0

    • if ext=2 or ‘raise’, raise a ValueError

    • if ext=3 or ‘const’, return the boundary value.

    The default value is 0, passed from the initialization of UnivariateSpline.

antiderivative(n=1)[source]#

Construct a new spline representing the antiderivative of this spline.

Parameters:

n (int, optional) – Order of antiderivative to evaluate. Default: 1

Returns:

spline – Spline of order k2=k+n representing the antiderivative of this spline.

Return type:

UnivariateSpline

derivative(n=1)[source]#

Construct a new spline representing the derivative of this spline.

Parameters:

n (int, optional) – Order of derivative to evaluate. Default: 1

Returns:

spline – Spline of order k2=k-n representing the derivative of this spline.

Return type:

UnivariateSpline

derivatives(x)[source]#

Return all derivatives of the spline at the point x.

Parameters:

x (float) – The point to evaluate the derivatives at.

Returns:

der – Derivatives of the orders 0 to k.

Return type:

ndarray, shape(k+1,)

get_coeffs()[source]#

Return spline coefficients.

get_knots()[source]#

Return positions of interior knots of the spline.

Internally, the knot vector contains 2*k additional boundary knots.

get_residual()[source]#

Return weighted sum of squared residuals of the spline approx.

This is equivalent to:

sum((w[i] * (y[i]-spl(x[i])))**2, axis=0
integral(a, b)[source]#

Return definite integral of the spline between two given points.

Parameters:
  • a (float) – Lower limit of integration.

  • b (float) – Upper limit of integration.

Returns:

integral – The value of the definite integral of the spline between limits.

Return type:

float

set_smoothing_factor(s, t=None)[source]#

Continue spline computation with the given smoothing factor s and with the knots found at the last call.

This routine modifies the spline in place.

__eq__(value, /)#

Return self==value.

__ne__(value, /)#

Return self!=value.

__lt__(value, /)#

Return self<value.

__le__(value, /)#

Return self<=value.

__gt__(value, /)#

Return self>value.

__ge__(value, /)#

Return self>=value.