cupyx.scipy.interpolate.BPoly#

class cupyx.scipy.interpolate.BPoly(c, x, extrapolate=None, axis=0)[source]#

Piecewise polynomial in terms of coefficients and breakpoints.

The polynomial between x[i] and x[i + 1] is written in the

Bernstein polynomial basis:

S = sum(c[a, i] * b(a, k; x) for a in range(k+1)),

where k is the degree of the polynomial, and:

b(a, k; x) = binom(k, a) * t**a * (1 - t)**(k - a),

with t = (x - x[i]) / (x[i+1] - x[i]) and binom is the binomial coefficient.

Parameters:
  • c (ndarray, shape (k, m, ...)) – Polynomial coefficients, order k and m intervals

  • x (ndarray, shape (m+1,)) – Polynomial breakpoints. Must be sorted in either increasing or decreasing order.

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

  • axis (int, optional) – Interpolation axis. Default is zero.

Variables:
  • x (ndarray) – Breakpoints.

  • c (ndarray) – Coefficients of the polynomials. They are reshaped to a 3-D array with the last dimension representing the trailing dimensions of the original coefficient array.

  • axis (int) – Interpolation axis.

See also

PPoly

piecewise polynomials in the power basis

Notes

Properties of Bernstein polynomials are well documented in the literature, see for example [1] [2] [3].

References

Examples

>>> from cupyx.scipy.interpolate import BPoly
>>> x = [0, 1]
>>> c = [[1], [2], [3]]
>>> bp = BPoly(c, x)

This creates a 2nd order polynomial

\[\begin{split}B(x) = 1 \times b_{0, 2}(x) + 2 \times b_{1, 2}(x) + 3 \times b_{2, 2}(x) \\ = 1 \times (1-x)^2 + 2 \times 2 x (1 - x) + 3 \times x^2\end{split}\]

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.

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:

bp – Piecewise polynomial of order k + nu representing the antiderivative of this polynomial.

Return type:

BPoly

Notes

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:

bp – Piecewise polynomial of order k - nu representing the derivative of this polynomial.

Return type:

BPoly

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_derivatives(xi, yi, orders=None, extrapolate=None)[source]#

Construct a piecewise polynomial in the Bernstein basis, compatible with the specified values and derivatives at breakpoints.

Parameters:
  • xi (array_like) – sorted 1-D array of x-coordinates

  • yi (array_like or list of array_likes) – yi[i][j] is the j th derivative known at xi[i]

  • orders (None or int or array_like of ints. Default: None.) – Specifies the degree of local polynomials. If not None, some derivatives are ignored.

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

Notes

If k derivatives are specified at a breakpoint x, the constructed polynomial is exactly k times continuously differentiable at x, unless the order is provided explicitly. In the latter case, the smoothness of the polynomial at the breakpoint is controlled by the order.

Deduces the number of derivatives to match at each end from order and the number of derivatives available. If possible it uses the same number of derivatives from each end; if the number is odd it tries to take the extra one from y2. In any case if not enough derivatives are available at one end or another it draws enough to make up the total from the other end.

If the order is too high and not enough derivatives are available, an exception is raised.

Examples

>>> from cupyx.scipy.interpolate import BPoly
>>> BPoly.from_derivatives([0, 1], [[1, 2], [3, 4]])

Creates a polynomial f(x) of degree 3, defined on [0, 1] such that f(0) = 1, df/dx(0) = 2, f(1) = 3, df/dx(1) = 4

>>> BPoly.from_derivatives([0, 1, 2], [[0, 1], [0], [2]])

Creates a piecewise polynomial f(x), such that f(0) = f(1) = 0, f(2) = 2, and df/dx(0) = 1. Based on the number of derivatives provided, the order of the local polynomials is 2 on [0, 1] and 1 on [1, 2]. Notice that no restriction is imposed on the derivatives at x = 1 and x = 2.

Indeed, the explicit form of the polynomial is:

f(x) = | x * (1 - x),  0 <= x < 1
       | 2 * (x - 1),  1 <= x <= 2

So that f’(1-0) = -1 and f’(1+0) = 2

classmethod from_power_basis(pp, extrapolate=None)[source]#

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

Parameters:
  • pp (PPoly) – A piecewise polynomial in the power basis

  • 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) – 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:

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

Return type:

array_like

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