scipy.interpolate.BPoly#

class 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:
cndarray, shape (k, m, …)

Polynomial coefficients, order k and m intervals

xndarray, shape (m+1,)

Polynomial breakpoints. Must be sorted in either increasing or decreasing order.

extrapolatebool, 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.

axisint, optional

Interpolation axis. Default is zero.

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

[3]

E. H. Doha, A. H. Bhrawy, and M. A. Saker, Boundary Value Problems, vol 2011, article ID 829546, DOI:10.1155/2011/829543.

Examples

>>> from 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}\]
Attributes:
xndarray

Breakpoints.

cndarray

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.

axisint

Interpolation axis.

Methods

__call__(x[, nu, extrapolate])

Evaluate the piecewise polynomial or its derivative.

extend(c, x)

Add additional breakpoints and coefficients to the polynomial.

derivative([nu])

Construct a new piecewise polynomial representing the derivative.

antiderivative([nu])

Construct a new piecewise polynomial representing the antiderivative.

integrate(a, b[, extrapolate])

Compute a definite integral over a piecewise polynomial.

construct_fast(c, x[, extrapolate, axis])

Construct the piecewise polynomial without making checks.

from_power_basis(pp[, extrapolate])

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

from_derivatives(xi, yi[, orders, extrapolate])

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