def _fix_real_abs_gt_1(x): """Convert `x` to complex if it has real components x_i with abs(x_i)>1. Otherwise, output is just the array version of the input (via asarray). Parameters ---------- x : array_like Returns ------- array Examples -------- >>> np.lib.scimath._fix_real_abs_gt_1([0,1]) array([0, 1]) >>> np.lib.scimath._fix_real_abs_gt_1([0,2]) array([ 0.+0.j, 2.+0.j]) """ x = asarray(x) if any(isreal(x) & (abs(x) > 1)): x = _tocomplex(x) return x
def _fix_real_lt_zero(x): """Convert `x` to complex if it has real, negative components. Otherwise, output is just the array version of the input (via asarray). Parameters ---------- x : array_like Returns ------- array Examples -------- >>> np.lib.scimath._fix_real_lt_zero([1,2]) array([1, 2]) >>> np.lib.scimath._fix_real_lt_zero([-1,2]) array([-1.+0.j, 2.+0.j]) """ x = asarray(x) if any(isreal(x) & (x < 0)): x = _tocomplex(x) return x
def asfarray(a, dtype=_nx.float_): """ Return an array converted to a float type. Parameters ---------- a : array_like The input array. dtype : str or dtype object, optional Float type code to coerce input array `a`. If `dtype` is one of the 'int' dtypes, it is replaced with float64. Returns ------- out : ndarray The input `a` as a float ndarray. Examples -------- >>> np.asfarray([2, 3]) array([ 2., 3.]) >>> np.asfarray([2, 3], dtype='float') array([ 2., 3.]) >>> np.asfarray([2, 3], dtype='int8') array([ 2., 3.]) """ if not _nx.issubdtype(dtype, _nx.inexact): dtype = _nx.float_ return asarray(a, dtype=dtype)
def mintypecode(typechars, typeset='GDFgdf', default='d'): """ Return the character for the minimum-size type to which given types can be safely cast. The returned type character must represent the smallest size dtype such that an array of the returned type can handle the data from an array of all types in `typechars` (or if `typechars` is an array, then its dtype.char). Parameters ---------- typechars : list of str or array_like If a list of strings, each string should represent a dtype. If array_like, the character representation of the array dtype is used. typeset : str or list of str, optional The set of characters that the returned character is chosen from. The default set is 'GDFgdf'. default : str, optional The default character, this is returned if none of the characters in `typechars` matches a character in `typeset`. Returns ------- typechar : str The character representing the minimum-size type that was found. See Also -------- dtype, sctype2char, maximum_sctype Examples -------- >>> np.mintypecode(['d', 'f', 'S']) 'd' >>> x = np.array([1.1, 2-3.j]) >>> np.mintypecode(x) 'D' >>> np.mintypecode('abceh', default='G') 'G' """ typecodes = [(isinstance(t, str) and t) or asarray(t).dtype.char for t in typechars] intersection = [t for t in typecodes if t in typeset] if not intersection: return default if 'F' in intersection and 'd' in intersection: return 'D' l = [] for t in intersection: i = _typecodes_by_elsize.index(t) l.append((i, t)) l.sort() return l[0][1]
def iscomplexobj(x): """ Check for a complex type or an array of complex numbers. The type of the input is checked, not the value. Even if the input has an imaginary part equal to zero, `iscomplexobj` evaluates to True. Parameters ---------- x : any The input can be of any type and shape. Returns ------- iscomplexobj : bool The return value, True if `x` is of a complex type or has at least one complex element. See Also -------- isrealobj, iscomplex Examples -------- >>> np.iscomplexobj(1) False >>> np.iscomplexobj(1+0j) True >>> np.iscomplexobj([3, 1+0j, True]) True """ try: dtype = x.dtype type_ = dtype.type except AttributeError: type_ = asarray(x).dtype.type return issubclass(type_, _nx.complexfloating)
def expand_dims(a, axis): """ Expand the shape of an array. Insert a new axis that will appear at the `axis` position in the expanded array shape. .. note:: Previous to NumPy 1.13.0, neither ``axis < -a.ndim - 1`` nor ``axis > a.ndim`` raised errors or put the new axis where documented. Those axis values are now deprecated and will raise an AxisError in the future. Parameters ---------- a : array_like Input array. axis : int Position in the expanded axes where the new axis is placed. Returns ------- res : ndarray Output array. The number of dimensions is one greater than that of the input array. See Also -------- squeeze : The inverse operation, removing singleton dimensions reshape : Insert, remove, and combine dimensions, and resize existing ones doc.indexing, atleast_1d, atleast_2d, atleast_3d Examples -------- >>> x = np.array([1,2]) >>> x.shape (2,) The following is equivalent to ``x[np.newaxis,:]`` or ``x[np.newaxis]``: >>> y = np.expand_dims(x, axis=0) >>> y array([[1, 2]]) >>> y.shape (1, 2) >>> y = np.expand_dims(x, axis=1) # Equivalent to x[:,np.newaxis] >>> y array([[1], [2]]) >>> y.shape (2, 1) Note that some examples may use ``None`` instead of ``np.newaxis``. These are the same objects: >>> np.newaxis is None True """ a = asarray(a) shape = a.shape if axis > a.ndim or axis < -a.ndim - 1: # 2017-05-17, 1.13.0 warnings.warn( "Both axis > a.ndim and axis < -a.ndim - 1 are " "deprecated and will raise an AxisError in the future.", DeprecationWarning, stacklevel=2) # When the deprecation period expires, delete this if block, if axis < 0: axis = axis + a.ndim + 1 # and uncomment the following line. # axis = normalize_axis_index(axis, a.ndim + 1) return a.reshape(shape[:axis] + (1, ) + shape[axis:])
def apply_over_axes(func, a, axes): """ Apply a function repeatedly over multiple axes. `func` is called as `res = func(a, axis)`, where `axis` is the first element of `axes`. The result `res` of the function call must have either the same dimensions as `a` or one less dimension. If `res` has one less dimension than `a`, a dimension is inserted before `axis`. The call to `func` is then repeated for each axis in `axes`, with `res` as the first argument. Parameters ---------- func : function This function must take two arguments, `func(a, axis)`. a : array_like Input array. axes : array_like Axes over which `func` is applied; the elements must be integers. Returns ------- apply_over_axis : ndarray The output array. The number of dimensions is the same as `a`, but the shape can be different. This depends on whether `func` changes the shape of its output with respect to its input. See Also -------- apply_along_axis : Apply a function to 1-D slices of an array along the given axis. Notes ------ This function is equivalent to tuple axis arguments to reorderable ufuncs with keepdims=True. Tuple axis arguments to ufuncs have been available since version 1.7.0. Examples -------- >>> a = np.arange(24).reshape(2,3,4) >>> a array([[[ 0, 1, 2, 3], [ 4, 5, 6, 7], [ 8, 9, 10, 11]], [[12, 13, 14, 15], [16, 17, 18, 19], [20, 21, 22, 23]]]) Sum over axes 0 and 2. The result has same number of dimensions as the original array: >>> np.apply_over_axes(np.sum, a, [0,2]) array([[[ 60], [ 92], [124]]]) Tuple axis arguments to ufuncs are equivalent: >>> np.sum(a, axis=(0,2), keepdims=True) array([[[ 60], [ 92], [124]]]) """ val = asarray(a) N = a.ndim if array(axes).ndim == 0: axes = (axes, ) for axis in axes: if axis < 0: axis = N + axis args = (val, axis) res = func(*args) if res.ndim == val.ndim: val = res else: res = expand_dims(res, axis) if res.ndim == val.ndim: val = res else: raise ValueError("function is not returning " "an array of the correct shape") return val
def __init__(self, arr): self.iter = asarray(arr).flat
def ix_(*args): """ Construct an open mesh from multiple sequences. This function takes N 1-D sequences and returns N outputs with N dimensions each, such that the shape is 1 in all but one dimension and the dimension with the non-unit shape value cycles through all N dimensions. Using `ix_` one can quickly construct index arrays that will index the cross product. ``a[np.ix_([1,3],[2,5])]`` returns the array ``[[a[1,2] a[1,5]], [a[3,2] a[3,5]]]``. Parameters ---------- args : 1-D sequences Each sequence should be of integer or boolean type. Boolean sequences will be interpreted as boolean masks for the corresponding dimension (equivalent to passing in ``np.nonzero(boolean_sequence)``). Returns ------- out : tuple of ndarrays N arrays with N dimensions each, with N the number of input sequences. Together these arrays form an open mesh. See Also -------- ogrid, mgrid, meshgrid Examples -------- >>> a = np.arange(10).reshape(2, 5) >>> a array([[0, 1, 2, 3, 4], [5, 6, 7, 8, 9]]) >>> ixgrid = np.ix_([0, 1], [2, 4]) >>> ixgrid (array([[0], [1]]), array([[2, 4]])) >>> ixgrid[0].shape, ixgrid[1].shape ((2, 1), (1, 2)) >>> a[ixgrid] array([[2, 4], [7, 9]]) >>> ixgrid = np.ix_([True, True], [2, 4]) >>> a[ixgrid] array([[2, 4], [7, 9]]) >>> ixgrid = np.ix_([True, True], [False, False, True, False, True]) >>> a[ixgrid] array([[2, 4], [7, 9]]) """ out = [] nd = len(args) for k, new in enumerate(args): new = asarray(new) if new.ndim != 1: raise ValueError("Cross index must be 1 dimensional") if new.size == 0: # Explicitly type empty arrays to avoid float default new = new.astype(_nx.intp) if issubdtype(new.dtype, _nx.bool_): new, = new.nonzero() new = new.reshape((1, ) * k + (new.size, ) + (1, ) * (nd - k - 1)) out.append(new) return tuple(out)
def polyval(p, x): """ Evaluate a polynomial at specific values. If `p` is of length N, this function returns the value: ``p[0]*x**(N-1) + p[1]*x**(N-2) + ... + p[N-2]*x + p[N-1]`` If `x` is a sequence, then `p(x)` is returned for each element of `x`. If `x` is another polynomial then the composite polynomial `p(x(t))` is returned. Parameters ---------- p : array_like or poly1d object 1D array of polynomial coefficients (including coefficients equal to zero) from highest degree to the constant term, or an instance of poly1d. x : array_like or poly1d object A number, an array of numbers, or an instance of poly1d, at which to evaluate `p`. Returns ------- values : ndarray or poly1d If `x` is a poly1d instance, the result is the composition of the two polynomials, i.e., `x` is "substituted" in `p` and the simplified result is returned. In addition, the type of `x` - array_like or poly1d - governs the type of the output: `x` array_like => `values` array_like, `x` a poly1d object => `values` is also. See Also -------- poly1d: A polynomial class. Notes ----- Horner's scheme [1]_ is used to evaluate the polynomial. Even so, for polynomials of high degree the values may be inaccurate due to rounding errors. Use carefully. References ---------- .. [1] I. N. Bronshtein, K. A. Semendyayev, and K. A. Hirsch (Eng. trans. Ed.), *Handbook of Mathematics*, New York, Van Nostrand Reinhold Co., 1985, pg. 720. Examples -------- >>> np.polyval([3,0,1], 5) # 3 * 5**2 + 0 * 5**1 + 1 76 >>> np.polyval([3,0,1], np.poly1d(5)) poly1d([ 76.]) >>> np.polyval(np.poly1d([3,0,1]), 5) 76 >>> np.polyval(np.poly1d([3,0,1]), np.poly1d(5)) poly1d([ 76.]) """ p = NX.asarray(p) if isinstance(x, poly1d): y = 0 else: x = NX.asarray(x) y = NX.zeros_like(x) for i in range(len(p)): y = y * x + p[i] return y
def polyfit(x, y, deg, rcond=None, full=False, w=None, cov=False): """ Least squares polynomial fit. Fit a polynomial ``p(x) = p[0] * x**deg + ... + p[deg]`` of degree `deg` to points `(x, y)`. Returns a vector of coefficients `p` that minimises the squared error. 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 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 to apply to the y-coordinates of the sample points. For gaussian uncertainties, use 1/sigma (not 1/sigma**2). cov : bool, optional Return the estimate and the covariance matrix of the estimate If full is True, then cov is not returned. Returns ------- p : ndarray, shape (deg + 1,) or (deg + 1, K) Polynomial coefficients, highest power first. If `y` was 2-D, the coefficients for `k`-th data set are in ``p[:,k]``. residuals, rank, singular_values, rcond Present only if `full` = True. Residuals of the least-squares fit, the effective rank of the scaled Vandermonde coefficient matrix, its singular values, and the specified value of `rcond`. For more details, see `linalg.lstsq`. V : ndarray, shape (M,M) or (M,M,K) Present only if `full` = False and `cov`=True. The covariance matrix of the polynomial coefficient estimates. The diagonal of this matrix are the variance estimates for each coefficient. If y is a 2-D array, then the covariance matrix for the `k`-th data set are in ``V[:,:,k]`` 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) See Also -------- polyval : Compute polynomial values. linalg.lstsq : Computes a least-squares fit. scipy.interpolate.UnivariateSpline : Computes spline fits. Notes ----- The solution minimizes the squared error .. math :: E = \\sum_{j=0}^k |p(x_j) - y_j|^2 in the equations:: x[0]**n * p[0] + ... + x[0] * p[n-1] + p[n] = y[0] x[1]**n * p[0] + ... + x[1] * p[n-1] + p[n] = y[1] ... x[k]**n * p[0] + ... + x[k] * p[n-1] + p[n] = y[k] The coefficient matrix of the coefficients `p` is a Vandermonde matrix. `polyfit` issues a `RankWarning` when the least-squares fit is badly conditioned. This implies that the best fit is not well-defined due to numerical error. The results may be improved by lowering the polynomial degree or by replacing `x` by `x` - `x`.mean(). The `rcond` parameter can also be set to a value smaller than its default, but the resulting fit may be spurious: including contributions from the small singular values can add numerical noise to the result. Note that fitting polynomial coefficients is inherently badly conditioned when the degree of the polynomial is large or the interval of sample points is badly centered. The quality of the fit should always be checked in these cases. When polynomial fits are not satisfactory, splines may be a good alternative. References ---------- .. [1] Wikipedia, "Curve fitting", http://en.wikipedia.org/wiki/Curve_fitting .. [2] Wikipedia, "Polynomial interpolation", http://en.wikipedia.org/wiki/Polynomial_interpolation Examples -------- >>> x = np.array([0.0, 1.0, 2.0, 3.0, 4.0, 5.0]) >>> y = np.array([0.0, 0.8, 0.9, 0.1, -0.8, -1.0]) >>> z = np.polyfit(x, y, 3) >>> z array([ 0.08703704, -0.81349206, 1.69312169, -0.03968254]) It is convenient to use `poly1d` objects for dealing with polynomials: >>> p = np.poly1d(z) >>> p(0.5) 0.6143849206349179 >>> p(3.5) -0.34732142857143039 >>> p(10) 22.579365079365115 High-order polynomials may oscillate wildly: >>> p30 = np.poly1d(np.polyfit(x, y, 30)) /... RankWarning: Polyfit may be poorly conditioned... >>> p30(4) -0.80000000000000204 >>> p30(5) -0.99999999999999445 >>> p30(4.5) -0.10547061179440398 Illustration: >>> import matplotlib.pyplot as plt >>> xp = np.linspace(-2, 6, 100) >>> _ = plt.plot(x, y, '.', xp, p(xp), '-', xp, p30(xp), '--') >>> plt.ylim(-2,2) (-2, 2) >>> plt.show() """ order = int(deg) + 1 x = NX.asarray(x) + 0.0 y = NX.asarray(y) + 0.0 # check arguments. if deg < 0: raise ValueError("expected deg >= 0") if x.ndim != 1: raise TypeError("expected 1D vector for x") if x.size == 0: raise TypeError("expected non-empty vector for x") if y.ndim < 1 or y.ndim > 2: raise TypeError("expected 1D or 2D array for y") if x.shape[0] != y.shape[0]: raise TypeError("expected x and y to have same length") # set rcond if rcond is None: rcond = len(x) * finfo(x.dtype).eps # set up least squares equation for powers of x lhs = vander(x, order) rhs = y # apply weighting if w is not None: w = NX.asarray(w) + 0.0 if w.ndim != 1: raise TypeError("expected a 1-d array for weights") if w.shape[0] != y.shape[0]: raise TypeError("expected w and y to have the same length") lhs *= w[:, NX.newaxis] if rhs.ndim == 2: rhs *= w[:, NX.newaxis] else: rhs *= w # scale lhs to improve condition number and solve scale = NX.sqrt((lhs * lhs).sum(axis=0)) lhs /= scale c, resids, rank, s = lstsq(lhs, rhs, rcond) c = (c.T / scale).T # broadcast scale coefficients # warn on rank reduction, which indicates an ill conditioned matrix if rank != order and not full: msg = "Polyfit may be poorly conditioned" warnings.warn(msg, RankWarning, stacklevel=2) if full: return c, resids, rank, s, rcond elif cov: Vbase = inv(dot(lhs.T, lhs)) Vbase /= NX.outer(scale, scale) # Some literature ignores the extra -2.0 factor in the denominator, but # it is included here because the covariance of Multivariate Student-T # (which is implied by a Bayesian uncertainty analysis) includes it. # Plus, it gives a slightly more conservative estimate of uncertainty. if len(x) <= order + 2: raise ValueError("the number of data points must exceed order + 2 " "for Bayesian estimate the covariance matrix") fac = resids / (len(x) - order - 2.0) if y.ndim == 1: return c, Vbase * fac else: return c, Vbase[:, :, NX.newaxis] * fac else: return c
def poly(seq_of_zeros): """ Find the coefficients of a polynomial with the given sequence of roots. Returns the coefficients of the polynomial whose leading coefficient is one for the given sequence of zeros (multiple roots must be included in the sequence as many times as their multiplicity; see Examples). A square matrix (or array, which will be treated as a matrix) can also be given, in which case the coefficients of the characteristic polynomial of the matrix are returned. Parameters ---------- seq_of_zeros : array_like, shape (N,) or (N, N) A sequence of polynomial roots, or a square array or matrix object. Returns ------- c : ndarray 1D array of polynomial coefficients from highest to lowest degree: ``c[0] * x**(N) + c[1] * x**(N-1) + ... + c[N-1] * x + c[N]`` where c[0] always equals 1. Raises ------ ValueError If input is the wrong shape (the input must be a 1-D or square 2-D array). See Also -------- polyval : Compute polynomial values. roots : Return the roots of a polynomial. polyfit : Least squares polynomial fit. poly1d : A one-dimensional polynomial class. Notes ----- Specifying the roots of a polynomial still leaves one degree of freedom, typically represented by an undetermined leading coefficient. [1]_ In the case of this function, that coefficient - the first one in the returned array - is always taken as one. (If for some reason you have one other point, the only automatic way presently to leverage that information is to use ``polyfit``.) The characteristic polynomial, :math:`p_a(t)`, of an `n`-by-`n` matrix **A** is given by :math:`p_a(t) = \\mathrm{det}(t\\, \\mathbf{I} - \\mathbf{A})`, where **I** is the `n`-by-`n` identity matrix. [2]_ References ---------- .. [1] M. Sullivan and M. Sullivan, III, "Algebra and Trignometry, Enhanced With Graphing Utilities," Prentice-Hall, pg. 318, 1996. .. [2] G. Strang, "Linear Algebra and Its Applications, 2nd Edition," Academic Press, pg. 182, 1980. Examples -------- Given a sequence of a polynomial's zeros: >>> np.poly((0, 0, 0)) # Multiple root example array([1, 0, 0, 0]) The line above represents z**3 + 0*z**2 + 0*z + 0. >>> np.poly((-1./2, 0, 1./2)) array([ 1. , 0. , -0.25, 0. ]) The line above represents z**3 - z/4 >>> np.poly((np.random.random(1.)[0], 0, np.random.random(1.)[0])) array([ 1. , -0.77086955, 0.08618131, 0. ]) #random Given a square array object: >>> P = np.array([[0, 1./3], [-1./2, 0]]) >>> np.poly(P) array([ 1. , 0. , 0.16666667]) Note how in all cases the leading coefficient is always 1. """ seq_of_zeros = atleast_1d(seq_of_zeros) sh = seq_of_zeros.shape if len(sh) == 2 and sh[0] == sh[1] and sh[0] != 0: seq_of_zeros = eigvals(seq_of_zeros) elif len(sh) == 1: dt = seq_of_zeros.dtype # Let object arrays slip through, e.g. for arbitrary precision if dt != object: seq_of_zeros = seq_of_zeros.astype(mintypecode(dt.char)) else: raise ValueError("input must be 1d or non-empty square 2d array.") if len(seq_of_zeros) == 0: return 1.0 dt = seq_of_zeros.dtype a = ones((1, ), dtype=dt) for k in range(len(seq_of_zeros)): a = NX.convolve(a, array([1, -seq_of_zeros[k]], dtype=dt), mode='full') if issubclass(a.dtype.type, NX.complexfloating): # if complex roots are all complex conjugates, the roots are real. roots = NX.asarray(seq_of_zeros, complex) if NX.all(NX.sort(roots) == NX.sort(roots.conjugate())): a = a.real.copy() return a
def polyder(p, m=1): """ Return the derivative of the specified order of a polynomial. Parameters ---------- p : poly1d or sequence Polynomial to differentiate. A sequence is interpreted as polynomial coefficients, see `poly1d`. m : int, optional Order of differentiation (default: 1) Returns ------- der : poly1d A new polynomial representing the derivative. See Also -------- polyint : Anti-derivative of a polynomial. poly1d : Class for one-dimensional polynomials. Examples -------- The derivative of the polynomial :math:`x^3 + x^2 + x^1 + 1` is: >>> p = np.poly1d([1,1,1,1]) >>> p2 = np.polyder(p) >>> p2 poly1d([3, 2, 1]) which evaluates to: >>> p2(2.) 17.0 We can verify this, approximating the derivative with ``(f(x + h) - f(x))/h``: >>> (p(2. + 0.001) - p(2.)) / 0.001 17.007000999997857 The fourth-order derivative of a 3rd-order polynomial is zero: >>> np.polyder(p, 2) poly1d([6, 2]) >>> np.polyder(p, 3) poly1d([6]) >>> np.polyder(p, 4) poly1d([ 0.]) """ m = int(m) if m < 0: raise ValueError("Order of derivative must be positive (see polyint)") truepoly = isinstance(p, poly1d) p = NX.asarray(p) n = len(p) - 1 y = p[:-1] * NX.arange(n, 0, -1) if m == 0: val = p else: val = polyder(y, m - 1) if truepoly: val = poly1d(val) return val
def polyint(p, m=1, k=None): """ Return an antiderivative (indefinite integral) of a polynomial. The returned order `m` antiderivative `P` of polynomial `p` satisfies :math:`\\frac{d^m}{dx^m}P(x) = p(x)` and is defined up to `m - 1` integration constants `k`. The constants determine the low-order polynomial part .. math:: \\frac{k_{m-1}}{0!} x^0 + \\ldots + \\frac{k_0}{(m-1)!}x^{m-1} of `P` so that :math:`P^{(j)}(0) = k_{m-j-1}`. Parameters ---------- p : array_like or poly1d Polynomial to differentiate. A sequence is interpreted as polynomial coefficients, see `poly1d`. m : int, optional Order of the antiderivative. (Default: 1) k : list of `m` scalars or scalar, optional Integration constants. They are given in the order of integration: those corresponding to highest-order terms come first. If ``None`` (default), all constants are assumed to be zero. If `m = 1`, a single scalar can be given instead of a list. See Also -------- polyder : derivative of a polynomial poly1d.integ : equivalent method Examples -------- The defining property of the antiderivative: >>> p = np.poly1d([1,1,1]) >>> P = np.polyint(p) >>> P poly1d([ 0.33333333, 0.5 , 1. , 0. ]) >>> np.polyder(P) == p True The integration constants default to zero, but can be specified: >>> P = np.polyint(p, 3) >>> P(0) 0.0 >>> np.polyder(P)(0) 0.0 >>> np.polyder(P, 2)(0) 0.0 >>> P = np.polyint(p, 3, k=[6,5,3]) >>> P poly1d([ 0.01666667, 0.04166667, 0.16666667, 3. , 5. , 3. ]) Note that 3 = 6 / 2!, and that the constants are given in the order of integrations. Constant of the highest-order polynomial term comes first: >>> np.polyder(P, 2)(0) 6.0 >>> np.polyder(P, 1)(0) 5.0 >>> P(0) 3.0 """ m = int(m) if m < 0: raise ValueError("Order of integral must be positive (see polyder)") if k is None: k = NX.zeros(m, float) k = atleast_1d(k) if len(k) == 1 and m > 1: k = k[0] * NX.ones(m, float) if len(k) < m: raise ValueError( "k must be a scalar or a rank-1 array of length 1 or >m.") truepoly = isinstance(p, poly1d) p = NX.asarray(p) if m == 0: if truepoly: return poly1d(p) return p else: # Note: this must work also with object and integer arrays y = NX.concatenate((p.__truediv__(NX.arange(len(p), 0, -1)), [k[0]])) val = polyint(y, m - 1, k=k[1:]) if truepoly: return poly1d(val) return val
def __array__(self, t=None): if t: return NX.asarray(self.coeffs, t) else: return NX.asarray(self.coeffs)