summaryrefslogtreecommitdiff
path: root/numpy/polynomial/laguerre.py
diff options
context:
space:
mode:
authorCharles Harris <charlesr.harris@gmail.com>2011-12-28 18:43:17 -0700
committerCharles Harris <charlesr.harris@gmail.com>2012-01-09 11:09:36 -0700
commitdc7719f66452288d7c0192f93c07c8b18d870b75 (patch)
tree54d6102e9dab5896fa402afa9e22807647173a59 /numpy/polynomial/laguerre.py
parentc462637f9b398600d25ca449aef8586d8d9d6210 (diff)
downloadnumpy-dc7719f66452288d7c0192f93c07c8b18d870b75.tar.gz
DOC: Finish documenting new functions in the polynomial package.
The old functions could use a review, but that isn't pressing.
Diffstat (limited to 'numpy/polynomial/laguerre.py')
-rw-r--r--numpy/polynomial/laguerre.py212
1 files changed, 140 insertions, 72 deletions
diff --git a/numpy/polynomial/laguerre.py b/numpy/polynomial/laguerre.py
index 489ecb8a2..710b480d9 100644
--- a/numpy/polynomial/laguerre.py
+++ b/numpy/polynomial/laguerre.py
@@ -650,6 +650,8 @@ def lagder(c, m=1, scl=1, axis=0) :
axis : int, optional
Axis over which the derivative is taken. (Default: 0).
+ .. versionadded:: 1.7.0
+
Returns
-------
der : ndarray
@@ -751,6 +753,8 @@ def lagint(c, m=1, k=[], lbnd=0, scl=1, axis=0):
axis : int, optional
Axis over which the derivative is taken. (Default: 0).
+ .. versionadded:: 1.7.0
+
Returns
-------
S : ndarray
@@ -950,8 +954,6 @@ def lagval2d(x, y, 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.
- .. versionadded::1.7.0
-
Parameters
----------
x, y : array_like, compatible objects
@@ -975,6 +977,11 @@ def lagval2d(x, y, c):
--------
lagval, laggrid2d, lagval3d, laggrid3d
+ Notes
+ -----
+
+ .. versionadded::1.7.0
+
"""
try:
x, y = np.array((x, y), copy=0)
@@ -1007,8 +1014,6 @@ def laggrid2d(x, y, c):
its shape to make it 2-D. The shape of the result will be c.shape[2:] +
x.shape + y.shape.
- .. versionadded:: 1.7.0
-
Parameters
----------
x, y : array_like, compatible objects
@@ -1032,6 +1037,11 @@ def laggrid2d(x, y, c):
--------
lagval, lagval2d, lagval3d, laggrid3d
+ Notes
+ -----
+
+ .. versionadded::1.7.0
+
"""
c = lagval(x, c)
c = lagval(y, c)
@@ -1056,8 +1066,6 @@ def lagval3d(x, y, z, c):
shape to make it 3-D. The shape of the result will be c.shape[3:] +
x.shape.
- .. versionadded::1.7.0
-
Parameters
----------
x, y, z : array_like, compatible object
@@ -1082,6 +1090,11 @@ def lagval3d(x, y, z, c):
--------
lagval, lagval2d, laggrid2d, laggrid3d
+ Notes
+ -----
+
+ .. versionadded::1.7.0
+
"""
try:
x, y, z = np.array((x, y, z), copy=0)
@@ -1117,8 +1130,6 @@ def laggrid3d(x, y, z, c):
its shape to make it 3-D. The shape of the result will be c.shape[3:] +
x.shape + yshape + z.shape.
- .. versionadded:: 1.7.0
-
Parameters
----------
x, y, z : array_like, compatible objects
@@ -1143,6 +1154,11 @@ def laggrid3d(x, y, z, c):
--------
lagval, lagval2d, laggrid2d, lagval3d
+ Notes
+ -----
+
+ .. versionadded::1.7.0
+
"""
c = lagval(x, c)
c = lagval(y, c)
@@ -1151,28 +1167,38 @@ def laggrid3d(x, y, z, c):
def lagvander(x, deg) :
- """Vandermonde matrix of given degree.
+ """Pseudo-Vandermonde matrix of given degree.
- Returns the Vandermonde matrix of degree `deg` and sample points `x`.
- This isn't a true Vandermonde matrix because `x` can be an arbitrary
- ndarray and the Laguerre polynomials aren't powers. If ``V`` is the
- returned matrix and `x` is a 2d array, then the elements of ``V`` are
- ``V[i,j,k] = P_k(x[i,j])``, where ``P_k`` is the Laguerre polynomial
- of degree ``k``.
+ Returns the pseudo-Vandermonde matrix of degree `deg` and sample points
+ `x`. The pseudo-Vandermonde matrix is defined by
+
+ .. math:: V[..., i] = L_i(x)
+
+ where `0 <= i <= deg`. The leading indices of `V` index the elements of
+ `x` and the last index is the degree of the Laguerre polynomial.
+
+ If `c` is a 1-D array of coefficients of length `n + 1` and `V` is the
+ array ``V = lagvander(x, n)``, then ``np.dot(V, c)`` and
+ ``lagval(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 Laguerre series of the same degree and sample points.
Parameters
----------
x : array_like
- Array of points. The values are converted to double or complex
- doubles. If x is scalar it is converted to a 1D array.
- deg : integer
+ Array 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 : int
Degree of the resulting matrix.
Returns
-------
- vander : Vandermonde matrix.
- The shape of the returned matrix is ``x.shape + (deg+1,)``. The last
- index is the degree.
+ vander: ndarray
+ 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 Laguerre polynomial. The dtype will be the same as
+ the converted `x`.
Examples
--------
@@ -1201,36 +1227,50 @@ def lagvander(x, deg) :
def lagvander2d(x, y, deg) :
- """Pseudo Vandermonde matrix of given degree.
-
- Returns the pseudo Vandermonde matrix for 2D Laguerre series in `x` and
- `y`. The sample point coordinates must all have the same shape after
- conversion to arrays and the dtype will be converted to either float64
- or complex128 depending on whether any of `x` or 'y' are complex. The
- maximum degrees of the 2D Laguerre series in each variable are specified in
- the list `deg` in the form ``[xdeg, ydeg]``. The return array has the
- shape ``x.shape + (order,)`` if `x`, and `y` are arrays or ``(1, order)
- if they are scalars. Here order is the number of elements in a
- flattened coefficient array of original shape ``(xdeg + 1, ydeg + 1)``.
- The flattening is done so that the resulting pseudo Vandermonde array
- can be easily used in least squares fits.
+ """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
+
+ .. math:: V[..., deg[1]*i + j] = L_i(x) * L_j(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 Laguerre polynomials.
+
+ If `c` is a 2-D array of coefficients of shape `(m + 1, n + 1)` and `V`
+ is the matrix ``V = lagvander2d(x, y, [m, n])``, then
+ ``np.dot(V, c.flat)`` and ``lagval2d(x, y, 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 2-D Laguerre series of the same
+ degrees and sample points.
Parameters
----------
- x,y : array_like
- Arrays of point coordinates, each of the same shape.
- deg : list
+ x, y : array_like
+ Arrays 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 ints
List of maximum degrees of the form [x_deg, y_deg].
Returns
-------
vander2d : ndarray
- The shape of the returned matrix is described above.
+ The shape of the returned matrix is ``x.shape + (order,)``, where
+ :math:`order = (deg[0]+1)*(deg([1]+1)`. The dtype will be the same
+ as the converted `x` and `y`.
See Also
--------
lagvander, lagvander3d. lagval2d, lagval3d
+ Notes
+ -----
+
+ .. versionadded::1.7.0
+
"""
ideg = [int(d) for d in deg]
is_valid = [id == d and id >= 0 for id, d in zip(ideg, deg)]
@@ -1246,37 +1286,51 @@ def lagvander2d(x, y, deg) :
def lagvander3d(x, y, z, deg) :
- """Pseudo Vandermonde matrix of given degree.
-
- Returns the pseudo Vandermonde matrix for 3D Laguerre series in `x`,
- `y`, or `z`. The sample point coordinates must all have the same shape
- after conversion to arrays and the dtype will be converted to either
- float64 or complex128 depending on whether any of `x`, `y`, or 'z' are
- complex. The maximum degrees of the 3D Laguerre series in each
- variable are specified in the list `deg` in the form ``[xdeg, ydeg,
- zdeg]``. The return array has the shape ``x.shape + (order,)`` if `x`,
- `y`, and `z` are arrays or ``(1, order) if they are scalars. Here order
- is the number of elements in a flattened coefficient array of original
- shape ``(xdeg + 1, ydeg + 1, zdeg + 1)``. The flattening is done so
- that the resulting pseudo Vandermonde array can be easily used in least
- squares fits.
+ """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
+
+ .. math:: V[..., (m+1)(n+1)i + (n+1)j + k] = L_i(x)*L_j(y)*L_k(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 Laguerre polynomials.
+
+ If `c` is a 3-D array of coefficients of shape `(l + 1, m + 1, n + 1)`
+ and `V` is the matrix ``V = lagvander3d(x, y, z, [l, m, n])``, then
+ ``np.dot(V, c.flat)`` and ``lagval3d(x, y, z, 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 3-D Laguerre series of the
+ same degrees and sample points.
Parameters
----------
- x,y,z : array_like
- Arrays of point coordinates, each of the same shape.
- deg : list
+ x, y, z : array_like
+ Arrays 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 ints
List of maximum degrees of the form [x_deg, y_deg, z_deg].
Returns
-------
vander3d : ndarray
- The shape of the returned matrix is described above.
+ The shape of the returned matrix is ``x.shape + (order,)``, where
+ :math:`order = (deg[0]+1)*(deg([1]+1)*(deg[2]+1)`. The dtype will
+ be the same as the converted `x`, `y`, and `z`.
See Also
--------
lagvander, lagvander3d. lagval2d, lagval3d
+ Notes
+ -----
+
+ .. versionadded::1.7.0
+
"""
ideg = [int(d) for d in deg]
is_valid = [id == d and id >= 0 for id, d in zip(ideg, deg)]
@@ -1462,11 +1516,12 @@ def lagfit(x, y, deg, rcond=None, full=False, w=None):
def lagcompanion(c):
- """Return the companion matrix of c.
+ """
+ Return the companion matrix of c.
- The unscaled companion matrix of the Laguerre polynomials is already
- symmetric when `c` represents a single Laguerre polynomial, so no
- further scaling is needed.
+ The usual companion matrix of the Laguerre polynomials is already
+ symmetric when `c` is a basis Laguerre polynomial, so no scaling is
+ applied.
Parameters
----------
@@ -1479,6 +1534,11 @@ def lagcompanion(c):
mat : ndarray
Companion matrix of dimensions (deg, deg).
+ Notes
+ -----
+
+ .. versionadded::1.7.0
+
"""
accprod = np.multiply.accumulate
# c is a trimmed copy
@@ -1560,12 +1620,13 @@ def lagroots(c):
def laggauss(deg):
- """Gauss Laguerre quadrature.
+ """
+ Gauss-Laguerre quadrature.
Computes the sample points and weights for Gauss-Laguerre quadrature.
These sample points and weights will correctly integrate polynomials of
- degree ``2*deg - 1`` or less over the interval ``[0, inf]`` with the
- weight function ``f(x) = exp(-x)``.
+ degree :math:`2*deg - 1` or less over the interval :math:`[0, \inf]` with the
+ weight function :math:`f(x) = \exp(-x)`.
Parameters
----------
@@ -1581,14 +1642,17 @@ def laggauss(deg):
Notes
-----
- The results have only been tested up to degree 100. Higher degrees may
+
+ .. versionadded::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
- w = c / (L'_n(x_k) * L_{n-1}(x_k))
+ .. math:: w_k = c / (L'_n(x_k) * L_{n-1}(x_k))
- where ``c`` is a constant independent of ``k`` and ``x_k`` is the k'th
- root of ``L_n``, and then scaling the results to get the right value
- when integrating 1.
+ where :math:`c` is a constant independent of :math:`k` and :math:`x_k`
+ is the k'th root of :math:`L_n`, and then scaling the results to get
+ the right value when integrating 1.
"""
ideg = int(deg)
@@ -1623,10 +1687,9 @@ def laggauss(deg):
def lagweight(x):
"""Weight function of the Laguerre polynomials.
- The weight function for which the Laguerre polynomials are orthogonal.
- In this case the weight function is ``exp(-x)``. Note that the Laguerre
- polynomials are not normalized, indeed, may be much greater than the
- normalized versions.
+ The weight function is :math:`exp(-x)` and the interval of integration
+ is :math:`[0, \inf]`. The Laguerre polynomials are orthogonal, but not
+ normalized, with respect to this weight function.
Parameters
----------
@@ -1638,6 +1701,11 @@ def lagweight(x):
w : ndarray
The weight function at `x`.
+ Notes
+ -----
+
+ .. versionadded::1.7.0
+
"""
w = np.exp(-x)
return w