summaryrefslogtreecommitdiff
path: root/numpy/polynomial/polyutils.py
diff options
context:
space:
mode:
authorJarrod Millman <millman@berkeley.edu>2010-02-17 23:53:04 +0000
committerJarrod Millman <millman@berkeley.edu>2010-02-17 23:53:04 +0000
commite2bb09430d90c73a7be6e47ea8c4528f094f693f (patch)
tree3ded297a6cbe634446d6a54afc4e95c8c71553e6 /numpy/polynomial/polyutils.py
parentdcc721a5bddde3afd4ce47d7a7b76ec6c7102b92 (diff)
downloadnumpy-e2bb09430d90c73a7be6e47ea8c4528f094f693f.tar.gz
more docstring updates from pydoc website (thanks to everyone who contributed!)
Diffstat (limited to 'numpy/polynomial/polyutils.py')
-rw-r--r--numpy/polynomial/polyutils.py237
1 files changed, 173 insertions, 64 deletions
diff --git a/numpy/polynomial/polyutils.py b/numpy/polynomial/polyutils.py
index 4d9516a03..5c65e03c2 100644
--- a/numpy/polynomial/polyutils.py
+++ b/numpy/polynomial/polyutils.py
@@ -1,30 +1,34 @@
-"""Utililty functions for polynomial modules.
+"""
+Utililty objects for the polynomial modules.
-This modules provides errors, warnings, and a polynomial base class along
-with some common routines that are used in both the polynomial and
-chebyshev modules.
+This module provides: error and warning objects; a polynomial base class;
+and some routines used in both the `polynomial` and `chebyshev` modules.
-Errors
-------
-- PolyError -- base class for errors
-- PolyDomainError -- mismatched domains
+Error objects
+-------------
+- `PolyError` -- base class for this sub-package's errors.
+- `PolyDomainError` -- raised when domains are "mismatched."
-Warnings
---------
-- RankWarning -- issued by least squares fits to warn of deficient rank
+Warning objects
+---------------
+- `RankWarning` -- raised by a least-squares fit when a rank-deficient
+ matrix is encountered.
-Base Class
+Base class
----------
-- PolyBase -- Base class for the Polynomial and Chebyshev classes.
+- `PolyBase` -- The base class for the `Polynomial` and `Chebyshev`
+ classes.
Functions
---------
-- as_series -- turns list of array_like into 1d arrays of common type
-- trimseq -- removes trailing zeros
-- trimcoef -- removes trailing coefficients less than given magnitude
-- getdomain -- finds appropriate domain for collection of points
-- mapdomain -- maps points between domains
-- mapparms -- parameters of the linear map between domains
+- `as_series` -- turns a list of array_likes into 1-D arrays of common
+ type.
+- `trimseq` -- removes trailing zeros.
+- `trimcoef` -- removes trailing coefficients that are less than a given
+ magnitude (thereby removing the corresponding terms).
+- `getdomain` -- returns a domain appropriate for a given set of abscissae.
+- `mapdomain` -- maps points between domains.
+- `mapparms` -- parameters of the linear map between domains.
"""
from __future__ import division
@@ -109,29 +113,44 @@ def trimseq(seq) :
def as_series(alist, trim=True) :
- """Return arguments as a list of 1d arrays.
+ """
+ Return argument as a list of 1-d arrays.
- The return type will always be an array of double, complex double. or
- object.
+ The returned list contains array(s) of dtype double, complex double, or
+ object. A 1-d argument of shape ``(N,)`` is parsed into ``N`` arrays of
+ size one; a 2-d argument of shape ``(M,N)`` is parsed into ``M`` arrays
+ of size ``N`` (i.e., is "parsed by row"); and a higher dimensional array
+ raises a Value Error if it is not first reshaped into either a 1-d or 2-d
+ array.
Parameters
----------
- [a1, a2,...] : list of array_like.
- The arrays must have no more than one dimension when converted.
- trim : boolean
+ a : array_like
+ A 1- or 2-d array_like
+ trim : boolean, optional
When True, trailing zeros are removed from the inputs.
When False, the inputs are passed through intact.
Returns
-------
[a1, a2,...] : list of 1d-arrays
- A copy of the input data as a 1d-arrays.
+ A copy of the input data as a list of 1-d arrays.
Raises
------
ValueError :
- Raised when an input can not be coverted to 1-d array or the
- resulting array is empty.
+ Raised when `as_series` cannot convert its input to 1-d arrays, or at
+ least one of the resulting arrays is empty.
+
+ Examples
+ --------
+ >>> from numpy import polynomial as P
+ >>> a = np.arange(4)
+ >>> P.as_series(a)
+ [array([ 0.]), array([ 1.]), array([ 2.]), array([ 3.])]
+ >>> b = np.arange(6).reshape((2,3))
+ >>> P.as_series(b)
+ [array([ 0., 1., 2.]), array([ 3., 4., 5.])]
"""
arrays = [np.array(a, ndmin=1, copy=0) for a in alist]
@@ -161,24 +180,47 @@ def as_series(alist, trim=True) :
def trimcoef(c, tol=0) :
- """Remove small trailing coefficients from a polynomial series.
+ """
+ Remove "small" "trailing" coefficients from a polynomial.
+
+ "Small" means "small in absolute value" and is controlled by the
+ parameter `tol`; "trailing" means highest order coefficient(s), e.g., in
+ ``[0, 1, 1, 0, 0]`` (which represents ``0 + x + x**2 + 0*x**3 + 0*x**4``)
+ both the 3-rd and 4-th order coefficients would be "trimmed."
Parameters
----------
c : array_like
- 1-d array of coefficients, ordered from low to high.
- tol : number
- Trailing elements with absolute value less than tol are removed.
+ 1-d array of coefficients, ordered from lowest order to highest.
+ tol : number, optional
+ Trailing (i.e., highest order) elements with absolute value less
+ than or equal to `tol` (default value is zero) are removed.
Returns
-------
trimmed : ndarray
- 1_d array with tailing zeros removed. If the resulting series would
- be empty, a series containing a singel zero is returned.
+ 1-d array with trailing zeros removed. If the resulting series
+ would be empty, a series containing a single zero is returned.
Raises
------
- ValueError : if tol < 0
+ ValueError
+ If `tol` < 0
+
+ See Also
+ --------
+ trimseq
+
+ Examples
+ --------
+ >>> from numpy import polynomial as P
+ >>> P.trimcoef((0,0,3,0,5,0,0))
+ array([ 0., 0., 3., 0., 5.])
+ >>> P.trimcoef((0,0,1e-3,0,1e-5,0,0),1e-3) # item == tol is trimmed
+ array([ 0.])
+ >>> i = complex(0,1) # works for complex
+ >>> P.trimcoef((3e-4,1e-3*(1-i),5e-4,2e-5*(1+i)), 1e-3)
+ array([ 0.0003+0.j , 0.0010-0.001j])
"""
if tol < 0 :
@@ -192,29 +234,42 @@ def trimcoef(c, tol=0) :
return c[:ind[-1] + 1].copy()
def getdomain(x) :
- """Determine suitable domain for given points.
+ """
+ Return a domain suitable for given abscissae.
+
+ Find a domain suitable for a polynomial or Chebyshev series
+ defined at the values supplied.
- Find a suitable domain in which to fit a function defined at the points
- `x` with a polynomial or Chebyshev series.
-
Parameters
----------
x : array_like
- 1D array of points whose domain will be determined.
+ 1-d array of abscissae whose domain will be determined.
Returns
-------
domain : ndarray
- 1D ndarray containing two values. If the inputs are complex, then
- the two points are the corners of the smallest rectangle alinged
- with the axes in the complex plane containing the points `x`. If
- the inputs are real, then the two points are the ends of the
- smallest interval containing the points `x`,
+ 1-d array containing two values. If the inputs are complex, then
+ the two returned points are the lower left and upper right corners
+ of the smallest rectangle (aligned with the axes) in the complex
+ plane containing the points `x`. If the inputs are real, then the
+ two points are the ends of the smallest interval containing the
+ points `x`.
See Also
--------
mapparms, mapdomain
+ Examples
+ --------
+ >>> from numpy.polynomial import polyutils as pu
+ >>> points = np.arange(4)**2 - 5; points
+ array([-5, -4, -1, 4])
+ >>> pu.getdomain(points)
+ array([-5., 4.])
+ >>> c = np.exp(complex(0,1)*np.pi*np.arange(12)/6) # unit circle
+ >>> pu.getdomain(c)
+ array([-1.-1.j, 1.+1.j])
+
"""
[x] = as_series([x], trim=False)
if x.dtype.char in np.typecodes['Complex'] :
@@ -225,27 +280,45 @@ def getdomain(x) :
return np.array((x.min(), x.max()))
def mapparms(old, new) :
- """Linear map between domains.
+ """
+ Linear map parameters between domains.
- Return the parameters of the linear map ``off + scl*x`` that maps the
- `old` domain to the `new` domain. The map is defined by the requirement
- that the left end of the old domain map to the left end of the new
- domain, and similarly for the right ends.
+ Return the parameters of the linear map ``offset + scale*x`` that maps
+ `old` to `new` such that ``old[i] -> new[i]``, ``i = 0, 1``.
Parameters
----------
old, new : array_like
- The two domains should convert as 1D arrays containing two values.
+ Each domain must (successfully) convert to a 1-d array containing
+ precisely two values.
Returns
-------
- off, scl : scalars
- The map `=``off + scl*x`` maps the first domain to the second.
+ offset, scale : scalars
+ The map ``L(x) = offset + scale*x`` maps the first domain to the
+ second.
See Also
--------
getdomain, mapdomain
+ Notes
+ -----
+ Also works for complex numbers, and thus can be used to calculate the
+ parameters required to map any line in the complex plane to any other
+ line therein.
+
+ Examples
+ --------
+ >>> from numpy import polynomial as P
+ >>> P.mapparms((-1,1),(-1,1))
+ (0.0, 1.0)
+ >>> P.mapparms((1,-1),(-1,1))
+ (0.0, -1.0)
+ >>> i = complex(0,1)
+ >>> P.mapparms((-i,-1),(1,i))
+ ((1+1j), (1+0j))
+
"""
oldlen = old[1] - old[0]
newlen = new[1] - new[0]
@@ -254,30 +327,66 @@ def mapparms(old, new) :
return off, scl
def mapdomain(x, old, new) :
- """Apply linear map to input points.
-
- The linear map of the form ``off + scl*x`` that takes the `old` domain
- to the `new` domain is applied to the points `x`.
+ """
+ Apply linear map to input points.
+
+ The linear map ``offset + scale*x`` that maps `old` to `new` is applied
+ to the points `x`.
Parameters
----------
x : array_like
- Points to be mapped
+ Points to be mapped.
old, new : array_like
- The two domains that determin the map. They should both convert as
- 1D arrays containing two values.
-
+ The two domains that determine the map. Each must (successfully)
+ convert to 1-d arrays containing precisely two values.
Returns
-------
- new_x : ndarray
- Array of points of the same shape as the input `x` after the linear
- map defined by the two domains is applied.
+ x_out : ndarray
+ Array of points of the same shape as `x`, after application of the
+ linear map between the two domains.
See Also
--------
getdomain, mapparms
+ Notes
+ -----
+ Effectively, this implements:
+
+ .. math ::
+ x\\_out = new[0] + m(x - old[0])
+
+ where
+
+ .. math ::
+ m = \\frac{new[1]-new[0]}{old[1]-old[0]}
+
+ Examples
+ --------
+ >>> from numpy import polynomial as P
+ >>> old_domain = (-1,1)
+ >>> new_domain = (0,2*np.pi)
+ >>> x = np.linspace(-1,1,6); x
+ array([-1. , -0.6, -0.2, 0.2, 0.6, 1. ])
+ >>> x_out = P.mapdomain(x, old_domain, new_domain); x_out
+ array([ 0. , 1.25663706, 2.51327412, 3.76991118, 5.02654825,
+ 6.28318531])
+ >>> x - P.mapdomain(x_out, new_domain, old_domain)
+ array([ 0., 0., 0., 0., 0., 0.])
+
+ Also works for complex numbers (and thus can be used to map any line in
+ the complex plane to any other line therein).
+
+ >>> i = complex(0,1)
+ >>> old = (-1 - i, 1 + i)
+ >>> new = (-1 + i, 1 - i)
+ >>> z = np.linspace(old[0], old[1], 6); z
+ array([-1.0-1.j , -0.6-0.6j, -0.2-0.2j, 0.2+0.2j, 0.6+0.6j, 1.0+1.j ])
+ >>> new_z = P.mapdomain(z, old, new); new_z
+ array([-1.0+1.j , -0.6+0.6j, -0.2+0.2j, 0.2-0.2j, 0.6-0.6j, 1.0-1.j ])
+
"""
[x] = as_series([x], trim=False)
off, scl = mapparms(old, new)