summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--doc/release/1.9.0-notes.rst14
-rw-r--r--numpy/polynomial/_polybase.py984
-rw-r--r--numpy/polynomial/chebyshev.py59
-rw-r--r--numpy/polynomial/hermite.py59
-rw-r--r--numpy/polynomial/hermite_e.py59
-rw-r--r--numpy/polynomial/laguerre.py59
-rw-r--r--numpy/polynomial/legendre.py59
-rw-r--r--numpy/polynomial/polynomial.py70
-rw-r--r--numpy/polynomial/polytemplate.py6
-rw-r--r--numpy/polynomial/polyutils.py57
10 files changed, 1360 insertions, 66 deletions
diff --git a/doc/release/1.9.0-notes.rst b/doc/release/1.9.0-notes.rst
index aa0431c6e..9a97dba09 100644
--- a/doc/release/1.9.0-notes.rst
+++ b/doc/release/1.9.0-notes.rst
@@ -18,6 +18,7 @@ Dropped Support
Future Changes
==============
+* The numpy/polynomial/polytemplate.py file will be removed in 1.10.0.
Compatibility notes
===================
@@ -70,6 +71,13 @@ The ``doc/swig`` directory moved
~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
The ``doc/swig`` directory has been moved to ``tools/swig``.
+Polynomial Classes no longer derived from PolyBase
+~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
+This may cause problems with folks who depended on the polynomial classes
+being derived from PolyBase. They are now all derived from the abstract
+base class ABCPolyBase. Strictly speaking, there should be a deprecation
+involved, but no external code making use of the old baseclass could be
+found.
New Features
============
@@ -151,6 +159,12 @@ Covariance check in ``np.random.multivariate_normal``
A ``RuntimeWarning`` warning is raised when the covariance matrix is not
positive-semidefinite.
+Polynomial Classes no longer template based
+~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
+The polynomial classes have been refactored to use an abstract base class
+rather than a template in order to implement a common interface. This makes
+importing the polynomial package faster as the classes do not need to be
+compiled on import.
Changes
=======
diff --git a/numpy/polynomial/_polybase.py b/numpy/polynomial/_polybase.py
new file mode 100644
index 000000000..4f174891a
--- /dev/null
+++ b/numpy/polynomial/_polybase.py
@@ -0,0 +1,984 @@
+"""
+Abstract base class for the various polynomial Classes.
+
+The ABCPolyBase class provides the methods needed to implement the common API
+for the various polynomial classes. It operates as a mixin, but uses the
+abc module from the stdlib, hence it is only available for Python >= 2.6.
+
+"""
+from __future__ import division, absolute_import, print_function
+
+from abc import ABCMeta, abstractmethod, abstractproperty
+
+import numpy as np
+from . import polyutils as pu
+
+__all__ = ['ABCPolyBase']
+
+class ABCPolyBase(object):
+ """An abstract base class for series classes.
+
+ ABCPolyBase provides the standard Python numerical methods
+ '+', '-', '*', '//', '%', 'divmod', '**', and '()' along with the
+ methods listed below.
+
+ .. versionadded:: 1.9.0
+
+ Parameters
+ ----------
+ coef : array_like
+ Series coefficients in order of increasing degree, i.e.,
+ ``(1, 2, 3)`` gives ``1*P_0(x) + 2*P_1(x) + 3*P_2(x)``, where
+ ``P_i`` is the basis polynomials of degree ``i``.
+ domain : (2,) array_like, optional
+ Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
+ to the interval ``[window[0], window[1]]`` by shifting and scaling.
+ The default value is the derived class domain.
+ window : (2,) array_like, optional
+ Window, see domain for its use. The default value is the
+ derived class window.
+
+ Attributes
+ ----------
+ coef : (N,) ndarray
+ Series coefficients in order of increasing degree.
+ domain : (2,) ndarray
+ Domain that is mapped to window.
+ window : (2,) ndarray
+ Window that domain is mapped to.
+
+ Class Attributes
+ ----------------
+ maxpower : int
+ Maximum power allowed, i.e., the largest number ``n`` such that
+ ``p(x)**n`` is allowed. This is to limit runaway polynomial size.
+ domain : (2,) ndarray
+ Default domain of the class.
+ window : (2,) ndarray
+ Default window of the class.
+
+ """
+ __metaclass__ = ABCMeta
+
+ # Not hashable
+ __hash__ = None
+
+ # Don't let participate in array operations. Value doesn't matter.
+ __array_priority__ = 1000
+
+ # Limit runaway size. T_n^m has degree n*m
+ maxpower = 100
+
+ @abstractproperty
+ def domain(self):
+ pass
+
+ @abstractproperty
+ def window(self):
+ pass
+
+ @abstractproperty
+ def nickname(self):
+ pass
+
+ @abstractmethod
+ def _add(self):
+ pass
+
+ @abstractmethod
+ def _sub(self):
+ pass
+
+ @abstractmethod
+ def _mul(self):
+ pass
+
+ @abstractmethod
+ def _div(self):
+ pass
+
+ @abstractmethod
+ def _pow(self):
+ pass
+
+ @abstractmethod
+ def _val(self):
+ pass
+
+ @abstractmethod
+ def _int(self):
+ pass
+
+ @abstractmethod
+ def _der(self):
+ pass
+
+ @abstractmethod
+ def _fit(self):
+ pass
+
+ @abstractmethod
+ def _line(self):
+ pass
+
+ @abstractmethod
+ def _roots(self):
+ pass
+
+ @abstractmethod
+ def _fromroots(self):
+ pass
+
+ def has_samecoef(self, other):
+ """Check if coefficients match.
+
+ .. versionadded:: 1.6.0
+
+ Parameters
+ ----------
+ other : class instance
+ The other class must have the ``coef`` attribute.
+
+ Returns
+ -------
+ bool : boolean
+ True if the coefficients are the same, False otherwise.
+
+ """
+ if len(self.coef) != len(other.coef):
+ return False
+ elif not np.all(self.coef == other.coef):
+ return False
+ else:
+ return True
+
+ def has_samedomain(self, other):
+ """Check if domains match.
+
+ .. versionadded:: 1.6.0
+
+ Parameters
+ ----------
+ other : class instance
+ The other class must have the ``domain`` attribute.
+
+ Returns
+ -------
+ bool : boolean
+ True if the domains are the same, False otherwise.
+
+ """
+ return np.all(self.domain == other.domain)
+
+ def has_samewindow(self, other):
+ """Check if windows match.
+
+ .. versionadded:: 1.6.0
+
+ Parameters
+ ----------
+ other : class instance
+ The other class must have the ``window`` attribute.
+
+ Returns
+ -------
+ bool : boolean
+ True if the windows are the same, False otherwise.
+
+ """
+ return np.all(self.window == other.window)
+
+ def has_sametype(self, other):
+ """Check if types match.
+
+ .. versionadded:: 1.7.0
+
+ Parameters
+ ----------
+ other : object
+ Class instance.
+
+ Returns
+ -------
+ bool : boolean
+ True if other is same class as self
+
+ """
+ return isinstance(other, self.__class__)
+
+ def __init__(self, coef, domain=None, window=None):
+ [coef] = pu.as_series([coef], trim=False)
+ self.coef = coef
+
+ if domain is not None:
+ [domain] = pu.as_series([domain], trim=False)
+ if len(domain) != 2:
+ raise ValueError("Domain has wrong number of elements.")
+ self.domain = domain
+
+ if window is not None:
+ [window] = pu.as_series([window], trim=False)
+ if len(window) != 2:
+ raise ValueError("Window has wrong number of elements.")
+ self.window = window
+
+ def __repr__(self):
+ format = "%s(%s, %s, %s)"
+ coef = repr(self.coef)[6:-1]
+ domain = repr(self.domain)[6:-1]
+ window = repr(self.window)[6:-1]
+ name = self.__class__.__name__
+ return format % (name, coef, domain, window)
+
+ def __str__(self):
+ format = "%s(%s)"
+ coef = str(self.coef)
+ name = self.nickname
+ return format % (name, coef)
+
+ # Pickle and copy
+
+ def __getstate__(self):
+ ret = self.__dict__.copy()
+ ret['coef'] = self.coef.copy()
+ ret['domain'] = self.domain.copy()
+ ret['window'] = self.window.copy()
+ return ret
+
+ def __setstate__(self, dict):
+ self.__dict__ = dict
+
+ # Call
+
+ def __call__(self, arg):
+ off, scl = pu.mapparms(self.domain, self.window)
+ arg = off + scl*arg
+ return self._val(arg, self.coef)
+
+ def __iter__(self):
+ return iter(self.coef)
+
+ def __len__(self):
+ return len(self.coef)
+
+ # Numeric properties.
+
+ def __neg__(self):
+ return self.__class__(-self.coef, self.domain, self.window)
+
+ def __pos__(self):
+ return self
+
+ def __add__(self, other):
+ if isinstance(other, ABCPolyBase):
+ if not self.has_sametype(other):
+ raise TypeError("Polynomial types differ")
+ elif not self.has_samedomain(other):
+ raise TypeError("Domains differ")
+ elif not self.has_samewindow(other):
+ raise TypeError("Windows differ")
+ else:
+ coef = self._add(self.coef, other.coef)
+ else:
+ try:
+ coef = self._add(self.coef, other)
+ except:
+ return NotImplemented
+ return self.__class__(coef, self.domain, self.window)
+
+ def __sub__(self, other):
+ if isinstance(other, ABCPolyBase):
+ if not self.has_sametype(other):
+ raise TypeError("Polynomial types differ")
+ elif not self.has_samedomain(other):
+ raise TypeError("Domains differ")
+ elif not self.has_samewindow(other):
+ raise TypeError("Windows differ")
+ else:
+ coef = self._sub(self.coef, other.coef)
+ else:
+ try:
+ coef = self._sub(self.coef, other)
+ except:
+ return NotImplemented
+ return self.__class__(coef, self.domain, self.window)
+
+ def __mul__(self, other):
+ if isinstance(other, ABCPolyBase):
+ if not self.has_sametype(other):
+ raise TypeError("Polynomial types differ")
+ elif not self.has_samedomain(other):
+ raise TypeError("Domains differ")
+ elif not self.has_samewindow(other):
+ raise TypeError("Windows differ")
+ else:
+ coef = self._mul(self.coef, other.coef)
+ else:
+ try:
+ coef = self._mul(self.coef, other)
+ except:
+ return NotImplemented
+ return self.__class__(coef, self.domain, self.window)
+
+ def __div__(self, other):
+ # set to __floordiv__, /, for now.
+ return self.__floordiv__(other)
+
+ def __truediv__(self, other):
+ # there is no true divide if the rhs is not a scalar, although it
+ # could return the first n elements of an infinite series.
+ # It is hard to see where n would come from, though.
+ if np.isscalar(other):
+ # this might be overly restrictive
+ coef = self.coef/other
+ return self.__class__(coef, self.domain, self.window)
+ else:
+ return NotImplemented
+
+ def __floordiv__(self, other):
+ if isinstance(other, ABCPolyBase):
+ if not self.has_sametype(other):
+ raise TypeError("Polynomial types differ")
+ elif not self.has_samedomain(other):
+ raise TypeError("Domains differ")
+ elif not self.has_samewindow(other):
+ raise TypeError("Windows differ")
+ else:
+ quo, rem = self._div(self.coef, other.coef)
+ else:
+ try:
+ quo, rem = self._div(self.coef, other)
+ except:
+ return NotImplemented
+ return self.__class__(quo, self.domain, self.window)
+
+ def __mod__(self, other):
+ if isinstance(other, ABCPolyBase):
+ if not self.has_sametype(other):
+ raise TypeError("Polynomial types differ")
+ elif not self.has_samedomain(other):
+ raise TypeError("Domains differ")
+ elif not self.has_samewindow(other):
+ raise TypeError("Windows differ")
+ else:
+ quo, rem = self._div(self.coef, other.coef)
+ else:
+ try:
+ quo, rem = self._div(self.coef, other)
+ except:
+ return NotImplemented
+ return self.__class__(rem, self.domain, self.window)
+
+ def __divmod__(self, other):
+ if isinstance(other, self.__class__):
+ if not self.has_samedomain(other):
+ raise TypeError("Domains are not equal")
+ elif not self.has_samewindow(other):
+ raise TypeError("Windows are not equal")
+ else:
+ quo, rem = self._div(self.coef, other.coef)
+ else:
+ try:
+ quo, rem = self._div(self.coef, other)
+ except:
+ return NotImplemented
+ quo = self.__class__(quo, self.domain, self.window)
+ rem = self.__class__(rem, self.domain, self.window)
+ return quo, rem
+
+ def __pow__(self, other):
+ try:
+ coef = self._pow(self.coef, other, maxpower = self.maxpower)
+ except:
+ raise
+ return self.__class__(coef, self.domain, self.window)
+
+ def __radd__(self, other):
+ try:
+ coef = self._add(other, self.coef)
+ except:
+ return NotImplemented
+ return self.__class__(coef, self.domain, self.window)
+
+ def __rsub__(self, other):
+ try:
+ coef = self._sub(other, self.coef)
+ except:
+ return NotImplemented
+ return self.__class__(coef, self.domain, self.window)
+
+ def __rmul__(self, other):
+ try:
+ coef = self._mul(other, self.coef)
+ except:
+ return NotImplemented
+ return self.__class__(coef, self.domain, self.window)
+
+ def __rdiv__(self, other):
+ # set to __floordiv__ /.
+ return self.__rfloordiv__(other)
+
+ def __rtruediv__(self, other):
+ # there is no true divide if the rhs is not a scalar, although it
+ # could return the first n elements of an infinite series.
+ # It is hard to see where n would come from, though.
+ if len(self.coef) == 1:
+ try:
+ quo, rem = self._div(other, self.coef[0])
+ except:
+ return NotImplemented
+ return self.__class__(quo, self.domain, self.window)
+
+ def __rfloordiv__(self, other):
+ try:
+ quo, rem = self._div(other, self.coef)
+ except:
+ return NotImplemented
+ return self.__class__(quo, self.domain, self.window)
+
+ def __rmod__(self, other):
+ try:
+ quo, rem = self._div(other, self.coef)
+ except:
+ return NotImplemented
+ return self.__class__(rem, self.domain, self.window)
+
+ def __rdivmod__(self, other):
+ try:
+ quo, rem = self._div(other, self.coef)
+ except:
+ return NotImplemented
+ quo = self.__class__(quo, self.domain, self.window)
+ rem = self.__class__(rem, self.domain, self.window)
+ return quo, rem
+
+ # Enhance me
+ # some augmented arithmetic operations could be added here
+
+ def __eq__(self, other):
+ res = isinstance(other, self.__class__) and\
+ self.has_samecoef(other) and \
+ self.has_samedomain(other) and\
+ self.has_samewindow(other)
+ return res
+
+ def __ne__(self, other):
+ return not self.__eq__(other)
+
+ #
+ # Extra methods.
+ #
+
+ def copy(self):
+ """Return a copy.
+
+ Returns
+ -------
+ new_series : series
+ Copy of self.
+
+ """
+ return self.__class__(self.coef, self.domain, self.window)
+
+ def degree(self):
+ """The degree of the series.
+
+ .. versionadded:: 1.5.0
+
+ Returns
+ -------
+ degree : int
+ Degree of the series, one less than the number of coefficients.
+
+ """
+ return len(self) - 1
+
+ def cutdeg(self, deg):
+ """Truncate series to the given degree.
+
+ Reduce the degree of the series to `deg` by discarding the
+ high order terms. If `deg` is greater than the current degree a
+ copy of the current series is returned. This can be useful in least
+ squares where the coefficients of the high degree terms may be very
+ small.
+
+ .. versionadded:: 1.5.0
+
+ Parameters
+ ----------
+ deg : non-negative int
+ The series is reduced to degree `deg` by discarding the high
+ order terms. The value of `deg` must be a non-negative integer.
+
+ Returns
+ -------
+ new_series : series
+ New instance of series with reduced degree.
+
+ """
+ return self.truncate(deg + 1)
+
+ def trim(self, tol=0):
+ """Remove trailing coefficients
+
+ Remove trailing coefficients until a coefficient is reached whose
+ absolute value greater than `tol` or the beginning of the series is
+ reached. If all the coefficients would be removed the series is set
+ to ``[0]``. A new series instance is returned with the new
+ coefficients. The current instance remains unchanged.
+
+ Parameters
+ ----------
+ tol : non-negative number.
+ All trailing coefficients less than `tol` will be removed.
+
+ Returns
+ -------
+ new_series : series
+ Contains the new set of coefficients.
+
+ """
+ coef = pu.trimcoef(self.coef, tol)
+ return self.__class__(coef, self.domain, self.window)
+
+ def truncate(self, size):
+ """Truncate series to length `size`.
+
+ Reduce the series to length `size` by discarding the high
+ degree terms. The value of `size` must be a positive integer. This
+ can be useful in least squares where the coefficients of the
+ high degree terms may be very small.
+
+ Parameters
+ ----------
+ size : positive int
+ The series is reduced to length `size` by discarding the high
+ degree terms. The value of `size` must be a positive integer.
+
+ Returns
+ -------
+ new_series : series
+ New instance of series with truncated coefficients.
+
+ """
+ isize = int(size)
+ if isize != size or isize < 1:
+ raise ValueError("size must be a positive integer")
+ if isize >= len(self.coef):
+ coef = self.coef
+ else:
+ coef = self.coef[:isize]
+ return self.__class__(coef, self.domain, self.window)
+
+ def convert(self, domain=None, kind=None, window=None):
+ """Convert series to a different kind and/or domain and/or window.
+
+ Parameters
+ ----------
+ domain : array_like, optional
+ The domain of the converted series. If the value is None,
+ the default domain of `kind` is used.
+ kind : class, optional
+ The polynomial series type class to which the current instance
+ should be converted. If kind is None, then the class of the
+ current instance is used.
+ window : array_like, optional
+ The window of the converted series. If the value is None,
+ the default window of `kind` is used.
+
+ Returns
+ -------
+ new_series : series
+ The returned class can be of different type than the current
+ instance and/or have a different domain and/or different
+ window.
+
+ Notes
+ -----
+ Conversion between domains and class types can result in
+ numerically ill defined series.
+
+ Examples
+ --------
+
+ """
+ if kind is None:
+ kind = self.__class__
+ if domain is None:
+ domain = kind.domain
+ if window is None:
+ window = kind.window
+ return self(kind.identity(domain, window=window))
+
+ def mapparms(self):
+ """Return the mapping parameters.
+
+ The returned values define a linear map ``off + scl*x`` that is
+ applied to the input arguments before the series is evaluated. The
+ map depends on the ``domain`` and ``window``; if the current
+ ``domain`` is equal to the ``window`` the resulting map is the
+ identity. If the coefficients of the series instance are to be
+ used by themselves outside this class, then the linear function
+ must be substituted for the ``x`` in the standard representation of
+ the base polynomials.
+
+ Returns
+ -------
+ off, scl : float or complex
+ The mapping function is defined by ``off + scl*x``.
+
+ Notes
+ -----
+ If the current domain is the interval ``[l1, r1]`` and the window
+ is ``[l2, r2]``, then the linear mapping function ``L`` is
+ defined by the equations::
+
+ L(l1) = l2
+ L(r1) = r2
+
+ """
+ return pu.mapparms(self.domain, self.window)
+
+ def integ(self, m=1, k=[], lbnd=None):
+ """Integrate.
+
+ Return a series instance that is the definite integral of the
+ current series.
+
+ Parameters
+ ----------
+ m : non-negative int
+ The number of integrations to perform.
+ k : array_like
+ Integration constants. The first constant is applied to the
+ first integration, the second to the second, and so on. The
+ list of values must less than or equal to `m` in length and any
+ missing values are set to zero.
+ lbnd : Scalar
+ The lower bound of the definite integral.
+
+ Returns
+ -------
+ new_series : series
+ A new series representing the integral. The domain is the same
+ as the domain of the integrated series.
+
+ """
+ off, scl = self.mapparms()
+ if lbnd is None:
+ lbnd = 0
+ else:
+ lbnd = off + scl*lbnd
+ coef = self._int(self.coef, m, k, lbnd, 1./scl)
+ return self.__class__(coef, self.domain, self.window)
+
+ def deriv(self, m=1):
+ """Differentiate.
+
+ Return a series instance of that is the derivative of the current
+ series.
+
+ Parameters
+ ----------
+ m : non-negative int
+ The number of integrations to perform.
+
+ Returns
+ -------
+ new_series : series
+ A new series representing the derivative. The domain is the same
+ as the domain of the differentiated series.
+
+ """
+ off, scl = self.mapparms()
+ coef = self._der(self.coef, m, scl)
+ return self.__class__(coef, self.domain, self.window)
+
+ def roots(self):
+ """Return the roots of the series polynomial.
+
+ Compute the roots for the series. Note that the accuracy of the
+ roots decrease the further outside the domain they lie.
+
+ Returns
+ -------
+ roots : ndarray
+ Array containing the roots of the series.
+
+ """
+ roots = self._roots(self.coef)
+ return pu.mapdomain(roots, self.window, self.domain)
+
+ def linspace(self, n=100, domain=None):
+ """Return x, y values at equally spaced points in domain.
+
+ Returns the x, y values at `n` linearly spaced points across the
+ domain. Here y is the value of the polynomial at the points x. By
+ default the domain is the same as that of the series instance.
+ This method is intended mostly as a plotting aid.
+
+ .. versionadded:: 1.5.0
+
+ Parameters
+ ----------
+ n : int, optional
+ Number of point pairs to return. The default value is 100.
+ domain : {None, array_like}, optional
+ If not None, the specified domain is used instead of that of
+ the calling instance. It should be of the form ``[beg,end]``.
+ The default is None which case the class domain is used.
+
+ Returns
+ -------
+ x, y : ndarray
+ x is equal to linspace(self.domain[0], self.domain[1], n) and
+ y is the series evaluated at element of x.
+
+ """
+ if domain is None:
+ domain = self.domain
+ x = np.linspace(domain[0], domain[1], n)
+ y = self(x)
+ return x, y
+
+
+
+ @classmethod
+ def fit(cls, x, y, deg, domain=None, rcond=None, full=False, w=None,
+ window=None):
+ """Least squares fit to data.
+
+ Return a series instance that is the least squares fit to the data
+ `y` sampled at `x`. The domain of the returned instance can be
+ specified and this will often result in a superior fit with less
+ chance of ill conditioning.
+
+ 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
+ Degree of the fitting polynomial.
+ domain : {None, [beg, end], []}, optional
+ Domain to use for the returned series. If ``None``,
+ then a minimal domain that covers the points `x` is chosen. If
+ ``[]`` the class domain is used. The default value was the
+ class domain in NumPy 1.4 and ``None`` in later versions.
+ The ``[]`` option was added in numpy 1.5.0.
+ rcond : float, optional
+ Relative 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, optional
+ Switch 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 contribution of each point
+ ``(x[i],y[i])`` to the fit is weighted by `w[i]`. Ideally the
+ weights are chosen so that the errors of the products
+ ``w[i]*y[i]`` all have the same variance. The default value is
+ None.
+
+ .. versionadded:: 1.5.0
+ window : {[beg, end]}, optional
+ Window to use for the returned series. The default
+ value is the default class domain
+
+ .. versionadded:: 1.6.0
+
+ Returns
+ -------
+ new_series : series
+ A series that represents the least squares fit to the data and
+ has the domain specified in the call.
+
+ [resid, rank, sv, rcond] : list
+ These values are only returned if `full` = True
+
+ resid -- sum of squared residuals of the least squares fit
+ rank -- the numerical rank of the scaled Vandermonde matrix
+ sv -- singular values of the scaled Vandermonde matrix
+ rcond -- value of `rcond`.
+
+ For more details, see `linalg.lstsq`.
+
+ """
+ if domain is None:
+ domain = pu.getdomain(x)
+ elif domain == []:
+ domain = cls.domain
+
+ if window is None:
+ window = cls.window
+
+ xnew = pu.mapdomain(x, domain, window)
+ res = cls._fit(xnew, y, deg, w=w, rcond=rcond, full=full)
+ if full:
+ [coef, status] = res
+ return cls(coef, domain=domain, window=window), status
+ else:
+ coef = res
+ return cls(coef, domain=domain, window=window)
+
+ @classmethod
+ def fromroots(cls, roots, domain=[], window=None):
+ """Return series instance that has the specified roots.
+
+ Returns a series representing the product
+ ``(x - r[0])*(x - r[1])*...*(x - r[n-1])``, where ``r`` is a
+ list of roots.
+
+ Parameters
+ ----------
+ roots : array_like
+ List of roots.
+ domain : {[], None, array_like}, optional
+ Domain for the resulting series. If None the domain is the
+ interval from the smallest root to the largest. If [] the
+ domain is the class domain. The default is [].
+ window : {None, array_like}, optional
+ Window for the returned series. If None the class window is
+ used. The default is None.
+
+ Returns
+ -------
+ new_series : series
+ Series with the specified roots.
+
+ """
+ [roots] = pu.as_series([roots], trim=False)
+ if domain is None:
+ domain = pu.getdomain(roots)
+ elif domain == []:
+ domain = cls.domain
+
+ if window is None:
+ window = cls.window
+
+ deg = len(roots)
+ off, scl = pu.mapparms(domain, window)
+ rnew = off + scl*roots
+ coef = cls._fromroots(rnew) / scl**deg
+ return cls(coef, domain=domain, window=window)
+
+ @classmethod
+ def identity(cls, domain=None, window=None):
+ """Identity function.
+
+ If ``p`` is the returned series, then ``p(x) == x`` for all
+ values of x.
+
+ Parameters
+ ----------
+ domain : {None, array_like}, optional
+ If given, the array must be of the form ``[beg, end]``, where
+ ``beg`` and ``end`` are the endpoints of the domain. If None is
+ given then the class domain is used. The default is None.
+ window : {None, array_like}, optional
+ If given, the resulting array must be if the form
+ ``[beg, end]``, where ``beg`` and ``end`` are the endpoints of
+ the window. If None is given then the class window is used. The
+ default is None.
+
+ Returns
+ -------
+ new_series : series
+ Series of representing the identity.
+
+ """
+ if domain is None:
+ domain = cls.domain
+ if window is None:
+ window = cls.window
+ off, scl = pu.mapparms(window, domain)
+ coef = cls._line(off, scl)
+ return cls(coef, domain, window)
+
+ @classmethod
+ def basis(cls, deg, domain=None, window=None):
+ """Series basis polynomial of degree `deg`.
+
+ Returns the series representing the basis polynomial of degree `deg`.
+
+ .. versionadded:: 1.7.0
+
+ Parameters
+ ----------
+ deg : int
+ Degree of the basis polynomial for the series. Must be >= 0.
+ domain : {None, array_like}, optional
+ If given, the array must be of the form ``[beg, end]``, where
+ ``beg`` and ``end`` are the endpoints of the domain. If None is
+ given then the class domain is used. The default is None.
+ window : {None, array_like}, optional
+ If given, the resulting array must be if the form
+ ``[beg, end]``, where ``beg`` and ``end`` are the endpoints of
+ the window. If None is given then the class window is used. The
+ default is None.
+
+ Returns
+ -------
+ new_series : series
+ A series with the coefficient of the `deg` term set to one and
+ all others zero.
+
+ """
+ if domain is None:
+ domain = cls.domain
+ if window is None:
+ window = cls.window
+ ideg = int(deg)
+
+ if ideg != deg or ideg < 0:
+ raise ValueError("deg must be non-negative integer")
+ return cls([0]*ideg + [1], domain, window)
+
+ @classmethod
+ def cast(cls, series, domain=None, window=None):
+ """Convert series to series of this class.
+
+ The `series` is expected to be an instance of some polynomial
+ series of one of the types supported by by the numpy.polynomial
+ module, but could be some other class that supports the convert
+ method.
+
+ .. versionadded:: 1.7.0
+
+ Parameters
+ ----------
+ series : series
+ The series instance to be converted.
+ domain : {None, array_like}, optional
+ If given, the array must be of the form ``[beg, end]``, where
+ ``beg`` and ``end`` are the endpoints of the domain. If None is
+ given then the class domain is used. The default is None.
+ window : {None, array_like}, optional
+ If given, the resulting array must be if the form
+ ``[beg, end]``, where ``beg`` and ``end`` are the endpoints of
+ the window. If None is given then the class window is used. The
+ default is None.
+
+ Returns
+ -------
+ new_series : series
+ A series of the same kind as the calling class and equal to
+ `series` when evaluated.
+
+ See Also
+ --------
+ convert : similar instance method
+
+ """
+ if domain is None:
+ domain = cls.domain
+ if window is None:
+ window = cls.window
+ return series.convert(domain, cls, window)
diff --git a/numpy/polynomial/chebyshev.py b/numpy/polynomial/chebyshev.py
index 6a2394382..ec3fde2f8 100644
--- a/numpy/polynomial/chebyshev.py
+++ b/numpy/polynomial/chebyshev.py
@@ -87,11 +87,12 @@ References
"""
from __future__ import division, absolute_import, print_function
+import warnings
import numpy as np
import numpy.linalg as la
+
from . import polyutils as pu
-import warnings
-from .polytemplate import polytemplate
+from ._polybase import ABCPolyBase
__all__ = ['chebzero', 'chebone', 'chebx', 'chebdomain', 'chebline',
'chebadd', 'chebsub', 'chebmulx', 'chebmul', 'chebdiv', 'chebpow',
@@ -1646,10 +1647,15 @@ def chebfit(x, y, deg, rcond=None, full=False, w=None):
the coefficients for the data in column k of `y` are in column
`k`.
- [residuals, rank, singular_values, rcond] : present when `full` = True
- Residuals of the least-squares fit, the effective rank of the
- scaled Vandermonde matrix and its singular values, and the
- specified value of `rcond`. For more details, see `linalg.lstsq`.
+ [residuals, rank, singular_values, rcond] : list
+ These values are only returned if `full` = True
+
+ resid -- sum of squared residuals of the least squares fit
+ rank -- the numerical rank of the scaled Vandermonde matrix
+ sv -- singular values of the scaled Vandermonde matrix
+ rcond -- value of `rcond`.
+
+ For more details, see `linalg.lstsq`.
Warns
-----
@@ -2012,4 +2018,43 @@ def chebpts2(npts):
# Chebyshev series class
#
-exec(polytemplate.substitute(name='Chebyshev', nick='cheb', domain='[-1,1]'))
+class Chebyshev(ABCPolyBase):
+ """A Chebyshev series class.
+
+ The Chebyshev class provides the standard Python numerical methods
+ '+', '-', '*', '//', '%', 'divmod', '**', and '()' as well as the
+ methods listed below.
+
+ Parameters
+ ----------
+ coef : array_like
+ Chebyshev coefficients in order of increasing degree, i.e.,
+ ``(1, 2, 3)`` gives ``1*T_0(x) + 2*T_1(x) + 3*T_2(x)``.
+ domain : (2,) array_like, optional
+ Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
+ to the interval ``[window[0], window[1]]`` by shifting and scaling.
+ The default value is [-1, 1].
+ window : (2,) array_like, optional
+ Window, see `domain` for its use. The default value is [-1, 1].
+
+ .. versionadded:: 1.6.0
+
+ """
+ # Virtual Functions
+ _add = staticmethod(chebadd)
+ _sub = staticmethod(chebsub)
+ _mul = staticmethod(chebmul)
+ _div = staticmethod(chebdiv)
+ _pow = staticmethod(chebpow)
+ _val = staticmethod(chebval)
+ _int = staticmethod(chebint)
+ _der = staticmethod(chebder)
+ _fit = staticmethod(chebfit)
+ _line = staticmethod(chebline)
+ _roots = staticmethod(chebroots)
+ _fromroots = staticmethod(chebfromroots)
+
+ # Virtual properties
+ nickname = 'cheb'
+ domain = np.array(chebdomain)
+ window = np.array(chebdomain)
diff --git a/numpy/polynomial/hermite.py b/numpy/polynomial/hermite.py
index 4140acfb7..7d2aa38a2 100644
--- a/numpy/polynomial/hermite.py
+++ b/numpy/polynomial/hermite.py
@@ -59,11 +59,12 @@ See also
"""
from __future__ import division, absolute_import, print_function
+import warnings
import numpy as np
import numpy.linalg as la
+
from . import polyutils as pu
-import warnings
-from .polytemplate import polytemplate
+from ._polybase import ABCPolyBase
__all__ = ['hermzero', 'hermone', 'hermx', 'hermdomain', 'hermline',
'hermadd', 'hermsub', 'hermmulx', 'hermmul', 'hermdiv', 'hermpow',
@@ -1416,10 +1417,15 @@ def hermfit(x, y, deg, rcond=None, full=False, w=None):
the coefficients for the data in column k of `y` are in column
`k`.
- [residuals, rank, singular_values, rcond] : present when `full` = True
- Residuals of the least-squares fit, the effective rank of the
- scaled Vandermonde matrix and its singular values, and the
- specified value of `rcond`. For more details, see `linalg.lstsq`.
+ [residuals, rank, singular_values, rcond] : list
+ These values are only returned if `full` = True
+
+ resid -- sum of squared residuals of the least squares fit
+ rank -- the numerical rank of the scaled Vandermonde matrix
+ sv -- singular values of the scaled Vandermonde matrix
+ rcond -- value of `rcond`.
+
+ For more details, see `linalg.lstsq`.
Warns
-----
@@ -1747,4 +1753,43 @@ def hermweight(x):
# Hermite series class
#
-exec(polytemplate.substitute(name='Hermite', nick='herm', domain='[-1,1]'))
+class Hermite(ABCPolyBase):
+ """An Hermite series class.
+
+ The Hermite class provides the standard Python numerical methods
+ '+', '-', '*', '//', '%', 'divmod', '**', and '()' as well as the
+ attributes and methods listed in the `ABCPolyBase` documentation.
+
+ Parameters
+ ----------
+ coef : array_like
+ Laguerre coefficients in order of increasing degree, i.e,
+ ``(1, 2, 3)`` gives ``1*H_0(x) + 2*H_1(X) + 3*H_2(x)``.
+ domain : (2,) array_like, optional
+ Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
+ to the interval ``[window[0], window[1]]`` by shifting and scaling.
+ The default value is [-1, 1].
+ window : (2,) array_like, optional
+ Window, see `domain` for its use. The default value is [-1, 1].
+
+ .. versionadded:: 1.6.0
+
+ """
+ # Virtual Functions
+ _add = staticmethod(hermadd)
+ _sub = staticmethod(hermsub)
+ _mul = staticmethod(hermmul)
+ _div = staticmethod(hermdiv)
+ _pow = staticmethod(hermpow)
+ _val = staticmethod(hermval)
+ _int = staticmethod(hermint)
+ _der = staticmethod(hermder)
+ _fit = staticmethod(hermfit)
+ _line = staticmethod(hermline)
+ _roots = staticmethod(hermroots)
+ _fromroots = staticmethod(hermfromroots)
+
+ # Virtual properties
+ nickname = 'herm'
+ domain = np.array(hermdomain)
+ window = np.array(hermdomain)
diff --git a/numpy/polynomial/hermite_e.py b/numpy/polynomial/hermite_e.py
index 735ca9470..c3a789809 100644
--- a/numpy/polynomial/hermite_e.py
+++ b/numpy/polynomial/hermite_e.py
@@ -59,11 +59,12 @@ See also
"""
from __future__ import division, absolute_import, print_function
+import warnings
import numpy as np
import numpy.linalg as la
+
from . import polyutils as pu
-import warnings
-from .polytemplate import polytemplate
+from ._polybase import ABCPolyBase
__all__ = ['hermezero', 'hermeone', 'hermex', 'hermedomain', 'hermeline',
'hermeadd', 'hermesub', 'hermemulx', 'hermemul', 'hermediv', 'hermpow',
@@ -1412,10 +1413,15 @@ def hermefit(x, y, deg, rcond=None, full=False, w=None):
the coefficients for the data in column k of `y` are in column
`k`.
- [residuals, rank, singular_values, rcond] : present when `full` = True
- Residuals of the least-squares fit, the effective rank of the
- scaled Vandermonde matrix and its singular values, and the
- specified value of `rcond`. For more details, see `linalg.lstsq`.
+ [residuals, rank, singular_values, rcond] : list
+ These values are only returned if `full` = True
+
+ resid -- sum of squared residuals of the least squares fit
+ rank -- the numerical rank of the scaled Vandermonde matrix
+ sv -- singular values of the scaled Vandermonde matrix
+ rcond -- value of `rcond`.
+
+ For more details, see `linalg.lstsq`.
Warns
-----
@@ -1743,4 +1749,43 @@ def hermeweight(x):
# HermiteE series class
#
-exec(polytemplate.substitute(name='HermiteE', nick='herme', domain='[-1,1]'))
+class HermiteE(ABCPolyBase):
+ """An HermiteE series class.
+
+ The HermiteE class provides the standard Python numerical methods
+ '+', '-', '*', '//', '%', 'divmod', '**', and '()' as well as the
+ attributes and methods listed in the `ABCPolyBase` documentation.
+
+ Parameters
+ ----------
+ coef : array_like
+ Laguerre coefficients in order of increasing degree, i.e,
+ ``(1, 2, 3)`` gives ``1*He_0(x) + 2*He_1(X) + 3*He_2(x)``.
+ domain : (2,) array_like, optional
+ Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
+ to the interval ``[window[0], window[1]]`` by shifting and scaling.
+ The default value is [-1, 1].
+ window : (2,) array_like, optional
+ Window, see `domain` for its use. The default value is [-1, 1].
+
+ .. versionadded:: 1.6.0
+
+ """
+ # Virtual Functions
+ _add = staticmethod(hermeadd)
+ _sub = staticmethod(hermesub)
+ _mul = staticmethod(hermemul)
+ _div = staticmethod(hermediv)
+ _pow = staticmethod(hermepow)
+ _val = staticmethod(hermeval)
+ _int = staticmethod(hermeint)
+ _der = staticmethod(hermeder)
+ _fit = staticmethod(hermefit)
+ _line = staticmethod(hermeline)
+ _roots = staticmethod(hermeroots)
+ _fromroots = staticmethod(hermefromroots)
+
+ # Virtual properties
+ nickname = 'herme'
+ domain = np.array(hermedomain)
+ window = np.array(hermedomain)
diff --git a/numpy/polynomial/laguerre.py b/numpy/polynomial/laguerre.py
index b7ffe9b0c..bf44dc5f4 100644
--- a/numpy/polynomial/laguerre.py
+++ b/numpy/polynomial/laguerre.py
@@ -59,11 +59,12 @@ See also
"""
from __future__ import division, absolute_import, print_function
+import warnings
import numpy as np
import numpy.linalg as la
+
from . import polyutils as pu
-import warnings
-from .polytemplate import polytemplate
+from ._polybase import ABCPolyBase
__all__ = ['lagzero', 'lagone', 'lagx', 'lagdomain', 'lagline',
'lagadd', 'lagsub', 'lagmulx', 'lagmul', 'lagdiv', 'lagpow',
@@ -1415,10 +1416,15 @@ def lagfit(x, y, deg, rcond=None, full=False, w=None):
the coefficients for the data in column k of `y` are in column
`k`.
- [residuals, rank, singular_values, rcond] : present when `full` = True
- Residuals of the least-squares fit, the effective rank of the
- scaled Vandermonde matrix and its singular values, and the
- specified value of `rcond`. For more details, see `linalg.lstsq`.
+ [residuals, rank, singular_values, rcond] : list
+ These values are only returned if `full` = True
+
+ resid -- sum of squared residuals of the least squares fit
+ rank -- the numerical rank of the scaled Vandermonde matrix
+ sv -- singular values of the scaled Vandermonde matrix
+ rcond -- value of `rcond`.
+
+ For more details, see `linalg.lstsq`.
Warns
-----
@@ -1739,4 +1745,43 @@ def lagweight(x):
# Laguerre series class
#
-exec(polytemplate.substitute(name='Laguerre', nick='lag', domain='[-1,1]'))
+class Laguerre(ABCPolyBase):
+ """A Laguerre series class.
+
+ The Laguerre class provides the standard Python numerical methods
+ '+', '-', '*', '//', '%', 'divmod', '**', and '()' as well as the
+ attributes and methods listed in the `ABCPolyBase` documentation.
+
+ Parameters
+ ----------
+ coef : array_like
+ Laguerre coefficients in order of increasing degree, i.e,
+ ``(1, 2, 3)`` gives ``1*L_0(x) + 2*L_1(X) + 3*L_2(x)``.
+ domain : (2,) array_like, optional
+ Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
+ to the interval ``[window[0], window[1]]`` by shifting and scaling.
+ The default value is [0, 1].
+ window : (2,) array_like, optional
+ Window, see `domain` for its use. The default value is [0, 1].
+
+ .. versionadded:: 1.6.0
+
+ """
+ # Virtual Functions
+ _add = staticmethod(lagadd)
+ _sub = staticmethod(lagsub)
+ _mul = staticmethod(lagmul)
+ _div = staticmethod(lagdiv)
+ _pow = staticmethod(lagpow)
+ _val = staticmethod(lagval)
+ _int = staticmethod(lagint)
+ _der = staticmethod(lagder)
+ _fit = staticmethod(lagfit)
+ _line = staticmethod(lagline)
+ _roots = staticmethod(lagroots)
+ _fromroots = staticmethod(lagfromroots)
+
+ # Virtual properties
+ nickname = 'lag'
+ domain = np.array(lagdomain)
+ window = np.array(lagdomain)
diff --git a/numpy/polynomial/legendre.py b/numpy/polynomial/legendre.py
index 8d89c8412..a54262823 100644
--- a/numpy/polynomial/legendre.py
+++ b/numpy/polynomial/legendre.py
@@ -83,11 +83,12 @@ numpy.polynomial.hermite_e
"""
from __future__ import division, absolute_import, print_function
+import warnings
import numpy as np
import numpy.linalg as la
+
from . import polyutils as pu
-import warnings
-from .polytemplate import polytemplate
+from ._polybase import ABCPolyBase
__all__ = ['legzero', 'legone', 'legx', 'legdomain', 'legline',
'legadd', 'legsub', 'legmulx', 'legmul', 'legdiv', 'legpow', 'legval',
@@ -1447,10 +1448,15 @@ def legfit(x, y, deg, rcond=None, full=False, w=None):
the coefficients for the data in column k of `y` are in column
`k`.
- [residuals, rank, singular_values, rcond] : present when `full` = True
- Residuals of the least-squares fit, the effective rank of the
- scaled Vandermonde matrix and its singular values, and the
- specified value of `rcond`. For more details, see `linalg.lstsq`.
+ [residuals, rank, singular_values, rcond] : list
+ These values are only returned if `full` = True
+
+ resid -- sum of squared residuals of the least squares fit
+ rank -- the numerical rank of the scaled Vandermonde matrix
+ sv -- singular values of the scaled Vandermonde matrix
+ rcond -- value of `rcond`.
+
+ For more details, see `linalg.lstsq`.
Warns
-----
@@ -1765,4 +1771,43 @@ def legweight(x):
# Legendre series class
#
-exec(polytemplate.substitute(name='Legendre', nick='leg', domain='[-1,1]'))
+class Legendre(ABCPolyBase):
+ """A Legendre series class.
+
+ The Legendre class provides the standard Python numerical methods
+ '+', '-', '*', '//', '%', 'divmod', '**', and '()' as well as the
+ attributes and methods listed in the `ABCPolyBase` documentation.
+
+ Parameters
+ ----------
+ coef : array_like
+ Legendre coefficients in order of increasing degree, i.e.,
+ ``(1, 2, 3)`` gives ``1*P_0(x) + 2*P_1(x) + 3*P_2(x)``.
+ domain : (2,) array_like, optional
+ Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
+ to the interval ``[window[0], window[1]]`` by shifting and scaling.
+ The default value is [-1, 1].
+ window : (2,) array_like, optional
+ Window, see `domain` for its use. The default value is [-1, 1].
+
+ .. versionadded:: 1.6.0
+
+ """
+ # Virtual Functions
+ _add = staticmethod(legadd)
+ _sub = staticmethod(legsub)
+ _mul = staticmethod(legmul)
+ _div = staticmethod(legdiv)
+ _pow = staticmethod(legpow)
+ _val = staticmethod(legval)
+ _int = staticmethod(legint)
+ _der = staticmethod(legder)
+ _fit = staticmethod(legfit)
+ _line = staticmethod(legline)
+ _roots = staticmethod(legroots)
+ _fromroots = staticmethod(legfromroots)
+
+ # Virtual properties
+ nickname = 'leg'
+ domain = np.array(legdomain)
+ window = np.array(legdomain)
diff --git a/numpy/polynomial/polynomial.py b/numpy/polynomial/polynomial.py
index c30fc6d0c..97f4fd148 100644
--- a/numpy/polynomial/polynomial.py
+++ b/numpy/polynomial/polynomial.py
@@ -61,11 +61,12 @@ __all__ = ['polyzero', 'polyone', 'polyx', 'polydomain', 'polyline',
'polyfit', 'polytrim', 'polyroots', 'Polynomial', 'polyval2d',
'polyval3d', 'polygrid2d', 'polygrid3d', 'polyvander2d', 'polyvander3d']
+import warnings
import numpy as np
import numpy.linalg as la
+
from . import polyutils as pu
-import warnings
-from .polytemplate import polytemplate
+from ._polybase import ABCPolyBase
polytrim = pu.trimcoef
@@ -823,7 +824,7 @@ def polyval2d(x, y, c):
Notes
-----
- .. versionadded::1.7.0
+ .. versionadded:: 1.7.0
"""
try:
@@ -883,7 +884,7 @@ def polygrid2d(x, y, c):
Notes
-----
- .. versionadded::1.7.0
+ .. versionadded:: 1.7.0
"""
c = polyval(x, c)
@@ -936,7 +937,7 @@ def polyval3d(x, y, z, c):
Notes
-----
- .. versionadded::1.7.0
+ .. versionadded:: 1.7.0
"""
try:
@@ -1000,7 +1001,7 @@ def polygrid3d(x, y, z, c):
Notes
-----
- .. versionadded::1.7.0
+ .. versionadded:: 1.7.0
"""
c = polyval(x, c)
@@ -1173,7 +1174,7 @@ def polyvander3d(x, y, z, deg) :
Notes
-----
- .. versionadded::1.7.0
+ .. versionadded:: 1.7.0
"""
ideg = [int(d) for d in deg]
@@ -1249,12 +1250,16 @@ def polyfit(x, y, deg, rcond=None, full=False, w=None):
the coefficients in column `k` of `coef` represent the polynomial
fit to the data in `y`'s `k`-th column.
- [residuals, rank, singular_values, rcond] : present when `full` == True
- Sum of the squared residuals (SSR) of the least-squares fit; the
- effective rank of the scaled Vandermonde matrix; its singular
- values; and the specified value of `rcond`. For more information,
- see `linalg.lstsq`.
+ [residuals, rank, singular_values, rcond] : list
+ These values are only returned if `full` = True
+ resid -- sum of squared residuals of the least squares fit
+ rank -- the numerical rank of the scaled Vandermonde matrix
+ sv -- singular values of the scaled Vandermonde matrix
+ rcond -- value of `rcond`.
+
+ For more details, see `linalg.lstsq`.
+
Raises
------
RankWarning
@@ -1490,4 +1495,43 @@ def polyroots(c):
# polynomial class
#
-exec(polytemplate.substitute(name='Polynomial', nick='poly', domain='[-1,1]'))
+class Polynomial(ABCPolyBase):
+ """A power series class.
+
+ The Polynomial class provides the standard Python numerical methods
+ '+', '-', '*', '//', '%', 'divmod', '**', and '()' as well as the
+ attributes and methods listed in the `ABCPolyBase` documentation.
+
+ Parameters
+ ----------
+ coef : array_like
+ Polynomial coefficients in order of increasing degree, i.e.,
+ ``(1, 2, 3)`` give ``1 + 2*x + 3*x**2``.
+ domain : (2,) array_like, optional
+ Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
+ to the interval ``[window[0], window[1]]`` by shifting and scaling.
+ The default value is [-1, 1].
+ window : (2,) array_like, optional
+ Window, see `domain` for its use. The default value is [-1, 1].
+
+ .. versionadded:: 1.6.0
+
+ """
+ # Virtual Functions
+ _add = staticmethod(polyadd)
+ _sub = staticmethod(polysub)
+ _mul = staticmethod(polymul)
+ _div = staticmethod(polydiv)
+ _pow = staticmethod(polypow)
+ _val = staticmethod(polyval)
+ _int = staticmethod(polyint)
+ _der = staticmethod(polyder)
+ _fit = staticmethod(polyfit)
+ _line = staticmethod(polyline)
+ _roots = staticmethod(polyroots)
+ _fromroots = staticmethod(polyfromroots)
+
+ # Virtual properties
+ nickname = 'poly'
+ domain = np.array(polydomain)
+ window = np.array(polydomain)
diff --git a/numpy/polynomial/polytemplate.py b/numpy/polynomial/polytemplate.py
index f315915d6..eeacf24fb 100644
--- a/numpy/polynomial/polytemplate.py
+++ b/numpy/polynomial/polytemplate.py
@@ -13,6 +13,12 @@ from __future__ import division, absolute_import, print_function
import string
import sys
+import warnings
+
+from numpy import ModuleDeprecationWarning
+
+warnings.warn("The polytemplate module will be removed in Numpy 1.10.0.",
+ ModuleDeprecationWarning)
polytemplate = string.Template('''
from __future__ import division, absolute_import, print_function
diff --git a/numpy/polynomial/polyutils.py b/numpy/polynomial/polyutils.py
index 63743bb40..99f508521 100644
--- a/numpy/polynomial/polyutils.py
+++ b/numpy/polynomial/polyutils.py
@@ -1,41 +1,53 @@
"""
-Utililty objects for the polynomial modules.
+Utililty classes and functions for the polynomial modules.
This module provides: error and warning objects; a polynomial base class;
and some routines used in both the `polynomial` and `chebyshev` modules.
Error objects
-------------
-- `PolyError` -- base class for this sub-package's errors.
-- `PolyDomainError` -- raised when domains are "mismatched."
+
+.. autosummary::
+ :toctree: generated/
+
+ PolyError base class for this sub-package's errors.
+ PolyDomainError raised when domains are mismatched.
Warning objects
---------------
-- `RankWarning` -- raised by a least-squares fit when a rank-deficient
- matrix is encountered.
+
+.. autosummary::
+ :toctree: generated/
+
+ RankWarning raised in least-squares fit for rank-deficient matrix.
Base class
----------
-- `PolyBase` -- The base class for the `Polynomial` and `Chebyshev`
- classes.
+
+.. autosummary::
+ :toctree: generated/
+
+ PolyBase Obsolete base class for the polynomial classes. Do not use.
Functions
---------
-- `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.
+
+.. autosummary::
+ :toctree: generated/
+
+ as_series convert list of array_likes into 1-D arrays of common type.
+ trimseq remove trailing zeros.
+ trimcoef remove small trailing coefficients.
+ getdomain return the 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, absolute_import, print_function
-__all__ = ['RankWarning', 'PolyError', 'PolyDomainError', 'PolyBase',
- 'as_series', 'trimseq', 'trimcoef', 'getdomain', 'mapdomain',
- 'mapparms']
+__all__ = ['RankWarning', 'PolyError', 'PolyDomainError', 'as_series',
+ 'trimseq', 'trimcoef', 'getdomain', 'mapdomain', 'mapparms',
+ 'PolyBase']
import warnings
import numpy as np
@@ -67,6 +79,15 @@ class PolyDomainError(PolyError) :
#
class PolyBase(object) :
+ """
+ Base class for all polynomial types.
+
+ Deprecated in numpy 1.9.0, use the abstract
+ ABCPolyBase class instead. Note that the latter
+ reguires a number of virtual functions to be
+ implemented.
+
+ """
pass
#