This module provides a number of objects (mostly functions) useful for
dealing with Chebyshev series, including a Chebyshev
class that
encapsulates the usual arithmetic operations. (General information
on how this module represents and works with such polynomials is in the
docstring for its "parent" sub-package, numpy.polynomial
).
The implementations of multiplication, division, integration, and differentiation use the algebraic identities [1]:
where
These identities allow a Chebyshev series to be expressed as a finite, symmetric Laurent series. In this module, this sort of Laurent series is referred to as a "z-series."
[1] | A. T. Benjamin, et al., "Combinatorial Trigonometry with Chebyshev Polynomials," Journal of Statistical Planning and Inference 14, 2008 (https://web.archive.org/web/20080221202153/https://www.math.hmc.edu/~benjamin/papers/CombTrig.pdf, pg. 4) |
Class | Chebyshev |
A Chebyshev series class. |
Function | cheb2poly |
Convert a Chebyshev series to a polynomial. |
Function | chebadd |
Add one Chebyshev series to another. |
Function | chebcompanion |
Return the scaled companion matrix of c. |
Function | chebder |
Differentiate a Chebyshev series. |
Function | chebdiv |
Divide one Chebyshev series by another. |
Function | chebfit |
Least squares fit of Chebyshev series to data. |
Function | chebfromroots |
Generate a Chebyshev series with given roots. |
Function | chebgauss |
Gauss-Chebyshev quadrature. |
Function | chebgrid2d |
Evaluate a 2-D Chebyshev series on the Cartesian product of x and y. |
Function | chebgrid3d |
Evaluate a 3-D Chebyshev series on the Cartesian product of x, y, and z. |
Function | chebint |
Integrate a Chebyshev series. |
Function | chebinterpolate |
Interpolate a function at the Chebyshev points of the first kind. |
Function | chebline |
Chebyshev series whose graph is a straight line. |
Function | chebmul |
Multiply one Chebyshev series by another. |
Function | chebmulx |
Multiply a Chebyshev series by x. |
Function | chebpow |
Raise a Chebyshev series to a power. |
Function | chebpts1 |
Chebyshev points of the first kind. |
Function | chebpts2 |
Chebyshev points of the second kind. |
Function | chebroots |
Compute the roots of a Chebyshev series. |
Function | chebsub |
Subtract one Chebyshev series from another. |
Function | chebval |
Evaluate a Chebyshev series at points x. |
Function | chebval2d |
Evaluate a 2-D Chebyshev series at points (x, y). |
Function | chebval3d |
Evaluate a 3-D Chebyshev series at points (x, y, z). |
Function | chebvander |
Pseudo-Vandermonde matrix of given degree. |
Function | chebvander2d |
Pseudo-Vandermonde matrix of given degrees. |
Function | chebvander3d |
Pseudo-Vandermonde matrix of given degrees. |
Function | chebweight |
The weight function of the Chebyshev polynomials. |
Function | poly2cheb |
Convert a polynomial to a Chebyshev series. |
Variable | chebdomain |
Undocumented |
Variable | chebone |
Undocumented |
Variable | chebx |
Undocumented |
Variable | chebzero |
Undocumented |
Function | _cseries_to_zseries |
Convert Chebyshev series to z-series. |
Function | _zseries_der |
Differentiate a z-series. |
Function | _zseries_div |
Divide the first z-series by the second. |
Function | _zseries_int |
Integrate a z-series. |
Function | _zseries_mul |
Multiply two z-series. |
Function | _zseries_to_cseries |
Convert z-series to a Chebyshev series. |
Convert a Chebyshev series to a polynomial.
Convert an array representing the coefficients of a Chebyshev series, ordered from lowest degree to highest, to an array of the coefficients of the equivalent polynomial (relative to the "standard" basis) ordered from lowest to highest degree.
poly2cheb
The easy way to do conversions between polynomial basis sets is to use the convert method of a class instance.
>>> from numpy import polynomial as P >>> c = P.Chebyshev(range(4)) >>> c Chebyshev([0., 1., 2., 3.], domain=[-1, 1], window=[-1, 1]) >>> p = c.convert(kind=P.Polynomial) >>> p Polynomial([-2., -8., 4., 12.], domain=[-1., 1.], window=[-1., 1.]) >>> P.chebyshev.cheb2poly(range(4)) array([-2., -8., 4., 12.])
Add one Chebyshev series to another.
Returns the sum of two Chebyshev series c1
+ c2
. The arguments
are sequences of coefficients ordered from lowest order term to
highest, i.e., [1,2,3] represents the series T_0 + 2*T_1 + 3*T_2.
chebsub, chebmulx, chebmul, chebdiv, chebpow
Unlike multiplication, division, etc., the sum of two Chebyshev series is a Chebyshev series (without having to "reproject" the result onto the basis set) so addition, just like that of "standard" polynomials, is simply "component-wise."
>>> from numpy.polynomial import chebyshev as C >>> c1 = (1,2,3) >>> c2 = (3,2,1) >>> C.chebadd(c1,c2) array([4., 4., 4.])
Return the scaled companion matrix of c.
The basis polynomials are scaled so that the companion matrix is
symmetric when c
is a Chebyshev basis polynomial. This provides
better eigenvalue estimates than the unscaled case and for basis
polynomials the eigenvalues are guaranteed to be real if
numpy.linalg.eigvalsh
is used to obtain them.
Differentiate a Chebyshev series.
Returns the Chebyshev series coefficients c
differentiated m
times
along axis
. At each iteration the result is multiplied by scl
(the
scaling factor is for use in a linear change of variable). The argument
c
is an array of coefficients from low to high degree along each
axis, e.g., [1,2,3] represents the series 1*T_0 + 2*T_1 + 3*T_2
while [[1,2],[1,2]] represents 1*T_0(x)*T_0(y) + 1*T_1(x)*T_0(y) +
2*T_0(x)*T_1(y) + 2*T_1(x)*T_1(y) if axis=0 is x and axis=1 is
y.
scl
. The end result is
multiplication by scl**m. This is for use in a linear change of
variable. (Default: 1)Axis over which the derivative is taken. (Default: 0).
chebint
In general, the result of differentiating a C-series needs to be "reprojected" onto the C-series basis set. Thus, typically, the result of this function is "unintuitive," albeit correct; see Examples section below.
>>> from numpy.polynomial import chebyshev as C >>> c = (1,2,3,4) >>> C.chebder(c) array([14., 12., 24.]) >>> C.chebder(c,3) array([96.]) >>> C.chebder(c,scl=-1) array([-14., -12., -24.]) >>> C.chebder(c,2,-1) array([12., 96.])
Divide one Chebyshev series by another.
Returns the quotient-with-remainder of two Chebyshev series
c1
/ c2
. The arguments are sequences of coefficients from lowest
order "term" to highest, e.g., [1,2,3] represents the series
T_0 + 2*T_1 + 3*T_2.
chebadd, chebsub, chebmulx, chebmul, chebpow
In general, the (polynomial) division of one C-series by another results in quotient and remainder terms that are not in the Chebyshev polynomial basis set. Thus, to express these results as C-series, it is typically necessary to "reproject" the results onto said basis set, which typically produces "unintuitive" (but correct) results; see Examples section below.
>>> from numpy.polynomial import chebyshev as C >>> c1 = (1,2,3) >>> c2 = (3,2,1) >>> C.chebdiv(c1,c2) # quotient "intuitive," remainder not (array([3.]), array([-8., -4.])) >>> c2 = (0,1,2,3) >>> C.chebdiv(c2,c1) # neither "intuitive" (array([0., 2.]), array([-2., -4.]))
Least squares fit of Chebyshev series to data.
Return the coefficients of a Chebyshev series of degree deg
that is the
least squares fit to the data values y
given at points x
. If y
is
1-D the returned coefficients will also be 1-D. If y
is 2-D multiple
fits are done, one for each column of y
, and the resulting
coefficients are stored in the corresponding columns of a 2-D return.
The fitted polynomial(s) are in the form
where n
is deg
.
deg
is a single integer,
all terms up to and including the deg
'th term are included in the
fit. For NumPy versions >= 1.11.0 a list of integers specifying the
degrees of the terms to include may be used instead.M
,), optionalWeights. If not None, the weight w[i] applies to the unsquared residual y[i] - y_hat[i] at x[i]. Ideally the weights are chosen so that the errors of the products w[i]*y[i] all have the same variance. When using inverse-variance weighting, use w[i] = 1/sigma(y[i]). The default value is None.
y
was 2-D,
the coefficients for the data in column k of y
are in column
k
.These values are only returned if full == True
rcond
.For more details, see numpy.linalg.lstsq
.
The rank of the coefficient matrix in the least-squares fit is deficient. The warning is only raised if full == False. The warnings can be turned off by
>>> import warnings >>> warnings.simplefilter('ignore', np.RankWarning)
numpy.polynomial.polynomial.polyfit numpy.polynomial.legendre.legfit numpy.polynomial.laguerre.lagfit numpy.polynomial.hermite.hermfit numpy.polynomial.hermite_e.hermefit chebval : Evaluates a Chebyshev series. chebvander : Vandermonde matrix of Chebyshev series. chebweight : Chebyshev weight function. numpy.linalg.lstsq : Computes a least-squares fit from the matrix. scipy.interpolate.UnivariateSpline : Computes spline fits.
The solution is the coefficients of the Chebyshev series p
that
minimizes the sum of the weighted squared errors
where wj are the weights. This problem is solved by setting up as the (typically) overdetermined matrix equation
where V
is the weighted pseudo Vandermonde matrix of x
, c
are the
coefficients to be solved for, w
are the weights, and y
are the
observed values. This equation is then solved using the singular value
decomposition of V
.
If some of the singular values of V
are so small that they are
neglected, then a RankWarning
will be issued. This means that the
coefficient values may be poorly determined. Using a lower order fit
will usually get rid of the warning. The rcond
parameter can also be
set to a value smaller than its default, but the resulting fit may be
spurious and have large contributions from roundoff error.
Fits using Chebyshev series are usually better conditioned than fits using power series, but much can depend on the distribution of the sample points and the smoothness of the data. If the quality of the fit is inadequate splines may be a good alternative.
[1] | Wikipedia, "Curve fitting", https://en.wikipedia.org/wiki/Curve_fitting |
Generate a Chebyshev series with given roots.
The function returns the coefficients of the polynomial
in Chebyshev form, where the r_n
are the roots specified in roots
.
If a zero has multiplicity n, then it must appear in roots
n times.
For instance, if 2 is a root of multiplicity three and 3 is a root of
multiplicity 2, then roots
looks something like [2, 2, 2, 3, 3]. The
roots can appear in any order.
If the returned coefficients are c
, then
The coefficient of the last term is not generally 1 for monic polynomials in Chebyshev form.
out
is a
real array, if some of the roots are complex, then out
is complex
even if all the coefficients in the result are real (see Examples
below).numpy.polynomial.polynomial.polyfromroots numpy.polynomial.legendre.legfromroots numpy.polynomial.laguerre.lagfromroots numpy.polynomial.hermite.hermfromroots numpy.polynomial.hermite_e.hermefromroots
>>> import numpy.polynomial.chebyshev as C >>> C.chebfromroots((-1,0,1)) # x^3 - x relative to the standard basis array([ 0. , -0.25, 0. , 0.25]) >>> j = complex(0,1) >>> C.chebfromroots((-j,j)) # x^2 + 1 relative to the standard basis array([1.5+0.j, 0. +0.j, 0.5+0.j])
Gauss-Chebyshev quadrature.
Computes the sample points and weights for Gauss-Chebyshev quadrature. These sample points and weights will correctly integrate polynomials of degree 2*deg − 1 or less over the interval [ − 1, 1] with the weight function f(x) = 1 ⁄ √(1 − x2).
The results have only been tested up to degree 100, higher degrees may
be problematic. For Gauss-Chebyshev there are closed form solutions for
the sample points and weights. If n = deg
, then
Evaluate a 2-D Chebyshev series on the Cartesian product of x and y.
This function returns the values:
where the points (a, b)
consist of all pairs formed by taking
a
from x
and b
from y
. The resulting points form a grid with
x
in the first dimension and y
in the second.
The parameters x
and y
are converted to arrays only if they are
tuples or a lists, otherwise they are treated as a scalars. In either
case, either x
and y
or their elements must support multiplication
and addition both with themselves and with the elements of c
.
If c
has fewer than two dimensions, ones are implicitly appended to
its shape to make it 2-D. The shape of the result will be c.shape[2:] +
x.shape + y.shape.
x
and y
. If x
or y
is a list or
tuple, it is first converted to an ndarray, otherwise it is left
unchanged and, if it isn't an ndarray, it is treated as a scalar.c[i,j]
. If c
has dimension
greater than two the remaining indices enumerate multiple sets of
coefficients.x
and y
.chebval, chebval2d, chebval3d, chebgrid3d
Evaluate a 3-D Chebyshev series on the Cartesian product of x, y, and z.
This function returns the values:
where the points (a, b, c)
consist of all triples formed by taking
a
from x
, b
from y
, and c
from z
. The resulting points form
a grid with x
in the first dimension, y
in the second, and z
in
the third.
The parameters x
, y
, and z
are converted to arrays only if they
are tuples or a lists, otherwise they are treated as a scalars. In
either case, either x
, y
, and z
or their elements must support
multiplication and addition both with themselves and with the elements
of c
.
If c
has fewer than three dimensions, ones are implicitly appended to
its shape to make it 3-D. The shape of the result will be c.shape[3:] +
x.shape + y.shape + z.shape.
x
, y
, and z
. If x
,`y`, or z
is a
list or tuple, it is first converted to an ndarray, otherwise it is
left unchanged and, if it isn't an ndarray, it is treated as a
scalar.c
has dimension
greater than two the remaining indices enumerate multiple sets of
coefficients.x
and y
.chebval, chebval2d, chebgrid2d, chebval3d
Integrate a Chebyshev series.
Returns the Chebyshev series coefficients c
integrated m
times from
lbnd
along axis
. At each iteration the resulting series is
multiplied by scl
and an integration constant, k
, is added.
The scaling factor is for use in a linear change of variable. ("Buyer
beware": note that, depending on what one is doing, one may want scl
to be the reciprocal of what one might expect; for more information,
see the Notes section below.) The argument c
is an array of
coefficients from low to high degree along each axis, e.g., [1,2,3]
represents the series T_0 + 2*T_1 + 3*T_2 while [[1,2],[1,2]]
represents 1*T_0(x)*T_0(y) + 1*T_1(x)*T_0(y) + 2*T_0(x)*T_1(y) +
2*T_1(x)*T_1(y) if axis=0 is x and axis=1 is y.
scl
before the integration constant is added. (Default: 1)Axis over which the integral is taken. (Default: 0).
chebder
Note that the result of each integration is multiplied by scl
.
Why is this important to note? Say one is making a linear change of
variable u = ax + b in an integral relative to x
. Then
dx = du ⁄ a, so one will need to set scl
equal to
1 ⁄ a- perhaps not what one would have first thought.
Also note that, in general, the result of integrating a C-series needs to be "reprojected" onto the C-series basis set. Thus, typically, the result of this function is "unintuitive," albeit correct; see Examples section below.
>>> from numpy.polynomial import chebyshev as C >>> c = (1,2,3) >>> C.chebint(c) array([ 0.5, -0.5, 0.5, 0.5]) >>> C.chebint(c,3) array([ 0.03125 , -0.1875 , 0.04166667, -0.05208333, 0.01041667, # may vary 0.00625 ]) >>> C.chebint(c, k=3) array([ 3.5, -0.5, 0.5, 0.5]) >>> C.chebint(c,lbnd=-2) array([ 8.5, -0.5, 0.5, 0.5]) >>> C.chebint(c,scl=-2) array([-1., 1., -1., -1.])
Interpolate a function at the Chebyshev points of the first kind.
Returns the Chebyshev series that interpolates func
at the Chebyshev
points of the first kind in the interval [-1, 1]. The interpolating
series tends to a minmax approximation to func
with increasing deg
if the function is continuous in the interval.
args
parameter.>>> import numpy.polynomial.chebyshev as C >>> C.chebfromfunction(lambda x: np.tanh(x) + 0.5, 8) array([ 5.00000000e-01, 8.11675684e-01, -9.86864911e-17, -5.42457905e-02, -2.71387850e-16, 4.51658839e-03, 2.46716228e-17, -3.79694221e-04, -3.26899002e-16])
The Chebyshev polynomials used in the interpolation are orthogonal when sampled at the Chebyshev points of the first kind. If it is desired to constrain some of the coefficients they can simply be set to the desired value after the interpolation, no new interpolation or fit is needed. This is especially useful if it is known apriori that some of coefficients are zero. For instance, if the function is even then the coefficients of the terms of odd degree in the result can be set to zero.
Chebyshev series whose graph is a straight line.
numpy.polynomial.polynomial.polyline numpy.polynomial.legendre.legline numpy.polynomial.laguerre.lagline numpy.polynomial.hermite.hermline numpy.polynomial.hermite_e.hermeline
>>> import numpy.polynomial.chebyshev as C >>> C.chebline(3,2) array([3, 2]) >>> C.chebval(-3, C.chebline(3,2)) # should be -3 -3.0
Multiply one Chebyshev series by another.
Returns the product of two Chebyshev series c1
* c2
. The arguments
are sequences of coefficients, from lowest order "term" to highest,
e.g., [1,2,3] represents the series T_0 + 2*T_1 + 3*T_2.
chebadd, chebsub, chebmulx, chebdiv, chebpow
In general, the (polynomial) product of two C-series results in terms that are not in the Chebyshev polynomial basis set. Thus, to express the product as a C-series, it is typically necessary to "reproject" the product onto said basis set, which typically produces "unintuitive live" (but correct) results; see Examples section below.
>>> from numpy.polynomial import chebyshev as C >>> c1 = (1,2,3) >>> c2 = (3,2,1) >>> C.chebmul(c1,c2) # multiplication requires "reprojection" array([ 6.5, 12. , 12. , 4. , 1.5])
Multiply a Chebyshev series by x.
Multiply the polynomial c
by x, where x is the independent
variable.
>>> from numpy.polynomial import chebyshev as C >>> C.chebmulx([1,2,3]) array([1. , 2.5, 1. , 1.5])
Raise a Chebyshev series to a power.
Returns the Chebyshev series c
raised to the power pow
. The
argument c
is a sequence of coefficients ordered from low to high.
i.e., [1,2,3] is the series T_0 + 2*T_1 + 3*T_2.
chebadd, chebsub, chebmulx, chebmul, chebdiv
>>> from numpy.polynomial import chebyshev as C >>> C.chebpow([1, 2, 3, 4], 2) array([15.5, 22. , 16. , ..., 12.5, 12. , 8. ])
Chebyshev points of the first kind.
The Chebyshev points of the first kind are the points cos(x), where x = [pi*(k + .5)/npts for k in range(npts)].
chebpts2
Chebyshev points of the second kind.
The Chebyshev points of the second kind are the points cos(x), where x = [pi*k/(npts - 1) for k in range(npts)].
Compute the roots of a Chebyshev series.
Return the roots (a.k.a. "zeros") of the polynomial
out
is also real, otherwise it is complex.numpy.polynomial.polynomial.polyroots numpy.polynomial.legendre.legroots numpy.polynomial.laguerre.lagroots numpy.polynomial.hermite.hermroots numpy.polynomial.hermite_e.hermeroots
The root estimates are obtained as the eigenvalues of the companion matrix, Roots far from the origin of the complex plane may have large errors due to the numerical instability of the series for such values. Roots with multiplicity greater than 1 will also show larger errors as the value of the series near such points is relatively insensitive to errors in the roots. Isolated roots near the origin can be improved by a few iterations of Newton's method.
The Chebyshev series basis polynomials aren't powers of x
so the
results of this function may seem unintuitive.
>>> import numpy.polynomial.chebyshev as cheb >>> cheb.chebroots((-1, 1,-1, 1)) # T3 - T2 + T1 - T0 has real roots array([ -5.00000000e-01, 2.60860684e-17, 1.00000000e+00]) # may vary
Subtract one Chebyshev series from another.
Returns the difference of two Chebyshev series c1
- c2
. The
sequences of coefficients are from lowest order term to highest, i.e.,
[1,2,3] represents the series T_0 + 2*T_1 + 3*T_2.
chebadd, chebmulx, chebmul, chebdiv, chebpow
Unlike multiplication, division, etc., the difference of two Chebyshev series is a Chebyshev series (without having to "reproject" the result onto the basis set) so subtraction, just like that of "standard" polynomials, is simply "component-wise."
>>> from numpy.polynomial import chebyshev as C >>> c1 = (1,2,3) >>> c2 = (3,2,1) >>> C.chebsub(c1,c2) array([-2., 0., 2.]) >>> C.chebsub(c2,c1) # -C.chebsub(c1,c2) array([ 2., 0., -2.])
Evaluate a Chebyshev series at points x.
If c
is of length n + 1
, this function returns the value:
The parameter x
is converted to an array only if it is a tuple or a
list, otherwise it is treated as a scalar. In either case, either x
or its elements must support multiplication and addition both with
themselves and with the elements of c
.
If c
is a 1-D array, then p(x)
will have the same shape as x
. If
c
is multidimensional, then the shape of the result depends on the
value of tensor
. If tensor
is true the shape will be c.shape[1:] +
x.shape. If tensor
is false the shape will be c.shape[1:]. Note that
scalars have shape (,).
Trailing zeros in the coefficients will be used in the evaluation, so they should be avoided if efficiency is a concern.
x
is a list or tuple, it is converted to an ndarray, otherwise
it is left unchanged and treated as a scalar. In either case, x
or its elements must support addition and multiplication with
with themselves and with the elements of c
.c
is multidimensional the
remaining indices enumerate multiple polynomials. In the two
dimensional case the coefficients may be thought of as stored in
the columns of c
.If True, the shape of the coefficient array is extended with ones
on the right, one for each dimension of x
. Scalars have dimension 0
for this action. The result is that every column of coefficients in
c
is evaluated for every element of x
. If False, x
is broadcast
over the columns of c
for the evaluation. This keyword is useful
when c
is multidimensional. The default value is True.
chebval2d, chebgrid2d, chebval3d, chebgrid3d
The evaluation uses Clenshaw recursion, aka synthetic division.
Evaluate a 2-D Chebyshev series at points (x, y).
This function returns the values:
The parameters x
and y
are converted to arrays only if they are
tuples or a lists, otherwise they are treated as a scalars and they
must have the same shape after conversion. In either case, either x
and y
or their elements must support multiplication and addition both
with themselves and with the elements of c
.
If c
is a 1-D array a one is implicitly appended to its shape to make
it 2-D. The shape of the result will be c.shape[2:] + x.shape.
(x, y)
,
where x
and y
must have the same shape. If x
or y
is a list
or tuple, it is first converted to an ndarray, otherwise it is left
unchanged and if it isn't an ndarray it is treated as a scalar.c
has
dimension greater than 2 the remaining indices enumerate multiple
sets of coefficients.x
and y
.chebval, chebgrid2d, chebval3d, chebgrid3d
Evaluate a 3-D Chebyshev series at points (x, y, z).
This function returns the values:
The parameters x
, y
, and z
are converted to arrays only if
they are tuples or a lists, otherwise they are treated as a scalars and
they must have the same shape after conversion. In either case, either
x
, y
, and z
or their elements must support multiplication and
addition both with themselves and with the elements of c
.
If c
has fewer than 3 dimensions, ones are implicitly appended to its
shape to make it 3-D. The shape of the result will be c.shape[3:] +
x.shape.
(x, y, z)
, where x
, y
, and z
must have the same shape. If
any of x
, y
, or z
is a list or tuple, it is first converted
to an ndarray, otherwise it is left unchanged and if it isn't an
ndarray it is treated as a scalar.c
has dimension
greater than 3 the remaining indices enumerate multiple sets of
coefficients.x
, y
, and z
.chebval, chebval2d, chebgrid2d, chebgrid3d
Pseudo-Vandermonde matrix of given degree.
Returns the pseudo-Vandermonde matrix of degree deg
and sample points
x
. The pseudo-Vandermonde matrix is defined by
where 0 <= i <= deg
. The leading indices of V
index the elements of
x
and the last index is the degree of the Chebyshev polynomial.
If c
is a 1-D array of coefficients of length n + 1
and V
is the
matrix V = chebvander(x, n), then np.dot(V, c) and
chebval(x, c) are the same up to roundoff. This equivalence is
useful both for least squares fitting and for the evaluation of a large
number of Chebyshev series of the same degree and sample points.
x
is
scalar it is converted to a 1-D array.x
.Pseudo-Vandermonde matrix of given degrees.
Returns the pseudo-Vandermonde matrix of degrees deg
and sample
points (x, y)
. The pseudo-Vandermonde matrix is defined by
where 0 <= i <= deg[0]
and 0 <= j <= deg[1]
. The leading indices of
V
index the points (x, y)
and the last index encodes the degrees of
the Chebyshev polynomials.
If V = chebvander2d(x, y, [xdeg, ydeg]), then the columns of V
correspond to the elements of a 2-D coefficient array c
of shape
(xdeg + 1, ydeg + 1) in the order
and np.dot(V, c.flat) and chebval2d(x, y, c) will be the same up to roundoff. This equivalence is useful both for least squares fitting and for the evaluation of a large number of 2-D Chebyshev series of the same degrees and sample points.
x
and y
.chebvander, chebvander3d, chebval2d, chebval3d
Pseudo-Vandermonde matrix of given degrees.
Returns the pseudo-Vandermonde matrix of degrees deg
and sample
points (x, y, z)
. If l, m, n
are the given degrees in x, y, z
,
then The pseudo-Vandermonde matrix is defined by
where 0 <= i <= l
, 0 <= j <= m
, and 0 <= j <= n
. The leading
indices of V
index the points (x, y, z)
and the last index encodes
the degrees of the Chebyshev polynomials.
If V = chebvander3d(x, y, z, [xdeg, ydeg, zdeg]), then the columns
of V
correspond to the elements of a 3-D coefficient array c
of
shape (xdeg + 1, ydeg + 1, zdeg + 1) in the order
and np.dot(V, c.flat) and chebval3d(x, y, z, c) will be the same up to roundoff. This equivalence is useful both for least squares fitting and for the evaluation of a large number of 3-D Chebyshev series of the same degrees and sample points.
x
, y
, and z
.chebvander, chebvander3d, chebval2d, chebval3d
The weight function of the Chebyshev polynomials.
The weight function is 1 ⁄ √(1 − x2) and the interval of integration is [ − 1, 1]. The Chebyshev polynomials are orthogonal, but not normalized, with respect to this weight function.
x
.Convert a polynomial to a Chebyshev series.
Convert an array representing the coefficients of a polynomial (relative to the "standard" basis) ordered from lowest degree to highest, to an array of the coefficients of the equivalent Chebyshev series, ordered from lowest to highest degree.
cheb2poly
The easy way to do conversions between polynomial basis sets is to use the convert method of a class instance.
>>> from numpy import polynomial as P >>> p = P.Polynomial(range(4)) >>> p Polynomial([0., 1., 2., 3.], domain=[-1, 1], window=[-1, 1]) >>> c = p.convert(kind=P.Chebyshev) >>> c Chebyshev([1. , 3.25, 1. , 0.75], domain=[-1., 1.], window=[-1., 1.]) >>> P.chebyshev.poly2cheb(range(4)) array([1. , 3.25, 1. , 0.75])
Convert Chebyshev series to z-series.
Convert a Chebyshev series to the equivalent z-series. The result is never an empty array. The dtype of the return is the same as that of the input. No checks are run on the arguments as this routine is for internal use.
Differentiate a z-series.
The derivative is with respect to x, not z. This is achieved using the chain rule and the value of dx/dz given in the module notes.
The zseries for x (ns) has been multiplied by two in order to avoid using floats that are incompatible with Decimal and likely other specialized scalar types. This scaling has been compensated by multiplying the value of zs by two also so that the two cancels in the division.
Divide the first z-series by the second.
Divide z1
by z2
and return the quotient and remainder as z-series.
Warning: this implementation only applies when both z1 and z2 have the
same symmetry, which is sufficient for present purposes.
This is not the same as polynomial division on account of the desired form of the remainder. If symmetric/anti-symmetric z-series are denoted by S/A then the following rules apply:
S/S -> S,S A/A -> S,A
The restriction to types of the same symmetry could be fixed but seems like unneeded generality. There is no natural form for the remainder in the case where there is no symmetry.
Integrate a z-series.
The integral is with respect to x, not z. This is achieved by a change of variable using dx/dz given in the module notes.
The zseries for x (ns) has been multiplied by two in order to avoid using floats that are incompatible with Decimal and likely other specialized scalar types. This scaling has been compensated by dividing the resulting zs by two.
Multiply two z-series.
Multiply two z-series to produce a z-series.
This is simply convolution. If symmetric/anti-symmetric z-series are denoted by S/A then the following rules apply:
S*S, A*A -> S S*A, A*S -> A
Convert z-series to a Chebyshev series.
Convert a z series to the equivalent Chebyshev series. The result is never an empty array. The dtype of the return is the same as that of the input. No checks are run on the arguments as this routine is for internal use.