cupyx.scipy.interpolate.CubicHermiteSpline#

class cupyx.scipy.interpolate.CubicHermiteSpline(x, y, dydx, axis=0, extrapolate=None)[source]#

Piecewise-cubic interpolator matching values and first derivatives.

The result is represented as a PPoly instance. [1]

Parameters:
  • x (array_like, shape (n,)) – 1-D array containing values of the independent variable. Values must be real, finite and in strictly increasing order.

  • y (array_like) – Array containing values of the dependent variable. It can have arbitrary number of dimensions, but the length along axis (see below) must match the length of x. Values must be finite.

  • dydx (array_like) – Array containing derivatives of the dependent variable. It can have arbitrary number of dimensions, but the length along axis (see below) must match the length of x. Values must be finite.

  • axis (int, optional) – Axis along which y is assumed to be varying. Meaning that for x[i] the corresponding values are cupy.take(y, i, axis=axis). Default is 0.

  • extrapolate ({bool, 'periodic', None}, 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. If None (default), it is set to True.

Variables:
  • x (ndarray, shape (n,)) – Breakpoints. The same x which was passed to the constructor.

  • c (ndarray, shape (4, n-1, ...)) – Coefficients of the polynomials on each segment. The trailing dimensions match the dimensions of y, excluding axis. For example, if y is 1-D, then c[k, i] is a coefficient for (x-x[i])**(3-k) on the segment between x[i] and x[i+1].

  • axis (int) – Interpolation axis. The same axis which was passed to the constructor.

See also

Akima1DInterpolator

Akima 1D interpolator.

PchipInterpolator

PCHIP 1-D monotonic cubic interpolator.

PPoly

Piecewise polynomial in terms of coefficients and breakpoints

Notes

If you want to create a higher-order spline matching higher-order derivatives, use BPoly.from_derivatives.

References

Methods

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

Evaluate the piecewise polynomial or its derivative.

Parameters:
  • x (array_like) – Points to evaluate the interpolant at.

  • nu (int, optional) – Order of derivative to evaluate. Must be non-negative.

  • extrapolate ({bool, 'periodic', None}, 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. If None (default), use self.extrapolate.

Returns:

y – Interpolated values. Shape is determined by replacing the interpolation axis in the original array with the shape of x.

Return type:

array_like

Notes

Derivatives are evaluated piecewise for each polynomial segment, even if the polynomial is not differentiable at the breakpoints. The polynomial intervals are considered half-open, [a, b), except for the last interval which is closed [a, b].

antiderivative(nu=1)[source]#

Construct a new piecewise polynomial representing the antiderivative. Antiderivative is also the indefinite integral of the function, and derivative is its inverse operation.

Parameters:

nu (int, optional) – Order of antiderivative to evaluate. Default is 1, i.e., compute the first integral. If negative, the derivative is returned.

Returns:

pp – Piecewise polynomial of order k2 = k + n representing the antiderivative of this polynomial.

Return type:

PPoly

Notes

The antiderivative returned by this function is continuous and continuously differentiable to order n-1, up to floating point rounding error.

If antiderivative is computed and self.extrapolate='periodic', it will be set to False for the returned instance. This is done because the antiderivative is no longer periodic and its correct evaluation outside of the initially given x interval is difficult.

classmethod construct_fast(c, x, extrapolate=None, axis=0)[source]#

Construct the piecewise polynomial without making checks. Takes the same parameters as the constructor. Input arguments c and x must be arrays of the correct shape and type. The c array can only be of dtypes float and complex, and x array must have dtype float.

derivative(nu=1)[source]#

Construct a new piecewise polynomial representing the derivative.

Parameters:

nu (int, optional) – Order of derivative to evaluate. Default is 1, i.e., compute the first derivative. If negative, the antiderivative is returned.

Returns:

pp – Piecewise polynomial of order k2 = k - n representing the derivative of this polynomial.

Return type:

PPoly

Notes

Derivatives are evaluated piecewise for each polynomial segment, even if the polynomial is not differentiable at the breakpoints. The polynomial intervals are considered half-open, [a, b), except for the last interval which is closed [a, b].

extend(c, x)[source]#

Add additional breakpoints and coefficients to the polynomial.

Parameters:
  • c (ndarray, size (k, m, ...)) – Additional coefficients for polynomials in intervals. Note that the first additional interval will be formed using one of the self.x end points.

  • x (ndarray, size (m,)) – Additional breakpoints. Must be sorted in the same order as self.x and either to the right or to the left of the current breakpoints.

classmethod from_bernstein_basis(bp, extrapolate=None)[source]#

Construct a piecewise polynomial in the power basis from a polynomial in Bernstein basis.

Parameters:
  • bp (BPoly) – A Bernstein basis polynomial, as created by BPoly

  • extrapolate (bool 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.

classmethod from_spline(tck, extrapolate=None)[source]#

Construct a piecewise polynomial from a spline

Parameters:
  • tck – A spline, as a (knots, coefficients, degree) tuple or a BSpline object.

  • extrapolate (bool 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.

integrate(a, b, extrapolate=None)[source]#

Compute a definite integral over a piecewise polynomial.

Parameters:
  • a (float) – Lower integration bound

  • b (float) – Upper integration bound

  • extrapolate ({bool, 'periodic', None}, 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. If None (default), use self.extrapolate.

Returns:

ig – Definite integral of the piecewise polynomial over [a, b]

Return type:

array_like

roots(discontinuity=True, extrapolate=None)[source]#

Find real roots of the piecewise polynomial.

Parameters:
  • discontinuity (bool, optional) – Whether to report sign changes across discontinuities at breakpoints as roots.

  • extrapolate ({bool, 'periodic', None}, optional) – If bool, determines whether to return roots from the polynomial extrapolated based on first and last intervals, ‘periodic’ works the same as False. If None (default), use self.extrapolate.

Returns:

roots – Roots of the polynomial(s). If the PPoly object describes multiple polynomials, the return value is an object array whose each element is an ndarray containing the roots.

Return type:

ndarray

See also

PPoly.solve

solve(y=0.0, discontinuity=True, extrapolate=None)[source]#

Find real solutions of the equation pp(x) == y.

Parameters:
  • y (float, optional) – Right-hand side. Default is zero.

  • discontinuity (bool, optional) – Whether to report sign changes across discontinuities at breakpoints as roots.

  • extrapolate ({bool, 'periodic', None}, optional) – If bool, determines whether to return roots from the polynomial extrapolated based on first and last intervals, ‘periodic’ works the same as False. If None (default), use self.extrapolate.

Returns:

roots – Roots of the polynomial(s). If the PPoly object describes multiple polynomials, the return value is an object array whose each element is an ndarray containing the roots.

Return type:

ndarray

Notes

This routine works only on real-valued polynomials. If the piecewise polynomial contains sections that are identically zero, the root list will contain the start point of the corresponding interval, followed by a nan value. If the polynomial is discontinuous across a breakpoint, and there is a sign change across the breakpoint, this is reported if the discont parameter is True.

At the moment, there is not an actual implementation.

__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.

Attributes

c#
x#
extrapolate#
axis#