module documentation

This module provides a number of objects (mostly functions) useful for dealing with Legendre series, including a Legendre 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).

Classes

Constants

Arithmetic

Calculus

Misc Functions

See Also

numpy.polynomial

Class Legendre A Legendre series class.
Function leg2poly Convert a Legendre series to a polynomial.
Function legadd Add one Legendre series to another.
Function legcompanion Return the scaled companion matrix of c.
Function legder Differentiate a Legendre series.
Function legdiv Divide one Legendre series by another.
Function legfit Least squares fit of Legendre series to data.
Function legfromroots Generate a Legendre series with given roots.
Function leggauss Gauss-Legendre quadrature.
Function leggrid2d Evaluate a 2-D Legendre series on the Cartesian product of x and y.
Function leggrid3d Evaluate a 3-D Legendre series on the Cartesian product of x, y, and z.
Function legint Integrate a Legendre series.
Function legline Legendre series whose graph is a straight line.
Function legmul Multiply one Legendre series by another.
Function legmulx Multiply a Legendre series by x.
Function legpow Raise a Legendre series to a power.
Function legroots Compute the roots of a Legendre series.
Function legsub Subtract one Legendre series from another.
Function legval Evaluate a Legendre series at points x.
Function legval2d Evaluate a 2-D Legendre series at points (x, y).
Function legval3d Evaluate a 3-D Legendre series at points (x, y, z).
Function legvander Pseudo-Vandermonde matrix of given degree.
Function legvander2d Pseudo-Vandermonde matrix of given degrees.
Function legvander3d Pseudo-Vandermonde matrix of given degrees.
Function legweight Weight function of the Legendre polynomials.
Function poly2leg Convert a polynomial to a Legendre series.
Variable legdomain Undocumented
Variable legone Undocumented
Variable legx Undocumented
Variable legzero Undocumented
def leg2poly(c): (source)

Convert a Legendre series to a polynomial.

Convert an array representing the coefficients of a Legendre 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.

See Also

poly2leg

Notes

The easy way to do conversions between polynomial basis sets is to use the convert method of a class instance.

Examples

>>> from numpy import polynomial as P
>>> c = P.Legendre(range(4))
>>> c
Legendre([0., 1., 2., 3.], domain=[-1,  1], window=[-1,  1])
>>> p = c.convert(kind=P.Polynomial)
>>> p
Polynomial([-1. , -3.5,  3. ,  7.5], domain=[-1.,  1.], window=[-1.,  1.])
>>> P.legendre.leg2poly(range(4))
array([-1. , -3.5,  3. ,  7.5])
Parameters
c:array_like1-D array containing the Legendre series coefficients, ordered from lowest order term to highest.
Returns
ndarraypol - 1-D array containing the coefficients of the equivalent polynomial (relative to the "standard" basis) ordered from lowest order term to highest.
def legadd(c1, c2): (source)

Add one Legendre series to another.

Returns the sum of two Legendre series c1 + c2. The arguments are sequences of coefficients ordered from lowest order term to highest, i.e., [1,2,3] represents the series P_0 + 2*P_1 + 3*P_2.

Notes

Unlike multiplication, division, etc., the sum of two Legendre series is a Legendre series (without having to "reproject" the result onto the basis set) so addition, just like that of "standard" polynomials, is simply "component-wise."

Examples

>>> from numpy.polynomial import legendre as L
>>> c1 = (1,2,3)
>>> c2 = (3,2,1)
>>> L.legadd(c1,c2)
array([4.,  4.,  4.])
Parameters
c1:array_like1-D arrays of Legendre series coefficients ordered from low to high.
c2:array_like1-D arrays of Legendre series coefficients ordered from low to high.
Returns
ndarrayout - Array representing the Legendre series of their sum.
def legcompanion(c): (source)

Return the scaled companion matrix of c.

The basis polynomials are scaled so that the companion matrix is symmetric when c is an Legendre 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.

Notes

New in version 1.7.0.
Parameters
c:array_like1-D array of Legendre series coefficients ordered from low to high degree.
Returns
ndarraymat - Scaled companion matrix of dimensions (deg, deg).
def legder(c, m=1, scl=1, axis=0): (source)

Differentiate a Legendre series.

Returns the Legendre 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*L_0 + 2*L_1 + 3*L_2 while [[1,2],[1,2]] represents 1*L_0(x)*L_0(y) + 1*L_1(x)*L_0(y) + 2*L_0(x)*L_1(y) + 2*L_1(x)*L_1(y) if axis=0 is x and axis=1 is y.

See Also

legint

Notes

In general, the result of differentiating a Legendre series does not resemble the same operation on a power series. Thus the result of this function may be "unintuitive," albeit correct; see Examples section below.

Examples

>>> from numpy.polynomial import legendre as L
>>> c = (1,2,3,4)
>>> L.legder(c)
array([  6.,   9.,  20.])
>>> L.legder(c, 3)
array([60.])
>>> L.legder(c, scl=-1)
array([ -6.,  -9., -20.])
>>> L.legder(c, 2,-1)
array([  9.,  60.])
Parameters
c:array_likeArray of Legendre series coefficients. If c is multidimensional the different axis correspond to different variables with the degree in each axis given by the corresponding index.
m:int, optionalNumber of derivatives taken, must be non-negative. (Default: 1)
scl:scalar, optionalEach differentiation is multiplied by scl. The end result is multiplication by scl**m. This is for use in a linear change of variable. (Default: 1)
axis:int, optional

Axis over which the derivative is taken. (Default: 0).

New in version 1.7.0.
Returns
ndarrayder - Legendre series of the derivative.
def legdiv(c1, c2): (source)

Divide one Legendre series by another.

Returns the quotient-with-remainder of two Legendre series c1 / c2. The arguments are sequences of coefficients from lowest order "term" to highest, e.g., [1,2,3] represents the series P_0 + 2*P_1 + 3*P_2.

Notes

In general, the (polynomial) division of one Legendre series by another results in quotient and remainder terms that are not in the Legendre polynomial basis set. Thus, to express these results as a Legendre series, it is necessary to "reproject" the results onto the Legendre basis set, which may produce "unintuitive" (but correct) results; see Examples section below.

Examples

>>> from numpy.polynomial import legendre as L
>>> c1 = (1,2,3)
>>> c2 = (3,2,1)
>>> L.legdiv(c1,c2) # quotient "intuitive," remainder not
(array([3.]), array([-8., -4.]))
>>> c2 = (0,1,2,3)
>>> L.legdiv(c2,c1) # neither "intuitive"
(array([-0.07407407,  1.66666667]), array([-1.03703704, -2.51851852])) # may vary
Parameters
c1:array_like1-D arrays of Legendre series coefficients ordered from low to high.
c2:array_like1-D arrays of Legendre series coefficients ordered from low to high.
Returns
ndarraysquo, rem - Of Legendre series coefficients representing the quotient and remainder.
def legfit(x, y, deg, rcond=None, full=False, w=None): (source)

Least squares fit of Legendre series to data.

Return the coefficients of a Legendre 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

p(x) = c0 + c1*L1(x) + ... + cn*Ln(x), 

where n is deg.

See Also

numpy.polynomial.polynomial.polyfit, numpy.polynomial.chebyshev.chebfit, numpy.polynomial.laguerre.lagfit, numpy.polynomial.hermite.hermfit, numpy.polynomial.hermite_e.hermefit

legval
Evaluates a Legendre series.
legvander
Vandermonde matrix of Legendre series.
legweight
Legendre weight function (= 1).
numpy.linalg.lstsq
Computes a least-squares fit from the matrix.
scipy.interpolate.UnivariateSpline
Computes spline fits.

Notes

The solution is the coefficients of the Legendre series p that minimizes the sum of the weighted squared errors

E = jw2j*|yj − p(xj)|2, 

where wj are the weights. This problem is solved by setting up as the (typically) overdetermined matrix equation

V(x)*c = w*y, 

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

References

[1]Wikipedia, "Curve fitting", https://en.wikipedia.org/wiki/Curve_fitting
Parameters
x:array_like, shape(M, )x-coordinates of the M sample points (x[i], y[i]).
y:array_like, shape(M, ) or (M, K)y-coordinates of the sample points. Several data sets of sample points sharing the same x-coordinates can be fitted at once by passing in a 2D-array that contains one dataset per column.
deg:int or 1-D array_likeDegree(s) of the fitting polynomials. If 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.
rcond:float, optionalRelative condition number of the fit. Singular values smaller than this relative to the largest singular value will be ignored. The default value is len(x)*eps, where eps is the relative precision of the float type, about 2e-16 in most cases.
full:bool, optionalSwitch determining nature of return value. When it is False (the default) just the coefficients are returned, when True diagnostic information from the singular value decomposition is also returned.
w:array_like, shape(M, ), optional

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

New in version 1.5.0.
Returns
  • coef: ndarray, shape (M, ) or (M, K) - Legendre coefficients ordered from low to high. If y was 2-D, the coefficients for the data in column k of y are in column k. If deg is specified as a list, coefficients for terms not included in the fit are set equal to zero in the returned coef.

  • [residuals, rank, singular_values, rcond]: list - These values are only returned if full == True

    • residuals -- sum of squared residuals of the least squares fit
    • rank -- the numerical rank of the scaled Vandermonde matrix
    • singular_values -- singular values of the scaled Vandermonde matrix
    • rcond -- value of rcond.

    For more details, see numpy.linalg.lstsq.

Warns
RankWarning

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)
def legfromroots(roots): (source)

Generate a Legendre series with given roots.

The function returns the coefficients of the polynomial

p(x) = (x − r0)*(x − r1)*...*(x − rn), 

in Legendre 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

p(x) = c0 + c1*L1(x) + ... + cn*Ln(x)

The coefficient of the last term is not generally 1 for monic polynomials in Legendre form.

Examples

>>> import numpy.polynomial.legendre as L
>>> L.legfromroots((-1,0,1)) # x^3 - x relative to the standard basis
array([ 0. , -0.4,  0. ,  0.4])
>>> j = complex(0,1)
>>> L.legfromroots((-j,j)) # x^2 + 1 relative to the standard basis
array([ 1.33333333+0.j,  0.00000000+0.j,  0.66666667+0.j]) # may vary
Parameters
roots:array_likeSequence containing the roots.
Returns
ndarrayout - 1-D array of coefficients. If all roots are real then 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).
def leggauss(deg): (source)

Gauss-Legendre quadrature.

Computes the sample points and weights for Gauss-Legendre 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.

Notes

New in version 1.7.0.

The results have only been tested up to degree 100, higher degrees may be problematic. The weights are determined by using the fact that

wk = c ⁄ (Ln(xk)*Ln − 1(xk))

where c is a constant independent of k and xk is the k'th root of Ln, and then scaling the results to get the right value when integrating 1.

Parameters
deg:intNumber of sample points and weights. It must be >= 1.
Returns
  • x: ndarray - 1-D ndarray containing the sample points.
  • y: ndarray - 1-D ndarray containing the weights.
def leggrid2d(x, y, c): (source)

Evaluate a 2-D Legendre series on the Cartesian product of x and y.

This function returns the values:

p(a, b) = i, jci, j*Li(a)*Lj(b)

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.

Notes

New in version 1.7.0.
Parameters
x:array_like, compatible objectsThe two dimensional series is evaluated at the points in the Cartesian product of 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.
y:array_like, compatible objectsThe two dimensional series is evaluated at the points in the Cartesian product of 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:array_likeArray of coefficients ordered so that the coefficient of the term of multi-degree i,j is contained in c[i,j]. If c has dimension greater than two the remaining indices enumerate multiple sets of coefficients.
Returns
ndarray, compatible objectvalues - The values of the two dimensional Chebyshev series at points in the Cartesian product of x and y.
def leggrid3d(x, y, z, c): (source)

Evaluate a 3-D Legendre series on the Cartesian product of x, y, and z.

This function returns the values:

p(a, b, c) = i, j, kci, j, k*Li(a)*Lj(b)*Lk(c)

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.

Notes

New in version 1.7.0.
Parameters
x:array_like, compatible objectsThe three dimensional series is evaluated at the points in the Cartesian product of 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.
y:array_like, compatible objectsThe three dimensional series is evaluated at the points in the Cartesian product of 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.
z:array_like, compatible objectsThe three dimensional series is evaluated at the points in the Cartesian product of 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:array_likeArray of coefficients ordered so that the coefficients for terms of degree i,j are contained in c[i,j]. If c has dimension greater than two the remaining indices enumerate multiple sets of coefficients.
Returns
ndarray, compatible objectvalues - The values of the two dimensional polynomial at points in the Cartesian product of x and y.
def legint(c, m=1, k=[], lbnd=0, scl=1, axis=0): (source)

Integrate a Legendre series.

Returns the Legendre 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 L_0 + 2*L_1 + 3*L_2 while [[1,2],[1,2]] represents 1*L_0(x)*L_0(y) + 1*L_1(x)*L_0(y) + 2*L_0(x)*L_1(y) + 2*L_1(x)*L_1(y) if axis=0 is x and axis=1 is y.

See Also

legder

Notes

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.

Examples

>>> from numpy.polynomial import legendre as L
>>> c = (1,2,3)
>>> L.legint(c)
array([ 0.33333333,  0.4       ,  0.66666667,  0.6       ]) # may vary
>>> L.legint(c, 3)
array([  1.66666667e-02,  -1.78571429e-02,   4.76190476e-02, # may vary
         -1.73472348e-18,   1.90476190e-02,   9.52380952e-03])
>>> L.legint(c, k=3)
 array([ 3.33333333,  0.4       ,  0.66666667,  0.6       ]) # may vary
>>> L.legint(c, lbnd=-2)
array([ 7.33333333,  0.4       ,  0.66666667,  0.6       ]) # may vary
>>> L.legint(c, scl=2)
array([ 0.66666667,  0.8       ,  1.33333333,  1.2       ]) # may vary
Parameters
c:array_likeArray of Legendre series coefficients. If c is multidimensional the different axis correspond to different variables with the degree in each axis given by the corresponding index.
m:int, optionalOrder of integration, must be positive. (Default: 1)
k:{[], list, scalar}, optionalIntegration constant(s). The value of the first integral at lbnd is the first value in the list, the value of the second integral at lbnd is the second value, etc. If k == [] (the default), all constants are set to zero. If m == 1, a single scalar can be given instead of a list.
lbnd:scalar, optionalThe lower bound of the integral. (Default: 0)
scl:scalar, optionalFollowing each integration the result is multiplied by scl before the integration constant is added. (Default: 1)
axis:int, optional

Axis over which the integral is taken. (Default: 0).

New in version 1.7.0.
Returns
ndarrayS - Legendre series coefficient array of the integral.
Raises
ValueErrorIf m < 0, len(k) > m, np.ndim(lbnd) != 0, or np.ndim(scl) != 0.
def legline(off, scl): (source)

Legendre series whose graph is a straight line.

Examples

>>> import numpy.polynomial.legendre as L
>>> L.legline(3,2)
array([3, 2])
>>> L.legval(-3, L.legline(3,2)) # should be -3
-3.0
Parameters
off:scalarsThe specified line is given by off + scl*x.
scl:scalarsThe specified line is given by off + scl*x.
Returns
ndarrayy - This module's representation of the Legendre series for off + scl*x.
def legmul(c1, c2): (source)

Multiply one Legendre series by another.

Returns the product of two Legendre series c1 * c2. The arguments are sequences of coefficients, from lowest order "term" to highest, e.g., [1,2,3] represents the series P_0 + 2*P_1 + 3*P_2.

Notes

In general, the (polynomial) product of two C-series results in terms that are not in the Legendre polynomial basis set. Thus, to express the product as a Legendre series, it is necessary to "reproject" the product onto said basis set, which may produce "unintuitive" (but correct) results; see Examples section below.

Examples

>>> from numpy.polynomial import legendre as L
>>> c1 = (1,2,3)
>>> c2 = (3,2)
>>> L.legmul(c1,c2) # multiplication requires "reprojection"
array([  4.33333333,  10.4       ,  11.66666667,   3.6       ]) # may vary
Parameters
c1:array_like1-D arrays of Legendre series coefficients ordered from low to high.
c2:array_like1-D arrays of Legendre series coefficients ordered from low to high.
Returns
ndarrayout - Of Legendre series coefficients representing their product.
def legmulx(c): (source)

Multiply a Legendre series by x.

Multiply the Legendre series c by x, where x is the independent variable.

See Also

legadd, legmul, legdiv, legpow

Notes

The multiplication uses the recursion relationship for Legendre polynomials in the form

xPi(x) = ((i + 1)*Pi + 1(x) + i*Pi − 1(x)) ⁄ (2i + 1)

Examples

>>> from numpy.polynomial import legendre as L
>>> L.legmulx([1,2,3])
array([ 0.66666667, 2.2, 1.33333333, 1.8]) # may vary
Parameters
c:array_like1-D array of Legendre series coefficients ordered from low to high.
Returns
ndarrayout - Array representing the result of the multiplication.
def legpow(c, pow, maxpower=16): (source)

Raise a Legendre series to a power.

Returns the Legendre 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 P_0 + 2*P_1 + 3*P_2.

Parameters
c:array_like1-D array of Legendre series coefficients ordered from low to high.
pow:integerPower to which the series will be raised
maxpower:integer, optionalMaximum power allowed. This is mainly to limit growth of the series to unmanageable size. Default is 16
Returns
ndarraycoef - Legendre series of power.
def legroots(c): (source)

Compute the roots of a Legendre series.

Return the roots (a.k.a. "zeros") of the polynomial

p(x) = ic[i]*Li(x).

Notes

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 Legendre series basis polynomials aren't powers of x so the results of this function may seem unintuitive.

Examples

>>> import numpy.polynomial.legendre as leg
>>> leg.legroots((1, 2, 3, 4)) # 4L_3 + 3L_2 + 2L_1 + 1L_0, all real roots
array([-0.85099543, -0.11407192,  0.51506735]) # may vary
Parameters
c:1-D array_like1-D array of coefficients.
Returns
ndarrayout - Array of the roots of the series. If all the roots are real, then out is also real, otherwise it is complex.
def legsub(c1, c2): (source)

Subtract one Legendre series from another.

Returns the difference of two Legendre series c1 - c2. The sequences of coefficients are from lowest order term to highest, i.e., [1,2,3] represents the series P_0 + 2*P_1 + 3*P_2.

Notes

Unlike multiplication, division, etc., the difference of two Legendre series is a Legendre series (without having to "reproject" the result onto the basis set) so subtraction, just like that of "standard" polynomials, is simply "component-wise."

Examples

>>> from numpy.polynomial import legendre as L
>>> c1 = (1,2,3)
>>> c2 = (3,2,1)
>>> L.legsub(c1,c2)
array([-2.,  0.,  2.])
>>> L.legsub(c2,c1) # -C.legsub(c1,c2)
array([ 2.,  0., -2.])
Parameters
c1:array_like1-D arrays of Legendre series coefficients ordered from low to high.
c2:array_like1-D arrays of Legendre series coefficients ordered from low to high.
Returns
ndarrayout - Of Legendre series coefficients representing their difference.
def legval(x, c, tensor=True): (source)

Evaluate a Legendre series at points x.

If c is of length n + 1, this function returns the value:

p(x) = c0*L0(x) + c1*L1(x) + ... + cn*Ln(x)

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.

Notes

The evaluation uses Clenshaw recursion, aka synthetic division.

Parameters
x:array_like, compatible objectIf 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 themselves and with the elements of c.
c:array_likeArray of coefficients ordered so that the coefficients for terms of degree n are contained in c[n]. If 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.
tensor:boolean, optional

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.

New in version 1.7.0.
Returns
ndarray, algebra_likevalues - The shape of the return value is described above.
def legval2d(x, y, c): (source)

Evaluate a 2-D Legendre series at points (x, y).

This function returns the values:

p(x, y) = i, jci, j*Li(x)*Lj(y)

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.

Notes

New in version 1.7.0.
Parameters
x:array_like, compatible objectsThe two dimensional series is evaluated at the points (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.
y:array_like, compatible objectsThe two dimensional series is evaluated at the points (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:array_likeArray of coefficients ordered so that the coefficient of the term of multi-degree i,j is contained in c[i,j]. If c has dimension greater than two the remaining indices enumerate multiple sets of coefficients.
Returns
ndarray, compatible objectvalues - The values of the two dimensional Legendre series at points formed from pairs of corresponding values from x and y.
def legval3d(x, y, z, c): (source)

Evaluate a 3-D Legendre series at points (x, y, z).

This function returns the values:

p(x, y, z) = i, j, kci, j, k*Li(x)*Lj(y)*Lk(z)

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.

Notes

New in version 1.7.0.
Parameters
x:array_like, compatible objectThe three dimensional series is evaluated at the points (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.
y:array_like, compatible objectThe three dimensional series is evaluated at the points (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.
z:array_like, compatible objectThe three dimensional series is evaluated at the points (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:array_likeArray of coefficients ordered so that the coefficient of the term of multi-degree i,j,k is contained in c[i,j,k]. If c has dimension greater than 3 the remaining indices enumerate multiple sets of coefficients.
Returns
ndarray, compatible objectvalues - The values of the multidimensional polynomial on points formed with triples of corresponding values from x, y, and z.
def legvander(x, deg): (source)

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

V[..., i] = Li(x)

where 0 <= i <= deg. The leading indices of V index the elements of x and the last index is the degree of the Legendre polynomial.

If c is a 1-D array of coefficients of length n + 1 and V is the array V = legvander(x, n), then np.dot(V, c) and legval(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 Legendre series of the same degree and sample points.

Parameters
x:array_likeArray of points. The dtype is converted to float64 or complex128 depending on whether any of the elements are complex. If x is scalar it is converted to a 1-D array.
deg:intDegree of the resulting matrix.
Returns
ndarrayvander - The pseudo-Vandermonde matrix. The shape of the returned matrix is x.shape + (deg + 1,), where The last index is the degree of the corresponding Legendre polynomial. The dtype will be the same as the converted x.
def legvander2d(x, y, deg): (source)

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

V[..., (deg[1] + 1)*i + j] = Li(x)*Lj(y), 

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 Legendre polynomials.

If V = legvander2d(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

c00, c01, c02..., c10, c11, c12...

and np.dot(V, c.flat) and legval2d(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 Legendre series of the same degrees and sample points.

Notes

New in version 1.7.0.
Parameters
x:array_likeArrays of point coordinates, all of the same shape. The dtypes will be converted to either float64 or complex128 depending on whether any of the elements are complex. Scalars are converted to 1-D arrays.
y:array_likeArrays of point coordinates, all of the same shape. The dtypes will be converted to either float64 or complex128 depending on whether any of the elements are complex. Scalars are converted to 1-D arrays.
deg:list of intsList of maximum degrees of the form [x_deg, y_deg].
Returns
ndarrayvander2d - The shape of the returned matrix is x.shape + (order,), where order = (deg[0] + 1)*(deg[1] + 1). The dtype will be the same as the converted x and y.
def legvander3d(x, y, z, deg): (source)

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

V[..., (m + 1)(n + 1)i + (n + 1)j + k] = Li(x)*Lj(y)*Lk(z), 

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 Legendre polynomials.

If V = legvander3d(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

c000, c001, c002, ..., c010, c011, c012, ...

and np.dot(V, c.flat) and legval3d(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 Legendre series of the same degrees and sample points.

Notes

New in version 1.7.0.
Parameters
x:array_likeArrays of point coordinates, all of the same shape. The dtypes will be converted to either float64 or complex128 depending on whether any of the elements are complex. Scalars are converted to 1-D arrays.
y:array_likeArrays of point coordinates, all of the same shape. The dtypes will be converted to either float64 or complex128 depending on whether any of the elements are complex. Scalars are converted to 1-D arrays.
z:array_likeArrays of point coordinates, all of the same shape. The dtypes will be converted to either float64 or complex128 depending on whether any of the elements are complex. Scalars are converted to 1-D arrays.
deg:list of intsList of maximum degrees of the form [x_deg, y_deg, z_deg].
Returns
ndarrayvander3d - The shape of the returned matrix is x.shape + (order,), where order = (deg[0] + 1)*(deg[1] + 1)*(deg[2] + 1). The dtype will be the same as the converted x, y, and z.
def legweight(x): (source)

Weight function of the Legendre polynomials.

The weight function is 1 and the interval of integration is [ − 1, 1]. The Legendre polynomials are orthogonal, but not normalized, with respect to this weight function.

Notes

New in version 1.7.0.
Parameters
x:array_likeValues at which the weight function will be computed.
Returns
ndarrayw - The weight function at x.
def poly2leg(pol): (source)

Convert a polynomial to a Legendre 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 Legendre series, ordered from lowest to highest degree.

See Also

leg2poly

Notes

The easy way to do conversions between polynomial basis sets is to use the convert method of a class instance.

Examples

>>> from numpy import polynomial as P
>>> p = P.Polynomial(np.arange(4))
>>> p
Polynomial([0.,  1.,  2.,  3.], domain=[-1,  1], window=[-1,  1])
>>> c = P.Legendre(P.legendre.poly2leg(p.coef))
>>> c
Legendre([ 1.  ,  3.25,  1.  ,  0.75], domain=[-1,  1], window=[-1,  1]) # may vary
Parameters
pol:array_like1-D array containing the polynomial coefficients
Returns
ndarrayc - 1-D array containing the coefficients of the equivalent Legendre series.
legdomain = (source)

Undocumented

Undocumented

Undocumented

Undocumented