SciPy

scipy.interpolate.BPoly

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

Piecewise polynomial in terms of coefficients and breakpoints.

The polynomial in the i-th interval x[i] <= xp < 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**k * (1 - t)**(k - a),

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

Parameters:

c : ndarray, shape (k, m, ...)

Polynomial coefficients, order k and m intervals

x : ndarray, shape (m+1,)

Polynomial breakpoints. These must be sorted in increasing 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.

See also

PPoly
piecewise polynomials in the power basis

Notes

Properties of Bernstein polynomials are well documented in the literature. Here’s a non-exhaustive list:

[R49]http://en.wikipedia.org/wiki/Bernstein_polynomial
[R50]Kenneth I. Joy, Bernstein polynomials, http://www.idav.ucdavis.edu/education/CAGDNotes/Bernstein-Polynomials.pdf
[R51]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

x (ndarray) Breakpoints.
c (ndarray) Coefficients of the polynomials. They are reshaped to a 3-dimensional array with the last dimension representing the trailing dimensions of the original coefficient array.
axis (int) Interpolation axis.

Methods

__call__(x[, nu, extrapolate]) Evaluate the piecewise polynomial or its derivative.
extend(c, x[, right]) 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.